Changeset d77b51d in git for modules/BTagging.cc
- Timestamp:
- Sep 29, 2015, 2:08:10 PM (9 years ago)
- Branches:
- ImprovedOutputFile, Timing, dual_readout, llp, master
- Children:
- a98c7ef
- Parents:
- d870fc5 (diff), 06ec139 (diff)
Note: this is a merge changeset, the changes displayed below correspond to the merge itself.
Use the(diff)
links above to see all the changes relative to each parent. - File:
-
- 1 edited
Legend:
- Unmodified
- Added
- Removed
-
modules/BTagging.cc
rd870fc5 rd77b51d 22 22 * Determines origin of jet, 23 23 * applies b-tagging efficiency (miss identification rate) formulas 24 * and sets b-tagging flags 24 * and sets b-tagging flags 25 25 * 26 26 * \author P. Demin - UCL, Louvain-la-Neuve … … 34 34 #include "classes/DelphesFormula.h" 35 35 36 #include "ExRootAnalysis/ExRootResult.h"37 #include "ExRootAnalysis/ExRootFilter.h"38 #include "ExRootAnalysis/ExRootClassifier.h"39 40 36 #include "TMath.h" 41 37 #include "TString.h" … … 46 42 #include "TLorentzVector.h" 47 43 48 #include <algorithm> 44 #include <algorithm> 49 45 #include <stdexcept> 50 46 #include <iostream> … … 55 51 //------------------------------------------------------------------------------ 56 52 57 class BTaggingPartonClassifier : public ExRootClassifier 53 BTagging::BTagging() : 54 fItJetInputArray(0) 58 55 { 59 public:60 61 BTaggingPartonClassifier() {}62 63 Int_t GetCategory(TObject *object);64 65 Double_t fEtaMax, fPTMin;66 };67 68 //------------------------------------------------------------------------------69 70 Int_t BTaggingPartonClassifier::GetCategory(TObject *object)71 {72 Candidate *parton = static_cast<Candidate*>(object);73 const TLorentzVector &momentum = parton->Momentum;74 Int_t pdgCode;75 76 if(momentum.Pt() <= fPTMin || TMath::Abs(momentum.Eta()) > fEtaMax) return -1;77 78 pdgCode = TMath::Abs(parton->PID);79 if(pdgCode != 21 && pdgCode > 5) return -1;80 81 return 0;82 }83 84 //------------------------------------------------------------------------------85 86 BTagging::BTagging() :87 fClassifier(0), fFilter(0),88 fItPartonInputArray(0), fItJetInputArray(0)89 {90 fClassifier = new BTaggingPartonClassifier;91 56 } 92 57 … … 95 60 BTagging::~BTagging() 96 61 { 97 if(fClassifier) delete fClassifier;98 62 } 99 63 … … 109 73 fBitNumber = GetInt("BitNumber", 0); 110 74 111 fDeltaR = GetDouble("DeltaR", 0.5);112 113 fClassifier->fPTMin = GetDouble("PartonPTMin", 1.0);114 fClassifier->fEtaMax = GetDouble("PartonEtaMax", 2.5);115 116 75 // read efficiency formulas 117 76 param = GetParam("EfficiencyFormula"); 118 77 size = param.GetSize(); 119 78 120 79 fEfficiencyMap.clear(); 121 80 for(i = 0; i < size/2; ++i) … … 139 98 // import input array(s) 140 99 141 fPartonInputArray = ImportArray(GetString("PartonInputArray", "Delphes/partons"));142 fItPartonInputArray = fPartonInputArray->MakeIterator();143 144 fFilter = new ExRootFilter(fPartonInputArray);145 146 100 fJetInputArray = ImportArray(GetString("JetInputArray", "FastJetFinder/jets")); 147 101 fItJetInputArray = fJetInputArray->MakeIterator(); … … 155 109 DelphesFormula *formula; 156 110 157 if(fFilter) delete fFilter;158 111 if(fItJetInputArray) delete fItJetInputArray; 159 if(fItPartonInputArray) delete fItPartonInputArray;160 112 161 113 for(itEfficiencyMap = fEfficiencyMap.begin(); itEfficiencyMap != fEfficiencyMap.end(); ++itEfficiencyMap) … … 170 122 void BTagging::Process() 171 123 { 172 Candidate *jet, *parton; 173 Double_t pt, eta, phi; 174 TObjArray *partonArray; 124 Candidate *jet; 125 Double_t pt, eta, phi, e; 175 126 map< Int_t, DelphesFormula * >::iterator itEfficiencyMap; 176 127 DelphesFormula *formula; 177 Int_t pdgCode, pdgCodeMax;178 128 179 // select quark and gluons180 fFilter->Reset();181 partonArray = fFilter->GetSubArray(fClassifier, 0);182 183 if(partonArray == 0) return;184 185 TIter itPartonArray(partonArray);186 187 129 // loop over all input jets 188 130 fItJetInputArray->Reset(); … … 190 132 { 191 133 const TLorentzVector &jetMomentum = jet->Momentum; 192 pdgCodeMax = -1;193 134 eta = jetMomentum.Eta(); 194 135 phi = jetMomentum.Phi(); 195 136 pt = jetMomentum.Pt(); 137 e = jetMomentum.E(); 196 138 197 // loop over all input partons 198 itPartonArray.Reset(); 199 while((parton = static_cast<Candidate*>(itPartonArray.Next()))) 200 { 201 pdgCode = TMath::Abs(parton->PID); 202 if(pdgCode == 21) pdgCode = 0; 203 if(jetMomentum.DeltaR(parton->Momentum) <= fDeltaR) 204 { 205 if(pdgCodeMax < pdgCode) pdgCodeMax = pdgCode; 206 } 207 } 208 if(pdgCodeMax == 0) pdgCodeMax = 21; 209 if(pdgCodeMax == -1) pdgCodeMax = 0; 210 211 // find an efficency formula 212 itEfficiencyMap = fEfficiencyMap.find(pdgCodeMax); 139 // find an efficiency formula 140 itEfficiencyMap = fEfficiencyMap.find(jet->Flavor); 213 141 if(itEfficiencyMap == fEfficiencyMap.end()) 214 142 { … … 217 145 formula = itEfficiencyMap->second; 218 146 219 // apply an efficency formula 220 jet->BTag |= (gRandom->Uniform() <= formula->Eval(pt, eta)) << fBitNumber; 147 // apply an efficiency formula 148 jet->BTag |= (gRandom->Uniform() <= formula->Eval(pt, eta, phi, e)) << fBitNumber; 149 150 // find an efficiency formula for algo flavor definition 151 itEfficiencyMap = fEfficiencyMap.find(jet->FlavorAlgo); 152 if(itEfficiencyMap == fEfficiencyMap.end()) 153 { 154 itEfficiencyMap = fEfficiencyMap.find(0); 155 } 156 formula = itEfficiencyMap->second; 157 158 // apply an efficiency formula 159 jet->BTagAlgo |= (gRandom->Uniform() <= formula->Eval(pt, eta, phi, e)) << fBitNumber; 160 161 // find an efficiency formula for phys flavor definition 162 itEfficiencyMap = fEfficiencyMap.find(jet->FlavorPhys); 163 if(itEfficiencyMap == fEfficiencyMap.end()) 164 { 165 itEfficiencyMap = fEfficiencyMap.find(0); 166 } 167 formula = itEfficiencyMap->second; 168 169 // apply an efficiency formula 170 jet->BTagPhys |= (gRandom->Uniform() <= formula->Eval(pt, eta, phi, e)) << fBitNumber; 221 171 } 222 172 }
Note:
See TracChangeset
for help on using the changeset viewer.