Fork me on GitHub

source: git/modules/TrackSmearing.cc@ befa255

ImprovedOutputFile Timing dual_readout llp
Last change on this file since befa255 was cde9f31, checked in by Kevin Pedro <kpedro88@…>, 7 years ago

propagate track smearing to Xd,Yd,Zd variables

  • Property mode set to 100644
File size: 12.7 KB
Line 
1/** \class TrackSmearing
2 *
3 * Performs d0, dZ, p, Theta, Phi smearing of tracks.
4 *
5 * \authors A. Hart, M. Selvaggi
6 *
7*/
8
9#include "modules/TrackSmearing.h"
10
11#include "classes/DelphesClasses.h"
12#include "classes/DelphesFactory.h"
13#include "classes/DelphesFormula.h"
14
15#include "ExRootAnalysis/ExRootResult.h"
16#include "ExRootAnalysis/ExRootFilter.h"
17#include "ExRootAnalysis/ExRootClassifier.h"
18
19#include "TMath.h"
20#include "TString.h"
21#include "TFormula.h"
22#include "TRandom3.h"
23#include "TObjArray.h"
24#include "TDatabasePDG.h"
25#include "TLorentzVector.h"
26#include "TFile.h"
27#include "TProfile2D.h"
28
29#include <algorithm>
30#include <stdexcept>
31#include <iostream>
32#include <sstream>
33
34using namespace std;
35
36//------------------------------------------------------------------------------
37
38TrackSmearing::TrackSmearing() :
39 fD0Formula(0), fDZFormula(0), fPFormula(0), fCtgThetaFormula(0), fPhiFormula(0), fItInputArray(0)
40{
41 fD0Formula = new DelphesFormula;
42 fDZFormula = new DelphesFormula;
43 fPFormula = new DelphesFormula;
44 fCtgThetaFormula = new DelphesFormula;
45 fPhiFormula = new DelphesFormula;
46}
47
48//------------------------------------------------------------------------------
49
50TrackSmearing::~TrackSmearing()
51{
52 if(fD0Formula) delete fD0Formula;
53 if(fDZFormula) delete fDZFormula;
54 if(fPFormula) delete fPFormula;
55 if(fCtgThetaFormula) delete fCtgThetaFormula;
56 if(fPhiFormula) delete fPhiFormula;
57}
58
59//------------------------------------------------------------------------------
60
61void TrackSmearing::Init()
62{
63 fBz = GetDouble("Bz", 0.0);
64
65 // read resolution formula
66
67 // !!! IF WE WANT TO KEEP ROOT INPUT !!!
68 if (string (GetString("D0ResolutionFormula", "0.0")) != "0.0")
69 {
70 fD0Formula->Compile(GetString("D0ResolutionFormula", "0.0"));
71 fUseD0Formula = true;
72 }
73 else
74 {
75 fD0ResolutionFile = GetString("D0ResolutionFile", "errors.root");
76 fD0ResolutionHist = GetString("D0ResolutionHist", "d0");
77 fUseD0Formula = false;
78 }
79 if (string (GetString("DZResolutionFormula", "0.0")) != "0.0")
80 {
81 fDZFormula->Compile(GetString("DZResolutionFormula", "0.0"));
82 fUseDZFormula = true;
83 }
84 else
85 {
86 fDZResolutionFile = GetString("DZResolutionFile", "errors.root");
87 fDZResolutionHist = GetString("DZResolutionHist", "dz");
88 fUseDZFormula = false;
89 }
90 if (string (GetString("PResolutionFormula", "0.0")) != "0.0")
91 {
92 fPFormula->Compile(GetString("PResolutionFormula", "0.0"));
93 fUsePFormula = true;
94 }
95 else
96 {
97 fPResolutionFile = GetString("PResolutionFile", "errors.root");
98 fPResolutionHist = GetString("PResolutionHist", "p");
99 fUsePFormula = false;
100 }
101 if (string (GetString("CtgThetaResolutionFormula", "0.0")) != "0.0")
102 {
103 fCtgThetaFormula->Compile(GetString("CtgThetaResolutionFormula", "0.0"));
104 fUseCtgThetaFormula = true;
105 }
106 else
107 {
108 fCtgThetaResolutionFile = GetString("CtgThetaResolutionFile", "errors.root");
109 fCtgThetaResolutionHist = GetString("CtgThetaResolutionHist", "ctgTheta");
110 fUseCtgThetaFormula = false;
111 }
112 if (string (GetString("PhiResolutionFormula", "0.0")) != "0.0")
113 {
114 fPhiFormula->Compile(GetString("PhiResolutionFormula", "0.0"));
115 fUsePhiFormula = true;
116 }
117 else
118 {
119 fPhiResolutionFile = GetString("PhiResolutionFile", "errors.root");
120 fPhiResolutionHist = GetString("PhiResolutionHist", "phi");
121 fUsePhiFormula = false;
122 }
123
124 fApplyToPileUp = GetBool("ApplyToPileUp", true);
125
126 // import input array
127
128 fInputArray = ImportArray(GetString("InputArray", "ParticlePropagator/stableParticles"));
129 fItInputArray = fInputArray->MakeIterator();
130
131 // import beamspot
132 try
133 {
134 fBeamSpotInputArray = ImportArray(GetString("BeamSpotInputArray", "BeamSpotFilter/beamSpotParticle"));
135 }
136 catch(runtime_error &e)
137 {
138 fBeamSpotInputArray = 0;
139 }
140
141 // create output array
142
143 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
144}
145
146//------------------------------------------------------------------------------
147
148void TrackSmearing::Finish()
149{
150 if(fItInputArray) delete fItInputArray;
151}
152
153//------------------------------------------------------------------------------
154
155void TrackSmearing::Process()
156{
157 Int_t iCandidate = 0;
158 TLorentzVector beamSpotPosition;
159 Candidate *candidate, *mother;
160 Double_t pt, eta, d0, d0Error, trueD0, dz, dzError, trueDZ, p, pError, trueP, ctgTheta, ctgThetaError, trueCtgTheta, phi, phiError, truePhi;
161 Double_t x, y, z, t, px, py, pz, theta;
162 Double_t q, r;
163 Double_t x_c, y_c, r_c, phi_0;
164 Double_t rcu, rc2, xd, yd, zd;
165 const Double_t c_light = 2.99792458E8;
166 TProfile2D *d0ErrorHist = NULL,
167 *dzErrorHist = NULL,
168 *pErrorHist = NULL,
169 *ctgThetaErrorHist = NULL,
170 *phiErrorHist = NULL;
171
172 if (!fBeamSpotInputArray || fBeamSpotInputArray->GetSize () == 0)
173 beamSpotPosition.SetXYZT(0.0, 0.0, 0.0, 0.0);
174 else
175 {
176 Candidate &beamSpotCandidate = *((Candidate *) fBeamSpotInputArray->At (0));
177 beamSpotPosition = beamSpotCandidate.Position;
178 }
179
180 if (!fUseD0Formula)
181 {
182 TFile *fin = TFile::Open (fD0ResolutionFile.c_str ());
183 d0ErrorHist = (TProfile2D *) fin->Get (fD0ResolutionHist.c_str ());
184 d0ErrorHist->SetDirectory (0);
185 fin->Close ();
186 }
187 if (!fUseDZFormula)
188 {
189 TFile *fin = TFile::Open (fDZResolutionFile.c_str ());
190 dzErrorHist = (TProfile2D *) fin->Get (fDZResolutionHist.c_str ());
191 dzErrorHist->SetDirectory (0);
192 fin->Close ();
193 }
194 if (!fUsePFormula)
195 {
196 TFile *fin = TFile::Open (fPResolutionFile.c_str ());
197 pErrorHist = (TProfile2D *) fin->Get (fPResolutionHist.c_str ());
198 pErrorHist->SetDirectory (0);
199 fin->Close ();
200 }
201 if (!fUseCtgThetaFormula)
202 {
203 TFile *fin = TFile::Open (fCtgThetaResolutionFile.c_str ());
204 ctgThetaErrorHist = (TProfile2D *) fin->Get (fCtgThetaResolutionHist.c_str ());
205 ctgThetaErrorHist->SetDirectory (0);
206 fin->Close ();
207 }
208 if (!fUsePhiFormula)
209 {
210 TFile *fin = TFile::Open (fPhiResolutionFile.c_str ());
211 phiErrorHist = (TProfile2D *) fin->Get (fPhiResolutionHist.c_str ());
212 phiErrorHist->SetDirectory (0);
213 fin->Close ();
214 }
215
216 fItInputArray->Reset();
217 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
218 {
219
220 const TLorentzVector &momentum = candidate->Momentum;
221 const TLorentzVector &position = candidate->InitialPosition;
222
223 pt = momentum.Pt();
224 eta = momentum.Eta();
225
226 d0 = trueD0 = candidate->D0;
227 dz = trueDZ = candidate->DZ;
228
229 p = trueP = candidate->P;
230 ctgTheta = trueCtgTheta = candidate->CtgTheta;
231 phi = truePhi = candidate->Phi;
232
233 if (fUseD0Formula)
234 d0Error = fD0Formula->Eval(pt, eta);
235 else
236 {
237 Int_t xbin, ybin;
238
239 xbin = pt < d0ErrorHist->GetXaxis ()->GetXmax () ? d0ErrorHist->GetXaxis ()->FindBin (pt)
240 : d0ErrorHist->GetXaxis ()->GetBinCenter (d0ErrorHist->GetXaxis ()->GetNbins ());
241 ybin = d0ErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
242 d0Error = d0ErrorHist->GetBinContent (xbin, ybin);
243 if (!d0Error)
244 d0Error = -1.0;
245 }
246 if (d0Error < 0.0)
247 continue;
248
249 if (fUseDZFormula)
250 dzError = fDZFormula->Eval(pt, eta);
251 else
252 {
253 Int_t xbin, ybin;
254
255 xbin = pt < dzErrorHist->GetXaxis ()->GetXmax () ? dzErrorHist->GetXaxis ()->FindBin (pt)
256 : dzErrorHist->GetXaxis ()->GetBinCenter (dzErrorHist->GetXaxis ()->GetNbins ());
257 ybin = dzErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
258 dzError = dzErrorHist->GetBinContent (xbin, ybin);
259 if (!dzError)
260 dzError = -1.0;
261 }
262 if (dzError < 0.0)
263 continue;
264
265 if (fUsePFormula)
266 pError = fPFormula->Eval(pt, eta) * p;
267 else
268 {
269 Int_t xbin, ybin;
270
271 xbin = pt < pErrorHist->GetXaxis ()->GetXmax () ? pErrorHist->GetXaxis ()->FindBin (pt)
272 : pErrorHist->GetXaxis ()->GetBinCenter (pErrorHist->GetXaxis ()->GetNbins ());
273 ybin = pErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
274 pError = pErrorHist->GetBinContent (xbin, ybin) * p;
275 if (!pError)
276 pError = -1.0;
277 }
278 if (pError < 0.0)
279 continue;
280
281 if (fUseCtgThetaFormula)
282 ctgThetaError = fCtgThetaFormula->Eval(pt, eta);
283 else
284 {
285 Int_t xbin, ybin;
286
287 xbin = pt < ctgThetaErrorHist->GetXaxis ()->GetXmax () ? ctgThetaErrorHist->GetXaxis ()->FindBin (pt)
288 : ctgThetaErrorHist->GetXaxis ()->GetBinCenter (ctgThetaErrorHist->GetXaxis ()->GetNbins ());
289 ybin = ctgThetaErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
290 ctgThetaError = ctgThetaErrorHist->GetBinContent (xbin, ybin);
291 if (!ctgThetaError)
292 ctgThetaError = -1.0;
293 }
294 if (ctgThetaError < 0.0)
295 continue;
296
297 if (fUsePhiFormula)
298 phiError = fPhiFormula->Eval(pt, eta);
299 else
300 {
301 Int_t xbin, ybin;
302
303 xbin = pt < phiErrorHist->GetXaxis ()->GetXmax () ? phiErrorHist->GetXaxis ()->FindBin (pt)
304 : phiErrorHist->GetXaxis ()->GetBinCenter (phiErrorHist->GetXaxis ()->GetNbins ());
305 ybin = phiErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
306 phiError = phiErrorHist->GetBinContent (xbin, ybin);
307 if (!phiError)
308 phiError = -1.0;
309 }
310 if (phiError < 0.0)
311 continue;
312
313 if (fApplyToPileUp || !candidate->IsPU)
314 {
315 d0 = gRandom->Gaus(d0, d0Error);
316 dz = gRandom->Gaus(dz, dzError);
317 p = gRandom->Gaus(p, pError);
318 ctgTheta = gRandom->Gaus(ctgTheta, ctgThetaError);
319 phi = gRandom->Gaus(phi, phiError);
320 }
321
322 if(p < 0.0) continue;
323 while (phi > TMath::Pi ()) phi -= TMath::TwoPi ();
324 while (phi <= -TMath::Pi ()) phi += TMath::TwoPi ();
325
326 mother = candidate;
327 candidate = static_cast<Candidate*>(candidate->Clone());
328 candidate->D0 = d0;
329 candidate->DZ = dz;
330 candidate->P = p;
331 candidate->CtgTheta = ctgTheta;
332 candidate->Phi = phi;
333
334 theta = TMath::ACos(ctgTheta / TMath::Sqrt (1.0 + ctgTheta * ctgTheta));
335 candidate->Momentum.SetPx (p * TMath::Cos (phi) * TMath::Sin (theta));
336 candidate->Momentum.SetPy (p * TMath::Sin (phi) * TMath::Sin (theta));
337 candidate->Momentum.SetPz (p * TMath::Cos (theta));
338 candidate->Momentum.SetE (candidate->Momentum.Pt () * TMath::CosH (eta));
339 candidate->PT = candidate->Momentum.Pt ();
340
341 x = position.X ();
342 y = position.Y ();
343 z = position.Z ();
344 t = position.T ();
345 px = candidate->Momentum.Px ();
346 py = candidate->Momentum.Py ();
347 pz = candidate->Momentum.Pz ();
348 pt = candidate->Momentum.Pt ();
349
350 // -- solve for delta: d0' = ( (x+delta)*py' - (y+delta)*px' )/pt'
351
352 candidate->InitialPosition.SetX (x + ((px * y - py * x + d0 * pt) / (py - px)));
353 candidate->InitialPosition.SetY (y + ((px * y - py * x + d0 * pt) / (py - px)));
354 x = candidate->InitialPosition.X ();
355 y = candidate->InitialPosition.Y ();
356 candidate->InitialPosition.SetZ (z + ((pz * (px * (x - beamSpotPosition.X ()) + py * (y - beamSpotPosition.Y ())) + pt * pt * (dz - z)) / (pt * pt)));
357 z = candidate->InitialPosition.Z ();
358
359 candidate->InitialPosition.SetT(t);
360
361 // update closest approach
362 x *= 1.0E-3;
363 y *= 1.0E-3;
364 z *= 1.0E-3;
365
366 q = candidate->Charge;
367
368 r = pt / (q * fBz) * 1.0E9/c_light; // in [m]
369 phi_0 = TMath::ATan2(py, px); // [rad] in [-pi, pi]
370
371 // 2. helix axis coordinates
372 x_c = x + r*TMath::Sin(phi_0);
373 y_c = y - r*TMath::Cos(phi_0);
374 r_c = TMath::Hypot(x_c, y_c);
375
376 rcu = TMath::Abs(r);
377 rc2 = r_c*r_c;
378
379 // calculate coordinates of closest approach to track circle in transverse plane xd, yd, zd
380 xd = x_c*x_c*x_c - x_c*rcu*r_c + x_c*y_c*y_c;
381 xd = (rc2 > 0.0) ? xd / rc2 : -999;
382 yd = y_c*(-rcu*r_c + rc2);
383 yd = (rc2 > 0.0) ? yd / rc2 : -999;
384 zd = z + (TMath::Sqrt(xd*xd + yd*yd) - TMath::Sqrt(x*x + y*y))*pz/pt;
385
386 candidate->Xd = xd*1.0E3;
387 candidate->Yd = yd*1.0E3;
388 candidate->Zd = zd*1.0E3;
389
390 if (fApplyToPileUp || !candidate->IsPU)
391 {
392 candidate->ErrorD0 = d0Error;
393 candidate->ErrorDZ = dzError;
394 candidate->ErrorP = pError;
395 candidate->ErrorCtgTheta = ctgThetaError;
396 candidate->ErrorPhi = phiError;
397 candidate->ErrorPT = ptError (p, ctgTheta, pError, ctgThetaError);
398 candidate->TrackResolution = pError/p;
399 }
400
401 candidate->AddCandidate(mother);
402 fOutputArray->Add(candidate);
403
404 iCandidate++;
405 }
406}
407
408Double_t TrackSmearing::ptError (const Double_t p, const Double_t ctgTheta, const Double_t dP, const Double_t dCtgTheta)
409{
410 Double_t a, b;
411 a = (p * p * ctgTheta * ctgTheta * dCtgTheta * dCtgTheta) / ((ctgTheta * ctgTheta + 1) * (ctgTheta * ctgTheta + 1) * (ctgTheta * ctgTheta + 1));
412 b = (dP * dP) / (ctgTheta * ctgTheta + 1);
413 return sqrt (a + b);
414}
415
416//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.