[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 |
|
---|
[24d005f] | 19 | /** \class PileUpJetID
|
---|
| 20 | *
|
---|
| 21 | * CMS PileUp Jet ID Variables, based on http://cds.cern.ch/record/1581583
|
---|
| 22 | *
|
---|
| 23 | * \author S. Zenz, December 2013
|
---|
| 24 | *
|
---|
| 25 | */
|
---|
| 26 |
|
---|
| 27 | #include "modules/PileUpJetID.h"
|
---|
| 28 |
|
---|
| 29 | #include "classes/DelphesClasses.h"
|
---|
| 30 | #include "classes/DelphesFactory.h"
|
---|
| 31 | #include "classes/DelphesFormula.h"
|
---|
| 32 |
|
---|
| 33 | #include "ExRootAnalysis/ExRootResult.h"
|
---|
| 34 | #include "ExRootAnalysis/ExRootFilter.h"
|
---|
| 35 | #include "ExRootAnalysis/ExRootClassifier.h"
|
---|
| 36 |
|
---|
| 37 | #include "TMath.h"
|
---|
| 38 | #include "TString.h"
|
---|
| 39 | #include "TFormula.h"
|
---|
| 40 | #include "TRandom3.h"
|
---|
| 41 | #include "TObjArray.h"
|
---|
| 42 | #include "TDatabasePDG.h"
|
---|
| 43 | #include "TLorentzVector.h"
|
---|
| 44 |
|
---|
| 45 | #include <algorithm>
|
---|
| 46 | #include <stdexcept>
|
---|
| 47 | #include <iostream>
|
---|
| 48 | #include <sstream>
|
---|
| 49 |
|
---|
| 50 | using namespace std;
|
---|
| 51 |
|
---|
| 52 | //------------------------------------------------------------------------------
|
---|
| 53 |
|
---|
| 54 | PileUpJetID::PileUpJetID() :
|
---|
[54b6dfc] | 55 | fItJetInputArray(0),fTrackInputArray(0),fNeutralInputArray(0),fItVertexInputArray(0)
|
---|
[24d005f] | 56 | {
|
---|
| 57 |
|
---|
| 58 | }
|
---|
| 59 |
|
---|
| 60 | //------------------------------------------------------------------------------
|
---|
| 61 |
|
---|
| 62 | PileUpJetID::~PileUpJetID()
|
---|
| 63 | {
|
---|
| 64 |
|
---|
| 65 | }
|
---|
| 66 |
|
---|
| 67 | //------------------------------------------------------------------------------
|
---|
| 68 |
|
---|
| 69 | void PileUpJetID::Init()
|
---|
| 70 | {
|
---|
| 71 | fJetPTMin = GetDouble("JetPTMin", 20.0);
|
---|
| 72 | fParameterR = GetDouble("ParameterR", 0.5);
|
---|
| 73 | fUseConstituents = GetInt("UseConstituents", 0);
|
---|
| 74 |
|
---|
| 75 | fAverageEachTower = false; // for timing
|
---|
| 76 |
|
---|
| 77 | // import input array(s)
|
---|
| 78 |
|
---|
| 79 | fJetInputArray = ImportArray(GetString("JetInputArray", "FastJetFinder/jets"));
|
---|
| 80 | fItJetInputArray = fJetInputArray->MakeIterator();
|
---|
| 81 |
|
---|
| 82 | fTrackInputArray = ImportArray(GetString("TrackInputArray", "Calorimeter/eflowTracks"));
|
---|
| 83 | fItTrackInputArray = fTrackInputArray->MakeIterator();
|
---|
| 84 |
|
---|
| 85 | fNeutralInputArray = ImportArray(GetString("NeutralInputArray", "Calorimeter/eflowTowers"));
|
---|
| 86 | fItNeutralInputArray = fNeutralInputArray->MakeIterator();
|
---|
[54b6dfc] | 87 |
|
---|
| 88 | fVertexInputArray = ImportArray(GetString("VertexInputArray", "PileUpMerger/vertices"));
|
---|
| 89 | fItVertexInputArray = fVertexInputArray->MakeIterator();
|
---|
| 90 |
|
---|
| 91 | fZVertexResolution = GetDouble("ZVertexResolution", 0.005)*1.0E3;
|
---|
| 92 | // create output array(s)
|
---|
[24d005f] | 93 |
|
---|
| 94 | fOutputArray = ExportArray(GetString("OutputArray", "jets"));
|
---|
| 95 |
|
---|
| 96 | }
|
---|
| 97 |
|
---|
| 98 | //------------------------------------------------------------------------------
|
---|
| 99 |
|
---|
| 100 | void PileUpJetID::Finish()
|
---|
| 101 | {
|
---|
| 102 |
|
---|
| 103 | if(fItJetInputArray) delete fItJetInputArray;
|
---|
| 104 | if(fItTrackInputArray) delete fItTrackInputArray;
|
---|
| 105 | if(fItNeutralInputArray) delete fItNeutralInputArray;
|
---|
[54b6dfc] | 106 | if(fItVertexInputArray) delete fItVertexInputArray;
|
---|
[24d005f] | 107 |
|
---|
| 108 | }
|
---|
| 109 |
|
---|
| 110 | //------------------------------------------------------------------------------
|
---|
| 111 |
|
---|
| 112 | void PileUpJetID::Process()
|
---|
| 113 | {
|
---|
| 114 | Candidate *candidate, *constituent;
|
---|
| 115 | TLorentzVector momentum, area;
|
---|
[54b6dfc] | 116 | Double_t zvtx=0;
|
---|
[24d005f] | 117 |
|
---|
| 118 | Candidate *trk;
|
---|
| 119 |
|
---|
[54b6dfc] | 120 | // find z position of primary vertex
|
---|
| 121 |
|
---|
| 122 | fItVertexInputArray->Reset();
|
---|
| 123 | while((candidate = static_cast<Candidate*>(fItVertexInputArray->Next())))
|
---|
| 124 | {
|
---|
| 125 | if(!candidate->IsPU)
|
---|
| 126 | {
|
---|
| 127 | zvtx = candidate->Position.Z();
|
---|
| 128 | break;
|
---|
| 129 | }
|
---|
| 130 | }
|
---|
| 131 |
|
---|
[24d005f] | 132 | // loop over all input candidates
|
---|
| 133 | fItJetInputArray->Reset();
|
---|
| 134 | while((candidate = static_cast<Candidate*>(fItJetInputArray->Next())))
|
---|
| 135 | {
|
---|
| 136 | momentum = candidate->Momentum;
|
---|
| 137 | area = candidate->Area;
|
---|
| 138 |
|
---|
| 139 | float sumpt = 0.;
|
---|
| 140 | float sumptch = 0.;
|
---|
| 141 | float sumptchpv = 0.;
|
---|
| 142 | float sumptchpu = 0.;
|
---|
| 143 | float sumdrsqptsq = 0.;
|
---|
| 144 | float sumptsq = 0.;
|
---|
| 145 | int nc = 0;
|
---|
| 146 | int nn = 0;
|
---|
| 147 | float pt_ann[5];
|
---|
| 148 |
|
---|
| 149 | for (int i = 0 ; i < 5 ; i++) {
|
---|
| 150 | pt_ann[i] = 0.;
|
---|
| 151 | }
|
---|
| 152 |
|
---|
| 153 | if (fUseConstituents) {
|
---|
| 154 | TIter itConstituents(candidate->GetCandidates());
|
---|
| 155 | while((constituent = static_cast<Candidate*>(itConstituents.Next()))) {
|
---|
| 156 | float pt = constituent->Momentum.Pt();
|
---|
| 157 | float dr = candidate->Momentum.DeltaR(constituent->Momentum);
|
---|
| 158 | sumpt += pt;
|
---|
| 159 | sumdrsqptsq += dr*dr*pt*pt;
|
---|
| 160 | sumptsq += pt*pt;
|
---|
| 161 | if (constituent->Charge == 0) {
|
---|
| 162 | // neutrals
|
---|
| 163 | nn++;
|
---|
| 164 | } else {
|
---|
| 165 | // charged
|
---|
[54b6dfc] | 166 | if (constituent->IsPU && TMath::Abs(constituent->Position.Z()-zvtx) > fZVertexResolution) {
|
---|
[24d005f] | 167 | sumptchpu += pt;
|
---|
| 168 | } else {
|
---|
| 169 | sumptchpv += pt;
|
---|
| 170 | }
|
---|
| 171 | sumptch += pt;
|
---|
| 172 | nc++;
|
---|
| 173 | }
|
---|
| 174 | for (int i = 0 ; i < 5 ; i++) {
|
---|
| 175 | if (dr > 0.1*i && dr < 0.1*(i+1)) {
|
---|
| 176 | pt_ann[i] += pt;
|
---|
| 177 | }
|
---|
| 178 | }
|
---|
| 179 | }
|
---|
| 180 | } else {
|
---|
| 181 | // Not using constituents, using dr
|
---|
| 182 | fItTrackInputArray->Reset();
|
---|
[54b6dfc] | 183 | while ((trk = static_cast<Candidate*>(fItTrackInputArray->Next()))) {
|
---|
[24d005f] | 184 | if (trk->Momentum.DeltaR(candidate->Momentum) < fParameterR) {
|
---|
| 185 | float pt = trk->Momentum.Pt();
|
---|
| 186 | sumpt += pt;
|
---|
| 187 | sumptch += pt;
|
---|
[54b6dfc] | 188 | if (trk->IsPU && TMath::Abs(trk->Position.Z()-zvtx) > fZVertexResolution) {
|
---|
[24d005f] | 189 | sumptchpu += pt;
|
---|
| 190 | } else {
|
---|
| 191 | sumptchpv += pt;
|
---|
| 192 | }
|
---|
| 193 | float dr = candidate->Momentum.DeltaR(trk->Momentum);
|
---|
| 194 | sumdrsqptsq += dr*dr*pt*pt;
|
---|
| 195 | sumptsq += pt*pt;
|
---|
| 196 | nc++;
|
---|
| 197 | for (int i = 0 ; i < 5 ; i++) {
|
---|
| 198 | if (dr > 0.1*i && dr < 0.1*(i+1)) {
|
---|
| 199 | pt_ann[i] += pt;
|
---|
| 200 | }
|
---|
| 201 | }
|
---|
| 202 | }
|
---|
| 203 | }
|
---|
| 204 | fItNeutralInputArray->Reset();
|
---|
| 205 | while ((constituent = static_cast<Candidate*>(fItNeutralInputArray->Next()))) {
|
---|
| 206 | if (constituent->Momentum.DeltaR(candidate->Momentum) < fParameterR) {
|
---|
| 207 | float pt = constituent->Momentum.Pt();
|
---|
| 208 | sumpt += pt;
|
---|
| 209 | float dr = candidate->Momentum.DeltaR(constituent->Momentum);
|
---|
| 210 | sumdrsqptsq += dr*dr*pt*pt;
|
---|
| 211 | sumptsq += pt*pt;
|
---|
| 212 | nn++;
|
---|
| 213 | for (int i = 0 ; i < 5 ; i++) {
|
---|
| 214 | if (dr > 0.1*i && dr < 0.1*(i+1)) {
|
---|
| 215 | pt_ann[i] += pt;
|
---|
| 216 | }
|
---|
| 217 | }
|
---|
| 218 | }
|
---|
| 219 | }
|
---|
| 220 | }
|
---|
| 221 |
|
---|
| 222 | if (sumptch > 0.) {
|
---|
| 223 | candidate->Beta = sumptchpu/sumptch;
|
---|
| 224 | candidate->BetaStar = sumptchpv/sumptch;
|
---|
| 225 | } else {
|
---|
| 226 | candidate->Beta = -999.;
|
---|
| 227 | candidate->BetaStar = -999.;
|
---|
| 228 | }
|
---|
| 229 | if (sumptsq > 0.) {
|
---|
| 230 | candidate->MeanSqDeltaR = sumdrsqptsq/sumptsq;
|
---|
| 231 | } else {
|
---|
| 232 | candidate->MeanSqDeltaR = -999.;
|
---|
| 233 | }
|
---|
| 234 | candidate->NCharged = nc;
|
---|
| 235 | candidate->NNeutrals = nn;
|
---|
| 236 | if (sumpt > 0.) {
|
---|
| 237 | candidate->PTD = TMath::Sqrt(sumptsq) / sumpt;
|
---|
| 238 | for (int i = 0 ; i < 5 ; i++) {
|
---|
| 239 | candidate->FracPt[i] = pt_ann[i]/sumpt;
|
---|
| 240 | }
|
---|
| 241 | } else {
|
---|
| 242 | candidate->PTD = -999.;
|
---|
| 243 | for (int i = 0 ; i < 5 ; i++) {
|
---|
| 244 | candidate->FracPt[i] = -999.;
|
---|
| 245 | }
|
---|
| 246 | }
|
---|
| 247 |
|
---|
| 248 | fOutputArray->Add(candidate);
|
---|
| 249 | }
|
---|
| 250 | }
|
---|
| 251 |
|
---|
| 252 | //------------------------------------------------------------------------------
|
---|
| 253 |
|
---|