Fork me on GitHub

source: svn/trunk/modules/ImpactParameterSmearing.cc

Last change on this file was 1372, checked in by Pavel Demin, 11 years ago

several minor corrections

File size: 3.6 KB
Line 
1/** \class ImpactParameterSmearing
2 *
3 * Performs transverse impact parameter smearing.
4 *
5 * $Date: 2014-16-03 14:57:44 +0100
6 *
7 *
8 * \author M. Selvaggi - UCL, Louvain-la-Neuve
9 *
10 */
11
12
13#include "modules/ImpactParameterSmearing.h"
14
15#include "classes/DelphesClasses.h"
16#include "classes/DelphesFactory.h"
17#include "classes/DelphesFormula.h"
18
19#include "ExRootAnalysis/ExRootResult.h"
20#include "ExRootAnalysis/ExRootFilter.h"
21#include "ExRootAnalysis/ExRootClassifier.h"
22
23#include "TMath.h"
24#include "TString.h"
25#include "TFormula.h"
26#include "TRandom3.h"
27#include "TObjArray.h"
28#include "TDatabasePDG.h"
29#include "TLorentzVector.h"
30
31#include <algorithm>
32#include <stdexcept>
33#include <iostream>
34#include <sstream>
35
36using namespace std;
37
38//------------------------------------------------------------------------------
39
40ImpactParameterSmearing::ImpactParameterSmearing() :
41 fFormula(0), fItInputArray(0)
42{
43 fFormula = new DelphesFormula;
44}
45
46//------------------------------------------------------------------------------
47
48ImpactParameterSmearing::~ImpactParameterSmearing()
49{
50 if(fFormula) delete fFormula;
51}
52
53//------------------------------------------------------------------------------
54
55void ImpactParameterSmearing::Init()
56{
57 // read resolution formula
58
59 fFormula->Compile(GetString("ResolutionFormula", "0.0"));
60
61 // import input array
62
63 fInputArray = ImportArray(GetString("InputArray", "TrackMerger/tracks"));
64 fItInputArray = fInputArray->MakeIterator();
65
66 // create output array
67
68 fOutputArray = ExportArray(GetString("OutputArray", "tracks"));
69}
70
71//------------------------------------------------------------------------------
72
73void ImpactParameterSmearing::Finish()
74{
75 if(fItInputArray) delete fItInputArray;
76}
77
78//------------------------------------------------------------------------------
79
80void ImpactParameterSmearing::Process()
81{
82 Candidate *candidate, *particle, *mother;
83 Double_t xd, yd, zd, dxy, dz, sx, sy, sz, ddxy, ddz;
84 Double_t pt, eta, px, py, ang_mom;
85
86 // cout<<"New event"<<endl;
87
88 fItInputArray->Reset();
89 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
90 {
91
92 //take momentum before smearing (otherwise apply double smearing on dxy)
93 particle = static_cast<Candidate*>(candidate->GetCandidates()->At(0));
94
95 const TLorentzVector &candidateMomentum = particle->Momentum;
96 // const TLorentzVector &candidateMomentum = candidate->Momentum;
97
98 eta = candidateMomentum.Eta();
99 pt = candidateMomentum.Pt();
100 px = candidateMomentum.Px();
101 py = candidateMomentum.Py();
102
103 // calculate coordinates of closest approach to track circle in transverse plane xd, yd, zd
104 xd = candidate->Xd;
105 yd = candidate->Yd;
106 zd = candidate->Zd;
107
108 // calculate smeared values
109 sx = gRandom->Gaus(0,fFormula->Eval(pt, eta));
110 sy = gRandom->Gaus(0,fFormula->Eval(pt, eta));
111 sz = gRandom->Gaus(0,fFormula->Eval(pt, eta));
112
113 xd += sx;
114 yd += sy;
115 zd += sz;
116
117 // calculate impact paramater (after-smearing)
118 ang_mom = (xd*py - yd*px);
119 dxy = ang_mom/pt;
120 dz = zd;
121
122 ddxy = gRandom->Gaus(0,fFormula->Eval(pt, eta));
123 ddz = gRandom->Gaus(0,fFormula->Eval(pt, eta));
124
125 //fill smeared values in candidate
126 mother = candidate;
127
128 candidate = static_cast<Candidate*>(candidate->Clone());
129 candidate->Xd = xd;
130 candidate->Yd = yd;
131 candidate->Zd = zd;
132
133 candidate->Dxy = dxy;
134 candidate->SDxy = ddxy;
135
136 candidate->AddCandidate(mother);
137 fOutputArray->Add(candidate);
138 }
139}
140
141//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.