[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;
|
---|
[ae93700] | 128 | Double_t x_c, y_c, r_c, phi_0;
|
---|
| 129 | Double_t x_t, y_t, z_t, r_t, phi_t;
|
---|
| 130 | Double_t t_r, t_z;
|
---|
| 131 | Double_t tmp;
|
---|
[38b4e15] | 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;
|
---|
[ae93700] | 136 | Double_t td, pio, phid, vz;
|
---|
[38b4e15] | 137 |
|
---|
[d7d2da3] | 138 | const Double_t c_light = 2.99792458E8;
|
---|
[e55f5b0] | 139 |
|
---|
[341014c] | 140 | if(!fBeamSpotInputArray || fBeamSpotInputArray->GetSize() == 0)
|
---|
[ae93700] | 141 | {
|
---|
[187fc41] | 142 | beamSpotPosition.SetXYZT(0.0, 0.0, 0.0, 0.0);
|
---|
[ae93700] | 143 | }
|
---|
[187fc41] | 144 | else
|
---|
| 145 | {
|
---|
[341014c] | 146 | Candidate &beamSpotCandidate = *((Candidate *)fBeamSpotInputArray->At(0));
|
---|
[187fc41] | 147 | beamSpotPosition = beamSpotCandidate.Position;
|
---|
| 148 | }
|
---|
[e55f5b0] | 149 |
|
---|
[d7d2da3] | 150 | fItInputArray->Reset();
|
---|
[341014c] | 151 | while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
|
---|
[d7d2da3] | 152 | {
|
---|
[60e1de6] | 153 | if(candidate->GetCandidates()->GetEntriesFast() == 0)
|
---|
| 154 | {
|
---|
| 155 | particle = candidate;
|
---|
| 156 | }
|
---|
| 157 | else
|
---|
| 158 | {
|
---|
[341014c] | 159 | particle = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
|
---|
[60e1de6] | 160 | }
|
---|
| 161 |
|
---|
| 162 | particlePosition = particle->Position;
|
---|
| 163 | particleMomentum = particle->Momentum;
|
---|
[38b4e15] | 164 |
|
---|
[341014c] | 165 | x = particlePosition.X() * 1.0E-3;
|
---|
| 166 | y = particlePosition.Y() * 1.0E-3;
|
---|
| 167 | z = particlePosition.Z() * 1.0E-3;
|
---|
[e55f5b0] | 168 |
|
---|
[341014c] | 169 | bsx = beamSpotPosition.X() * 1.0E-3;
|
---|
| 170 | bsy = beamSpotPosition.Y() * 1.0E-3;
|
---|
| 171 | bsz = beamSpotPosition.Z() * 1.0E-3;
|
---|
[e55f5b0] | 172 |
|
---|
[60e1de6] | 173 | q = particle->Charge;
|
---|
[d7d2da3] | 174 |
|
---|
| 175 | // check that particle position is inside the cylinder
|
---|
[9330b7b] | 176 | if(TMath::Hypot(x, y) > fRadiusMax || TMath::Abs(z) > fHalfLengthMax)
|
---|
[d7d2da3] | 177 | {
|
---|
| 178 | continue;
|
---|
| 179 | }
|
---|
| 180 |
|
---|
[60e1de6] | 181 | px = particleMomentum.Px();
|
---|
| 182 | py = particleMomentum.Py();
|
---|
| 183 | pz = particleMomentum.Pz();
|
---|
| 184 | pt = particleMomentum.Pt();
|
---|
| 185 | pt2 = particleMomentum.Perp2();
|
---|
| 186 | e = particleMomentum.E();
|
---|
[d7d2da3] | 187 |
|
---|
| 188 | if(pt2 < 1.0E-9)
|
---|
| 189 | {
|
---|
| 190 | continue;
|
---|
| 191 | }
|
---|
| 192 |
|
---|
[9330b7b] | 193 | if(TMath::Hypot(x, y) > fRadius || TMath::Abs(z) > fHalfLength)
|
---|
| 194 | {
|
---|
| 195 | mother = candidate;
|
---|
[341014c] | 196 | candidate = static_cast<Candidate *>(candidate->Clone());
|
---|
[9330b7b] | 197 |
|
---|
[60e1de6] | 198 | candidate->InitialPosition = particlePosition;
|
---|
| 199 | candidate->Position = particlePosition;
|
---|
[9330b7b] | 200 | candidate->L = 0.0;
|
---|
| 201 |
|
---|
[60e1de6] | 202 | candidate->Momentum = particleMomentum;
|
---|
[9330b7b] | 203 | candidate->AddCandidate(mother);
|
---|
| 204 |
|
---|
| 205 | fOutputArray->Add(candidate);
|
---|
| 206 | }
|
---|
| 207 | else if(TMath::Abs(q) < 1.0E-9 || TMath::Abs(fBz) < 1.0E-9)
|
---|
[d7d2da3] | 208 | {
|
---|
[ae93700] | 209 | // solve pt2*t^2 + 2*(px*x + py*y)*t - (fRadius2 - x*x - y*y) = 0
|
---|
| 210 | tmp = px * y - py * x;
|
---|
| 211 | t_r = (TMath::Sqrt(pt2 * fRadius2 - tmp * tmp) - px * x - py * y) / pt2;
|
---|
[d41ba4a] | 212 |
|
---|
[ae93700] | 213 | t_z = (TMath::Sign(fHalfLength, pz) - z) / pz;
|
---|
[d7d2da3] | 214 |
|
---|
[ae93700] | 215 | t = TMath::Min(t_r, t_z);
|
---|
[d7d2da3] | 216 |
|
---|
[ae93700] | 217 | x_t = x + px * t;
|
---|
| 218 | y_t = y + py * t;
|
---|
| 219 | z_t = z + pz * t;
|
---|
[38b4e15] | 220 |
|
---|
[ae93700] | 221 | l = TMath::Sqrt((x_t - x) * (x_t - x) + (y_t - y) * (y_t - y) + (z_t - z) * (z_t - z));
|
---|
[d7d2da3] | 222 |
|
---|
| 223 | mother = candidate;
|
---|
[ae93700] | 224 | candidate = static_cast<Candidate *>(candidate->Clone());
|
---|
[d7d2da3] | 225 |
|
---|
[60e1de6] | 226 | candidate->InitialPosition = particlePosition;
|
---|
[ae93700] | 227 | candidate->Position.SetXYZT(x_t * 1.0E3, y_t * 1.0E3, z_t * 1.0E3, particlePosition.T() + t * e * 1.0E3);
|
---|
| 228 | candidate->L = l * 1.0E3;
|
---|
[e55f5b0] | 229 |
|
---|
[60e1de6] | 230 | candidate->Momentum = particleMomentum;
|
---|
[d7d2da3] | 231 | candidate->AddCandidate(mother);
|
---|
[d41ba4a] | 232 |
|
---|
[d7d2da3] | 233 | fOutputArray->Add(candidate);
|
---|
[38b4e15] | 234 |
|
---|
[d41ba4a] | 235 | if(TMath::Abs(q) > 1.0E-9)
|
---|
[d7d2da3] | 236 | {
|
---|
| 237 | switch(TMath::Abs(candidate->PID))
|
---|
| 238 | {
|
---|
[341014c] | 239 | case 11:
|
---|
| 240 | fElectronOutputArray->Add(candidate);
|
---|
| 241 | break;
|
---|
| 242 | case 13:
|
---|
| 243 | fMuonOutputArray->Add(candidate);
|
---|
| 244 | break;
|
---|
| 245 | default:
|
---|
| 246 | fChargedHadronOutputArray->Add(candidate);
|
---|
[d7d2da3] | 247 | }
|
---|
| 248 | }
|
---|
[751cb9c] | 249 | else
|
---|
| 250 | {
|
---|
[341014c] | 251 | fNeutralOutputArray->Add(candidate);
|
---|
[751cb9c] | 252 | }
|
---|
[d7d2da3] | 253 | }
|
---|
| 254 | else
|
---|
| 255 | {
|
---|
| 256 |
|
---|
[ae93700] | 257 | // 1. initial transverse momentum p_{T0}: Part->pt
|
---|
| 258 | // initial transverse momentum direction phi_0 = -atan(p_{X0} / p_{Y0})
|
---|
| 259 | // relativistic gamma: gamma = E / mc^2; gammam = gamma * m
|
---|
| 260 | // gyration frequency omega = q * Bz / (gammam)
|
---|
| 261 | // helix radius r = p_{T0} / (omega * gammam)
|
---|
[d7d2da3] | 262 |
|
---|
[ae93700] | 263 | gammam = e * 1.0E9 / (c_light * c_light); // gammam in [eV/c^2]
|
---|
| 264 | omega = q * fBz / gammam; // omega is here in [89875518/s]
|
---|
| 265 | r = pt / (q * fBz) * 1.0E9 / c_light; // in [m]
|
---|
[d7d2da3] | 266 |
|
---|
[b594101] | 267 | phi_0 = TMath::ATan2(py, px); // [rad] in [-pi, pi]
|
---|
[d7d2da3] | 268 |
|
---|
| 269 | // 2. helix axis coordinates
|
---|
[ae93700] | 270 | x_c = x + r * TMath::Sin(phi_0);
|
---|
| 271 | y_c = y - r * TMath::Cos(phi_0);
|
---|
[d7d2da3] | 272 | r_c = TMath::Hypot(x_c, y_c);
|
---|
[e55f5b0] | 273 |
|
---|
[ae93700] | 274 | // time of closest approach
|
---|
| 275 | td = (phi_0 + TMath::ATan2(x_c, y_c)) / omega;
|
---|
[e55f5b0] | 276 |
|
---|
[ae93700] | 277 | // remove all the modulo pi that might have come from the atan
|
---|
| 278 | pio = TMath::Abs(TMath::Pi() / omega);
|
---|
| 279 | while(TMath::Abs(td) > 0.5 * pio)
|
---|
[38b4e15] | 280 | {
|
---|
[ae93700] | 281 | td -= TMath::Sign(1.0, td) * pio;
|
---|
[38b4e15] | 282 | }
|
---|
| 283 |
|
---|
[ae93700] | 284 | vz = pz * c_light / e;
|
---|
[38b4e15] | 285 |
|
---|
[ae93700] | 286 | // calculate coordinates of closest approach to z axis
|
---|
| 287 | phid = phi_0 - omega * td;
|
---|
| 288 | xd = x_c - r * TMath::Sin(phid);
|
---|
| 289 | yd = y_c + r * TMath::Cos(phid);
|
---|
| 290 | zd = z + vz * td;
|
---|
| 291 |
|
---|
| 292 | // momentum at closest approach
|
---|
| 293 | px = pt * TMath::Cos(phid);
|
---|
| 294 | py = pt * TMath::Sin(phid);
|
---|
[38b4e15] | 295 |
|
---|
| 296 | particleMomentum.SetPtEtaPhiE(pt, particleMomentum.Eta(), phid, particleMomentum.E());
|
---|
| 297 |
|
---|
| 298 | // calculate additional track parameters (correct for beamspot position)
|
---|
| 299 | d0 = ((xd - bsx) * py - (yd - bsy) * px) / pt;
|
---|
| 300 | dz = zd - bsz;
|
---|
[ae93700] | 301 | ctgTheta = 1.0 / TMath::Tan(particleMomentum.Theta());
|
---|
[e55f5b0] | 302 |
|
---|
[d7d2da3] | 303 | // 3. time evaluation t = TMath::Min(t_r, t_z)
|
---|
| 304 | // t_r : time to exit from the sides
|
---|
| 305 | // t_z : time to exit from the front or the back
|
---|
[ae93700] | 306 | t_z = (vz == 0.0) ? 1.0E99 : (TMath::Sign(fHalfLength, pz) - z) / vz;
|
---|
[d7d2da3] | 307 |
|
---|
[ae93700] | 308 | if(r_c + TMath::Abs(r) < fRadius)
|
---|
[d7d2da3] | 309 | {
|
---|
[ae93700] | 310 | // helix does not cross the cylinder sides
|
---|
[d7d2da3] | 311 | t = t_z;
|
---|
| 312 | }
|
---|
| 313 | else
|
---|
| 314 | {
|
---|
[ae93700] | 315 | alpha = TMath::ACos((r * r + r_c * r_c - fRadius * fRadius) / (2 * TMath::Abs(r) * r_c));
|
---|
| 316 | t_r = td + TMath::Abs(alpha / omega);
|
---|
[38b4e15] | 317 |
|
---|
[d41ba4a] | 318 | t = TMath::Min(t_r, t_z);
|
---|
[d7d2da3] | 319 | }
|
---|
| 320 |
|
---|
[ae93700] | 321 | // 4. position in terms of x(t), y(t), z(t)
|
---|
| 322 | phi_t = phi_0 - omega * t;
|
---|
| 323 | x_t = x_c - r * TMath::Sin(phi_t);
|
---|
| 324 | y_t = y_c + r * TMath::Cos(phi_t);
|
---|
| 325 | z_t = z + vz * t;
|
---|
| 326 | r_t = TMath::Hypot(x_t, y_t);
|
---|
[d7d2da3] | 327 |
|
---|
[ae93700] | 328 | // lenght of the path from production to tracker
|
---|
| 329 | l = t * TMath::Hypot(vz, r * omega);
|
---|
[e55f5b0] | 330 |
|
---|
[d7d2da3] | 331 | if(r_t > 0.0)
|
---|
| 332 | {
|
---|
[acd0621] | 333 | // store these variables before cloning
|
---|
[5b51d33] | 334 | if(particle == candidate)
|
---|
| 335 | {
|
---|
[341014c] | 336 | particle->D0 = d0 * 1.0E3;
|
---|
| 337 | particle->DZ = dz * 1.0E3;
|
---|
[38b4e15] | 338 | particle->P = particleMomentum.P();
|
---|
[5b51d33] | 339 | particle->PT = pt;
|
---|
| 340 | particle->CtgTheta = ctgTheta;
|
---|
[38b4e15] | 341 | particle->Phi = particleMomentum.Phi();
|
---|
[5b51d33] | 342 | }
|
---|
[e55f5b0] | 343 |
|
---|
[d7d2da3] | 344 | mother = candidate;
|
---|
[341014c] | 345 | candidate = static_cast<Candidate *>(candidate->Clone());
|
---|
[d7d2da3] | 346 |
|
---|
[60e1de6] | 347 | candidate->InitialPosition = particlePosition;
|
---|
[341014c] | 348 | candidate->Position.SetXYZT(x_t * 1.0E3, y_t * 1.0E3, z_t * 1.0E3, particlePosition.T() + t * c_light * 1.0E3);
|
---|
[d7d2da3] | 349 |
|
---|
[60e1de6] | 350 | candidate->Momentum = particleMomentum;
|
---|
[e55f5b0] | 351 |
|
---|
[341014c] | 352 | candidate->L = l * 1.0E3;
|
---|
[e55f5b0] | 353 |
|
---|
[341014c] | 354 | candidate->Xd = xd * 1.0E3;
|
---|
| 355 | candidate->Yd = yd * 1.0E3;
|
---|
| 356 | candidate->Zd = zd * 1.0E3;
|
---|
[b594101] | 357 |
|
---|
| 358 | candidate->AddCandidate(mother);
|
---|
[d7d2da3] | 359 |
|
---|
| 360 | fOutputArray->Add(candidate);
|
---|
| 361 | switch(TMath::Abs(candidate->PID))
|
---|
| 362 | {
|
---|
[341014c] | 363 | case 11:
|
---|
| 364 | fElectronOutputArray->Add(candidate);
|
---|
| 365 | break;
|
---|
| 366 | case 13:
|
---|
| 367 | fMuonOutputArray->Add(candidate);
|
---|
| 368 | break;
|
---|
| 369 | default:
|
---|
| 370 | fChargedHadronOutputArray->Add(candidate);
|
---|
[d7d2da3] | 371 | }
|
---|
| 372 | }
|
---|
| 373 | }
|
---|
| 374 | }
|
---|
| 375 | }
|
---|
| 376 |
|
---|
| 377 | //------------------------------------------------------------------------------
|
---|