Fork me on GitHub

source: git/modules/Isolation.cc@ 9c491e1

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 9c491e1 was df35033, checked in by Pavel Demin <pavel.demin@…>, 10 years ago

replace call to Overlaps with UniqueID comparison (close #362)

  • Property mode set to 100644
File size: 5.9 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 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 - 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
54using namespace std;
55
56//------------------------------------------------------------------------------
57
58class IsolationClassifier : public ExRootClassifier
59{
60public:
61
62 IsolationClassifier() {}
63
64 Int_t GetCategory(TObject *object);
65
66 Double_t fPTMin;
67};
68
69//------------------------------------------------------------------------------
70
71Int_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
83Isolation::Isolation() :
84 fClassifier(0), fFilter(0),
[e9023b9]85 fItIsolationInputArray(0), fItCandidateInputArray(0),
86 fItRhoInputArray(0)
[d7d2da3]87{
88 fClassifier = new IsolationClassifier;
89}
90
91//------------------------------------------------------------------------------
92
93Isolation::~Isolation()
94{
95}
96
97//------------------------------------------------------------------------------
98
99void Isolation::Init()
100{
101 const char *rhoInputArrayName;
102
103 fDeltaRMax = GetDouble("DeltaRMax", 0.5);
104
105 fPTRatioMax = GetDouble("PTRatioMax", 0.1);
106
[b286067]107 fPTSumMax = GetDouble("PTSumMax", 5.0);
108
109 fUsePTSum = GetBool("UsePTSum", false);
110
[d7d2da3]111 fClassifier->fPTMin = GetDouble("PTMin", 0.5);
112
113 // import input array(s)
114
115 fIsolationInputArray = ImportArray(GetString("IsolationInputArray", "Delphes/partons"));
116 fItIsolationInputArray = fIsolationInputArray->MakeIterator();
117
118 fFilter = new ExRootFilter(fIsolationInputArray);
119
120 fCandidateInputArray = ImportArray(GetString("CandidateInputArray", "Calorimeter/electrons"));
121 fItCandidateInputArray = fCandidateInputArray->MakeIterator();
122
123 rhoInputArrayName = GetString("RhoInputArray", "");
124 if(rhoInputArrayName[0] != '\0')
125 {
126 fRhoInputArray = ImportArray(rhoInputArrayName);
[e9023b9]127 fItRhoInputArray = fRhoInputArray->MakeIterator();
[d7d2da3]128 }
129 else
130 {
131 fRhoInputArray = 0;
132 }
133
134 // create output array
135
136 fOutputArray = ExportArray(GetString("OutputArray", "electrons"));
137}
138
139//------------------------------------------------------------------------------
140
141void Isolation::Finish()
142{
[e9023b9]143 if(fItRhoInputArray) delete fItRhoInputArray;
[d7d2da3]144 if(fFilter) delete fFilter;
145 if(fItCandidateInputArray) delete fItCandidateInputArray;
146 if(fItIsolationInputArray) delete fItIsolationInputArray;
147}
148
149//------------------------------------------------------------------------------
150
151void Isolation::Process()
152{
[e9023b9]153 Candidate *candidate, *isolation, *object;
[d7d2da3]154 TObjArray *isolationArray;
[b286067]155 Double_t sum, ratio;
[d7d2da3]156 Int_t counter;
[e9023b9]157 Double_t eta = 0.0;
[d7d2da3]158 Double_t rho = 0.0;
159
160 if(fRhoInputArray && fRhoInputArray->GetEntriesFast() > 0)
161 {
162 candidate = static_cast<Candidate*>(fRhoInputArray->At(0));
163 rho = candidate->Momentum.Pt();
164 }
165
166 // select isolation objects
167 fFilter->Reset();
168 isolationArray = fFilter->GetSubArray(fClassifier, 0);
169
170 if(isolationArray == 0) return;
171
172 TIter itIsolationArray(isolationArray);
173
174 // loop over all input jets
175 fItCandidateInputArray->Reset();
176 while((candidate = static_cast<Candidate*>(fItCandidateInputArray->Next())))
177 {
178 const TLorentzVector &candidateMomentum = candidate->Momentum;
[e9023b9]179 eta = TMath::Abs(candidateMomentum.Eta());
[d7d2da3]180
181 // loop over all input tracks
[b286067]182 sum = 0.0;
[d7d2da3]183 counter = 0;
184 itIsolationArray.Reset();
185 while((isolation = static_cast<Candidate*>(itIsolationArray.Next())))
186 {
187 const TLorentzVector &isolationMomentum = isolation->Momentum;
188
189 if(candidateMomentum.DeltaR(isolationMomentum) <= fDeltaRMax &&
[df35033]190 candidate->GetUniqueID() != isolation->GetUniqueID())
[d7d2da3]191 {
[b286067]192 sum += isolationMomentum.Pt();
[d7d2da3]193 ++counter;
194 }
195 }
196
[e9023b9]197 // find rho
198 rho = 0.0;
199 if(fRhoInputArray)
200 {
201 fItRhoInputArray->Reset();
202 while((object = static_cast<Candidate*>(fItRhoInputArray->Next())))
203 {
204 if(eta >= object->Edges[0] && eta < object->Edges[1])
205 {
206 rho = object->Momentum.Pt();
207 }
208 }
209 }
210
[b286067]211 // correct sum for pile-up contamination
[e9023b9]212 sum = sum - rho*fDeltaRMax*fDeltaRMax*TMath::Pi();
[d7d2da3]213
[b286067]214 ratio = sum/candidateMomentum.Pt();
215 if((fUsePTSum && sum > fPTSumMax) || ratio > fPTRatioMax) continue;
[d7d2da3]216
217 fOutputArray->Add(candidate);
218 }
219}
220
221//------------------------------------------------------------------------------
Note: See TracBrowser for help on using the repository browser.