Fork me on GitHub

source: git/modules/TrackSmearing.cc@ 6e91f77

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

added comment on smeared initial position

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