Fork me on GitHub

source: git/modules/MomentumSmearing.cc@ 4e8e72b

Last change on this file since 4e8e72b was 27197df, checked in by Michele Selvaggi <michele.selvaggi@…>, 4 years ago

added UseMomentumVector option to momentum smearing

  • Property mode set to 100644
File size: 4.3 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 MomentumSmearing
20 *
21 * Performs transverse momentum resolution smearing.
22 *
23 * \author P. Demin - UCL, Louvain-la-Neuve
24 *
25 */
26
27#include "modules/MomentumSmearing.h"
28
29#include "classes/DelphesClasses.h"
30#include "classes/DelphesFactory.h"
31#include "classes/DelphesFormula.h"
32
33#include "ExRootAnalysis/ExRootClassifier.h"
34#include "ExRootAnalysis/ExRootFilter.h"
35#include "ExRootAnalysis/ExRootResult.h"
36
37#include "TDatabasePDG.h"
38#include "TFormula.h"
39#include "TLorentzVector.h"
40#include "TMath.h"
41#include "TObjArray.h"
42#include "TRandom3.h"
43#include "TString.h"
44
45#include <algorithm>
46#include <iostream>
47#include <sstream>
48#include <stdexcept>
49
50using namespace std;
51
52//------------------------------------------------------------------------------
53
54MomentumSmearing::MomentumSmearing() :
55 fFormula(0), fItInputArray(0)
56{
57 fFormula = new DelphesFormula;
58}
59
60//------------------------------------------------------------------------------
61
62MomentumSmearing::~MomentumSmearing()
63{
64 if(fFormula) delete fFormula;
65}
66
67//------------------------------------------------------------------------------
68
69void MomentumSmearing::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", "ParticlePropagator/stableParticles"));
78 fItInputArray = fInputArray->MakeIterator();
79
80 // switch to compute momentum smearing based on momentum vector eta, phi
81 fUseMomentumVector = GetBool("UseMomentumVector", false);
82
83 // create output array
84
85 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
86}
87
88//------------------------------------------------------------------------------
89
90void MomentumSmearing::Finish()
91{
92 if(fItInputArray) delete fItInputArray;
93}
94
95//------------------------------------------------------------------------------
96
97void MomentumSmearing::Process()
98{
99 Candidate *candidate, *mother;
100 Double_t pt, eta, phi, e, m, res;
101
102 fItInputArray->Reset();
103 while((candidate = static_cast<Candidate *>(fItInputArray->Next())))
104 {
105 const TLorentzVector &candidatePosition = candidate->Position;
106 const TLorentzVector &candidateMomentum = candidate->Momentum;
107 eta = candidatePosition.Eta();
108 phi = candidatePosition.Phi();
109
110 if (fUseMomentumVector){
111 eta = candidateMomentum.Eta();
112 phi = candidateMomentum.Phi();
113 }
114
115 pt = candidateMomentum.Pt();
116 e = candidateMomentum.E();
117 m = candidateMomentum.M();
118 res = fFormula->Eval(pt, eta, phi, e, candidate);
119
120 // apply smearing formula
121 //pt = gRandom->Gaus(pt, fFormula->Eval(pt, eta, phi, e) * pt);
122
123 res = (res > 1.0) ? 1.0 : res;
124
125 pt = LogNormal(pt, res * pt);
126
127 //if(pt <= 0.0) continue;
128
129 mother = candidate;
130 candidate = static_cast<Candidate *>(candidate->Clone());
131 eta = candidateMomentum.Eta();
132 phi = candidateMomentum.Phi();
133 candidate->Momentum.SetPtEtaPhiM(pt, eta, phi, m);
134 //candidate->TrackResolution = fFormula->Eval(pt, eta, phi, e);
135 candidate->TrackResolution = res;
136 candidate->AddCandidate(mother);
137
138 fOutputArray->Add(candidate);
139 }
140}
141//----------------------------------------------------------------
142
143Double_t MomentumSmearing::LogNormal(Double_t mean, Double_t sigma)
144{
145 Double_t a, b;
146
147 if(mean > 0.0)
148 {
149 b = TMath::Sqrt(TMath::Log((1.0 + (sigma * sigma) / (mean * mean))));
150 a = TMath::Log(mean) - 0.5 * b * b;
151
152 return TMath::Exp(a + b * gRandom->Gaus(0.0, 1.0));
153 }
154 else
155 {
156 return 0.0;
157 }
158}
159
160//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.