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 |
|
---|
54 | using namespace std;
|
---|
55 |
|
---|
56 | //------------------------------------------------------------------------------
|
---|
57 |
|
---|
58 | ParticlePropagator::ParticlePropagator() :
|
---|
59 | fItInputArray(0)
|
---|
60 | {
|
---|
61 | }
|
---|
62 |
|
---|
63 | //------------------------------------------------------------------------------
|
---|
64 |
|
---|
65 | ParticlePropagator::~ParticlePropagator()
|
---|
66 | {
|
---|
67 | }
|
---|
68 |
|
---|
69 |
|
---|
70 | //------------------------------------------------------------------------------
|
---|
71 |
|
---|
72 | void 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 | try
|
---|
96 | {
|
---|
97 | fBeamSpotInputArray = ImportArray(GetString("BeamSpotInputArray", "BeamSpotFilter/beamSpotParticle"));
|
---|
98 | }
|
---|
99 | catch(runtime_error &e)
|
---|
100 | {
|
---|
101 | fBeamSpotInputArray = 0;
|
---|
102 | }
|
---|
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 |
|
---|
113 | void ParticlePropagator::Finish()
|
---|
114 | {
|
---|
115 | if(fItInputArray) delete fItInputArray;
|
---|
116 | }
|
---|
117 |
|
---|
118 | //------------------------------------------------------------------------------
|
---|
119 |
|
---|
120 | void ParticlePropagator::Process()
|
---|
121 | {
|
---|
122 | Candidate *candidate, *mother;
|
---|
123 | TLorentzVector candidatePosition, candidateMomentum, beamSpotPosition;
|
---|
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;
|
---|
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 |
|
---|
136 | const Double_t c_light = 2.99792458E8;
|
---|
137 |
|
---|
138 | if (!fBeamSpotInputArray || fBeamSpotInputArray->GetSize () == 0)
|
---|
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 |
|
---|
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;
|
---|
154 |
|
---|
155 | bsx = beamSpotPosition.X()*1.0E-3;
|
---|
156 | bsy = beamSpotPosition.Y()*1.0E-3;
|
---|
157 | bsz = beamSpotPosition.Z()*1.0E-3;
|
---|
158 |
|
---|
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 | {
|
---|
181 | // solve pt2*t^2 + 2*(px*x + py*y)*t - (fRadius2 - x*x - y*y) = 0
|
---|
182 | tmp = px*y - py*x;
|
---|
183 | discr2 = pt2*fRadius2 - tmp*tmp;
|
---|
184 |
|
---|
185 | if(discr2 < 0.0)
|
---|
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;
|
---|
195 | t = (t1 < 0.0) ? t2 : t1;
|
---|
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;
|
---|
202 | t = (t3 < 0.0) ? t4 : t3;
|
---|
203 | }
|
---|
204 |
|
---|
205 | x_t = x + px*t;
|
---|
206 | y_t = y + py*t;
|
---|
207 | z_t = z + pz*t;
|
---|
208 |
|
---|
209 | l = TMath::Sqrt( (x_t - x)*(x_t - x) + (y_t - y)*(y_t - y) + (z_t - z)*(z_t - z));
|
---|
210 |
|
---|
211 | mother = candidate;
|
---|
212 | candidate = static_cast<Candidate*>(candidate->Clone());
|
---|
213 |
|
---|
214 | candidate->InitialPosition = candidatePosition;
|
---|
215 | candidate->Position.SetXYZT(x_t*1.0E3, y_t*1.0E3, z_t*1.0E3, candidatePosition.T() + t*e*1.0E3);
|
---|
216 | candidate->L = l*1.0E3;
|
---|
217 |
|
---|
218 | candidate->Momentum = candidateMomentum;
|
---|
219 | candidate->AddCandidate(mother);
|
---|
220 |
|
---|
221 | fOutputArray->Add(candidate);
|
---|
222 | if(TMath::Abs(q) > 1.0E-9)
|
---|
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 |
|
---|
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)
|
---|
245 |
|
---|
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]
|
---|
248 | r = pt / (q * fBz) * 1.0E9/c_light; // in [m]
|
---|
249 |
|
---|
250 | phi_0 = TMath::ATan2(py, px); // [rad] in [-pi, pi]
|
---|
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 |
|
---|
260 | rcu = TMath::Abs(r);
|
---|
261 | rc2 = r_c*r_c;
|
---|
262 |
|
---|
263 | // calculate coordinates of closest approach to track circle in transverse plane xd, yd, zd
|
---|
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;
|
---|
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;
|
---|
283 | dz = z - ((x - bsx) * px + (y - bsy) * py) / pt * (pz / pt);
|
---|
284 | p = candidateMomentum.P();
|
---|
285 | ctgTheta = 1.0 / TMath::Tan (candidateMomentum.Theta ());
|
---|
286 |
|
---|
287 |
|
---|
288 | // 3. time evaluation t = TMath::Min(t_r, t_z)
|
---|
289 | // t_r : time to exit from the sides
|
---|
290 | // t_z : time to exit from the front or the back
|
---|
291 | t_r = 0.0; // in [ns]
|
---|
292 | int sign_pz = (pz > 0.0) ? 1 : -1;
|
---|
293 | if(pz == 0.0) t_z = 1.0E99;
|
---|
294 | else t_z = gammam / (pz*1.0E9/c_light) * (-z + fHalfLength*sign_pz);
|
---|
295 |
|
---|
296 | if(r_c + TMath::Abs(r) < fRadius)
|
---|
297 | {
|
---|
298 | // helix does not cross the cylinder sides
|
---|
299 | t = t_z;
|
---|
300 | }
|
---|
301 | else
|
---|
302 | {
|
---|
303 | asinrho = TMath::ASin( (fRadius*fRadius - r_c*r_c - r*r) / (2*TMath::Abs(r)*r_c) );
|
---|
304 | delta = phi_0 - phi;
|
---|
305 | if(delta <-TMath::Pi()) delta += 2*TMath::Pi();
|
---|
306 | if(delta > TMath::Pi()) delta -= 2*TMath::Pi();
|
---|
307 | t1 = (delta + asinrho) / omega;
|
---|
308 | t2 = (delta + TMath::Pi() - asinrho) / omega;
|
---|
309 | t3 = (delta + TMath::Pi() + asinrho) / omega;
|
---|
310 | t4 = (delta - asinrho) / omega;
|
---|
311 | t5 = (delta - TMath::Pi() - asinrho) / omega;
|
---|
312 | t6 = (delta - TMath::Pi() + asinrho) / omega;
|
---|
313 |
|
---|
314 | if(t1 < 0.0) t1 = 1.0E99;
|
---|
315 | if(t2 < 0.0) t2 = 1.0E99;
|
---|
316 | if(t3 < 0.0) t3 = 1.0E99;
|
---|
317 | if(t4 < 0.0) t4 = 1.0E99;
|
---|
318 | if(t5 < 0.0) t5 = 1.0E99;
|
---|
319 | if(t6 < 0.0) t6 = 1.0E99;
|
---|
320 |
|
---|
321 | t_ra = TMath::Min(t1, TMath::Min(t2, t3));
|
---|
322 | t_rb = TMath::Min(t4, TMath::Min(t5, t6));
|
---|
323 | t_r = TMath::Min(t_ra, t_rb);
|
---|
324 | t = TMath::Min(t_r, t_z);
|
---|
325 | }
|
---|
326 |
|
---|
327 | // 4. position in terms of x(t), y(t), z(t)
|
---|
328 | x_t = x_c + r * TMath::Sin(omega * t - phi_0);
|
---|
329 | y_t = y_c + r * TMath::Cos(omega * t - phi_0);
|
---|
330 | z_t = z + pz*1.0E9 / c_light / gammam * t;
|
---|
331 | r_t = TMath::Hypot(x_t, y_t);
|
---|
332 |
|
---|
333 |
|
---|
334 | // compute path length for an helix
|
---|
335 |
|
---|
336 | alpha = pz*1.0E9 / c_light / gammam;
|
---|
337 | l = t * TMath::Sqrt(alpha*alpha + r*r*omega*omega);
|
---|
338 |
|
---|
339 | if(r_t > 0.0)
|
---|
340 | {
|
---|
341 |
|
---|
342 | // store these variables before cloning
|
---|
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 |
|
---|
350 | mother = candidate;
|
---|
351 | candidate = static_cast<Candidate*>(candidate->Clone());
|
---|
352 |
|
---|
353 | candidate->InitialPosition = candidatePosition;
|
---|
354 | candidate->Position.SetXYZT(x_t*1.0E3, y_t*1.0E3, z_t*1.0E3, candidatePosition.T() + t*c_light*1.0E3);
|
---|
355 |
|
---|
356 | candidate->Momentum = candidateMomentum;
|
---|
357 |
|
---|
358 | candidate->L = l*1.0E3;
|
---|
359 |
|
---|
360 | candidate->Xd = xd*1.0E3;
|
---|
361 | candidate->Yd = yd*1.0E3;
|
---|
362 | candidate->Zd = zd*1.0E3;
|
---|
363 |
|
---|
364 | candidate->AddCandidate(mother);
|
---|
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 | //------------------------------------------------------------------------------
|
---|
384 |
|
---|