Fork me on GitHub

source: git/modules/ParticlePropagator.cc@ cb46568

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

added BeamSpot array in propagator, and calculation of track parameters

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