Fork me on GitHub

source: git/modules/MomentumSmearing.cc@ e5ea42e

Last change on this file since e5ea42e was 4827699, checked in by Michele Selvaggi <michele.selvaggi@…>, 8 years ago

store correct TrackResolution

  • Property mode set to 100644
File size: 4.1 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
[d7d2da3]19
20/** \class MomentumSmearing
21 *
22 * Performs transverse momentum resolution smearing.
23 *
24 * \author P. Demin - UCL, Louvain-la-Neuve
25 *
26 */
27
28#include "modules/MomentumSmearing.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
55MomentumSmearing::MomentumSmearing() :
56 fFormula(0), fItInputArray(0)
57{
58 fFormula = new DelphesFormula;
59}
60
61//------------------------------------------------------------------------------
62
63MomentumSmearing::~MomentumSmearing()
64{
65 if(fFormula) delete fFormula;
66}
67
68//------------------------------------------------------------------------------
69
70void MomentumSmearing::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", "ParticlePropagator/stableParticles"));
79 fItInputArray = fInputArray->MakeIterator();
80
81 // create output array
82
83 fOutputArray = ExportArray(GetString("OutputArray", "stableParticles"));
84}
85
86//------------------------------------------------------------------------------
87
88void MomentumSmearing::Finish()
89{
90 if(fItInputArray) delete fItInputArray;
91}
92
93//------------------------------------------------------------------------------
94
95void MomentumSmearing::Process()
96{
97 Candidate *candidate, *mother;
[79d4123]98 Double_t pt, eta, phi, e, res;
[d7d2da3]99
100 fItInputArray->Reset();
101 while((candidate = static_cast<Candidate*>(fItInputArray->Next())))
102 {
103 const TLorentzVector &candidatePosition = candidate->Position;
104 const TLorentzVector &candidateMomentum = candidate->Momentum;
105 eta = candidatePosition.Eta();
106 phi = candidatePosition.Phi();
107 pt = candidateMomentum.Pt();
[95aa610]108 e = candidateMomentum.E();
[79d4123]109 res = fFormula->Eval(pt, eta, phi, e);
110
[d7d2da3]111 // apply smearing formula
[bc58cf5]112 //pt = gRandom->Gaus(pt, fFormula->Eval(pt, eta, phi, e) * pt);
[79d4123]113
114 res = ( res > 1.0 ) ? 1.0 : res;
115
116 pt = LogNormal(pt, res * pt );
[d7d2da3]117
[bc58cf5]118 //if(pt <= 0.0) continue;
[d7d2da3]119
120 mother = candidate;
121 candidate = static_cast<Candidate*>(candidate->Clone());
122 eta = candidateMomentum.Eta();
123 phi = candidateMomentum.Phi();
124 candidate->Momentum.SetPtEtaPhiE(pt, eta, phi, pt*TMath::CosH(eta));
[4827699]125 //candidate->TrackResolution = fFormula->Eval(pt, eta, phi, e);
126 candidate->TrackResolution = res;
[d7d2da3]127 candidate->AddCandidate(mother);
128
129 fOutputArray->Add(candidate);
130 }
131}
[bc58cf5]132//----------------------------------------------------------------
133
134Double_t MomentumSmearing::LogNormal(Double_t mean, Double_t sigma)
135{
136 Double_t a, b;
137
138 if(mean > 0.0)
139 {
140 b = TMath::Sqrt(TMath::Log((1.0 + (sigma*sigma)/(mean*mean))));
141 a = TMath::Log(mean) - 0.5*b*b;
142
143 return TMath::Exp(a + b*gRandom->Gaus(0.0, 1.0));
144 }
145 else
146 {
147 return 0.0;
148 }
149}
150
[d7d2da3]151
152//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.