Fork me on GitHub

source: git/modules/ImpactParameterSmearing.cc@ 505a779

ImprovedOutputFile Timing
Last change on this file since 505a779 was 341014c, checked in by Pavel Demin <pavel-demin@…>, 6 years ago

apply .clang-format to all .h, .cc and .cpp files

  • Property mode set to 100644
File size: 4.2 KB
RevLine 
[b443089]1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
[1fa50c2]4 *
[b443089]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.
[1fa50c2]9 *
[b443089]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.
[1fa50c2]14 *
[b443089]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
[a0431dc]19/** \class ImpactParameterSmearing
20 *
21 * Performs transverse impact parameter smearing.
22 *
23 * \author M. Selvaggi - UCL, Louvain-la-Neuve
24 *
25 */
[4d753a1]26
[a0431dc]27#include "modules/ImpactParameterSmearing.h"
28
29#include "classes/DelphesClasses.h"
30#include "classes/DelphesFactory.h"
31#include "classes/DelphesFormula.h"
32
33#include "ExRootAnalysis/ExRootClassifier.h"
[341014c]34#include "ExRootAnalysis/ExRootFilter.h"
35#include "ExRootAnalysis/ExRootResult.h"
[a0431dc]36
37#include "TDatabasePDG.h"
[341014c]38#include "TFormula.h"
[a0431dc]39#include "TLorentzVector.h"
[341014c]40#include "TMath.h"
41#include "TObjArray.h"
42#include "TRandom3.h"
43#include "TString.h"
[a0431dc]44
[4d753a1]45#include <algorithm>
[a0431dc]46#include <iostream>
47#include <sstream>
[341014c]48#include <stdexcept>
[a0431dc]49
50using namespace std;
51
52//------------------------------------------------------------------------------
53
54ImpactParameterSmearing::ImpactParameterSmearing() :
55 fFormula(0), fItInputArray(0)
56{
57 fFormula = new DelphesFormula;
58}
59
60//------------------------------------------------------------------------------
61
62ImpactParameterSmearing::~ImpactParameterSmearing()
63{
64 if(fFormula) delete fFormula;
65}
66
67//------------------------------------------------------------------------------
68
69void ImpactParameterSmearing::Init()
70{
71 // read resolution formula
72
73 fFormula->Compile(GetString("ResolutionFormula", "0.0"));
74
75 // import input array
76
77 fInputArray = ImportArray(GetString("InputArray", "TrackMerger/tracks"));
78 fItInputArray = fInputArray->MakeIterator();
79
80 // create output array
81
82 fOutputArray = ExportArray(GetString("OutputArray", "tracks"));
83}
84
85//------------------------------------------------------------------------------
86
87void ImpactParameterSmearing::Finish()
88{
89 if(fItInputArray) delete fItInputArray;
90}
91
92//------------------------------------------------------------------------------
93
94void ImpactParameterSmearing::Process()
95{
96 Candidate *candidate, *particle, *mother;
[632da30]97 Double_t xd, yd, zd, d0, sx, sy, sz, dd0;
[95aa610]98 Double_t pt, eta, px, py, phi, e;
[a0431dc]99
100 fItInputArray->Reset();
[341014c]101 while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
[a0431dc]102 {
[4d753a1]103
[632da30]104 // take momentum before smearing (otherwise apply double smearing on d0)
[341014c]105 particle = static_cast<Candidate *>(candidate->GetCandidates()->At(0));
[4d753a1]106
[a0431dc]107 const TLorentzVector &candidateMomentum = particle->Momentum;
[4d753a1]108
[a0431dc]109 eta = candidateMomentum.Eta();
110 pt = candidateMomentum.Pt();
[95aa610]111 phi = candidateMomentum.Phi();
112 e = candidateMomentum.E();
[341014c]113
[a0431dc]114 px = candidateMomentum.Px();
115 py = candidateMomentum.Py();
[4d753a1]116
[a0431dc]117 // calculate coordinates of closest approach to track circle in transverse plane xd, yd, zd
[341014c]118 xd = candidate->Xd;
119 yd = candidate->Yd;
120 zd = candidate->Zd;
[4d753a1]121
122 // calculate smeared values
[95aa610]123 sx = gRandom->Gaus(0.0, fFormula->Eval(pt, eta, phi, e));
124 sy = gRandom->Gaus(0.0, fFormula->Eval(pt, eta, phi, e));
125 sz = gRandom->Gaus(0.0, fFormula->Eval(pt, eta, phi, e));
[4d753a1]126
[a0431dc]127 xd += sx;
128 yd += sy;
[4d753a1]129 zd += sz;
130
131 // calculate impact parameter (after-smearing)
[341014c]132 d0 = (xd * py - yd * px) / pt;
[4d753a1]133
[632da30]134 dd0 = gRandom->Gaus(0.0, fFormula->Eval(pt, eta, phi, e));
[4d753a1]135
136 // fill smeared values in candidate
[a0431dc]137 mother = candidate;
[4d753a1]138
[341014c]139 candidate = static_cast<Candidate *>(candidate->Clone());
[e4c3fef]140 candidate->Xd = xd;
141 candidate->Yd = yd;
142 candidate->Zd = zd;
[4d753a1]143
[632da30]144 candidate->D0 = d0;
145 candidate->ErrorD0 = dd0;
[4d753a1]146
[a0431dc]147 candidate->AddCandidate(mother);
148 fOutputArray->Add(candidate);
149 }
150}
151
152//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.