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 |
|
---|
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 |
|
---|
51 | using namespace std;
|
---|
52 |
|
---|
53 | //------------------------------------------------------------------------------
|
---|
54 |
|
---|
55 | MomentumSmearing::MomentumSmearing() :
|
---|
56 | fFormula(0), fItInputArray(0)
|
---|
57 | {
|
---|
58 | fFormula = new DelphesFormula;
|
---|
59 | }
|
---|
60 |
|
---|
61 | //------------------------------------------------------------------------------
|
---|
62 |
|
---|
63 | MomentumSmearing::~MomentumSmearing()
|
---|
64 | {
|
---|
65 | if(fFormula) delete fFormula;
|
---|
66 | }
|
---|
67 |
|
---|
68 | //------------------------------------------------------------------------------
|
---|
69 |
|
---|
70 | void 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 |
|
---|
88 | void MomentumSmearing::Finish()
|
---|
89 | {
|
---|
90 | if(fItInputArray) delete fItInputArray;
|
---|
91 | }
|
---|
92 |
|
---|
93 | //------------------------------------------------------------------------------
|
---|
94 |
|
---|
95 | void MomentumSmearing::Process()
|
---|
96 | {
|
---|
97 | Candidate *candidate, *mother;
|
---|
98 | Double_t pt, eta, phi, e;
|
---|
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();
|
---|
108 | e = candidateMomentum.E();
|
---|
109 |
|
---|
110 | // apply smearing formula
|
---|
111 | //pt = gRandom->Gaus(pt, fFormula->Eval(pt, eta, phi, e) * pt);
|
---|
112 | pt = LogNormal(pt, fFormula->Eval(pt, eta, phi, e) * pt );
|
---|
113 |
|
---|
114 | //if(pt <= 0.0) continue;
|
---|
115 |
|
---|
116 | mother = candidate;
|
---|
117 | candidate = static_cast<Candidate*>(candidate->Clone());
|
---|
118 | eta = candidateMomentum.Eta();
|
---|
119 | phi = candidateMomentum.Phi();
|
---|
120 | candidate->Momentum.SetPtEtaPhiE(pt, eta, phi, pt*TMath::CosH(eta));
|
---|
121 | candidate->TrackResolution = fFormula->Eval(pt, eta, phi, e);
|
---|
122 | candidate->AddCandidate(mother);
|
---|
123 |
|
---|
124 | fOutputArray->Add(candidate);
|
---|
125 | }
|
---|
126 | }
|
---|
127 | //----------------------------------------------------------------
|
---|
128 |
|
---|
129 | Double_t MomentumSmearing::LogNormal(Double_t mean, Double_t sigma)
|
---|
130 | {
|
---|
131 | Double_t a, b;
|
---|
132 |
|
---|
133 | if(mean > 0.0)
|
---|
134 | {
|
---|
135 | b = TMath::Sqrt(TMath::Log((1.0 + (sigma*sigma)/(mean*mean))));
|
---|
136 | a = TMath::Log(mean) - 0.5*b*b;
|
---|
137 |
|
---|
138 | return TMath::Exp(a + b*gRandom->Gaus(0.0, 1.0));
|
---|
139 | }
|
---|
140 | else
|
---|
141 | {
|
---|
142 | return 0.0;
|
---|
143 | }
|
---|
144 | }
|
---|
145 |
|
---|
146 |
|
---|
147 | //------------------------------------------------------------------------------
|
---|