1 | /*
|
---|
2 | * Delphes: a framework for fast simulation of a generic collider experiment
|
---|
3 | * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
|
---|
4 | *
|
---|
5 | * This program is free software: you can redistribute it and/or modify
|
---|
6 | * it under the terms of the GNU General Public License as published by
|
---|
7 | * the Free Software Foundation, either version 3 of the License, or
|
---|
8 | * (at your option) any later version.
|
---|
9 | *
|
---|
10 | * This program is distributed in the hope that it will be useful,
|
---|
11 | * but WITHOUT ANY WARRANTY; without even the implied warranty of
|
---|
12 | * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
|
---|
13 | * GNU General Public License for more details.
|
---|
14 | *
|
---|
15 | * You should have received a copy of the GNU General Public License
|
---|
16 | * along with this program. If not, see <http://www.gnu.org/licenses/>.
|
---|
17 | */
|
---|
18 |
|
---|
19 | /** \class ImpactParameterSmearing
|
---|
20 | *
|
---|
21 | * Performs transverse impact parameter smearing.
|
---|
22 | *
|
---|
23 | * \author M. Selvaggi - UCL, Louvain-la-Neuve
|
---|
24 | *
|
---|
25 | */
|
---|
26 |
|
---|
27 |
|
---|
28 | #include "modules/ImpactParameterSmearing.h"
|
---|
29 |
|
---|
30 | #include "classes/DelphesClasses.h"
|
---|
31 | #include "classes/DelphesFactory.h"
|
---|
32 | #include "classes/DelphesFormula.h"
|
---|
33 |
|
---|
34 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
35 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
36 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
37 |
|
---|
38 | #include "TMath.h"
|
---|
39 | #include "TString.h"
|
---|
40 | #include "TFormula.h"
|
---|
41 | #include "TRandom3.h"
|
---|
42 | #include "TObjArray.h"
|
---|
43 | #include "TDatabasePDG.h"
|
---|
44 | #include "TLorentzVector.h"
|
---|
45 |
|
---|
46 | #include <algorithm>
|
---|
47 | #include <stdexcept>
|
---|
48 | #include <iostream>
|
---|
49 | #include <sstream>
|
---|
50 |
|
---|
51 | using namespace std;
|
---|
52 |
|
---|
53 | //------------------------------------------------------------------------------
|
---|
54 |
|
---|
55 | ImpactParameterSmearing::ImpactParameterSmearing() :
|
---|
56 | fFormula(0), fItInputArray(0)
|
---|
57 | {
|
---|
58 | fFormula = new DelphesFormula;
|
---|
59 | }
|
---|
60 |
|
---|
61 | //------------------------------------------------------------------------------
|
---|
62 |
|
---|
63 | ImpactParameterSmearing::~ImpactParameterSmearing()
|
---|
64 | {
|
---|
65 | if(fFormula) delete fFormula;
|
---|
66 | }
|
---|
67 |
|
---|
68 | //------------------------------------------------------------------------------
|
---|
69 |
|
---|
70 | void ImpactParameterSmearing::Init()
|
---|
71 | {
|
---|
72 | // read resolution formula
|
---|
73 |
|
---|
74 | fFormula->Compile(GetString("ResolutionFormula", "0.0"));
|
---|
75 |
|
---|
76 | // import input array
|
---|
77 |
|
---|
78 | fInputArray = ImportArray(GetString("InputArray", "TrackMerger/tracks"));
|
---|
79 | fItInputArray = fInputArray->MakeIterator();
|
---|
80 |
|
---|
81 | // create output array
|
---|
82 |
|
---|
83 | fOutputArray = ExportArray(GetString("OutputArray", "tracks"));
|
---|
84 | }
|
---|
85 |
|
---|
86 | //------------------------------------------------------------------------------
|
---|
87 |
|
---|
88 | void ImpactParameterSmearing::Finish()
|
---|
89 | {
|
---|
90 | if(fItInputArray) delete fItInputArray;
|
---|
91 | }
|
---|
92 |
|
---|
93 | //------------------------------------------------------------------------------
|
---|
94 |
|
---|
95 | void ImpactParameterSmearing::Process()
|
---|
96 | {
|
---|
97 | Candidate *candidate, *particle, *mother;
|
---|
98 | Double_t xd, yd, zd, dxy, sx, sy, sz, ddxy;
|
---|
99 | Double_t pt, eta, px, py;
|
---|
100 |
|
---|
101 | fItInputArray->Reset();
|
---|
102 | while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
|
---|
103 | {
|
---|
104 |
|
---|
105 | // take momentum before smearing (otherwise apply double smearing on dxy)
|
---|
106 | particle = static_cast<Candidate*>(candidate->GetCandidates()->At(0));
|
---|
107 |
|
---|
108 | const TLorentzVector &candidateMomentum = particle->Momentum;
|
---|
109 |
|
---|
110 | eta = candidateMomentum.Eta();
|
---|
111 | pt = candidateMomentum.Pt();
|
---|
112 | px = candidateMomentum.Px();
|
---|
113 | py = candidateMomentum.Py();
|
---|
114 |
|
---|
115 | // calculate coordinates of closest approach to track circle in transverse plane xd, yd, zd
|
---|
116 | xd = candidate->Xd;
|
---|
117 | yd = candidate->Yd;
|
---|
118 | zd = candidate->Zd;
|
---|
119 |
|
---|
120 | // calculate smeared values
|
---|
121 | sx = gRandom->Gaus(0.0, fFormula->Eval(pt, eta));
|
---|
122 | sy = gRandom->Gaus(0.0, fFormula->Eval(pt, eta));
|
---|
123 | sz = gRandom->Gaus(0.0, fFormula->Eval(pt, eta));
|
---|
124 |
|
---|
125 | xd += sx;
|
---|
126 | yd += sy;
|
---|
127 | zd += sz;
|
---|
128 |
|
---|
129 | // calculate impact parameter (after-smearing)
|
---|
130 | dxy = (xd*py - yd*px)/pt;
|
---|
131 |
|
---|
132 | ddxy = gRandom->Gaus(0.0, fFormula->Eval(pt, eta));
|
---|
133 |
|
---|
134 | // fill smeared values in candidate
|
---|
135 | mother = candidate;
|
---|
136 |
|
---|
137 | candidate = static_cast<Candidate*>(candidate->Clone());
|
---|
138 | candidate->Xd = xd;
|
---|
139 | candidate->Yd = yd;
|
---|
140 | candidate->Zd = zd;
|
---|
141 |
|
---|
142 | candidate->Dxy = dxy;
|
---|
143 | candidate->SDxy = ddxy;
|
---|
144 |
|
---|
145 | candidate->AddCandidate(mother);
|
---|
146 | fOutputArray->Add(candidate);
|
---|
147 | }
|
---|
148 | }
|
---|
149 |
|
---|
150 | //------------------------------------------------------------------------------
|
---|