Fork me on GitHub

source: git/modules/ParticlePropagator.cc@ 17826f2

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 17826f2 was fcdb8bc, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

fix comment in ParticlePropagator

  • Property mode set to 100644
File size: 9.5 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
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
[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 * $Date$
28 * $Revision$
29 *
30 *
31 * \author P. Demin - UCL, Louvain-la-Neuve
32 *
33 */
34
35#include "modules/ParticlePropagator.h"
36
37#include "classes/DelphesClasses.h"
38#include "classes/DelphesFactory.h"
39#include "classes/DelphesFormula.h"
40
41#include "ExRootAnalysis/ExRootResult.h"
42#include "ExRootAnalysis/ExRootFilter.h"
43#include "ExRootAnalysis/ExRootClassifier.h"
44
45#include "TMath.h"
46#include "TString.h"
47#include "TFormula.h"
48#include "TRandom3.h"
49#include "TObjArray.h"
50#include "TDatabasePDG.h"
51#include "TLorentzVector.h"
52
[d41ba4a]53#include <algorithm>
[d7d2da3]54#include <stdexcept>
55#include <iostream>
56#include <sstream>
57
58using namespace std;
59
60//------------------------------------------------------------------------------
61
62ParticlePropagator::ParticlePropagator() :
63 fItInputArray(0)
64{
65}
66
67//------------------------------------------------------------------------------
68
69ParticlePropagator::~ParticlePropagator()
70{
71}
72
73//------------------------------------------------------------------------------
74
75void ParticlePropagator::Init()
76{
77 fRadius = GetDouble("Radius", 1.0);
78 fRadius2 = fRadius*fRadius;
79 fHalfLength = GetDouble("HalfLength", 3.0);
80 fBz = GetDouble("Bz", 0.0);
81 if(fRadius < 1.0E-2)
[d41ba4a]82 {
[d7d2da3]83 cout << "ERROR: magnetic field radius is too low\n";
84 return;
85 }
86 if(fHalfLength < 1.0E-2)
87 {
88 cout << "ERROR: magnetic field length is too low\n";
89 return;
90 }
91
92 // import array with output from filter/classifier module
93
94 fInputArray = ImportArray(GetString("InputArray", "Delphes/stableParticles"));
95 fItInputArray = fInputArray->MakeIterator();
96
97 // create output arrays
98
99 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
100 fChargedHadronOutputArray = ExportArray(GetString("ChargedHadronOutputArray", "chargedHadrons"));
101 fElectronOutputArray = ExportArray(GetString("ElectronOutputArray", "electrons"));
102 fMuonOutputArray = ExportArray(GetString("MuonOutputArray", "muons"));
103}
104
105//------------------------------------------------------------------------------
106
107void ParticlePropagator::Finish()
108{
109 if(fItInputArray) delete fItInputArray;
110}
111
112//------------------------------------------------------------------------------
113
114void ParticlePropagator::Process()
115{
116 Candidate *candidate, *mother;
117 TLorentzVector candidatePosition, candidateMomentum;
118 Double_t px, py, pz, pt, pt2, e, q;
119 Double_t x, y, z, t, r, phi;
120 Double_t x_c, y_c, r_c, phi_c, phi_0;
121 Double_t x_t, y_t, z_t, r_t;
122 Double_t t1, t2, t3, t4, t5, t6;
123 Double_t t_z, t_r, t_ra, t_rb;
124 Double_t tmp, discr, discr2;
125 Double_t delta, gammam, omega, asinrho;
[a0431dc]126 Double_t ang_mom, rcu, rc2, dxy, xd, yd, zd;
127
[d7d2da3]128 const Double_t c_light = 2.99792458E8;
[d41ba4a]129
[d7d2da3]130 fItInputArray->Reset();
131 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
132 {
133 candidatePosition = candidate->Position;
134 candidateMomentum = candidate->Momentum;
135 x = candidatePosition.X()*1.0E-3;
136 y = candidatePosition.Y()*1.0E-3;
137 z = candidatePosition.Z()*1.0E-3;
138 q = candidate->Charge;
139
140 // check that particle position is inside the cylinder
141 if(TMath::Hypot(x, y) > fRadius || TMath::Abs(z) > fHalfLength)
142 {
143 continue;
144 }
145
146 px = candidateMomentum.Px();
147 py = candidateMomentum.Py();
148 pz = candidateMomentum.Pz();
149 pt = candidateMomentum.Pt();
150 pt2 = candidateMomentum.Perp2();
151 e = candidateMomentum.E();
152
153 if(pt2 < 1.0E-9)
154 {
155 continue;
156 }
157
158 if(TMath::Abs(q) < 1.0E-9 || TMath::Abs(fBz) < 1.0E-9)
159 {
[fcdb8bc]160 // solve pt2*t^2 + 2*(px*x + py*y)*t - (fRadius2 - x*x - y*y) = 0
[d7d2da3]161 tmp = px*y - py*x;
162 discr2 = pt2*fRadius2 - tmp*tmp;
[d41ba4a]163
[d7d2da3]164 if(discr2 < 0)
165 {
166 // no solutions
167 continue;
168 }
169
170 tmp = px*x + py*y;
171 discr = TMath::Sqrt(discr2);
172 t1 = (-tmp + discr)/pt2;
173 t2 = (-tmp - discr)/pt2;
[d41ba4a]174 t = (t1 < 0) ? t2 : t1;
[d7d2da3]175
176 z_t = z + pz*t;
177 if(TMath::Abs(z_t) > fHalfLength)
178 {
179 t3 = (+fHalfLength - z) / pz;
180 t4 = (-fHalfLength - z) / pz;
[d41ba4a]181 t = (t3 < 0) ? t4 : t3;
[d7d2da3]182 }
183
184 x_t = x + px*t;
185 y_t = y + py*t;
186 z_t = z + pz*t;
187
188 mother = candidate;
189 candidate = static_cast<Candidate*>(candidate->Clone());
190
[d41ba4a]191 candidate->Position.SetXYZT(x_t*1.0E3, y_t*1.0E3, z_t*1.0E3, candidatePosition.T() + t*e*1.0E3);
[d7d2da3]192
193 candidate->Momentum = candidateMomentum;
194 candidate->AddCandidate(mother);
[d41ba4a]195
[d7d2da3]196 fOutputArray->Add(candidate);
[d41ba4a]197 if(TMath::Abs(q) > 1.0E-9)
[d7d2da3]198 {
199 switch(TMath::Abs(candidate->PID))
200 {
201 case 11:
202 fElectronOutputArray->Add(candidate);
203 break;
204 case 13:
205 fMuonOutputArray->Add(candidate);
206 break;
207 default:
208 fChargedHadronOutputArray->Add(candidate);
209 }
210 }
211 }
212 else
213 {
214
215 // 1. initial transverse momentum p_{T0} : Part->pt
[d41ba4a]216 // initial transverse momentum direction \phi_0 = -atan(p_X0/p_Y0)
[d7d2da3]217 // relativistic gamma : gamma = E/mc² ; gammam = gamma \times m
218 // giration frequency \omega = q/(gamma m) fBz
219 // helix radius r = p_T0 / (omega gamma m)
220
221 gammam = e*1.0E9 / (c_light*c_light); // gammam in [eV/c²]
[d41ba4a]222 omega = q * fBz / (gammam); // omega is here in [ 89875518 / s]
223 r = pt / (q * fBz) * 1.0E9/c_light; // in [m]
[d7d2da3]224
225 phi_0 = TMath::ATan2(py, px); // [rad] in [-pi; pi]
226
227 // 2. helix axis coordinates
228 x_c = x + r*TMath::Sin(phi_0);
229 y_c = y - r*TMath::Cos(phi_0);
230 r_c = TMath::Hypot(x_c, y_c);
231 phi_c = TMath::ATan2(y_c, x_c);
232 phi = phi_c;
233 if(x_c < 0.0) phi += TMath::Pi();
234
[a0431dc]235 rcu = TMath::Abs(r);
236 rc2 = r_c*r_c;
237
238 // calculate coordinates of closest approach to track circle in transverse plane xd, yd, zd
239 xd = x_c*x_c*x_c - x_c*rcu*r_c + x_c*y_c*y_c;
240 xd = ( rc2 > 0.0 ) ? xd / rc2 : -999;
241 yd = y_c*(-rcu*r_c + rc2);
242 yd = ( rc2 > 0.0 ) ? yd / rc2 : -999;
243 zd = z + (TMath::Sqrt(xd*xd+yd*yd) - TMath::Sqrt(x*x+y*y))*pz/pt;
244
245 // calculate impact paramater
246 ang_mom = (xd*py - yd*px);
247 dxy = ang_mom/pt;
248
249
[d7d2da3]250 // 3. time evaluation t = TMath::Min(t_r, t_z)
251 // t_r : time to exit from the sides
252 // t_z : time to exit from the front or the back
253 t_r = 0.0; // in [ns]
254 int sign_pz = (pz > 0.0) ? 1 : -1;
255 if(pz == 0.0) t_z = 1.0E99;
256 else t_z = gammam / (pz*1.0E9/c_light) * (-z + fHalfLength*sign_pz);
257
258 if(r_c + TMath::Abs(r) < fRadius)
259 {
260 // helix does not cross the cylinder sides
261 t = t_z;
262 }
263 else
264 {
265 asinrho = TMath::ASin( (fRadius*fRadius - r_c*r_c - r*r) / (2*TMath::Abs(r)*r_c) );
266 delta = phi_0 - phi;
267 if(delta <-TMath::Pi()) delta += 2*TMath::Pi();
268 if(delta > TMath::Pi()) delta -= 2*TMath::Pi();
269 t1 = (delta + asinrho) / omega;
270 t2 = (delta + TMath::Pi() - asinrho) / omega;
271 t3 = (delta + TMath::Pi() + asinrho) / omega;
272 t4 = (delta - asinrho) / omega;
273 t5 = (delta - TMath::Pi() - asinrho) / omega;
274 t6 = (delta - TMath::Pi() + asinrho) / omega;
275
276 if(t1 < 0) t1 = 1.0E99;
277 if(t2 < 0) t2 = 1.0E99;
278 if(t3 < 0) t3 = 1.0E99;
279 if(t4 < 0) t4 = 1.0E99;
280 if(t5 < 0) t5 = 1.0E99;
281 if(t6 < 0) t6 = 1.0E99;
282
283 t_ra = TMath::Min(t1, TMath::Min(t2, t3));
284 t_rb = TMath::Min(t4, TMath::Min(t5, t6));
285 t_r = TMath::Min(t_ra, t_rb);
[d41ba4a]286 t = TMath::Min(t_r, t_z);
[d7d2da3]287 }
288
289 // 4. position in terms of x(t), y(t), z(t)
290 x_t = x_c + r * TMath::Sin(omega * t - phi_0);
291 y_t = y_c + r * TMath::Cos(omega * t - phi_0);
292 z_t = z + pz*1.0E9 / c_light / gammam * t;
293 r_t = TMath::Hypot(x_t, y_t);
294
295 if(r_t > 0.0)
296 {
297 mother = candidate;
298 candidate = static_cast<Candidate*>(candidate->Clone());
299
[d41ba4a]300 candidate->Position.SetXYZT(x_t*1.0E3, y_t*1.0E3, z_t*1.0E3, candidatePosition.T() + t*c_light*1.0E3);
[d7d2da3]301
302 candidate->Momentum = candidateMomentum;
[a0431dc]303 candidate->Xd = xd*1.0E3;
304 candidate->Yd = yd*1.0E3;
305 candidate->Zd = zd*1.0E3;
306
307 candidate->AddCandidate(mother);
[d7d2da3]308
309 fOutputArray->Add(candidate);
310 switch(TMath::Abs(candidate->PID))
311 {
312 case 11:
313 fElectronOutputArray->Add(candidate);
314 break;
315 case 13:
316 fMuonOutputArray->Add(candidate);
317 break;
318 default:
319 fChargedHadronOutputArray->Add(candidate);
320 }
321 }
322 }
323 }
324}
325
326//------------------------------------------------------------------------------
[a0431dc]327
Note: See TracBrowser for help on using the repository browser.