[b443089] | 1 | /*
|
---|
| 2 | * Delphes: a framework for fast simulation of a generic collider experiment
|
---|
| 3 | * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
|
---|
[1fa50c2] | 4 | *
|
---|
[b443089] | 5 | * This program is free software: you can redistribute it and/or modify
|
---|
| 6 | * it under the terms of the GNU General Public License as published by
|
---|
| 7 | * the Free Software Foundation, either version 3 of the License, or
|
---|
| 8 | * (at your option) any later version.
|
---|
[1fa50c2] | 9 | *
|
---|
[b443089] | 10 | * This program is distributed in the hope that it will be useful,
|
---|
| 11 | * but WITHOUT ANY WARRANTY; without even the implied warranty of
|
---|
| 12 | * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
---|
| 13 | * GNU General Public License for more details.
|
---|
[1fa50c2] | 14 | *
|
---|
[b443089] | 15 | * You should have received a copy of the GNU General Public License
|
---|
| 16 | * along with this program. If not, see <http://www.gnu.org/licenses/>.
|
---|
| 17 | */
|
---|
| 18 |
|
---|
[d7d2da3] | 19 | /** \class ParticlePropagator
|
---|
| 20 | *
|
---|
| 21 | * Propagates charged and neutral particles
|
---|
[d41ba4a] | 22 | * from a given vertex to a cylinder defined by its radius,
|
---|
[d7d2da3] | 23 | * its half-length, centered at (0,0,0) and with its axis
|
---|
| 24 | * oriented along the z-axis.
|
---|
| 25 | *
|
---|
| 26 | * \author P. Demin - UCL, Louvain-la-Neuve
|
---|
| 27 | *
|
---|
| 28 | */
|
---|
| 29 |
|
---|
| 30 | #include "modules/ParticlePropagator.h"
|
---|
| 31 |
|
---|
| 32 | #include "classes/DelphesClasses.h"
|
---|
| 33 | #include "classes/DelphesFactory.h"
|
---|
| 34 | #include "classes/DelphesFormula.h"
|
---|
| 35 |
|
---|
| 36 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
[341014c] | 37 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
| 38 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
[d7d2da3] | 39 |
|
---|
| 40 | #include "TDatabasePDG.h"
|
---|
[341014c] | 41 | #include "TFormula.h"
|
---|
[d7d2da3] | 42 | #include "TLorentzVector.h"
|
---|
[341014c] | 43 | #include "TMath.h"
|
---|
| 44 | #include "TObjArray.h"
|
---|
| 45 | #include "TRandom3.h"
|
---|
| 46 | #include "TString.h"
|
---|
[d7d2da3] | 47 |
|
---|
[d41ba4a] | 48 | #include <algorithm>
|
---|
[d7d2da3] | 49 | #include <iostream>
|
---|
| 50 | #include <sstream>
|
---|
[341014c] | 51 | #include <stdexcept>
|
---|
[d7d2da3] | 52 |
|
---|
| 53 | using namespace std;
|
---|
| 54 |
|
---|
| 55 | //------------------------------------------------------------------------------
|
---|
| 56 |
|
---|
| 57 | ParticlePropagator::ParticlePropagator() :
|
---|
| 58 | fItInputArray(0)
|
---|
| 59 | {
|
---|
| 60 | }
|
---|
| 61 |
|
---|
| 62 | //------------------------------------------------------------------------------
|
---|
| 63 |
|
---|
| 64 | ParticlePropagator::~ParticlePropagator()
|
---|
| 65 | {
|
---|
| 66 | }
|
---|
| 67 |
|
---|
| 68 | //------------------------------------------------------------------------------
|
---|
| 69 |
|
---|
| 70 | void ParticlePropagator::Init()
|
---|
| 71 | {
|
---|
| 72 | fRadius = GetDouble("Radius", 1.0);
|
---|
[341014c] | 73 | fRadius2 = fRadius * fRadius;
|
---|
[d7d2da3] | 74 | fHalfLength = GetDouble("HalfLength", 3.0);
|
---|
| 75 | fBz = GetDouble("Bz", 0.0);
|
---|
| 76 | if(fRadius < 1.0E-2)
|
---|
[d41ba4a] | 77 | {
|
---|
[d7d2da3] | 78 | cout << "ERROR: magnetic field radius is too low\n";
|
---|
| 79 | return;
|
---|
| 80 | }
|
---|
| 81 | if(fHalfLength < 1.0E-2)
|
---|
| 82 | {
|
---|
| 83 | cout << "ERROR: magnetic field length is too low\n";
|
---|
| 84 | return;
|
---|
| 85 | }
|
---|
| 86 |
|
---|
[9330b7b] | 87 | fRadiusMax = GetDouble("RadiusMax", fRadius);
|
---|
| 88 | fHalfLengthMax = GetDouble("HalfLengthMax", fHalfLength);
|
---|
| 89 |
|
---|
[d7d2da3] | 90 | // import array with output from filter/classifier module
|
---|
| 91 |
|
---|
| 92 | fInputArray = ImportArray(GetString("InputArray", "Delphes/stableParticles"));
|
---|
| 93 | fItInputArray = fInputArray->MakeIterator();
|
---|
| 94 |
|
---|
[187fc41] | 95 | // import beamspot
|
---|
[bff2e33] | 96 | try
|
---|
| 97 | {
|
---|
| 98 | fBeamSpotInputArray = ImportArray(GetString("BeamSpotInputArray", "BeamSpotFilter/beamSpotParticle"));
|
---|
| 99 | }
|
---|
| 100 | catch(runtime_error &e)
|
---|
| 101 | {
|
---|
| 102 | fBeamSpotInputArray = 0;
|
---|
[e55f5b0] | 103 | }
|
---|
[d7d2da3] | 104 | // create output arrays
|
---|
| 105 |
|
---|
| 106 | fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
|
---|
[751cb9c] | 107 | fNeutralOutputArray = ExportArray(GetString("NeutralOutputArray", "neutralParticles"));
|
---|
[d7d2da3] | 108 | fChargedHadronOutputArray = ExportArray(GetString("ChargedHadronOutputArray", "chargedHadrons"));
|
---|
| 109 | fElectronOutputArray = ExportArray(GetString("ElectronOutputArray", "electrons"));
|
---|
| 110 | fMuonOutputArray = ExportArray(GetString("MuonOutputArray", "muons"));
|
---|
| 111 | }
|
---|
| 112 |
|
---|
| 113 | //------------------------------------------------------------------------------
|
---|
| 114 |
|
---|
| 115 | void ParticlePropagator::Finish()
|
---|
| 116 | {
|
---|
| 117 | if(fItInputArray) delete fItInputArray;
|
---|
| 118 | }
|
---|
| 119 |
|
---|
| 120 | //------------------------------------------------------------------------------
|
---|
| 121 |
|
---|
| 122 | void ParticlePropagator::Process()
|
---|
| 123 | {
|
---|
[60e1de6] | 124 | Candidate *candidate, *mother, *particle;
|
---|
| 125 | TLorentzVector particlePosition, particleMomentum, beamSpotPosition;
|
---|
[d7d2da3] | 126 | Double_t px, py, pz, pt, pt2, e, q;
|
---|
[38b4e15] | 127 | Double_t x, y, z, t, r;
|
---|
[d7d2da3] | 128 | Double_t x_c, y_c, r_c, phi_c, phi_0;
|
---|
| 129 | Double_t x_t, y_t, z_t, r_t;
|
---|
[38b4e15] | 130 | Double_t t_z, t_r;
|
---|
| 131 | Double_t discr;
|
---|
| 132 | Double_t gammam, omega;
|
---|
| 133 | Double_t xd, yd, zd;
|
---|
| 134 | Double_t l, d0, dz, ctgTheta, alpha;
|
---|
[187fc41] | 135 | Double_t bsx, bsy, bsz;
|
---|
[38b4e15] | 136 | Double_t rxp, rdp, t_R;
|
---|
| 137 | Double_t td, pio, phid, sign_pz, vz;
|
---|
| 138 |
|
---|
[d7d2da3] | 139 | const Double_t c_light = 2.99792458E8;
|
---|
[e55f5b0] | 140 |
|
---|
[341014c] | 141 | if(!fBeamSpotInputArray || fBeamSpotInputArray->GetSize() == 0)
|
---|
[187fc41] | 142 | beamSpotPosition.SetXYZT(0.0, 0.0, 0.0, 0.0);
|
---|
| 143 | else
|
---|
| 144 | {
|
---|
[341014c] | 145 | Candidate &beamSpotCandidate = *((Candidate *)fBeamSpotInputArray->At(0));
|
---|
[187fc41] | 146 | beamSpotPosition = beamSpotCandidate.Position;
|
---|
| 147 | }
|
---|
[e55f5b0] | 148 |
|
---|
[d7d2da3] | 149 | fItInputArray->Reset();
|
---|
[341014c] | 150 | while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
|
---|
[d7d2da3] | 151 | {
|
---|
[60e1de6] | 152 | if(candidate->GetCandidates()->GetEntriesFast() == 0)
|
---|
| 153 | {
|
---|
| 154 | particle = candidate;
|
---|
| 155 | }
|
---|
| 156 | else
|
---|
| 157 | {
|
---|
[341014c] | 158 | particle = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
|
---|
[60e1de6] | 159 | }
|
---|
| 160 |
|
---|
| 161 | particlePosition = particle->Position;
|
---|
| 162 | particleMomentum = particle->Momentum;
|
---|
[38b4e15] | 163 |
|
---|
| 164 | // Constants
|
---|
| 165 |
|
---|
[341014c] | 166 | x = particlePosition.X() * 1.0E-3;
|
---|
| 167 | y = particlePosition.Y() * 1.0E-3;
|
---|
| 168 | z = particlePosition.Z() * 1.0E-3;
|
---|
[e55f5b0] | 169 |
|
---|
[341014c] | 170 | bsx = beamSpotPosition.X() * 1.0E-3;
|
---|
| 171 | bsy = beamSpotPosition.Y() * 1.0E-3;
|
---|
| 172 | bsz = beamSpotPosition.Z() * 1.0E-3;
|
---|
[e55f5b0] | 173 |
|
---|
[60e1de6] | 174 | q = particle->Charge;
|
---|
[d7d2da3] | 175 |
|
---|
| 176 | // check that particle position is inside the cylinder
|
---|
[9330b7b] | 177 | if(TMath::Hypot(x, y) > fRadiusMax || TMath::Abs(z) > fHalfLengthMax)
|
---|
[d7d2da3] | 178 | {
|
---|
| 179 | continue;
|
---|
| 180 | }
|
---|
| 181 |
|
---|
[60e1de6] | 182 | px = particleMomentum.Px();
|
---|
| 183 | py = particleMomentum.Py();
|
---|
| 184 | pz = particleMomentum.Pz();
|
---|
| 185 | pt = particleMomentum.Pt();
|
---|
| 186 | pt2 = particleMomentum.Perp2();
|
---|
| 187 | e = particleMomentum.E();
|
---|
[d7d2da3] | 188 |
|
---|
| 189 | if(pt2 < 1.0E-9)
|
---|
| 190 | {
|
---|
| 191 | continue;
|
---|
| 192 | }
|
---|
| 193 |
|
---|
[9330b7b] | 194 | if(TMath::Hypot(x, y) > fRadius || TMath::Abs(z) > fHalfLength)
|
---|
| 195 | {
|
---|
| 196 | mother = candidate;
|
---|
[341014c] | 197 | candidate = static_cast<Candidate *>(candidate->Clone());
|
---|
[9330b7b] | 198 |
|
---|
[60e1de6] | 199 | candidate->InitialPosition = particlePosition;
|
---|
| 200 | candidate->Position = particlePosition;
|
---|
[9330b7b] | 201 | candidate->L = 0.0;
|
---|
| 202 |
|
---|
[60e1de6] | 203 | candidate->Momentum = particleMomentum;
|
---|
[9330b7b] | 204 | candidate->AddCandidate(mother);
|
---|
| 205 |
|
---|
| 206 | fOutputArray->Add(candidate);
|
---|
| 207 | }
|
---|
| 208 | else if(TMath::Abs(q) < 1.0E-9 || TMath::Abs(fBz) < 1.0E-9)
|
---|
[d7d2da3] | 209 | {
|
---|
[d41ba4a] | 210 |
|
---|
[38b4e15] | 211 | rxp = x*py - y*px;
|
---|
| 212 | rdp = x*px + y*py;
|
---|
[d7d2da3] | 213 |
|
---|
[38b4e15] | 214 | discr = fRadius*fRadius*pt*pt - rxp*rxp;
|
---|
[d7d2da3] | 215 |
|
---|
[38b4e15] | 216 | t_R = e * (sqrt(discr) - rdp) / (c_light * pt * pt);
|
---|
| 217 | t_z = e * (TMath::Sign(fHalfLengthMax, pz) - z) / ( c_light * pz);
|
---|
| 218 |
|
---|
| 219 | t = TMath::Min(t_R, t_z);
|
---|
[d7d2da3] | 220 |
|
---|
[38b4e15] | 221 | x_t = x + px*t*c_light/e;
|
---|
| 222 | y_t = y + py*t*c_light/e;
|
---|
| 223 | z_t = z + pz*t*c_light/e;
|
---|
| 224 | r_t = TMath::Hypot(x_t, y_t);
|
---|
[e55f5b0] | 225 |
|
---|
[38b4e15] | 226 | l = TMath::Sqrt( (x_t - x)*(x_t - x) + (y_t - y)*(y_t - y) + (z_t - z)*(z_t - z));
|
---|
[d7d2da3] | 227 |
|
---|
| 228 | mother = candidate;
|
---|
[38b4e15] | 229 | candidate = static_cast<Candidate*>(candidate->Clone());
|
---|
[d7d2da3] | 230 |
|
---|
[60e1de6] | 231 | candidate->InitialPosition = particlePosition;
|
---|
[38b4e15] | 232 | candidate->Position.SetXYZT(x_t*1.0E3, y_t*1.0E3, z_t*1.0E3, particlePosition.T() + t*c_light*1.0E3);
|
---|
| 233 | candidate->L = l*1.0E3;
|
---|
[e55f5b0] | 234 |
|
---|
[60e1de6] | 235 | candidate->Momentum = particleMomentum;
|
---|
[d7d2da3] | 236 | candidate->AddCandidate(mother);
|
---|
[d41ba4a] | 237 |
|
---|
[d7d2da3] | 238 | fOutputArray->Add(candidate);
|
---|
[38b4e15] | 239 |
|
---|
[d41ba4a] | 240 | if(TMath::Abs(q) > 1.0E-9)
|
---|
[d7d2da3] | 241 | {
|
---|
| 242 | switch(TMath::Abs(candidate->PID))
|
---|
| 243 | {
|
---|
[341014c] | 244 | case 11:
|
---|
| 245 | fElectronOutputArray->Add(candidate);
|
---|
| 246 | break;
|
---|
| 247 | case 13:
|
---|
| 248 | fMuonOutputArray->Add(candidate);
|
---|
| 249 | break;
|
---|
| 250 | default:
|
---|
| 251 | fChargedHadronOutputArray->Add(candidate);
|
---|
[d7d2da3] | 252 | }
|
---|
| 253 | }
|
---|
[751cb9c] | 254 | else
|
---|
| 255 | {
|
---|
[341014c] | 256 | fNeutralOutputArray->Add(candidate);
|
---|
[751cb9c] | 257 | }
|
---|
[d7d2da3] | 258 | }
|
---|
| 259 | else
|
---|
| 260 | {
|
---|
| 261 |
|
---|
[b594101] | 262 | // 1. initial transverse momentum p_{T0}: Part->pt
|
---|
| 263 | // initial transverse momentum direction phi_0 = -atan(p_X0/p_Y0)
|
---|
| 264 | // relativistic gamma: gamma = E/mc^2; gammam = gamma * m
|
---|
| 265 | // gyration frequency omega = q/(gamma m) fBz
|
---|
| 266 | // helix radius r = p_{T0} / (omega gamma m)
|
---|
[d7d2da3] | 267 |
|
---|
[38b4e15] | 268 | gammam = e*1.0E9 / (c_light*c_light); // gammam in [eV/c^2]
|
---|
| 269 | omega = q * fBz / (gammam); // omega is here in [89875518/s]
|
---|
| 270 | r = pt / (q * fBz) * 1.0E9/c_light; // in [m]
|
---|
[d7d2da3] | 271 |
|
---|
[b594101] | 272 | phi_0 = TMath::ATan2(py, px); // [rad] in [-pi, pi]
|
---|
[d7d2da3] | 273 |
|
---|
| 274 | // 2. helix axis coordinates
|
---|
[38b4e15] | 275 | x_c = x + r*TMath::Sin(phi_0);
|
---|
| 276 | y_c = y - r*TMath::Cos(phi_0);
|
---|
[d7d2da3] | 277 | r_c = TMath::Hypot(x_c, y_c);
|
---|
[38b4e15] | 278 | phi_c = TMath::ATan(y_c/x_c);
|
---|
| 279 | if(x_c < 0.0) phi_c -= TMath::Sign(1., phi_c)*TMath::Pi();
|
---|
[e55f5b0] | 280 |
|
---|
[38b4e15] | 281 | //Find the time of closest approach
|
---|
| 282 | td = (phi_0 - TMath::ATan(-x_c/y_c))/omega;
|
---|
[e55f5b0] | 283 |
|
---|
[38b4e15] | 284 | //Remove all the modulo pi that might have come from the atan
|
---|
| 285 | pio = fabs(TMath::Pi()/omega);
|
---|
| 286 | while(fabs(td) > 0.5*pio)
|
---|
| 287 | {
|
---|
| 288 | td -= TMath::Sign(1., td)*pio;
|
---|
| 289 | }
|
---|
| 290 |
|
---|
| 291 | //Compute the coordinate of closed approach to z axis
|
---|
| 292 | //if wants wtr beamline need to be changedto re-center with a traslation of the z axis
|
---|
| 293 | phid = phi_0 - omega*td;
|
---|
| 294 | xd = x_c - r*TMath::Sin(phid);
|
---|
| 295 | yd = y_c + r*TMath::Cos(phid);
|
---|
| 296 | zd = z + c_light*(pz/e)*td;
|
---|
| 297 |
|
---|
| 298 | //Compute momentum at closest approach (perigee??)
|
---|
| 299 | px = pt*TMath::Cos(phid);
|
---|
| 300 | py = pt*TMath::Sin(phid);
|
---|
| 301 |
|
---|
| 302 | particleMomentum.SetPtEtaPhiE(pt, particleMomentum.Eta(), phid, particleMomentum.E());
|
---|
| 303 |
|
---|
| 304 | // calculate additional track parameters (correct for beamspot position)
|
---|
| 305 | d0 = ((xd - bsx) * py - (yd - bsy) * px) / pt;
|
---|
| 306 | dz = zd - bsz;
|
---|
| 307 | ctgTheta = 1.0 / TMath::Tan (particleMomentum.Theta());
|
---|
[e55f5b0] | 308 |
|
---|
[d7d2da3] | 309 | // 3. time evaluation t = TMath::Min(t_r, t_z)
|
---|
| 310 | // t_r : time to exit from the sides
|
---|
| 311 | // t_z : time to exit from the front or the back
|
---|
[38b4e15] | 312 | t = 0;
|
---|
| 313 | t_z = 0;
|
---|
| 314 | sign_pz = (pz > 0.0) ? 1 : -1;
|
---|
| 315 | if(pz == 0.0) t_z = 1.0E99;
|
---|
| 316 | else t_z = gammam / (pz*1.0E9/c_light) * (-z + fHalfLength*sign_pz);
|
---|
[d7d2da3] | 317 |
|
---|
[38b4e15] | 318 | if(r_c + TMath::Abs(r) < fRadius) // helix does not cross the cylinder sides
|
---|
[d7d2da3] | 319 | {
|
---|
| 320 | t = t_z;
|
---|
| 321 | }
|
---|
| 322 | else
|
---|
| 323 | {
|
---|
[38b4e15] | 324 | alpha = -(fRadius*fRadius - r*r - r_c*r_c)/(2*fabs(r)*r_c);
|
---|
| 325 | alpha = fabs(TMath::ACos(alpha));
|
---|
| 326 | t_r = td + alpha/fabs(omega);
|
---|
| 327 |
|
---|
[d41ba4a] | 328 | t = TMath::Min(t_r, t_z);
|
---|
[d7d2da3] | 329 | }
|
---|
| 330 |
|
---|
[38b4e15] | 331 | x_t = x_c - r*TMath::Sin(phi_0 - omega*t);
|
---|
| 332 | y_t = y_c + r*TMath::Cos(phi_0 - omega*t);
|
---|
| 333 | z_t = z + c_light*t*pz/e;
|
---|
| 334 | r_t = TMath::Hypot(x_t, y_t);
|
---|
[d7d2da3] | 335 |
|
---|
[187fc41] | 336 | // compute path length for an helix
|
---|
[38b4e15] | 337 | vz = pz*1.0E9 / c_light / gammam;
|
---|
| 338 | //lenght of the path from production to tracker
|
---|
| 339 | l = t * TMath::Sqrt(vz*vz + r*r*omega*omega);
|
---|
[e55f5b0] | 340 |
|
---|
[d7d2da3] | 341 | if(r_t > 0.0)
|
---|
| 342 | {
|
---|
[acd0621] | 343 | // store these variables before cloning
|
---|
[5b51d33] | 344 | if(particle == candidate)
|
---|
| 345 | {
|
---|
[341014c] | 346 | particle->D0 = d0 * 1.0E3;
|
---|
| 347 | particle->DZ = dz * 1.0E3;
|
---|
[38b4e15] | 348 | particle->P = particleMomentum.P();
|
---|
[5b51d33] | 349 | particle->PT = pt;
|
---|
| 350 | particle->CtgTheta = ctgTheta;
|
---|
[38b4e15] | 351 | particle->Phi = particleMomentum.Phi();
|
---|
[5b51d33] | 352 | }
|
---|
[e55f5b0] | 353 |
|
---|
[d7d2da3] | 354 | mother = candidate;
|
---|
[341014c] | 355 | candidate = static_cast<Candidate *>(candidate->Clone());
|
---|
[d7d2da3] | 356 |
|
---|
[60e1de6] | 357 | candidate->InitialPosition = particlePosition;
|
---|
[341014c] | 358 | candidate->Position.SetXYZT(x_t * 1.0E3, y_t * 1.0E3, z_t * 1.0E3, particlePosition.T() + t * c_light * 1.0E3);
|
---|
[d7d2da3] | 359 |
|
---|
[60e1de6] | 360 | candidate->Momentum = particleMomentum;
|
---|
[e55f5b0] | 361 |
|
---|
[341014c] | 362 | candidate->L = l * 1.0E3;
|
---|
[e55f5b0] | 363 |
|
---|
[341014c] | 364 | candidate->Xd = xd * 1.0E3;
|
---|
| 365 | candidate->Yd = yd * 1.0E3;
|
---|
| 366 | candidate->Zd = zd * 1.0E3;
|
---|
[b594101] | 367 |
|
---|
| 368 | candidate->AddCandidate(mother);
|
---|
[d7d2da3] | 369 |
|
---|
| 370 | fOutputArray->Add(candidate);
|
---|
| 371 | switch(TMath::Abs(candidate->PID))
|
---|
| 372 | {
|
---|
[341014c] | 373 | case 11:
|
---|
| 374 | fElectronOutputArray->Add(candidate);
|
---|
| 375 | break;
|
---|
| 376 | case 13:
|
---|
| 377 | fMuonOutputArray->Add(candidate);
|
---|
| 378 | break;
|
---|
| 379 | default:
|
---|
| 380 | fChargedHadronOutputArray->Add(candidate);
|
---|
[d7d2da3] | 381 | }
|
---|
| 382 | }
|
---|
| 383 | }
|
---|
| 384 | }
|
---|
| 385 | }
|
---|
| 386 |
|
---|
| 387 | //------------------------------------------------------------------------------
|
---|