[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 | *
|
---|
[7f9ae0a] | 27 | * \author P. Demin, M. Selvaggi, R. Gerosa - UCL, Louvain-la-Neuve
|
---|
[d7d2da3] | 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),
|
---|
[e9023b9] | 85 | fItIsolationInputArray(0), fItCandidateInputArray(0),
|
---|
| 86 | fItRhoInputArray(0)
|
---|
[d7d2da3] | 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 |
|
---|
[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 |
|
---|
| 141 | void 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 |
|
---|
| 151 | void Isolation::Process()
|
---|
| 152 | {
|
---|
[e9023b9] | 153 | Candidate *candidate, *isolation, *object;
|
---|
[d7d2da3] | 154 | TObjArray *isolationArray;
|
---|
[7f9ae0a] | 155 | Double_t sumCharged, sumNeutral, sumAllParticles, sumChargedPU, sumDBeta, ratioDBeta, sumRhoCorr, ratioRhoCorr;
|
---|
[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
|
---|
[7f9ae0a] | 182 |
|
---|
[3db5282] | 183 | sumNeutral = 0.0;
|
---|
[7f9ae0a] | 184 | sumCharged = 0.0;
|
---|
| 185 | sumChargedPU = 0.0;
|
---|
| 186 | sumAllParticles = 0.0;
|
---|
| 187 |
|
---|
[d7d2da3] | 188 | counter = 0;
|
---|
| 189 | itIsolationArray.Reset();
|
---|
[7f9ae0a] | 190 |
|
---|
[d7d2da3] | 191 | while((isolation = static_cast<Candidate*>(itIsolationArray.Next())))
|
---|
| 192 | {
|
---|
| 193 | const TLorentzVector &isolationMomentum = isolation->Momentum;
|
---|
| 194 |
|
---|
| 195 | if(candidateMomentum.DeltaR(isolationMomentum) <= fDeltaRMax &&
|
---|
[e8abecd] | 196 | candidate->GetUniqueID() != isolation->GetUniqueID())
|
---|
[d7d2da3] | 197 | {
|
---|
[7f9ae0a] | 198 | sumAllParticles += isolationMomentum.Pt();
|
---|
| 199 | if(isolation->Charge !=0)
|
---|
| 200 | {
|
---|
| 201 | sumCharged += isolationMomentum.Pt();
|
---|
| 202 | if(isolation->IsRecoPU != 0) sumChargedPU += isolationMomentum.Pt();
|
---|
| 203 | }
|
---|
[e8abecd] | 204 | else
|
---|
| 205 | {
|
---|
| 206 | sumNeutral += isolationMomentum.Pt();
|
---|
| 207 | }
|
---|
[d7d2da3] | 208 | ++counter;
|
---|
| 209 | }
|
---|
| 210 | }
|
---|
| 211 |
|
---|
[e9023b9] | 212 | // find rho
|
---|
| 213 | rho = 0.0;
|
---|
| 214 | if(fRhoInputArray)
|
---|
| 215 | {
|
---|
| 216 | fItRhoInputArray->Reset();
|
---|
| 217 | while((object = static_cast<Candidate*>(fItRhoInputArray->Next())))
|
---|
| 218 | {
|
---|
| 219 | if(eta >= object->Edges[0] && eta < object->Edges[1])
|
---|
| 220 | {
|
---|
| 221 | rho = object->Momentum.Pt();
|
---|
| 222 | }
|
---|
| 223 | }
|
---|
| 224 | }
|
---|
| 225 |
|
---|
[7f9ae0a] | 226 | // correct sum for pile-up contamination
|
---|
| 227 | sumDBeta = sumCharged + TMath::Max(sumNeutral-0.5*sumChargedPU,0.0);
|
---|
| 228 | sumRhoCorr = sumCharged + TMath::Max(sumNeutral-TMath::Max(rho,0.0)*fDeltaRMax*fDeltaRMax*TMath::Pi(),0.0);
|
---|
| 229 | ratioDBeta = sumDBeta/candidateMomentum.Pt();
|
---|
| 230 | ratioRhoCorr = sumRhoCorr/candidateMomentum.Pt();
|
---|
| 231 |
|
---|
[b62c2da] | 232 | candidate->IsolationVar = ratioDBeta;
|
---|
[7f9ae0a] | 233 | candidate->IsolationVarRhoCorr = ratioRhoCorr;
|
---|
[b62c2da] | 234 | candidate->SumPtCharged = sumCharged;
|
---|
| 235 | candidate->SumPtNeutral = sumNeutral;
|
---|
| 236 | candidate->SumPtChargedPU = sumChargedPU;
|
---|
| 237 | candidate->SumPt = sumAllParticles;
|
---|
[7f9ae0a] | 238 |
|
---|
[e8abecd] | 239 | if((fUsePTSum && sumDBeta > fPTSumMax) || (!fUsePTSum && ratioDBeta > fPTRatioMax)) continue;
|
---|
[d7d2da3] | 240 | fOutputArray->Add(candidate);
|
---|
| 241 | }
|
---|
| 242 | }
|
---|
| 243 |
|
---|
| 244 | //------------------------------------------------------------------------------
|
---|