Fork me on GitHub

source: git/modules/ImpactParameterSmearing.cc@ 4170544

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

change dxy variable into d0

  • Property mode set to 100644
File size: 4.2 KB
Line 
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
51using namespace std;
52
53//------------------------------------------------------------------------------
54
55ImpactParameterSmearing::ImpactParameterSmearing() :
56 fFormula(0), fItInputArray(0)
57{
58 fFormula = new DelphesFormula;
59}
60
61//------------------------------------------------------------------------------
62
63ImpactParameterSmearing::~ImpactParameterSmearing()
64{
65 if(fFormula) delete fFormula;
66}
67
68//------------------------------------------------------------------------------
69
70void 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
88void ImpactParameterSmearing::Finish()
89{
90 if(fItInputArray) delete fItInputArray;
91}
92
93//------------------------------------------------------------------------------
94
95void ImpactParameterSmearing::Process()
96{
97 Candidate *candidate, *particle, *mother;
98 Double_t xd, yd, zd, d0, sx, sy, sz, dd0;
99 Double_t pt, eta, px, py, phi, e;
100
101 fItInputArray->Reset();
102 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
103 {
104
105 // take momentum before smearing (otherwise apply double smearing on d0)
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 phi = candidateMomentum.Phi();
113 e = candidateMomentum.E();
114
115 px = candidateMomentum.Px();
116 py = candidateMomentum.Py();
117
118 // calculate coordinates of closest approach to track circle in transverse plane xd, yd, zd
119 xd = candidate->Xd;
120 yd = candidate->Yd;
121 zd = candidate->Zd;
122
123 // calculate smeared values
124 sx = gRandom->Gaus(0.0, fFormula->Eval(pt, eta, phi, e));
125 sy = gRandom->Gaus(0.0, fFormula->Eval(pt, eta, phi, e));
126 sz = gRandom->Gaus(0.0, fFormula->Eval(pt, eta, phi, e));
127
128 xd += sx;
129 yd += sy;
130 zd += sz;
131
132 // calculate impact parameter (after-smearing)
133 d0 = (xd*py - yd*px)/pt;
134
135 dd0 = gRandom->Gaus(0.0, fFormula->Eval(pt, eta, phi, e));
136
137 // fill smeared values in candidate
138 mother = candidate;
139
140 candidate = static_cast<Candidate*>(candidate->Clone());
141 candidate->Xd = xd;
142 candidate->Yd = yd;
143 candidate->Zd = zd;
144
145 candidate->D0 = d0;
146 candidate->ErrorD0 = dd0;
147
148 candidate->AddCandidate(mother);
149 fOutputArray->Add(candidate);
150 }
151}
152
153//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.