Fork me on GitHub

source: git/modules/ImpactParameterSmearing.cc@ 493aa9b

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 493aa9b was 7c0fcd5, checked in by Pavel Demin <demin@…>, 10 years ago

delete duplicate license file and prepend GPLv3 header to all source code files

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