Fork me on GitHub

source: git/modules/ParticlePropagator.cc@ 952bbbc

Last change on this file since 952bbbc was a07b54c, checked in by Roberto Preghenella <preghenella@…>, 5 years ago

Proper calculation of the DCAz coordinate

  • Property mode set to 100644
File size: 12.7 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/** \class ParticlePropagator
20 *
21 * Propagates charged and neutral particles
[d41ba4a]22 * from a given vertex to a cylinder defined by its radius,
[d7d2da3]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"
[341014c]37#include "ExRootAnalysis/ExRootFilter.h"
38#include "ExRootAnalysis/ExRootResult.h"
[d7d2da3]39
40#include "TDatabasePDG.h"
[341014c]41#include "TFormula.h"
[d7d2da3]42#include "TLorentzVector.h"
[341014c]43#include "TMath.h"
44#include "TObjArray.h"
45#include "TRandom3.h"
46#include "TString.h"
[d7d2da3]47
[d41ba4a]48#include <algorithm>
[d7d2da3]49#include <iostream>
50#include <sstream>
[341014c]51#include <stdexcept>
[d7d2da3]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);
[341014c]73 fRadius2 = fRadius * fRadius;
[d7d2da3]74 fHalfLength = GetDouble("HalfLength", 3.0);
75 fBz = GetDouble("Bz", 0.0);
76 if(fRadius < 1.0E-2)
[d41ba4a]77 {
[d7d2da3]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
[9330b7b]87 fRadiusMax = GetDouble("RadiusMax", fRadius);
88 fHalfLengthMax = GetDouble("HalfLengthMax", fHalfLength);
89
[d7d2da3]90 // import array with output from filter/classifier module
91
92 fInputArray = ImportArray(GetString("InputArray", "Delphes/stableParticles"));
93 fItInputArray = fInputArray->MakeIterator();
94
[187fc41]95 // import beamspot
[bff2e33]96 try
97 {
98 fBeamSpotInputArray = ImportArray(GetString("BeamSpotInputArray", "BeamSpotFilter/beamSpotParticle"));
99 }
100 catch(runtime_error &e)
101 {
102 fBeamSpotInputArray = 0;
[e55f5b0]103 }
[d7d2da3]104 // create output arrays
105
106 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
[751cb9c]107 fNeutralOutputArray = ExportArray(GetString("NeutralOutputArray", "neutralParticles"));
[d7d2da3]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{
[60e1de6]124 Candidate *candidate, *mother, *particle;
125 TLorentzVector particlePosition, particleMomentum, beamSpotPosition;
[d7d2da3]126 Double_t px, py, pz, pt, pt2, e, q;
127 Double_t x, y, z, t, r, phi;
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 t1, t2, t3, t4, t5, t6;
131 Double_t t_z, t_r, t_ra, t_rb;
132 Double_t tmp, discr, discr2;
133 Double_t delta, gammam, omega, asinrho;
[187fc41]134 Double_t rcu, rc2, xd, yd, zd;
135 Double_t l, d0, dz, p, ctgTheta, phip, etap, alpha;
136 Double_t bsx, bsy, bsz;
[a07b54c]137 Double_t s0, s1, sd;
138
[d7d2da3]139 const Double_t c_light = 2.99792458E8;
[e55f5b0]140
[341014c]141 if(!fBeamSpotInputArray || fBeamSpotInputArray->GetSize() == 0)
[187fc41]142 beamSpotPosition.SetXYZT(0.0, 0.0, 0.0, 0.0);
143 else
144 {
[341014c]145 Candidate &beamSpotCandidate = *((Candidate *)fBeamSpotInputArray->At(0));
[187fc41]146 beamSpotPosition = beamSpotCandidate.Position;
147 }
[e55f5b0]148
[d7d2da3]149 fItInputArray->Reset();
[341014c]150 while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
[d7d2da3]151 {
[60e1de6]152 if(candidate->GetCandidates()->GetEntriesFast() == 0)
153 {
154 particle = candidate;
155 }
156 else
157 {
[341014c]158 particle = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
[60e1de6]159 }
160
161 particlePosition = particle->Position;
162 particleMomentum = particle->Momentum;
[341014c]163 x = particlePosition.X() * 1.0E-3;
164 y = particlePosition.Y() * 1.0E-3;
165 z = particlePosition.Z() * 1.0E-3;
[e55f5b0]166
[341014c]167 bsx = beamSpotPosition.X() * 1.0E-3;
168 bsy = beamSpotPosition.Y() * 1.0E-3;
169 bsz = beamSpotPosition.Z() * 1.0E-3;
[e55f5b0]170
[60e1de6]171 q = particle->Charge;
[d7d2da3]172
173 // check that particle position is inside the cylinder
[9330b7b]174 if(TMath::Hypot(x, y) > fRadiusMax || TMath::Abs(z) > fHalfLengthMax)
[d7d2da3]175 {
176 continue;
177 }
178
[60e1de6]179 px = particleMomentum.Px();
180 py = particleMomentum.Py();
181 pz = particleMomentum.Pz();
182 pt = particleMomentum.Pt();
183 pt2 = particleMomentum.Perp2();
184 e = particleMomentum.E();
[d7d2da3]185
186 if(pt2 < 1.0E-9)
187 {
188 continue;
189 }
190
[9330b7b]191 if(TMath::Hypot(x, y) > fRadius || TMath::Abs(z) > fHalfLength)
192 {
193 mother = candidate;
[341014c]194 candidate = static_cast<Candidate *>(candidate->Clone());
[9330b7b]195
[60e1de6]196 candidate->InitialPosition = particlePosition;
197 candidate->Position = particlePosition;
[9330b7b]198 candidate->L = 0.0;
199
[60e1de6]200 candidate->Momentum = particleMomentum;
[9330b7b]201 candidate->AddCandidate(mother);
202
203 fOutputArray->Add(candidate);
204 }
205 else if(TMath::Abs(q) < 1.0E-9 || TMath::Abs(fBz) < 1.0E-9)
[d7d2da3]206 {
[fcdb8bc]207 // solve pt2*t^2 + 2*(px*x + py*y)*t - (fRadius2 - x*x - y*y) = 0
[341014c]208 tmp = px * y - py * x;
209 discr2 = pt2 * fRadius2 - tmp * tmp;
[d41ba4a]210
[b594101]211 if(discr2 < 0.0)
[d7d2da3]212 {
213 // no solutions
214 continue;
215 }
216
[341014c]217 tmp = px * x + py * y;
[d7d2da3]218 discr = TMath::Sqrt(discr2);
[341014c]219 t1 = (-tmp + discr) / pt2;
220 t2 = (-tmp - discr) / pt2;
[b594101]221 t = (t1 < 0.0) ? t2 : t1;
[d7d2da3]222
[341014c]223 z_t = z + pz * t;
[d7d2da3]224 if(TMath::Abs(z_t) > fHalfLength)
225 {
226 t3 = (+fHalfLength - z) / pz;
227 t4 = (-fHalfLength - z) / pz;
[b594101]228 t = (t3 < 0.0) ? t4 : t3;
[d7d2da3]229 }
230
[341014c]231 x_t = x + px * t;
232 y_t = y + py * t;
233 z_t = z + pz * t;
[e55f5b0]234
[341014c]235 l = TMath::Sqrt((x_t - x) * (x_t - x) + (y_t - y) * (y_t - y) + (z_t - z) * (z_t - z));
[d7d2da3]236
237 mother = candidate;
[341014c]238 candidate = static_cast<Candidate *>(candidate->Clone());
[d7d2da3]239
[60e1de6]240 candidate->InitialPosition = particlePosition;
[341014c]241 candidate->Position.SetXYZT(x_t * 1.0E3, y_t * 1.0E3, z_t * 1.0E3, particlePosition.T() + t * e * 1.0E3);
242 candidate->L = l * 1.0E3;
[e55f5b0]243
[60e1de6]244 candidate->Momentum = particleMomentum;
[d7d2da3]245 candidate->AddCandidate(mother);
[d41ba4a]246
[d7d2da3]247 fOutputArray->Add(candidate);
[d41ba4a]248 if(TMath::Abs(q) > 1.0E-9)
[d7d2da3]249 {
250 switch(TMath::Abs(candidate->PID))
251 {
[341014c]252 case 11:
253 fElectronOutputArray->Add(candidate);
254 break;
255 case 13:
256 fMuonOutputArray->Add(candidate);
257 break;
258 default:
259 fChargedHadronOutputArray->Add(candidate);
[d7d2da3]260 }
261 }
[751cb9c]262 else
263 {
[341014c]264 fNeutralOutputArray->Add(candidate);
[751cb9c]265 }
[d7d2da3]266 }
267 else
268 {
269
[b594101]270 // 1. initial transverse momentum p_{T0}: Part->pt
271 // initial transverse momentum direction phi_0 = -atan(p_X0/p_Y0)
272 // relativistic gamma: gamma = E/mc^2; gammam = gamma * m
273 // gyration frequency omega = q/(gamma m) fBz
274 // helix radius r = p_{T0} / (omega gamma m)
[d7d2da3]275
[341014c]276 gammam = e * 1.0E9 / (c_light * c_light); // gammam in [eV/c^2]
277 omega = q * fBz / (gammam); // omega is here in [89875518/s]
278 r = pt / (q * fBz) * 1.0E9 / c_light; // in [m]
[d7d2da3]279
[b594101]280 phi_0 = TMath::ATan2(py, px); // [rad] in [-pi, pi]
[d7d2da3]281
282 // 2. helix axis coordinates
[341014c]283 x_c = x + r * TMath::Sin(phi_0);
284 y_c = y - r * TMath::Cos(phi_0);
[d7d2da3]285 r_c = TMath::Hypot(x_c, y_c);
286 phi_c = TMath::ATan2(y_c, x_c);
287 phi = phi_c;
288 if(x_c < 0.0) phi += TMath::Pi();
289
[a0431dc]290 rcu = TMath::Abs(r);
[341014c]291 rc2 = r_c * r_c;
[b594101]292
[a0431dc]293 // calculate coordinates of closest approach to track circle in transverse plane xd, yd, zd
[341014c]294 xd = x_c * x_c * x_c - x_c * rcu * r_c + x_c * y_c * y_c;
[b594101]295 xd = (rc2 > 0.0) ? xd / rc2 : -999;
[341014c]296 yd = y_c * (-rcu * r_c + rc2);
[b594101]297 yd = (rc2 > 0.0) ? yd / rc2 : -999;
[341014c]298 zd = z + (TMath::Sqrt(xd * xd + yd * yd) - TMath::Sqrt(x * x + y * y)) * pz / pt;
[e55f5b0]299
[a07b54c]300 // proper calculation of the DCAz coordinate
301 // s0: track circle parameter at the track origin
302 // s1: track circle parameter at the closest approach to beam pipe
303 // sd: s1-s0 signed angular difference
304 s0 = atan2(y - y_c, x - x_c);
305 s1 = atan2(yd - y_c, xd - x_c);
306 sd = atan2(sin(s1 - s0), cos(s1 - s0));
307 zd = z - r * pz / pt * sd;
308
[187fc41]309 // use perigee momentum rather than original particle
[e55f5b0]310 // momentum, since the orignal particle momentum isn't known
311
[9330b7b]312 px = TMath::Sign(1.0, r) * pt * (-y_c / r_c);
313 py = TMath::Sign(1.0, r) * pt * (x_c / r_c);
[60e1de6]314 etap = particleMomentum.Eta();
[187fc41]315 phip = TMath::ATan2(py, px);
[e55f5b0]316
[60e1de6]317 particleMomentum.SetPtEtaPhiE(pt, etap, phip, particleMomentum.E());
[e55f5b0]318
[187fc41]319 // calculate additional track parameters (correct for beamspot position)
[e55f5b0]320
[341014c]321 d0 = ((x - bsx) * py - (y - bsy) * px) / pt;
322 dz = z - ((x - bsx) * px + (y - bsy) * py) / pt * (pz / pt);
323 p = particleMomentum.P();
324 ctgTheta = 1.0 / TMath::Tan(particleMomentum.Theta());
[e55f5b0]325
[d7d2da3]326 // 3. time evaluation t = TMath::Min(t_r, t_z)
327 // t_r : time to exit from the sides
328 // t_z : time to exit from the front or the back
329 t_r = 0.0; // in [ns]
330 int sign_pz = (pz > 0.0) ? 1 : -1;
[341014c]331 if(pz == 0.0)
332 t_z = 1.0E99;
333 else
334 t_z = gammam / (pz * 1.0E9 / c_light) * (-z + fHalfLength * sign_pz);
[d7d2da3]335
[341014c]336 if(r_c + TMath::Abs(r) < fRadius)
[d7d2da3]337 {
338 // helix does not cross the cylinder sides
339 t = t_z;
340 }
341 else
342 {
[341014c]343 asinrho = TMath::ASin((fRadius * fRadius - r_c * r_c - r * r) / (2 * TMath::Abs(r) * r_c));
[d7d2da3]344 delta = phi_0 - phi;
[341014c]345 if(delta < -TMath::Pi()) delta += 2 * TMath::Pi();
346 if(delta > TMath::Pi()) delta -= 2 * TMath::Pi();
[d7d2da3]347 t1 = (delta + asinrho) / omega;
348 t2 = (delta + TMath::Pi() - asinrho) / omega;
349 t3 = (delta + TMath::Pi() + asinrho) / omega;
350 t4 = (delta - asinrho) / omega;
351 t5 = (delta - TMath::Pi() - asinrho) / omega;
352 t6 = (delta - TMath::Pi() + asinrho) / omega;
353
[b594101]354 if(t1 < 0.0) t1 = 1.0E99;
355 if(t2 < 0.0) t2 = 1.0E99;
356 if(t3 < 0.0) t3 = 1.0E99;
357 if(t4 < 0.0) t4 = 1.0E99;
358 if(t5 < 0.0) t5 = 1.0E99;
359 if(t6 < 0.0) t6 = 1.0E99;
[d7d2da3]360
361 t_ra = TMath::Min(t1, TMath::Min(t2, t3));
362 t_rb = TMath::Min(t4, TMath::Min(t5, t6));
363 t_r = TMath::Min(t_ra, t_rb);
[d41ba4a]364 t = TMath::Min(t_r, t_z);
[d7d2da3]365 }
366
367 // 4. position in terms of x(t), y(t), z(t)
368 x_t = x_c + r * TMath::Sin(omega * t - phi_0);
369 y_t = y_c + r * TMath::Cos(omega * t - phi_0);
[341014c]370 z_t = z + pz * 1.0E9 / c_light / gammam * t;
[d7d2da3]371 r_t = TMath::Hypot(x_t, y_t);
372
[187fc41]373 // compute path length for an helix
[e55f5b0]374
[341014c]375 alpha = pz * 1.0E9 / c_light / gammam;
376 l = t * TMath::Sqrt(alpha * alpha + r * r * omega * omega);
[e55f5b0]377
[d7d2da3]378 if(r_t > 0.0)
379 {
[e55f5b0]380
[acd0621]381 // store these variables before cloning
[5b51d33]382 if(particle == candidate)
383 {
[341014c]384 particle->D0 = d0 * 1.0E3;
385 particle->DZ = dz * 1.0E3;
[5b51d33]386 particle->P = p;
387 particle->PT = pt;
388 particle->CtgTheta = ctgTheta;
389 particle->Phi = phip;
390 }
[e55f5b0]391
[d7d2da3]392 mother = candidate;
[341014c]393 candidate = static_cast<Candidate *>(candidate->Clone());
[d7d2da3]394
[60e1de6]395 candidate->InitialPosition = particlePosition;
[341014c]396 candidate->Position.SetXYZT(x_t * 1.0E3, y_t * 1.0E3, z_t * 1.0E3, particlePosition.T() + t * c_light * 1.0E3);
[d7d2da3]397
[60e1de6]398 candidate->Momentum = particleMomentum;
[e55f5b0]399
[341014c]400 candidate->L = l * 1.0E3;
[e55f5b0]401
[341014c]402 candidate->Xd = xd * 1.0E3;
403 candidate->Yd = yd * 1.0E3;
404 candidate->Zd = zd * 1.0E3;
[b594101]405
406 candidate->AddCandidate(mother);
[d7d2da3]407
408 fOutputArray->Add(candidate);
409 switch(TMath::Abs(candidate->PID))
410 {
[341014c]411 case 11:
412 fElectronOutputArray->Add(candidate);
413 break;
414 case 13:
415 fMuonOutputArray->Add(candidate);
416 break;
417 default:
418 fChargedHadronOutputArray->Add(candidate);
[d7d2da3]419 }
420 }
421 }
422 }
423}
424
425//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.