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