Fork me on GitHub

source: git/modules/TrackSmearing.cc@ 5a697dde

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 5a697dde was 9fccf4f, checked in by Michele Selvaggi <michele.selvaggi@…>, 9 years ago

remove empty spaces

  • 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 * \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 // read resolution formula
64
65 // !!! IF WE WANT TO KEEP ROOT INPUT !!!
66 if (string (GetString("D0ResolutionFormula", "0.0")) != "0.0")
67 {
68 fD0Formula->Compile(GetString("D0ResolutionFormula", "0.0"));
69 fUseD0Formula = true;
70 }
71 else
72 {
73 fD0ResolutionFile = GetString("D0ResolutionFile", "errors.root");
74 fD0ResolutionHist = GetString("D0ResolutionHist", "d0");
75 fUseD0Formula = false;
76 }
77 if (string (GetString("DZResolutionFormula", "0.0")) != "0.0")
78 {
79 fDZFormula->Compile(GetString("DZResolutionFormula", "0.0"));
80 fUseDZFormula = true;
81 }
82 else
83 {
84 fDZResolutionFile = GetString("DZResolutionFile", "errors.root");
85 fDZResolutionHist = GetString("DZResolutionHist", "dz");
86 fUseDZFormula = false;
87 }
88 if (string (GetString("PResolutionFormula", "0.0")) != "0.0")
89 {
90 fPFormula->Compile(GetString("PResolutionFormula", "0.0"));
91 fUsePFormula = true;
92 }
93 else
94 {
95 fPResolutionFile = GetString("PResolutionFile", "errors.root");
96 fPResolutionHist = GetString("PResolutionHist", "p");
97 fUsePFormula = false;
98 }
99 if (string (GetString("CtgThetaResolutionFormula", "0.0")) != "0.0")
100 {
101 fCtgThetaFormula->Compile(GetString("CtgThetaResolutionFormula", "0.0"));
102 fUseCtgThetaFormula = true;
103 }
104 else
105 {
106 fCtgThetaResolutionFile = GetString("CtgThetaResolutionFile", "errors.root");
107 fCtgThetaResolutionHist = GetString("CtgThetaResolutionHist", "ctgTheta");
108 fUseCtgThetaFormula = false;
109 }
110 if (string (GetString("PhiResolutionFormula", "0.0")) != "0.0")
111 {
112 fPhiFormula->Compile(GetString("PhiResolutionFormula", "0.0"));
113 fUsePhiFormula = true;
114 }
115 else
116 {
117 fPhiResolutionFile = GetString("PhiResolutionFile", "errors.root");
118 fPhiResolutionHist = GetString("PhiResolutionHist", "phi");
119 fUsePhiFormula = false;
120 }
121
122 fApplyToPileUp = GetBool("ApplyToPileUp", true);
123
124 // import input array
125
126 fInputArray = ImportArray(GetString("InputArray", "ParticlePropagator/stableParticles"));
127 fItInputArray = fInputArray->MakeIterator();
128
129 // import beamspot
130 try
131 {
132 fBeamSpotInputArray = ImportArray(GetString("BeamSpotInputArray", "BeamSpotFilter/beamSpotParticle"));
133 }
134 catch(runtime_error &e)
135 {
136 fBeamSpotInputArray = 0;
137 }
138
139 // create output array
140
141 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
142}
143
144//------------------------------------------------------------------------------
145
146void TrackSmearing::Finish()
147{
148 if(fItInputArray) delete fItInputArray;
149}
150
151//------------------------------------------------------------------------------
152
153void TrackSmearing::Process()
154{
155 Int_t iCandidate = 0;
156 TLorentzVector beamSpotPosition;
157 Candidate *candidate, *mother;
158 Double_t pt, eta, d0, d0Error, trueD0, dz, dzError, trueDZ, p, pError, trueP, ctgTheta, ctgThetaError, trueCtgTheta, phi, phiError, truePhi;
159 Double_t x, y, z, t, px, py, pz, theta;
160 TProfile2D *d0ErrorHist = NULL,
161 *dzErrorHist = NULL,
162 *pErrorHist = NULL,
163 *ctgThetaErrorHist = NULL,
164 *phiErrorHist = NULL;
165
166 if (!fBeamSpotInputArray || fBeamSpotInputArray->GetSize () == 0)
167 beamSpotPosition.SetXYZT(0.0, 0.0, 0.0, 0.0);
168 else
169 {
170 Candidate &beamSpotCandidate = *((Candidate *) fBeamSpotInputArray->At (0));
171 beamSpotPosition = beamSpotCandidate.Position;
172 }
173
174 if (!fUseD0Formula)
175 {
176 TFile *fin = TFile::Open (fD0ResolutionFile.c_str ());
177 d0ErrorHist = (TProfile2D *) fin->Get (fD0ResolutionHist.c_str ());
178 d0ErrorHist->SetDirectory (0);
179 fin->Close ();
180 }
181 if (!fUseDZFormula)
182 {
183 TFile *fin = TFile::Open (fDZResolutionFile.c_str ());
184 dzErrorHist = (TProfile2D *) fin->Get (fDZResolutionHist.c_str ());
185 dzErrorHist->SetDirectory (0);
186 fin->Close ();
187 }
188 if (!fUsePFormula)
189 {
190 TFile *fin = TFile::Open (fPResolutionFile.c_str ());
191 pErrorHist = (TProfile2D *) fin->Get (fPResolutionHist.c_str ());
192 pErrorHist->SetDirectory (0);
193 fin->Close ();
194 }
195 if (!fUseCtgThetaFormula)
196 {
197 TFile *fin = TFile::Open (fCtgThetaResolutionFile.c_str ());
198 ctgThetaErrorHist = (TProfile2D *) fin->Get (fCtgThetaResolutionHist.c_str ());
199 ctgThetaErrorHist->SetDirectory (0);
200 fin->Close ();
201 }
202 if (!fUsePhiFormula)
203 {
204 TFile *fin = TFile::Open (fPhiResolutionFile.c_str ());
205 phiErrorHist = (TProfile2D *) fin->Get (fPhiResolutionHist.c_str ());
206 phiErrorHist->SetDirectory (0);
207 fin->Close ();
208 }
209
210 fItInputArray->Reset();
211 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
212 {
213
214 const TLorentzVector &momentum = candidate->Momentum;
215 const TLorentzVector &position = candidate->InitialPosition;
216
217 pt = momentum.Pt();
218 eta = momentum.Eta();
219
220 d0 = trueD0 = candidate->D0;
221 dz = trueDZ = candidate->DZ;
222
223 p = trueP = candidate->P;
224 ctgTheta = trueCtgTheta = candidate->CtgTheta;
225 phi = truePhi = candidate->Phi;
226
227 if (fUseD0Formula)
228 d0Error = fD0Formula->Eval(pt, eta);
229 else
230 {
231 Int_t xbin, ybin;
232
233 xbin = pt < d0ErrorHist->GetXaxis ()->GetXmax () ? d0ErrorHist->GetXaxis ()->FindBin (pt)
234 : d0ErrorHist->GetXaxis ()->GetBinCenter (d0ErrorHist->GetXaxis ()->GetNbins ());
235 ybin = d0ErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
236 d0Error = d0ErrorHist->GetBinContent (xbin, ybin);
237 if (!d0Error)
238 d0Error = -1.0;
239 }
240 if (d0Error < 0.0)
241 continue;
242
243 if (fUseDZFormula)
244 dzError = fDZFormula->Eval(pt, eta);
245 else
246 {
247 Int_t xbin, ybin;
248
249 xbin = pt < dzErrorHist->GetXaxis ()->GetXmax () ? dzErrorHist->GetXaxis ()->FindBin (pt)
250 : dzErrorHist->GetXaxis ()->GetBinCenter (dzErrorHist->GetXaxis ()->GetNbins ());
251 ybin = dzErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
252 dzError = dzErrorHist->GetBinContent (xbin, ybin);
253 if (!dzError)
254 dzError = -1.0;
255 }
256 if (dzError < 0.0)
257 continue;
258
259 if (fUsePFormula)
260 pError = fPFormula->Eval(pt, eta) * p;
261 else
262 {
263 Int_t xbin, ybin;
264
265 xbin = pt < pErrorHist->GetXaxis ()->GetXmax () ? pErrorHist->GetXaxis ()->FindBin (pt)
266 : pErrorHist->GetXaxis ()->GetBinCenter (pErrorHist->GetXaxis ()->GetNbins ());
267 ybin = pErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
268 pError = pErrorHist->GetBinContent (xbin, ybin) * p;
269 if (!pError)
270 pError = -1.0;
271 }
272 if (pError < 0.0)
273 continue;
274
275 if (fUseCtgThetaFormula)
276 ctgThetaError = fCtgThetaFormula->Eval(pt, eta);
277 else
278 {
279 Int_t xbin, ybin;
280
281 xbin = pt < ctgThetaErrorHist->GetXaxis ()->GetXmax () ? ctgThetaErrorHist->GetXaxis ()->FindBin (pt)
282 : ctgThetaErrorHist->GetXaxis ()->GetBinCenter (ctgThetaErrorHist->GetXaxis ()->GetNbins ());
283 ybin = ctgThetaErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
284 ctgThetaError = ctgThetaErrorHist->GetBinContent (xbin, ybin);
285 if (!ctgThetaError)
286 ctgThetaError = -1.0;
287 }
288 if (ctgThetaError < 0.0)
289 continue;
290
291 if (fUsePhiFormula)
292 phiError = fPhiFormula->Eval(pt, eta);
293 else
294 {
295 Int_t xbin, ybin;
296
297 xbin = pt < phiErrorHist->GetXaxis ()->GetXmax () ? phiErrorHist->GetXaxis ()->FindBin (pt)
298 : phiErrorHist->GetXaxis ()->GetBinCenter (phiErrorHist->GetXaxis ()->GetNbins ());
299 ybin = phiErrorHist->GetYaxis ()->FindBin (TMath::Abs (eta));
300 phiError = phiErrorHist->GetBinContent (xbin, ybin);
301 if (!phiError)
302 phiError = -1.0;
303 }
304 if (phiError < 0.0)
305 continue;
306
307 if (fApplyToPileUp || !candidate->IsPU)
308 {
309 d0 = gRandom->Gaus(d0, d0Error);
310 dz = gRandom->Gaus(dz, dzError);
311 p = gRandom->Gaus(p, pError);
312 ctgTheta = gRandom->Gaus(ctgTheta, ctgThetaError);
313 phi = gRandom->Gaus(phi, phiError);
314 }
315
316 if(p < 0.0) continue;
317 while (phi > TMath::Pi ()) phi -= TMath::TwoPi ();
318 while (phi <= -TMath::Pi ()) phi += TMath::TwoPi ();
319
320 mother = candidate;
321 candidate = static_cast<Candidate*>(candidate->Clone());
322 candidate->D0 = d0;
323 candidate->DZ = dz;
324 candidate->P = p;
325 candidate->CtgTheta = ctgTheta;
326 candidate->Phi = phi;
327
328 theta = TMath::ACos(ctgTheta / TMath::Sqrt (1.0 + ctgTheta * ctgTheta));
329 candidate->Momentum.SetPx (p * TMath::Cos (phi) * TMath::Sin (theta));
330 candidate->Momentum.SetPy (p * TMath::Sin (phi) * TMath::Sin (theta));
331 candidate->Momentum.SetPz (p * TMath::Cos (theta));
332 candidate->Momentum.SetE (candidate->Momentum.Pt () * TMath::CosH (eta));
333 candidate->PT = candidate->Momentum.Pt ();
334
335 x = position.X ();
336 y = position.Y ();
337 z = position.Z ();
338 t = position.T ();
339 px = candidate->Momentum.Px ();
340 py = candidate->Momentum.Py ();
341 pz = candidate->Momentum.Pz ();
342 pt = candidate->Momentum.Pt ();
343
344 // -- solve for delta: d0' = ( (x+delta)*py' - (y+delta)*px' )/pt'
345
346 candidate->InitialPosition.SetX (x + ((px * y - py * x + d0 * pt) / (py - px)));
347 candidate->InitialPosition.SetY (y + ((px * y - py * x + d0 * pt) / (py - px)));
348 x = candidate->InitialPosition.X ();
349 y = candidate->InitialPosition.Y ();
350 candidate->InitialPosition.SetZ (z + ((pz * (px * (x - beamSpotPosition.X ()) + py * (y - beamSpotPosition.Y ())) + pt * pt * (dz - z)) / (pt * pt)));
351
352 candidate->InitialPosition.SetT(t);
353
354 if (fApplyToPileUp || !candidate->IsPU)
355 {
356 candidate->ErrorD0 = d0Error;
357 candidate->ErrorDZ = dzError;
358 candidate->ErrorP = pError;
359 candidate->ErrorCtgTheta = ctgThetaError;
360 candidate->ErrorPhi = phiError;
361 candidate->ErrorPT = ptError (p, ctgTheta, pError, ctgThetaError);
362 candidate->TrackResolution = pError/p;
363 }
364
365 candidate->AddCandidate(mother);
366 fOutputArray->Add(candidate);
367
368 iCandidate++;
369 }
370}
371
372Double_t TrackSmearing::ptError (const Double_t p, const Double_t ctgTheta, const Double_t dP, const Double_t dCtgTheta)
373{
374 Double_t a, b;
375 a = (p * p * ctgTheta * ctgTheta * dCtgTheta * dCtgTheta) / ((ctgTheta * ctgTheta + 1) * (ctgTheta * ctgTheta + 1) * (ctgTheta * ctgTheta + 1));
376 b = (dP * dP) / (ctgTheta * ctgTheta + 1);
377 return sqrt (a + b);
378}
379
380//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.