Fork me on GitHub

source: git/modules/ParticlePropagator.cc@ c6f9311

Last change on this file since c6f9311 was 38b4e15, checked in by michele <michele.selvaggi@…>, 4 years ago

fixed bug in propagation at low pT

  • Property mode set to 100644
File size: 11.2 KB
Line 
1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
4 *
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.
9 *
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.
14 *
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
19/** \class ParticlePropagator
20 *
21 * Propagates charged and neutral particles
22 * from a given vertex to a cylinder defined by its radius,
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"
37#include "ExRootAnalysis/ExRootFilter.h"
38#include "ExRootAnalysis/ExRootResult.h"
39
40#include "TDatabasePDG.h"
41#include "TFormula.h"
42#include "TLorentzVector.h"
43#include "TMath.h"
44#include "TObjArray.h"
45#include "TRandom3.h"
46#include "TString.h"
47
48#include <algorithm>
49#include <iostream>
50#include <sstream>
51#include <stdexcept>
52
53using namespace std;
54
55//------------------------------------------------------------------------------
56
57ParticlePropagator::ParticlePropagator() :
58 fItInputArray(0)
59{
60}
61
62//------------------------------------------------------------------------------
63
64ParticlePropagator::~ParticlePropagator()
65{
66}
67
68//------------------------------------------------------------------------------
69
70void ParticlePropagator::Init()
71{
72 fRadius = GetDouble("Radius", 1.0);
73 fRadius2 = fRadius * fRadius;
74 fHalfLength = GetDouble("HalfLength", 3.0);
75 fBz = GetDouble("Bz", 0.0);
76 if(fRadius < 1.0E-2)
77 {
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
87 fRadiusMax = GetDouble("RadiusMax", fRadius);
88 fHalfLengthMax = GetDouble("HalfLengthMax", fHalfLength);
89
90 // import array with output from filter/classifier module
91
92 fInputArray = ImportArray(GetString("InputArray", "Delphes/stableParticles"));
93 fItInputArray = fInputArray->MakeIterator();
94
95 // import beamspot
96 try
97 {
98 fBeamSpotInputArray = ImportArray(GetString("BeamSpotInputArray", "BeamSpotFilter/beamSpotParticle"));
99 }
100 catch(runtime_error &e)
101 {
102 fBeamSpotInputArray = 0;
103 }
104 // create output arrays
105
106 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
107 fNeutralOutputArray = ExportArray(GetString("NeutralOutputArray", "neutralParticles"));
108 fChargedHadronOutputArray = ExportArray(GetString("ChargedHadronOutputArray", "chargedHadrons"));
109 fElectronOutputArray = ExportArray(GetString("ElectronOutputArray", "electrons"));
110 fMuonOutputArray = ExportArray(GetString("MuonOutputArray", "muons"));
111}
112
113//------------------------------------------------------------------------------
114
115void ParticlePropagator::Finish()
116{
117 if(fItInputArray) delete fItInputArray;
118}
119
120//------------------------------------------------------------------------------
121
122void ParticlePropagator::Process()
123{
124 Candidate *candidate, *mother, *particle;
125 TLorentzVector particlePosition, particleMomentum, beamSpotPosition;
126 Double_t px, py, pz, pt, pt2, e, q;
127 Double_t x, y, z, t, r;
128 Double_t x_c, y_c, r_c, phi_c, phi_0;
129 Double_t x_t, y_t, z_t, r_t;
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;
135 Double_t bsx, bsy, bsz;
136 Double_t rxp, rdp, t_R;
137 Double_t td, pio, phid, sign_pz, vz;
138
139 const Double_t c_light = 2.99792458E8;
140
141 if(!fBeamSpotInputArray || fBeamSpotInputArray->GetSize() == 0)
142 beamSpotPosition.SetXYZT(0.0, 0.0, 0.0, 0.0);
143 else
144 {
145 Candidate &beamSpotCandidate = *((Candidate *)fBeamSpotInputArray->At(0));
146 beamSpotPosition = beamSpotCandidate.Position;
147 }
148
149 fItInputArray->Reset();
150 while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
151 {
152 if(candidate->GetCandidates()->GetEntriesFast() == 0)
153 {
154 particle = candidate;
155 }
156 else
157 {
158 particle = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
159 }
160
161 particlePosition = particle->Position;
162 particleMomentum = particle->Momentum;
163
164 // Constants
165
166 x = particlePosition.X() * 1.0E-3;
167 y = particlePosition.Y() * 1.0E-3;
168 z = particlePosition.Z() * 1.0E-3;
169
170 bsx = beamSpotPosition.X() * 1.0E-3;
171 bsy = beamSpotPosition.Y() * 1.0E-3;
172 bsz = beamSpotPosition.Z() * 1.0E-3;
173
174 q = particle->Charge;
175
176 // check that particle position is inside the cylinder
177 if(TMath::Hypot(x, y) > fRadiusMax || TMath::Abs(z) > fHalfLengthMax)
178 {
179 continue;
180 }
181
182 px = particleMomentum.Px();
183 py = particleMomentum.Py();
184 pz = particleMomentum.Pz();
185 pt = particleMomentum.Pt();
186 pt2 = particleMomentum.Perp2();
187 e = particleMomentum.E();
188
189 if(pt2 < 1.0E-9)
190 {
191 continue;
192 }
193
194 if(TMath::Hypot(x, y) > fRadius || TMath::Abs(z) > fHalfLength)
195 {
196 mother = candidate;
197 candidate = static_cast<Candidate *>(candidate->Clone());
198
199 candidate->InitialPosition = particlePosition;
200 candidate->Position = particlePosition;
201 candidate->L = 0.0;
202
203 candidate->Momentum = particleMomentum;
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)
209 {
210
211 rxp = x*py - y*px;
212 rdp = x*px + y*py;
213
214 discr = fRadius*fRadius*pt*pt - rxp*rxp;
215
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);
220
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);
225
226 l = TMath::Sqrt( (x_t - x)*(x_t - x) + (y_t - y)*(y_t - y) + (z_t - z)*(z_t - z));
227
228 mother = candidate;
229 candidate = static_cast<Candidate*>(candidate->Clone());
230
231 candidate->InitialPosition = particlePosition;
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;
234
235 candidate->Momentum = particleMomentum;
236 candidate->AddCandidate(mother);
237
238 fOutputArray->Add(candidate);
239
240 if(TMath::Abs(q) > 1.0E-9)
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 }
254 else
255 {
256 fNeutralOutputArray->Add(candidate);
257 }
258 }
259 else
260 {
261
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)
267
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]
271
272 phi_0 = TMath::ATan2(py, px); // [rad] in [-pi, pi]
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::ATan(y_c/x_c);
279 if(x_c < 0.0) phi_c -= TMath::Sign(1., phi_c)*TMath::Pi();
280
281 //Find the time of closest approach
282 td = (phi_0 - TMath::ATan(-x_c/y_c))/omega;
283
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());
308
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
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);
317
318 if(r_c + TMath::Abs(r) < fRadius) // helix does not cross the cylinder sides
319 {
320 t = t_z;
321 }
322 else
323 {
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
328 t = TMath::Min(t_r, t_z);
329 }
330
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);
335
336 // compute path length for an helix
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);
340
341 if(r_t > 0.0)
342 {
343 // store these variables before cloning
344 if(particle == candidate)
345 {
346 particle->D0 = d0 * 1.0E3;
347 particle->DZ = dz * 1.0E3;
348 particle->P = particleMomentum.P();
349 particle->PT = pt;
350 particle->CtgTheta = ctgTheta;
351 particle->Phi = particleMomentum.Phi();
352 }
353
354 mother = candidate;
355 candidate = static_cast<Candidate *>(candidate->Clone());
356
357 candidate->InitialPosition = particlePosition;
358 candidate->Position.SetXYZT(x_t * 1.0E3, y_t * 1.0E3, z_t * 1.0E3, particlePosition.T() + t * c_light * 1.0E3);
359
360 candidate->Momentum = particleMomentum;
361
362 candidate->L = l * 1.0E3;
363
364 candidate->Xd = xd * 1.0E3;
365 candidate->Yd = yd * 1.0E3;
366 candidate->Zd = zd * 1.0E3;
367
368 candidate->AddCandidate(mother);
369
370 fOutputArray->Add(candidate);
371 switch(TMath::Abs(candidate->PID))
372 {
373 case 11:
374 fElectronOutputArray->Add(candidate);
375 break;
376 case 13:
377 fMuonOutputArray->Add(candidate);
378 break;
379 default:
380 fChargedHadronOutputArray->Add(candidate);
381 }
382 }
383 }
384 }
385}
386
387//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.