Fork me on GitHub

source: git/modules/TrackSmearing.cc@ 641cb3d

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 641cb3d was 786028a, checked in by Michele Selvaggi <michele.selvaggi@…>, 8 years ago

fixed typo in TrackResolution

  • Property mode set to 100644
File size: 11.7 KB
Line 
1/** \class TrackSmearing
2 *
3 * Performs d0, dZ, p, Theta, Phi smearing of tracks.
4 *
5 *
6 *
7 * \author A. Hart, M. Selvaggi
8 *
9*/
10
11#include "modules/TrackSmearing.h"
12
13#include "classes/DelphesClasses.h"
14#include "classes/DelphesFactory.h"
15#include "classes/DelphesFormula.h"
16
17#include "ExRootAnalysis/ExRootResult.h"
18#include "ExRootAnalysis/ExRootFilter.h"
19#include "ExRootAnalysis/ExRootClassifier.h"
20
21#include "TMath.h"
22#include "TString.h"
23#include "TFormula.h"
24#include "TRandom3.h"
25#include "TObjArray.h"
26#include "TDatabasePDG.h"
27#include "TLorentzVector.h"
28#include "TFile.h"
29#include "TProfile2D.h"
30
31#include <algorithm>
32#include <stdexcept>
33#include <iostream>
34#include <sstream>
35
36using namespace std;
37
38//------------------------------------------------------------------------------
39
40TrackSmearing::TrackSmearing() :
41 fD0Formula(0), fDZFormula(0), fPFormula(0), fCtgThetaFormula(0), fPhiFormula(0), fItInputArray(0)
42{
43 fD0Formula = new DelphesFormula;
44 fDZFormula = new DelphesFormula;
45 fPFormula = new DelphesFormula;
46 fCtgThetaFormula = new DelphesFormula;
47 fPhiFormula = new DelphesFormula;
48}
49
50//------------------------------------------------------------------------------
51
52TrackSmearing::~TrackSmearing()
53{
54 if(fD0Formula) delete fD0Formula;
55 if(fDZFormula) delete fDZFormula;
56 if(fPFormula) delete fPFormula;
57 if(fCtgThetaFormula) delete fCtgThetaFormula;
58 if(fPhiFormula) delete fPhiFormula;
59}
60
61//------------------------------------------------------------------------------
62
63void TrackSmearing::Init()
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 TProfile2D *d0ErrorHist = NULL,
163 *dzErrorHist = NULL,
164 *pErrorHist = NULL,
165 *ctgThetaErrorHist = NULL,
166 *phiErrorHist = NULL;
167
168 if (!fBeamSpotInputArray || fBeamSpotInputArray->GetSize () == 0)
169 beamSpotPosition.SetXYZT(0.0, 0.0, 0.0, 0.0);
170 else
171 {
172 Candidate &beamSpotCandidate = *((Candidate *) fBeamSpotInputArray->At (0));
173 beamSpotPosition = beamSpotCandidate.Position;
174 }
175
176 if (!fUseD0Formula)
177 {
178 TFile *fin = TFile::Open (fD0ResolutionFile.c_str ());
179 d0ErrorHist = (TProfile2D *) fin->Get (fD0ResolutionHist.c_str ());
180 d0ErrorHist->SetDirectory (0);
181 fin->Close ();
182 }
183 if (!fUseDZFormula)
184 {
185 TFile *fin = TFile::Open (fDZResolutionFile.c_str ());
186 dzErrorHist = (TProfile2D *) fin->Get (fDZResolutionHist.c_str ());
187 dzErrorHist->SetDirectory (0);
188 fin->Close ();
189 }
190 if (!fUsePFormula)
191 {
192 TFile *fin = TFile::Open (fPResolutionFile.c_str ());
193 pErrorHist = (TProfile2D *) fin->Get (fPResolutionHist.c_str ());
194 pErrorHist->SetDirectory (0);
195 fin->Close ();
196 }
197 if (!fUseCtgThetaFormula)
198 {
199 TFile *fin = TFile::Open (fCtgThetaResolutionFile.c_str ());
200 ctgThetaErrorHist = (TProfile2D *) fin->Get (fCtgThetaResolutionHist.c_str ());
201 ctgThetaErrorHist->SetDirectory (0);
202 fin->Close ();
203 }
204 if (!fUsePhiFormula)
205 {
206 TFile *fin = TFile::Open (fPhiResolutionFile.c_str ());
207 phiErrorHist = (TProfile2D *) fin->Get (fPhiResolutionHist.c_str ());
208 phiErrorHist->SetDirectory (0);
209 fin->Close ();
210 }
211
212 fItInputArray->Reset();
213 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
214 {
215
216 const TLorentzVector &momentum = candidate->Momentum;
217 const TLorentzVector &position = candidate->InitialPosition;
218
219 pt = momentum.Pt();
220 eta = momentum.Eta();
221
222 d0 = trueD0 = candidate->D0;
223 dz = trueDZ = candidate->DZ;
224
225 p = trueP = candidate->P;
226 ctgTheta = trueCtgTheta = candidate->CtgTheta;
227 phi = truePhi = candidate->Phi;
228
229 if (fUseD0Formula)
230 d0Error = fD0Formula->Eval(pt, eta);
231 else
232 {
233 Int_t xbin, ybin;
234
235 xbin = pt < d0ErrorHist->GetXaxis ()->GetXmax () ? d0ErrorHist->GetXaxis ()->FindBin (pt)
236 : d0ErrorHist->GetXaxis ()->GetBinCenter (d0ErrorHist->GetXaxis ()->GetNbins ());
237 ybin = d0ErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
238 d0Error = d0ErrorHist->GetBinContent (xbin, ybin);
239 if (!d0Error)
240 d0Error = -1.0;
241 }
242 if (d0Error < 0.0)
243 continue;
244
245 if (fUseDZFormula)
246 dzError = fDZFormula->Eval(pt, eta);
247 else
248 {
249 Int_t xbin, ybin;
250
251 xbin = pt < dzErrorHist->GetXaxis ()->GetXmax () ? dzErrorHist->GetXaxis ()->FindBin (pt)
252 : dzErrorHist->GetXaxis ()->GetBinCenter (dzErrorHist->GetXaxis ()->GetNbins ());
253 ybin = dzErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
254 dzError = dzErrorHist->GetBinContent (xbin, ybin);
255 if (!dzError)
256 dzError = -1.0;
257 }
258 if (dzError < 0.0)
259 continue;
260
261 if (fUsePFormula)
262 pError = fPFormula->Eval(pt, eta) * p;
263 else
264 {
265 Int_t xbin, ybin;
266
267 xbin = pt < pErrorHist->GetXaxis ()->GetXmax () ? pErrorHist->GetXaxis ()->FindBin (pt)
268 : pErrorHist->GetXaxis ()->GetBinCenter (pErrorHist->GetXaxis ()->GetNbins ());
269 ybin = pErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
270 pError = pErrorHist->GetBinContent (xbin, ybin) * p;
271 if (!pError)
272 pError = -1.0;
273 }
274 if (pError < 0.0)
275 continue;
276
277 if (fUseCtgThetaFormula)
278 ctgThetaError = fCtgThetaFormula->Eval(pt, eta);
279 else
280 {
281 Int_t xbin, ybin;
282
283 xbin = pt < ctgThetaErrorHist->GetXaxis ()->GetXmax () ? ctgThetaErrorHist->GetXaxis ()->FindBin (pt)
284 : ctgThetaErrorHist->GetXaxis ()->GetBinCenter (ctgThetaErrorHist->GetXaxis ()->GetNbins ());
285 ybin = ctgThetaErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
286 ctgThetaError = ctgThetaErrorHist->GetBinContent (xbin, ybin);
287 if (!ctgThetaError)
288 ctgThetaError = -1.0;
289 }
290 if (ctgThetaError < 0.0)
291 continue;
292
293 if (fUsePhiFormula)
294 phiError = fPhiFormula->Eval(pt, eta);
295 else
296 {
297 Int_t xbin, ybin;
298
299 xbin = pt < phiErrorHist->GetXaxis ()->GetXmax () ? phiErrorHist->GetXaxis ()->FindBin (pt)
300 : phiErrorHist->GetXaxis ()->GetBinCenter (phiErrorHist->GetXaxis ()->GetNbins ());
301 ybin = phiErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
302 phiError = phiErrorHist->GetBinContent (xbin, ybin);
303 if (!phiError)
304 phiError = -1.0;
305 }
306 if (phiError < 0.0)
307 continue;
308
309 if (fApplyToPileUp || !candidate->IsPU)
310 {
311 d0 = gRandom->Gaus(d0, d0Error);
312 dz = gRandom->Gaus(dz, dzError);
313 p = gRandom->Gaus(p, pError);
314 ctgTheta = gRandom->Gaus(ctgTheta, ctgThetaError);
315 phi = gRandom->Gaus(phi, phiError);
316 }
317
318 if(p < 0.0) continue;
319 while (phi > TMath::Pi ()) phi -= TMath::TwoPi ();
320 while (phi <= -TMath::Pi ()) phi += TMath::TwoPi ();
321
322 mother = candidate;
323 candidate = static_cast<Candidate*>(candidate->Clone());
324 candidate->D0 = d0;
325 candidate->DZ = dz;
326 candidate->P = p;
327 candidate->CtgTheta = ctgTheta;
328 candidate->Phi = phi;
329
330 theta = TMath::ACos(ctgTheta / TMath::Sqrt (1.0 + ctgTheta * ctgTheta));
331 candidate->Momentum.SetPx (p * TMath::Cos (phi) * TMath::Sin (theta));
332 candidate->Momentum.SetPy (p * TMath::Sin (phi) * TMath::Sin (theta));
333 candidate->Momentum.SetPz (p * TMath::Cos (theta));
334 candidate->Momentum.SetE (candidate->Momentum.Pt () * TMath::CosH (eta));
335 candidate->PT = candidate->Momentum.Pt ();
336
337 x = position.X ();
338 y = position.Y ();
339 z = position.Z ();
340 t = position.T ();
341 px = candidate->Momentum.Px ();
342 py = candidate->Momentum.Py ();
343 pz = candidate->Momentum.Pz ();
344 pt = candidate->Momentum.Pt ();
345
346 // -- solve for delta: d0' = ( (x+delta)*py' - (y+delta)*px' )/pt'
347
348 candidate->InitialPosition.SetX (x + ((px * y - py * x + d0 * pt) / (py - px)));
349 candidate->InitialPosition.SetY (y + ((px * y - py * x + d0 * pt) / (py - px)));
350 x = candidate->InitialPosition.X ();
351 y = candidate->InitialPosition.Y ();
352 candidate->InitialPosition.SetZ (z + ((pz * (px * (x - beamSpotPosition.X ()) + py * (y - beamSpotPosition.Y ())) + pt * pt * (dz - z)) / (pt * pt)));
353
354
355 candidate->InitialPosition.SetT (t);
356
357 if (fApplyToPileUp || !candidate->IsPU)
358 {
359 candidate->ErrorD0 = d0Error;
360 candidate->ErrorDZ = dzError;
361 candidate->ErrorP = pError;
362 candidate->ErrorCtgTheta = ctgThetaError;
363 candidate->ErrorPhi = phiError;
364 candidate->ErrorPT = ptError (p, ctgTheta, pError, ctgThetaError);
365 candidate->TrackResolution = pError/p;
366 }
367
368 candidate->AddCandidate(mother);
369 fOutputArray->Add(candidate);
370
371 iCandidate++;
372 }
373}
374
375Double_t TrackSmearing::ptError (const Double_t p, const Double_t ctgTheta, const Double_t dP, const Double_t dCtgTheta)
376{
377 Double_t a, b;
378 a = (p * p * ctgTheta * ctgTheta * dCtgTheta * dCtgTheta) / ((ctgTheta * ctgTheta + 1) * (ctgTheta * ctgTheta + 1) * (ctgTheta * ctgTheta + 1));
379 b = (dP * dP) / (ctgTheta * ctgTheta + 1);
380 return sqrt (a + b);
381}
382
383//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.