Fork me on GitHub

source: git/modules/ParticlePropagator.cc@ e0f8f99

ImprovedOutputFile Timing dual_readout llp
Last change on this file since e0f8f99 was 3f7ef54, checked in by Michele Selvaggi <michele.selvaggi@…>, 9 years ago

fixed typo in dZ formula

  • Property mode set to 100644
File size: 11.4 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
89 // import array with output from filter/classifier module
90
91 fInputArray = ImportArray(GetString("InputArray", "Delphes/stableParticles"));
92 fItInputArray = fInputArray->MakeIterator();
93
[187fc41]94 // import beamspot
[bff2e33]95 try
96 {
97 fBeamSpotInputArray = ImportArray(GetString("BeamSpotInputArray", "BeamSpotFilter/beamSpotParticle"));
98 }
99 catch(runtime_error &e)
100 {
101 fBeamSpotInputArray = 0;
102 }
[d7d2da3]103 // create output arrays
104
105 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
106 fChargedHadronOutputArray = ExportArray(GetString("ChargedHadronOutputArray", "chargedHadrons"));
107 fElectronOutputArray = ExportArray(GetString("ElectronOutputArray", "electrons"));
108 fMuonOutputArray = ExportArray(GetString("MuonOutputArray", "muons"));
109}
110
111//------------------------------------------------------------------------------
112
113void ParticlePropagator::Finish()
114{
115 if(fItInputArray) delete fItInputArray;
116}
117
118//------------------------------------------------------------------------------
119
120void ParticlePropagator::Process()
121{
122 Candidate *candidate, *mother;
[187fc41]123 TLorentzVector candidatePosition, candidateMomentum, beamSpotPosition;
[d7d2da3]124 Double_t px, py, pz, pt, pt2, e, q;
125 Double_t x, y, z, t, r, phi;
126 Double_t x_c, y_c, r_c, phi_c, phi_0;
127 Double_t x_t, y_t, z_t, r_t;
128 Double_t t1, t2, t3, t4, t5, t6;
129 Double_t t_z, t_r, t_ra, t_rb;
130 Double_t tmp, discr, discr2;
131 Double_t delta, gammam, omega, asinrho;
[187fc41]132 Double_t rcu, rc2, xd, yd, zd;
133 Double_t l, d0, dz, p, ctgTheta, phip, etap, alpha;
134 Double_t bsx, bsy, bsz;
135
[d7d2da3]136 const Double_t c_light = 2.99792458E8;
[187fc41]137
[bff2e33]138 if (!fBeamSpotInputArray || fBeamSpotInputArray->GetSize () == 0)
[187fc41]139 beamSpotPosition.SetXYZT(0.0, 0.0, 0.0, 0.0);
140 else
141 {
142 Candidate &beamSpotCandidate = *((Candidate *) fBeamSpotInputArray->At(0));
143 beamSpotPosition = beamSpotCandidate.Position;
144 }
145
[d7d2da3]146 fItInputArray->Reset();
147 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
148 {
149 candidatePosition = candidate->Position;
150 candidateMomentum = candidate->Momentum;
151 x = candidatePosition.X()*1.0E-3;
152 y = candidatePosition.Y()*1.0E-3;
153 z = candidatePosition.Z()*1.0E-3;
[187fc41]154
155 bsx = beamSpotPosition.X()*1.0E-3;
156 bsy = beamSpotPosition.Y()*1.0E-3;
157 bsz = beamSpotPosition.Z()*1.0E-3;
158
[d7d2da3]159 q = candidate->Charge;
160
161 // check that particle position is inside the cylinder
162 if(TMath::Hypot(x, y) > fRadius || TMath::Abs(z) > fHalfLength)
163 {
164 continue;
165 }
166
167 px = candidateMomentum.Px();
168 py = candidateMomentum.Py();
169 pz = candidateMomentum.Pz();
170 pt = candidateMomentum.Pt();
171 pt2 = candidateMomentum.Perp2();
172 e = candidateMomentum.E();
173
174 if(pt2 < 1.0E-9)
175 {
176 continue;
177 }
178
179 if(TMath::Abs(q) < 1.0E-9 || TMath::Abs(fBz) < 1.0E-9)
180 {
[fcdb8bc]181 // solve pt2*t^2 + 2*(px*x + py*y)*t - (fRadius2 - x*x - y*y) = 0
[d7d2da3]182 tmp = px*y - py*x;
183 discr2 = pt2*fRadius2 - tmp*tmp;
[d41ba4a]184
[b594101]185 if(discr2 < 0.0)
[d7d2da3]186 {
187 // no solutions
188 continue;
189 }
190
191 tmp = px*x + py*y;
192 discr = TMath::Sqrt(discr2);
193 t1 = (-tmp + discr)/pt2;
194 t2 = (-tmp - discr)/pt2;
[b594101]195 t = (t1 < 0.0) ? t2 : t1;
[d7d2da3]196
197 z_t = z + pz*t;
198 if(TMath::Abs(z_t) > fHalfLength)
199 {
200 t3 = (+fHalfLength - z) / pz;
201 t4 = (-fHalfLength - z) / pz;
[b594101]202 t = (t3 < 0.0) ? t4 : t3;
[d7d2da3]203 }
204
205 x_t = x + px*t;
206 y_t = y + py*t;
207 z_t = z + pz*t;
[187fc41]208
209 l = TMath::Sqrt( (x_t - x)*(x_t - x) + (y_t - y)*(y_t - y) + (z_t - z)*(z_t - z));
[d7d2da3]210
211 mother = candidate;
212 candidate = static_cast<Candidate*>(candidate->Clone());
213
[3f7ef54]214 candidate->InitialPosition = candidatePosition;
[d41ba4a]215 candidate->Position.SetXYZT(x_t*1.0E3, y_t*1.0E3, z_t*1.0E3, candidatePosition.T() + t*e*1.0E3);
[187fc41]216 candidate->L = l*1.0E3;
217
[d7d2da3]218 candidate->Momentum = candidateMomentum;
219 candidate->AddCandidate(mother);
[d41ba4a]220
[d7d2da3]221 fOutputArray->Add(candidate);
[d41ba4a]222 if(TMath::Abs(q) > 1.0E-9)
[d7d2da3]223 {
224 switch(TMath::Abs(candidate->PID))
225 {
226 case 11:
227 fElectronOutputArray->Add(candidate);
228 break;
229 case 13:
230 fMuonOutputArray->Add(candidate);
231 break;
232 default:
233 fChargedHadronOutputArray->Add(candidate);
234 }
235 }
236 }
237 else
238 {
239
[b594101]240 // 1. initial transverse momentum p_{T0}: Part->pt
241 // initial transverse momentum direction phi_0 = -atan(p_X0/p_Y0)
242 // relativistic gamma: gamma = E/mc^2; gammam = gamma * m
243 // gyration frequency omega = q/(gamma m) fBz
244 // helix radius r = p_{T0} / (omega gamma m)
[d7d2da3]245
[b594101]246 gammam = e*1.0E9 / (c_light*c_light); // gammam in [eV/c^2]
247 omega = q * fBz / (gammam); // omega is here in [89875518/s]
[d41ba4a]248 r = pt / (q * fBz) * 1.0E9/c_light; // in [m]
[d7d2da3]249
[b594101]250 phi_0 = TMath::ATan2(py, px); // [rad] in [-pi, pi]
[d7d2da3]251
252 // 2. helix axis coordinates
253 x_c = x + r*TMath::Sin(phi_0);
254 y_c = y - r*TMath::Cos(phi_0);
255 r_c = TMath::Hypot(x_c, y_c);
256 phi_c = TMath::ATan2(y_c, x_c);
257 phi = phi_c;
258 if(x_c < 0.0) phi += TMath::Pi();
259
[a0431dc]260 rcu = TMath::Abs(r);
261 rc2 = r_c*r_c;
[b594101]262
[a0431dc]263 // calculate coordinates of closest approach to track circle in transverse plane xd, yd, zd
[b594101]264 xd = x_c*x_c*x_c - x_c*rcu*r_c + x_c*y_c*y_c;
265 xd = (rc2 > 0.0) ? xd / rc2 : -999;
266 yd = y_c*(-rcu*r_c + rc2);
267 yd = (rc2 > 0.0) ? yd / rc2 : -999;
268 zd = z + (TMath::Sqrt(xd*xd + yd*yd) - TMath::Sqrt(x*x + y*y))*pz/pt;
[187fc41]269
270 // use perigee momentum rather than original particle
271 // momentum, since the orignal particle momentum isn't known
272
273 px = TMath::Sign(1.0,r) * pt * (-y_c / r_c);
274 py = TMath::Sign(1.0,r) * pt * (x_c / r_c);
275 etap = candidateMomentum.Eta();
276 phip = TMath::ATan2(py, px);
277
278 candidateMomentum.SetPtEtaPhiE(pt, etap, phip, candidateMomentum.E());
279
280 // calculate additional track parameters (correct for beamspot position)
281
282 d0 = ( (x - bsx) * py - (y - bsy) * px) / pt;
[3f7ef54]283 dz = z - ((x - bsx) * px + (y - bsy) * py) / pt * (pz / pt);
[187fc41]284 p = candidateMomentum.P();
285 ctgTheta = 1.0 / TMath::Tan (candidateMomentum.Theta ());
[acd0621]286
[d7d2da3]287 // 3. time evaluation t = TMath::Min(t_r, t_z)
288 // t_r : time to exit from the sides
289 // t_z : time to exit from the front or the back
290 t_r = 0.0; // in [ns]
291 int sign_pz = (pz > 0.0) ? 1 : -1;
292 if(pz == 0.0) t_z = 1.0E99;
293 else t_z = gammam / (pz*1.0E9/c_light) * (-z + fHalfLength*sign_pz);
294
295 if(r_c + TMath::Abs(r) < fRadius)
296 {
297 // helix does not cross the cylinder sides
298 t = t_z;
299 }
300 else
301 {
302 asinrho = TMath::ASin( (fRadius*fRadius - r_c*r_c - r*r) / (2*TMath::Abs(r)*r_c) );
303 delta = phi_0 - phi;
304 if(delta <-TMath::Pi()) delta += 2*TMath::Pi();
305 if(delta > TMath::Pi()) delta -= 2*TMath::Pi();
306 t1 = (delta + asinrho) / omega;
307 t2 = (delta + TMath::Pi() - asinrho) / omega;
308 t3 = (delta + TMath::Pi() + asinrho) / omega;
309 t4 = (delta - asinrho) / omega;
310 t5 = (delta - TMath::Pi() - asinrho) / omega;
311 t6 = (delta - TMath::Pi() + asinrho) / omega;
312
[b594101]313 if(t1 < 0.0) t1 = 1.0E99;
314 if(t2 < 0.0) t2 = 1.0E99;
315 if(t3 < 0.0) t3 = 1.0E99;
316 if(t4 < 0.0) t4 = 1.0E99;
317 if(t5 < 0.0) t5 = 1.0E99;
318 if(t6 < 0.0) t6 = 1.0E99;
[d7d2da3]319
320 t_ra = TMath::Min(t1, TMath::Min(t2, t3));
321 t_rb = TMath::Min(t4, TMath::Min(t5, t6));
322 t_r = TMath::Min(t_ra, t_rb);
[d41ba4a]323 t = TMath::Min(t_r, t_z);
[d7d2da3]324 }
325
326 // 4. position in terms of x(t), y(t), z(t)
327 x_t = x_c + r * TMath::Sin(omega * t - phi_0);
328 y_t = y_c + r * TMath::Cos(omega * t - phi_0);
329 z_t = z + pz*1.0E9 / c_light / gammam * t;
330 r_t = TMath::Hypot(x_t, y_t);
331
[187fc41]332
333 // compute path length for an helix
334
335 alpha = pz*1.0E9 / c_light / gammam;
336 l = t * TMath::Sqrt(alpha*alpha + r*r*omega*omega);
337
[d7d2da3]338 if(r_t > 0.0)
339 {
[acd0621]340
341 // store these variables before cloning
342
343 candidate->D0 = d0*1.0E3;
344 candidate->DZ = dz*1.0E3;
345 candidate->P = p;
346 candidate->PT = pt;
347 candidate->CtgTheta = ctgTheta;
348 candidate->Phi = phip;
349
[d7d2da3]350 mother = candidate;
351 candidate = static_cast<Candidate*>(candidate->Clone());
352
[3f7ef54]353 candidate->InitialPosition = candidatePosition;
[d41ba4a]354 candidate->Position.SetXYZT(x_t*1.0E3, y_t*1.0E3, z_t*1.0E3, candidatePosition.T() + t*c_light*1.0E3);
[d7d2da3]355
356 candidate->Momentum = candidateMomentum;
[187fc41]357
[acd0621]358 candidate->L = l*1.0E3;
359
360 candidate->Xd = xd*1.0E3;
[b594101]361 candidate->Yd = yd*1.0E3;
[a0431dc]362 candidate->Zd = zd*1.0E3;
[b594101]363
364 candidate->AddCandidate(mother);
[d7d2da3]365
366 fOutputArray->Add(candidate);
367 switch(TMath::Abs(candidate->PID))
368 {
369 case 11:
370 fElectronOutputArray->Add(candidate);
371 break;
372 case 13:
373 fMuonOutputArray->Add(candidate);
374 break;
375 default:
376 fChargedHadronOutputArray->Add(candidate);
377 }
378 }
379 }
380 }
381}
382
383//------------------------------------------------------------------------------
[a0431dc]384
Note: See TracBrowser for help on using the repository browser.