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 Isolation
|
---|
21 | *
|
---|
22 | * Sums transverse momenta of isolation objects (tracks, calorimeter towers, etc)
|
---|
23 | * within a DeltaR cone around a candidate and calculates fraction of this sum
|
---|
24 | * to the candidate's transverse momentum. outputs candidates that have
|
---|
25 | * the transverse momenta fraction within (PTRatioMin, PTRatioMax].
|
---|
26 | *
|
---|
27 | * \author P. Demin, M. Selvaggi, R. Gerosa - UCL, Louvain-la-Neuve
|
---|
28 | *
|
---|
29 | */
|
---|
30 |
|
---|
31 | #include "modules/Isolation.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 |
|
---|
54 | using namespace std;
|
---|
55 |
|
---|
56 | //------------------------------------------------------------------------------
|
---|
57 |
|
---|
58 | class IsolationClassifier : public ExRootClassifier
|
---|
59 | {
|
---|
60 | public:
|
---|
61 |
|
---|
62 | IsolationClassifier() {}
|
---|
63 |
|
---|
64 | Int_t GetCategory(TObject *object);
|
---|
65 |
|
---|
66 | Double_t fPTMin;
|
---|
67 | };
|
---|
68 |
|
---|
69 | //------------------------------------------------------------------------------
|
---|
70 |
|
---|
71 | Int_t IsolationClassifier::GetCategory(TObject *object)
|
---|
72 | {
|
---|
73 | Candidate *track = static_cast<Candidate*>(object);
|
---|
74 | const TLorentzVector &momentum = track->Momentum;
|
---|
75 |
|
---|
76 | if(momentum.Pt() < fPTMin) return -1;
|
---|
77 |
|
---|
78 | return 0;
|
---|
79 | }
|
---|
80 |
|
---|
81 | //------------------------------------------------------------------------------
|
---|
82 |
|
---|
83 | Isolation::Isolation() :
|
---|
84 | fClassifier(0), fFilter(0),
|
---|
85 | fItIsolationInputArray(0), fItCandidateInputArray(0),
|
---|
86 | fItRhoInputArray(0)
|
---|
87 | {
|
---|
88 | fClassifier = new IsolationClassifier;
|
---|
89 | }
|
---|
90 |
|
---|
91 | //------------------------------------------------------------------------------
|
---|
92 |
|
---|
93 | Isolation::~Isolation()
|
---|
94 | {
|
---|
95 | }
|
---|
96 |
|
---|
97 | //------------------------------------------------------------------------------
|
---|
98 |
|
---|
99 | void Isolation::Init()
|
---|
100 | {
|
---|
101 | const char *rhoInputArrayName;
|
---|
102 |
|
---|
103 | fDeltaRMax = GetDouble("DeltaRMax", 0.5);
|
---|
104 |
|
---|
105 | fPTRatioMax = GetDouble("PTRatioMax", 0.1);
|
---|
106 |
|
---|
107 | fPTSumMax = GetDouble("PTSumMax", 5.0);
|
---|
108 |
|
---|
109 | fUsePTSum = GetBool("UsePTSum", false);
|
---|
110 |
|
---|
111 | fVetoLeptons = GetBool("VetoLeptons", true);
|
---|
112 |
|
---|
113 | fClassifier->fPTMin = GetDouble("PTMin", 0.5);
|
---|
114 |
|
---|
115 |
|
---|
116 | // import input array(s)
|
---|
117 |
|
---|
118 | fIsolationInputArray = ImportArray(GetString("IsolationInputArray", "Delphes/partons"));
|
---|
119 | fItIsolationInputArray = fIsolationInputArray->MakeIterator();
|
---|
120 |
|
---|
121 | fFilter = new ExRootFilter(fIsolationInputArray);
|
---|
122 |
|
---|
123 | fCandidateInputArray = ImportArray(GetString("CandidateInputArray", "Calorimeter/electrons"));
|
---|
124 | fItCandidateInputArray = fCandidateInputArray->MakeIterator();
|
---|
125 |
|
---|
126 | rhoInputArrayName = GetString("RhoInputArray", "");
|
---|
127 | if(rhoInputArrayName[0] != '\0')
|
---|
128 | {
|
---|
129 | fRhoInputArray = ImportArray(rhoInputArrayName);
|
---|
130 | fItRhoInputArray = fRhoInputArray->MakeIterator();
|
---|
131 | }
|
---|
132 | else
|
---|
133 | {
|
---|
134 | fRhoInputArray = 0;
|
---|
135 | }
|
---|
136 |
|
---|
137 | // create output array
|
---|
138 |
|
---|
139 | fOutputArray = ExportArray(GetString("OutputArray", "electrons"));
|
---|
140 | }
|
---|
141 |
|
---|
142 | //------------------------------------------------------------------------------
|
---|
143 |
|
---|
144 | void Isolation::Finish()
|
---|
145 | {
|
---|
146 | if(fItRhoInputArray) delete fItRhoInputArray;
|
---|
147 | if(fFilter) delete fFilter;
|
---|
148 | if(fItCandidateInputArray) delete fItCandidateInputArray;
|
---|
149 | if(fItIsolationInputArray) delete fItIsolationInputArray;
|
---|
150 | }
|
---|
151 |
|
---|
152 | //------------------------------------------------------------------------------
|
---|
153 |
|
---|
154 | void Isolation::Process()
|
---|
155 | {
|
---|
156 | Candidate *candidate, *isolation, *object;
|
---|
157 | TObjArray *isolationArray;
|
---|
158 | Double_t sumCharged, sumNeutral, sumAllParticles, sumChargedPU, sumDBeta, ratioDBeta, sumRhoCorr, ratioRhoCorr;
|
---|
159 | Int_t counter;
|
---|
160 | Double_t eta = 0.0;
|
---|
161 | Double_t rho = 0.0;
|
---|
162 |
|
---|
163 | if(fRhoInputArray && fRhoInputArray->GetEntriesFast() > 0)
|
---|
164 | {
|
---|
165 | candidate = static_cast<Candidate*>(fRhoInputArray->At(0));
|
---|
166 | rho = candidate->Momentum.Pt();
|
---|
167 | }
|
---|
168 |
|
---|
169 | // select isolation objects
|
---|
170 | fFilter->Reset();
|
---|
171 | isolationArray = fFilter->GetSubArray(fClassifier, 0);
|
---|
172 |
|
---|
173 | if(isolationArray == 0) return;
|
---|
174 |
|
---|
175 | TIter itIsolationArray(isolationArray);
|
---|
176 |
|
---|
177 | // loop over all input jets
|
---|
178 | fItCandidateInputArray->Reset();
|
---|
179 | while((candidate = static_cast<Candidate*>(fItCandidateInputArray->Next())))
|
---|
180 | {
|
---|
181 | const TLorentzVector &candidateMomentum = candidate->Momentum;
|
---|
182 | eta = TMath::Abs(candidateMomentum.Eta());
|
---|
183 |
|
---|
184 | // loop over all input tracks
|
---|
185 |
|
---|
186 | sumNeutral = 0.0;
|
---|
187 | sumCharged = 0.0;
|
---|
188 | sumChargedPU = 0.0;
|
---|
189 | sumAllParticles = 0.0;
|
---|
190 |
|
---|
191 | counter = 0;
|
---|
192 | itIsolationArray.Reset();
|
---|
193 |
|
---|
194 | while((isolation = static_cast<Candidate*>(itIsolationArray.Next())))
|
---|
195 | {
|
---|
196 | const TLorentzVector &isolationMomentum = isolation->Momentum;
|
---|
197 |
|
---|
198 | if(candidateMomentum.DeltaR(isolationMomentum) <= fDeltaRMax &&
|
---|
199 | candidate->GetUniqueID() != isolation->GetUniqueID() &&
|
---|
200 | ( !fVetoLeptons || (TMath::Abs(candidate->PID) != 11 && (TMath::Abs(candidate->PID) != 13)) ) )
|
---|
201 | {
|
---|
202 |
|
---|
203 | sumAllParticles += isolationMomentum.Pt();
|
---|
204 | if(isolation->Charge !=0)
|
---|
205 | {
|
---|
206 | sumCharged += isolationMomentum.Pt();
|
---|
207 | if(isolation->IsRecoPU != 0) sumChargedPU += isolationMomentum.Pt();
|
---|
208 | }
|
---|
209 |
|
---|
210 | else sumNeutral += isolationMomentum.Pt();
|
---|
211 |
|
---|
212 | ++counter;
|
---|
213 | }
|
---|
214 | }
|
---|
215 |
|
---|
216 | // find rho
|
---|
217 | rho = 0.0;
|
---|
218 | if(fRhoInputArray)
|
---|
219 | {
|
---|
220 | fItRhoInputArray->Reset();
|
---|
221 | while((object = static_cast<Candidate*>(fItRhoInputArray->Next())))
|
---|
222 | {
|
---|
223 | if(eta >= object->Edges[0] && eta < object->Edges[1])
|
---|
224 | {
|
---|
225 | rho = object->Momentum.Pt();
|
---|
226 | }
|
---|
227 | }
|
---|
228 | }
|
---|
229 |
|
---|
230 |
|
---|
231 | // correct sum for pile-up contamination
|
---|
232 | sumDBeta = sumCharged + TMath::Max(sumNeutral-0.5*sumChargedPU,0.0);
|
---|
233 | sumRhoCorr = sumCharged + TMath::Max(sumNeutral-TMath::Max(rho,0.0)*fDeltaRMax*fDeltaRMax*TMath::Pi(),0.0);
|
---|
234 | ratioDBeta = sumDBeta/candidateMomentum.Pt();
|
---|
235 | ratioRhoCorr = sumRhoCorr/candidateMomentum.Pt();
|
---|
236 |
|
---|
237 | candidate->IsolationVar = ratioDBeta;
|
---|
238 | candidate->IsolationVarRhoCorr = ratioRhoCorr;
|
---|
239 | candidate->SumPtCharged = sumCharged;
|
---|
240 | candidate->SumPtNeutral = sumNeutral;
|
---|
241 | candidate->SumPtChargedPU = sumChargedPU;
|
---|
242 | candidate->SumPt = sumAllParticles;
|
---|
243 |
|
---|
244 | if((fUsePTSum && sumDBeta > fPTSumMax) || ratioDBeta > fPTRatioMax) continue;
|
---|
245 | fOutputArray->Add(candidate);
|
---|
246 | }
|
---|
247 | }
|
---|
248 |
|
---|
249 | //------------------------------------------------------------------------------
|
---|