Fork me on GitHub

source: git/modules/ParticlePropagator.cc@ 80306e6

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 80306e6 was 48f7a77, checked in by Michele Selvaggi <michele.selvaggi@…>, 8 years ago

added compute beamspot function in propagator

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