Fork me on GitHub

source: git/modules/TrackSmearing.cc@ c438758

Last change on this file since c438758 was fd4b326, checked in by michele <michele.selvaggi@…>, 4 years ago

tracks have now montecarlo mass

  • Property mode set to 100644
File size: 12.4 KB
RevLine 
[d1e6379]1/** \class TrackSmearing
2 *
3 * Performs d0, dZ, p, Theta, Phi smearing of tracks.
4 *
[5658083]5 * \authors A. Hart, M. Selvaggi
[d1e6379]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/ExRootClassifier.h"
[341014c]16#include "ExRootAnalysis/ExRootFilter.h"
17#include "ExRootAnalysis/ExRootResult.h"
[d1e6379]18
19#include "TDatabasePDG.h"
20#include "TFile.h"
[341014c]21#include "TFormula.h"
22#include "TLorentzVector.h"
23#include "TMath.h"
24#include "TObjArray.h"
[d1e6379]25#include "TProfile2D.h"
[341014c]26#include "TRandom3.h"
27#include "TString.h"
[d1e6379]28
[341014c]29#include <algorithm>
[d1e6379]30#include <iostream>
31#include <sstream>
[341014c]32#include <stdexcept>
[d1e6379]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{
[cde9f31]63 fBz = GetDouble("Bz", 0.0);
64
[d1e6379]65 // read resolution formula
66
67 // !!! IF WE WANT TO KEEP ROOT INPUT !!!
[341014c]68 if(string(GetString("D0ResolutionFormula", "0.0")) != "0.0")
69 {
70 fD0Formula->Compile(GetString("D0ResolutionFormula", "0.0"));
71 fUseD0Formula = true;
72 }
[d1e6379]73 else
[341014c]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 }
[d1e6379]84 else
[341014c]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 }
[d1e6379]95 else
[341014c]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 }
[d1e6379]106 else
[341014c]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 }
[d1e6379]117 else
[341014c]118 {
119 fPhiResolutionFile = GetString("PhiResolutionFile", "errors.root");
120 fPhiResolutionHist = GetString("PhiResolutionHist", "phi");
121 fUsePhiFormula = false;
122 }
[d1e6379]123
124 fApplyToPileUp = GetBool("ApplyToPileUp", true);
125
126 // import input array
127
128 fInputArray = ImportArray(GetString("InputArray", "ParticlePropagator/stableParticles"));
129 fItInputArray = fInputArray->MakeIterator();
[341014c]130
[bff2e33]131 // import beamspot
132 try
133 {
134 fBeamSpotInputArray = ImportArray(GetString("BeamSpotInputArray", "BeamSpotFilter/beamSpotParticle"));
135 }
136 catch(runtime_error &e)
137 {
138 fBeamSpotInputArray = 0;
[341014c]139 }
140
[d1e6379]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;
[fd4b326]160 Double_t pt, eta, e, m, d0, d0Error, trueD0, dz, dzError, trueDZ, p, pError, trueP, ctgTheta, ctgThetaError, trueCtgTheta, phi, phiError, truePhi;
[d1e6379]161 Double_t x, y, z, t, px, py, pz, theta;
[cde9f31]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;
[d1e6379]166 TProfile2D *d0ErrorHist = NULL,
167 *dzErrorHist = NULL,
168 *pErrorHist = NULL,
169 *ctgThetaErrorHist = NULL,
170 *phiErrorHist = NULL;
171
[341014c]172 if(!fBeamSpotInputArray || fBeamSpotInputArray->GetSize() == 0)
[d1e6379]173 beamSpotPosition.SetXYZT(0.0, 0.0, 0.0, 0.0);
174 else
175 {
[341014c]176 Candidate &beamSpotCandidate = *((Candidate *)fBeamSpotInputArray->At(0));
[d1e6379]177 beamSpotPosition = beamSpotCandidate.Position;
178 }
179
[341014c]180 if(!fUseD0Formula)
[d1e6379]181 {
[341014c]182 TFile *fin = TFile::Open(fD0ResolutionFile.c_str());
183 d0ErrorHist = (TProfile2D *)fin->Get(fD0ResolutionHist.c_str());
184 d0ErrorHist->SetDirectory(0);
185 fin->Close();
[d1e6379]186 }
[341014c]187 if(!fUseDZFormula)
[d1e6379]188 {
[341014c]189 TFile *fin = TFile::Open(fDZResolutionFile.c_str());
190 dzErrorHist = (TProfile2D *)fin->Get(fDZResolutionHist.c_str());
191 dzErrorHist->SetDirectory(0);
192 fin->Close();
[d1e6379]193 }
[341014c]194 if(!fUsePFormula)
[d1e6379]195 {
[341014c]196 TFile *fin = TFile::Open(fPResolutionFile.c_str());
197 pErrorHist = (TProfile2D *)fin->Get(fPResolutionHist.c_str());
198 pErrorHist->SetDirectory(0);
199 fin->Close();
[d1e6379]200 }
[341014c]201 if(!fUseCtgThetaFormula)
[d1e6379]202 {
[341014c]203 TFile *fin = TFile::Open(fCtgThetaResolutionFile.c_str());
204 ctgThetaErrorHist = (TProfile2D *)fin->Get(fCtgThetaResolutionHist.c_str());
205 ctgThetaErrorHist->SetDirectory(0);
206 fin->Close();
[d1e6379]207 }
[341014c]208 if(!fUsePhiFormula)
[d1e6379]209 {
[341014c]210 TFile *fin = TFile::Open(fPhiResolutionFile.c_str());
211 phiErrorHist = (TProfile2D *)fin->Get(fPhiResolutionHist.c_str());
212 phiErrorHist->SetDirectory(0);
213 fin->Close();
[d1e6379]214 }
215
216 fItInputArray->Reset();
[341014c]217 while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
[d1e6379]218 {
[341014c]219
[d1e6379]220 const TLorentzVector &momentum = candidate->Momentum;
221 const TLorentzVector &position = candidate->InitialPosition;
[341014c]222
[d1e6379]223 pt = momentum.Pt();
224 eta = momentum.Eta();
[f298e48]225 e = momentum.E();
[fd4b326]226 m = momentum.M();
227
[d1e6379]228 d0 = trueD0 = candidate->D0;
229 dz = trueDZ = candidate->DZ;
[341014c]230
[d1e6379]231 p = trueP = candidate->P;
232 ctgTheta = trueCtgTheta = candidate->CtgTheta;
233 phi = truePhi = candidate->Phi;
234
[341014c]235 if(fUseD0Formula)
[f298e48]236 d0Error = fD0Formula->Eval(pt, eta, phi, e, candidate);
[d1e6379]237 else
238 {
[341014c]239 Int_t xbin, ybin;
240
241 xbin = pt < d0ErrorHist->GetXaxis()->GetXmax() ? d0ErrorHist->GetXaxis()->FindBin(pt) : d0ErrorHist->GetXaxis()->GetBinCenter(d0ErrorHist->GetXaxis()->GetNbins());
242 ybin = d0ErrorHist->GetYaxis()->FindBin(TMath::Abs(eta));
243 d0Error = d0ErrorHist->GetBinContent(xbin, ybin);
244 if(!d0Error)
245 d0Error = -1.0;
[d1e6379]246 }
[341014c]247 if(d0Error < 0.0)
[d1e6379]248 continue;
249
[341014c]250 if(fUseDZFormula)
[f298e48]251 dzError = fDZFormula->Eval(pt, eta, phi, e, candidate);
[d1e6379]252 else
253 {
[341014c]254 Int_t xbin, ybin;
255
256 xbin = pt < dzErrorHist->GetXaxis()->GetXmax() ? dzErrorHist->GetXaxis()->FindBin(pt) : 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;
[d1e6379]261 }
[341014c]262 if(dzError < 0.0)
[d1e6379]263 continue;
264
[341014c]265 if(fUsePFormula)
[f298e48]266 pError = fPFormula->Eval(pt, eta, phi, e, candidate) * p;
[d1e6379]267 else
268 {
[341014c]269 Int_t xbin, ybin;
270
271 xbin = pt < pErrorHist->GetXaxis()->GetXmax() ? pErrorHist->GetXaxis()->FindBin(pt) : pErrorHist->GetXaxis()->GetBinCenter(pErrorHist->GetXaxis()->GetNbins());
272 ybin = pErrorHist->GetYaxis()->FindBin(TMath::Abs(eta));
273 pError = pErrorHist->GetBinContent(xbin, ybin) * p;
274 if(!pError)
275 pError = -1.0;
[d1e6379]276 }
[341014c]277 if(pError < 0.0)
[d1e6379]278 continue;
279
[341014c]280 if(fUseCtgThetaFormula)
[f298e48]281 ctgThetaError = fCtgThetaFormula->Eval(pt, eta, phi, e, candidate);
[d1e6379]282 else
283 {
[341014c]284 Int_t xbin, ybin;
285
286 xbin = pt < ctgThetaErrorHist->GetXaxis()->GetXmax() ? ctgThetaErrorHist->GetXaxis()->FindBin(pt) : ctgThetaErrorHist->GetXaxis()->GetBinCenter(ctgThetaErrorHist->GetXaxis()->GetNbins());
287 ybin = ctgThetaErrorHist->GetYaxis()->FindBin(TMath::Abs(eta));
288 ctgThetaError = ctgThetaErrorHist->GetBinContent(xbin, ybin);
289 if(!ctgThetaError)
290 ctgThetaError = -1.0;
[d1e6379]291 }
[341014c]292 if(ctgThetaError < 0.0)
[d1e6379]293 continue;
294
[341014c]295 if(fUsePhiFormula)
[f298e48]296 phiError = fPhiFormula->Eval(pt, eta, phi, e, candidate);
[d1e6379]297 else
298 {
[341014c]299 Int_t xbin, ybin;
300
301 xbin = pt < phiErrorHist->GetXaxis()->GetXmax() ? phiErrorHist->GetXaxis()->FindBin(pt) : phiErrorHist->GetXaxis()->GetBinCenter(phiErrorHist->GetXaxis()->GetNbins());
302 ybin = phiErrorHist->GetYaxis()->FindBin(TMath::Abs(eta));
303 phiError = phiErrorHist->GetBinContent(xbin, ybin);
304 if(!phiError)
305 phiError = -1.0;
[d1e6379]306 }
[341014c]307 if(phiError < 0.0)
[d1e6379]308 continue;
309
[341014c]310 if(fApplyToPileUp || !candidate->IsPU)
[d1e6379]311 {
[341014c]312 d0 = gRandom->Gaus(d0, d0Error);
313 dz = gRandom->Gaus(dz, dzError);
314 p = gRandom->Gaus(p, pError);
315 ctgTheta = gRandom->Gaus(ctgTheta, ctgThetaError);
316 phi = gRandom->Gaus(phi, phiError);
[d1e6379]317 }
318
319 if(p < 0.0) continue;
[341014c]320 while(phi > TMath::Pi()) phi -= TMath::TwoPi();
321 while(phi <= -TMath::Pi()) phi += TMath::TwoPi();
[d1e6379]322
323 mother = candidate;
[341014c]324 candidate = static_cast<Candidate *>(candidate->Clone());
[d1e6379]325 candidate->D0 = d0;
326 candidate->DZ = dz;
327 candidate->P = p;
328 candidate->CtgTheta = ctgTheta;
329 candidate->Phi = phi;
330
[341014c]331 theta = TMath::ACos(ctgTheta / TMath::Sqrt(1.0 + ctgTheta * ctgTheta));
332 candidate->Momentum.SetPx(p * TMath::Cos(phi) * TMath::Sin(theta));
333 candidate->Momentum.SetPy(p * TMath::Sin(phi) * TMath::Sin(theta));
334 candidate->Momentum.SetPz(p * TMath::Cos(theta));
[fd4b326]335 candidate->Momentum.SetE(TMath::Sqrt(p*p + m*m));
[341014c]336 candidate->PT = candidate->Momentum.Pt();
337
338 x = position.X();
339 y = position.Y();
340 z = position.Z();
341 t = position.T();
342 px = candidate->Momentum.Px();
343 py = candidate->Momentum.Py();
344 pz = candidate->Momentum.Pz();
345 pt = candidate->Momentum.Pt();
346
[6e91f77]347 // -- solve for delta: d0' = ( (x+delta)*py' - (y+delta)*px' )/pt'
[341014c]348
349 candidate->InitialPosition.SetX(x + ((px * y - py * x + d0 * pt) / (py - px)));
350 candidate->InitialPosition.SetY(y + ((px * y - py * x + d0 * pt) / (py - px)));
351 x = candidate->InitialPosition.X();
352 y = candidate->InitialPosition.Y();
353 candidate->InitialPosition.SetZ(z + ((pz * (px * (x - beamSpotPosition.X()) + py * (y - beamSpotPosition.Y())) + pt * pt * (dz - z)) / (pt * pt)));
354 z = candidate->InitialPosition.Z();
[cde9f31]355
[9fccf4f]356 candidate->InitialPosition.SetT(t);
[cde9f31]357
358 // update closest approach
359 x *= 1.0E-3;
360 y *= 1.0E-3;
361 z *= 1.0E-3;
362
363 q = candidate->Charge;
364
[341014c]365 r = pt / (q * fBz) * 1.0E9 / c_light; // in [m]
[cde9f31]366 phi_0 = TMath::ATan2(py, px); // [rad] in [-pi, pi]
367
368 // 2. helix axis coordinates
[341014c]369 x_c = x + r * TMath::Sin(phi_0);
370 y_c = y - r * TMath::Cos(phi_0);
[cde9f31]371 r_c = TMath::Hypot(x_c, y_c);
372
373 rcu = TMath::Abs(r);
[341014c]374 rc2 = r_c * r_c;
[cde9f31]375
376 // calculate coordinates of closest approach to track circle in transverse plane xd, yd, zd
[341014c]377 xd = x_c * x_c * x_c - x_c * rcu * r_c + x_c * y_c * y_c;
[cde9f31]378 xd = (rc2 > 0.0) ? xd / rc2 : -999;
[341014c]379 yd = y_c * (-rcu * r_c + rc2);
[cde9f31]380 yd = (rc2 > 0.0) ? yd / rc2 : -999;
[341014c]381 zd = z + (TMath::Sqrt(xd * xd + yd * yd) - TMath::Sqrt(x * x + y * y)) * pz / pt;
[cde9f31]382
[341014c]383 candidate->Xd = xd * 1.0E3;
384 candidate->Yd = yd * 1.0E3;
385 candidate->Zd = zd * 1.0E3;
[cde9f31]386
[341014c]387 if(fApplyToPileUp || !candidate->IsPU)
[d1e6379]388 {
[341014c]389 candidate->ErrorD0 = d0Error;
390 candidate->ErrorDZ = dzError;
391 candidate->ErrorP = pError;
392 candidate->ErrorCtgTheta = ctgThetaError;
393 candidate->ErrorPhi = phiError;
394 candidate->ErrorPT = ptError(p, ctgTheta, pError, ctgThetaError);
395 candidate->TrackResolution = pError / p;
[d1e6379]396 }
[341014c]397
[d1e6379]398 candidate->AddCandidate(mother);
399 fOutputArray->Add(candidate);
400
401 iCandidate++;
402 }
403}
404
[341014c]405Double_t TrackSmearing::ptError(const Double_t p, const Double_t ctgTheta, const Double_t dP, const Double_t dCtgTheta)
[d1e6379]406{
407 Double_t a, b;
408 a = (p * p * ctgTheta * ctgTheta * dCtgTheta * dCtgTheta) / ((ctgTheta * ctgTheta + 1) * (ctgTheta * ctgTheta + 1) * (ctgTheta * ctgTheta + 1));
409 b = (dP * dP) / (ctgTheta * ctgTheta + 1);
[341014c]410 return sqrt(a + b);
[d1e6379]411}
412
413//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.