Fork me on GitHub

source: git/modules/ParticlePropagator.cc@ 3874a9c

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 3874a9c was 751cb9c, checked in by Michele Selvaggi <michele.selvaggi@…>, 7 years ago

added NeutralOutputArray in ParticlePropagator

  • Property mode set to 100644
File size: 12.0 KB
RevLine 
[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
20/** \class ParticlePropagator
21 *
22 * Propagates charged and neutral particles
[d41ba4a]23 * from a given vertex to a cylinder defined by its radius,
[d7d2da3]24 * its half-length, centered at (0,0,0) and with its axis
25 * oriented along the z-axis.
26 *
27 * \author P. Demin - UCL, Louvain-la-Neuve
28 *
29 */
30
31#include "modules/ParticlePropagator.h"
32
33#include "classes/DelphesClasses.h"
34#include "classes/DelphesFactory.h"
35#include "classes/DelphesFormula.h"
36
37#include "ExRootAnalysis/ExRootResult.h"
38#include "ExRootAnalysis/ExRootFilter.h"
39#include "ExRootAnalysis/ExRootClassifier.h"
40
41#include "TMath.h"
42#include "TString.h"
43#include "TFormula.h"
44#include "TRandom3.h"
45#include "TObjArray.h"
46#include "TDatabasePDG.h"
47#include "TLorentzVector.h"
48
[d41ba4a]49#include <algorithm>
[d7d2da3]50#include <stdexcept>
51#include <iostream>
52#include <sstream>
53
54using namespace std;
55
56//------------------------------------------------------------------------------
57
58ParticlePropagator::ParticlePropagator() :
59 fItInputArray(0)
60{
61}
62
63//------------------------------------------------------------------------------
64
65ParticlePropagator::~ParticlePropagator()
66{
67}
68
[48f7a77]69
[d7d2da3]70//------------------------------------------------------------------------------
71
72void ParticlePropagator::Init()
73{
74 fRadius = GetDouble("Radius", 1.0);
75 fRadius2 = fRadius*fRadius;
76 fHalfLength = GetDouble("HalfLength", 3.0);
77 fBz = GetDouble("Bz", 0.0);
78 if(fRadius < 1.0E-2)
[d41ba4a]79 {
[d7d2da3]80 cout << "ERROR: magnetic field radius is too low\n";
81 return;
82 }
83 if(fHalfLength < 1.0E-2)
84 {
85 cout << "ERROR: magnetic field length is too low\n";
86 return;
87 }
88
[9330b7b]89 fRadiusMax = GetDouble("RadiusMax", fRadius);
90 fHalfLengthMax = GetDouble("HalfLengthMax", fHalfLength);
91
[d7d2da3]92 // import array with output from filter/classifier module
93
94 fInputArray = ImportArray(GetString("InputArray", "Delphes/stableParticles"));
95 fItInputArray = fInputArray->MakeIterator();
96
[187fc41]97 // import beamspot
[bff2e33]98 try
99 {
100 fBeamSpotInputArray = ImportArray(GetString("BeamSpotInputArray", "BeamSpotFilter/beamSpotParticle"));
101 }
102 catch(runtime_error &e)
103 {
104 fBeamSpotInputArray = 0;
[e55f5b0]105 }
[d7d2da3]106 // create output arrays
107
108 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
[751cb9c]109 fNeutralOutputArray = ExportArray(GetString("NeutralOutputArray", "neutralParticles"));
[d7d2da3]110 fChargedHadronOutputArray = ExportArray(GetString("ChargedHadronOutputArray", "chargedHadrons"));
111 fElectronOutputArray = ExportArray(GetString("ElectronOutputArray", "electrons"));
112 fMuonOutputArray = ExportArray(GetString("MuonOutputArray", "muons"));
113}
114
115//------------------------------------------------------------------------------
116
117void ParticlePropagator::Finish()
118{
119 if(fItInputArray) delete fItInputArray;
120}
121
122//------------------------------------------------------------------------------
123
124void ParticlePropagator::Process()
125{
126 Candidate *candidate, *mother;
[187fc41]127 TLorentzVector candidatePosition, candidateMomentum, beamSpotPosition;
[d7d2da3]128 Double_t px, py, pz, pt, pt2, e, q;
129 Double_t x, y, z, t, r, phi;
130 Double_t x_c, y_c, r_c, phi_c, phi_0;
131 Double_t x_t, y_t, z_t, r_t;
132 Double_t t1, t2, t3, t4, t5, t6;
133 Double_t t_z, t_r, t_ra, t_rb;
134 Double_t tmp, discr, discr2;
135 Double_t delta, gammam, omega, asinrho;
[187fc41]136 Double_t rcu, rc2, xd, yd, zd;
137 Double_t l, d0, dz, p, ctgTheta, phip, etap, alpha;
138 Double_t bsx, bsy, bsz;
[e55f5b0]139
[d7d2da3]140 const Double_t c_light = 2.99792458E8;
[e55f5b0]141
142 if (!fBeamSpotInputArray || fBeamSpotInputArray->GetSize () == 0)
[187fc41]143 beamSpotPosition.SetXYZT(0.0, 0.0, 0.0, 0.0);
144 else
145 {
146 Candidate &beamSpotCandidate = *((Candidate *) fBeamSpotInputArray->At(0));
147 beamSpotPosition = beamSpotCandidate.Position;
148 }
[e55f5b0]149
[d7d2da3]150 fItInputArray->Reset();
151 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
152 {
153 candidatePosition = candidate->Position;
154 candidateMomentum = candidate->Momentum;
155 x = candidatePosition.X()*1.0E-3;
156 y = candidatePosition.Y()*1.0E-3;
157 z = candidatePosition.Z()*1.0E-3;
[e55f5b0]158
[187fc41]159 bsx = beamSpotPosition.X()*1.0E-3;
160 bsy = beamSpotPosition.Y()*1.0E-3;
161 bsz = beamSpotPosition.Z()*1.0E-3;
[e55f5b0]162
[d7d2da3]163 q = candidate->Charge;
164
165 // check that particle position is inside the cylinder
[9330b7b]166 if(TMath::Hypot(x, y) > fRadiusMax || TMath::Abs(z) > fHalfLengthMax)
[d7d2da3]167 {
168 continue;
169 }
170
171 px = candidateMomentum.Px();
172 py = candidateMomentum.Py();
173 pz = candidateMomentum.Pz();
174 pt = candidateMomentum.Pt();
175 pt2 = candidateMomentum.Perp2();
176 e = candidateMomentum.E();
177
178 if(pt2 < 1.0E-9)
179 {
180 continue;
181 }
182
[9330b7b]183 if(TMath::Hypot(x, y) > fRadius || TMath::Abs(z) > fHalfLength)
184 {
185 mother = candidate;
186 candidate = static_cast<Candidate*>(candidate->Clone());
187
188 candidate->InitialPosition = candidatePosition;
189 candidate->Position = candidatePosition;
190 candidate->L = 0.0;
191
192 candidate->Momentum = candidateMomentum;
193 candidate->AddCandidate(mother);
194
195 fOutputArray->Add(candidate);
196 }
197 else if(TMath::Abs(q) < 1.0E-9 || TMath::Abs(fBz) < 1.0E-9)
[d7d2da3]198 {
[fcdb8bc]199 // solve pt2*t^2 + 2*(px*x + py*y)*t - (fRadius2 - x*x - y*y) = 0
[d7d2da3]200 tmp = px*y - py*x;
201 discr2 = pt2*fRadius2 - tmp*tmp;
[d41ba4a]202
[b594101]203 if(discr2 < 0.0)
[d7d2da3]204 {
205 // no solutions
206 continue;
207 }
208
209 tmp = px*x + py*y;
210 discr = TMath::Sqrt(discr2);
211 t1 = (-tmp + discr)/pt2;
212 t2 = (-tmp - discr)/pt2;
[b594101]213 t = (t1 < 0.0) ? t2 : t1;
[d7d2da3]214
215 z_t = z + pz*t;
216 if(TMath::Abs(z_t) > fHalfLength)
217 {
218 t3 = (+fHalfLength - z) / pz;
219 t4 = (-fHalfLength - z) / pz;
[b594101]220 t = (t3 < 0.0) ? t4 : t3;
[d7d2da3]221 }
222
223 x_t = x + px*t;
224 y_t = y + py*t;
225 z_t = z + pz*t;
[e55f5b0]226
[187fc41]227 l = TMath::Sqrt( (x_t - x)*(x_t - x) + (y_t - y)*(y_t - y) + (z_t - z)*(z_t - z));
[d7d2da3]228
229 mother = candidate;
230 candidate = static_cast<Candidate*>(candidate->Clone());
231
[3f7ef54]232 candidate->InitialPosition = candidatePosition;
[d41ba4a]233 candidate->Position.SetXYZT(x_t*1.0E3, y_t*1.0E3, z_t*1.0E3, candidatePosition.T() + t*e*1.0E3);
[187fc41]234 candidate->L = l*1.0E3;
[e55f5b0]235
[d7d2da3]236 candidate->Momentum = candidateMomentum;
237 candidate->AddCandidate(mother);
[d41ba4a]238
[d7d2da3]239 fOutputArray->Add(candidate);
[d41ba4a]240 if(TMath::Abs(q) > 1.0E-9)
[d7d2da3]241 {
242 switch(TMath::Abs(candidate->PID))
243 {
244 case 11:
245 fElectronOutputArray->Add(candidate);
246 break;
247 case 13:
248 fMuonOutputArray->Add(candidate);
249 break;
250 default:
251 fChargedHadronOutputArray->Add(candidate);
252 }
253 }
[751cb9c]254 else
255 {
256 fNeutralOutputArray->Add(candidate);
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
[b594101]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]
[d41ba4a]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
275 x_c = x + r*TMath::Sin(phi_0);
276 y_c = y - r*TMath::Cos(phi_0);
277 r_c = TMath::Hypot(x_c, y_c);
278 phi_c = TMath::ATan2(y_c, x_c);
279 phi = phi_c;
280 if(x_c < 0.0) phi += TMath::Pi();
281
[a0431dc]282 rcu = TMath::Abs(r);
283 rc2 = r_c*r_c;
[b594101]284
[a0431dc]285 // calculate coordinates of closest approach to track circle in transverse plane xd, yd, zd
[b594101]286 xd = x_c*x_c*x_c - x_c*rcu*r_c + x_c*y_c*y_c;
287 xd = (rc2 > 0.0) ? xd / rc2 : -999;
288 yd = y_c*(-rcu*r_c + rc2);
289 yd = (rc2 > 0.0) ? yd / rc2 : -999;
290 zd = z + (TMath::Sqrt(xd*xd + yd*yd) - TMath::Sqrt(x*x + y*y))*pz/pt;
[e55f5b0]291
[187fc41]292 // use perigee momentum rather than original particle
[e55f5b0]293 // momentum, since the orignal particle momentum isn't known
294
[9330b7b]295 px = TMath::Sign(1.0, r) * pt * (-y_c / r_c);
296 py = TMath::Sign(1.0, r) * pt * (x_c / r_c);
[187fc41]297 etap = candidateMomentum.Eta();
298 phip = TMath::ATan2(py, px);
[e55f5b0]299
[187fc41]300 candidateMomentum.SetPtEtaPhiE(pt, etap, phip, candidateMomentum.E());
[e55f5b0]301
[187fc41]302 // calculate additional track parameters (correct for beamspot position)
[e55f5b0]303
[9330b7b]304 d0 = ((x - bsx) * py - (y - bsy) * px) / pt;
[3f7ef54]305 dz = z - ((x - bsx) * px + (y - bsy) * py) / pt * (pz / pt);
[187fc41]306 p = candidateMomentum.P();
[9330b7b]307 ctgTheta = 1.0 / TMath::Tan (candidateMomentum.Theta());
[e55f5b0]308
309
[d7d2da3]310 // 3. time evaluation t = TMath::Min(t_r, t_z)
311 // t_r : time to exit from the sides
312 // t_z : time to exit from the front or the back
313 t_r = 0.0; // in [ns]
314 int 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);
317
318 if(r_c + TMath::Abs(r) < fRadius)
319 {
320 // helix does not cross the cylinder sides
321 t = t_z;
322 }
323 else
324 {
[9330b7b]325 asinrho = TMath::ASin((fRadius*fRadius - r_c*r_c - r*r) / (2*TMath::Abs(r)*r_c));
[d7d2da3]326 delta = phi_0 - phi;
327 if(delta <-TMath::Pi()) delta += 2*TMath::Pi();
328 if(delta > TMath::Pi()) delta -= 2*TMath::Pi();
329 t1 = (delta + asinrho) / omega;
330 t2 = (delta + TMath::Pi() - asinrho) / omega;
331 t3 = (delta + TMath::Pi() + asinrho) / omega;
332 t4 = (delta - asinrho) / omega;
333 t5 = (delta - TMath::Pi() - asinrho) / omega;
334 t6 = (delta - TMath::Pi() + asinrho) / omega;
335
[b594101]336 if(t1 < 0.0) t1 = 1.0E99;
337 if(t2 < 0.0) t2 = 1.0E99;
338 if(t3 < 0.0) t3 = 1.0E99;
339 if(t4 < 0.0) t4 = 1.0E99;
340 if(t5 < 0.0) t5 = 1.0E99;
341 if(t6 < 0.0) t6 = 1.0E99;
[d7d2da3]342
343 t_ra = TMath::Min(t1, TMath::Min(t2, t3));
344 t_rb = TMath::Min(t4, TMath::Min(t5, t6));
345 t_r = TMath::Min(t_ra, t_rb);
[d41ba4a]346 t = TMath::Min(t_r, t_z);
[d7d2da3]347 }
348
349 // 4. position in terms of x(t), y(t), z(t)
350 x_t = x_c + r * TMath::Sin(omega * t - phi_0);
351 y_t = y_c + r * TMath::Cos(omega * t - phi_0);
352 z_t = z + pz*1.0E9 / c_light / gammam * t;
353 r_t = TMath::Hypot(x_t, y_t);
354
[e55f5b0]355
[187fc41]356 // compute path length for an helix
[e55f5b0]357
[187fc41]358 alpha = pz*1.0E9 / c_light / gammam;
359 l = t * TMath::Sqrt(alpha*alpha + r*r*omega*omega);
[e55f5b0]360
[d7d2da3]361 if(r_t > 0.0)
362 {
[e55f5b0]363
[acd0621]364 // store these variables before cloning
365 candidate->D0 = d0*1.0E3;
366 candidate->DZ = dz*1.0E3;
367 candidate->P = p;
368 candidate->PT = pt;
369 candidate->CtgTheta = ctgTheta;
370 candidate->Phi = phip;
[e55f5b0]371
[d7d2da3]372 mother = candidate;
373 candidate = static_cast<Candidate*>(candidate->Clone());
374
[3f7ef54]375 candidate->InitialPosition = candidatePosition;
[d41ba4a]376 candidate->Position.SetXYZT(x_t*1.0E3, y_t*1.0E3, z_t*1.0E3, candidatePosition.T() + t*c_light*1.0E3);
[d7d2da3]377
378 candidate->Momentum = candidateMomentum;
[e55f5b0]379
[9330b7b]380 candidate->L = l*1.0E3;
[e55f5b0]381
[9330b7b]382 candidate->Xd = xd*1.0E3;
[b594101]383 candidate->Yd = yd*1.0E3;
[a0431dc]384 candidate->Zd = zd*1.0E3;
[b594101]385
386 candidate->AddCandidate(mother);
[d7d2da3]387
388 fOutputArray->Add(candidate);
389 switch(TMath::Abs(candidate->PID))
390 {
391 case 11:
392 fElectronOutputArray->Add(candidate);
393 break;
394 case 13:
395 fMuonOutputArray->Add(candidate);
396 break;
397 default:
398 fChargedHadronOutputArray->Add(candidate);
399 }
400 }
401 }
402 }
403}
404
405//------------------------------------------------------------------------------
[a0431dc]406
Note: See TracBrowser for help on using the repository browser.