CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/HLTrigger/HLTanalyzers/src/HLTBJet.cc

Go to the documentation of this file.
00001 #include <cmath>
00002 #include <algorithm>
00003 #include <utility>
00004 #include <cstring>
00005 #include <boost/foreach.hpp>
00006 
00007 #include <TTree.h>
00008 
00009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00010 #include "FWCore/Framework/interface/ESHandle.h" 
00011 #include "FWCore/Framework/interface/EventSetup.h"
00012 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00013 
00014 #include "DataFormats/Common/interface/View.h"
00015 #include "DataFormats/JetReco/interface/Jet.h"
00016 #include "DataFormats/BTauReco/interface/JetTag.h"
00017 
00018 #include "HLTrigger/HLTanalyzers/interface/HLTBJet.h"
00019 
00020 static const size_t kMaxBJets = 10;
00021 
00022 #include "HLTMessages.h"
00023 
00024 HLTBJet::HLTBJet() 
00025 {
00026   // set of variables for uncorrected L2 jets
00027   NohBJetL2                 = 0;
00028   ohBJetL2Energy            = new float[kMaxBJets];
00029   ohBJetL2Et                = new float[kMaxBJets];
00030   ohBJetL2Pt                = new float[kMaxBJets];
00031   ohBJetL2Eta               = new float[kMaxBJets];
00032   ohBJetL2Phi               = new float[kMaxBJets];
00033               
00034   // set of variables for corrected L2 jets
00035   NohBJetL2Corrected        = 0;
00036   ohBJetL2CorrectedEnergy   = new float[kMaxBJets];
00037   ohBJetL2CorrectedEt       = new float[kMaxBJets];
00038   ohBJetL2CorrectedPt       = new float[kMaxBJets];
00039   ohBJetL2CorrectedEta      = new float[kMaxBJets];
00040   ohBJetL2CorrectedPhi      = new float[kMaxBJets];
00041   
00042   // set of variables for lifetime-based b-tag
00043   ohBJetIPL25Tag            = new float[kMaxBJets];
00044   ohBJetIPL3Tag             = new float[kMaxBJets];
00045 
00046   // set of variables for lifetime-based b-tag Single Track
00047   ohBJetIPL25TagSingleTrack = new float[kMaxBJets];
00048   ohBJetIPL3TagSingleTrack  = new float[kMaxBJets];
00049   
00050   // set of variables for b-tagging performance measurements
00051   // SoftMuonbyDR
00052   ohBJetPerfL25Tag          = new int[kMaxBJets];          // do not optimize 
00053   ohBJetPerfL3Tag           = new int[kMaxBJets];          // do not optimize
00054 }
00055 
00056 void HLTBJet::clear() 
00057 {
00058   NohBJetL2          = 0;
00059   NohBJetL2Corrected = 0;
00060   std::memset(ohBJetL2Energy,            '\0', kMaxBJets * sizeof(float));
00061   std::memset(ohBJetL2Et,                '\0', kMaxBJets * sizeof(float));
00062   std::memset(ohBJetL2Et,                '\0', kMaxBJets * sizeof(float));
00063   std::memset(ohBJetL2Pt,                '\0', kMaxBJets * sizeof(float));
00064   std::memset(ohBJetL2Eta,               '\0', kMaxBJets * sizeof(float));
00065   std::memset(ohBJetL2Phi,               '\0', kMaxBJets * sizeof(float));
00066   std::memset(ohBJetL2CorrectedEnergy,   '\0', kMaxBJets * sizeof(float));
00067   std::memset(ohBJetL2CorrectedEt,       '\0', kMaxBJets * sizeof(float));
00068   std::memset(ohBJetL2CorrectedPt,       '\0', kMaxBJets * sizeof(float));
00069   std::memset(ohBJetL2CorrectedEta,      '\0', kMaxBJets * sizeof(float));
00070   std::memset(ohBJetL2CorrectedPhi,      '\0', kMaxBJets * sizeof(float));
00071 
00072   std::memset(ohBJetIPL25Tag,            '\0', kMaxBJets * sizeof(float));
00073   std::memset(ohBJetIPL3Tag,             '\0', kMaxBJets * sizeof(float));
00074 
00075   std::memset(ohBJetIPL25TagSingleTrack, '\0', kMaxBJets * sizeof(float));
00076   std::memset(ohBJetIPL3TagSingleTrack,  '\0', kMaxBJets * sizeof(float));
00077 
00078   std::memset(ohBJetPerfL25Tag,          '\0', kMaxBJets * sizeof(int));
00079   std::memset(ohBJetPerfL3Tag,           '\0', kMaxBJets * sizeof(int));
00080 }
00081 
00082 HLTBJet::~HLTBJet() 
00083 { 
00084   delete[] ohBJetL2Energy;
00085   delete[] ohBJetL2Et;
00086   delete[] ohBJetL2Pt;
00087   delete[] ohBJetL2Eta;
00088   delete[] ohBJetL2Phi;
00089   delete[] ohBJetL2CorrectedEnergy;
00090   delete[] ohBJetL2CorrectedEt;
00091   delete[] ohBJetL2CorrectedPt;
00092   delete[] ohBJetL2CorrectedEta;
00093   delete[] ohBJetL2CorrectedPhi;
00094   delete[] ohBJetIPL25Tag;
00095   delete[] ohBJetIPL3Tag;
00096   delete[] ohBJetIPL25TagSingleTrack;
00097   delete[] ohBJetIPL3TagSingleTrack;
00098   delete[] ohBJetPerfL25Tag;
00099   delete[] ohBJetPerfL3Tag;
00100 }
00101 
00102 void HLTBJet::setup(const edm::ParameterSet & config, TTree * tree)
00103 {
00104   // create the TTree branches
00105   if (tree) {
00106     tree->Branch("NohBJetL2",               & NohBJetL2,             "NohBJetL2/I");
00107     tree->Branch("ohBJetL2Energy",          ohBJetL2Energy,          "ohBJetL2Energy[NohBJetL2]/F");
00108     tree->Branch("ohBJetL2Et",              ohBJetL2Et,              "ohBJetL2Et[NohBJetL2]/F");
00109     tree->Branch("ohBJetL2Pt",              ohBJetL2Pt,              "ohBJetL2Pt[NohBJetL2]/F");
00110     tree->Branch("ohBJetL2Eta",             ohBJetL2Eta,             "ohBJetL2Eta[NohBJetL2]/F");
00111     tree->Branch("ohBJetL2Phi",             ohBJetL2Phi,             "ohBJetL2Phi[NohBJetL2]/F");
00112                 
00113     tree->Branch("NohBJetL2Corrected",      & NohBJetL2Corrected,    "NohBJetL2Corrected/I");
00114     tree->Branch("ohBJetL2CorrectedEnergy", ohBJetL2CorrectedEnergy, "ohBJetL2CorrectedEnergy[NohBJetL2Corrected]/F");
00115     tree->Branch("ohBJetL2CorrectedEt",     ohBJetL2CorrectedEt,     "ohBJetL2CorrectedEt[NohBJetL2Corrected]/F");
00116     tree->Branch("ohBJetL2CorrectedPt",     ohBJetL2CorrectedPt,     "ohBJetL2CorrectedPt[NohBJetL2Corrected]/F");
00117     tree->Branch("ohBJetL2CorrectedEta",    ohBJetL2CorrectedEta,    "ohBJetL2CorrectedEta[NohBJetL2Corrected]/F");
00118     tree->Branch("ohBJetL2CorrectedPhi",    ohBJetL2CorrectedPhi,    "ohBJetL2CorrectedPhi[NohBJetL2Corrected]/F");
00119     
00120     tree->Branch("ohBJetIPL25Tag",          ohBJetIPL25Tag,          "ohBJetIPL25Tag[NohBJetL2Corrected]/F");
00121     tree->Branch("ohBJetIPL3Tag",           ohBJetIPL3Tag,           "ohBJetIPL3Tag[NohBJetL2Corrected]/F");
00122 
00123     tree->Branch("ohBJetIPL25TagSingleTrack", ohBJetIPL25TagSingleTrack,      "ohBJetIPL25TagSingleTrack[NohBJetL2Corrected]/F");
00124     tree->Branch("ohBJetIPL3TagSingleTrack",  ohBJetIPL3TagSingleTrack,       "ohBJetIPL3TagSingleTrack[NohBJetL2Corrected]/F");
00125 
00126     tree->Branch("ohBJetPerfL25Tag",        ohBJetPerfL25Tag,        "ohBJetPerfL25Tag[NohBJetL2Corrected]/I");
00127     tree->Branch("ohBJetPerfL3Tag",         ohBJetPerfL3Tag,         "ohBJetPerfL3Tag[NohBJetL2Corrected]/I");
00128   }
00129 }
00130 
00131 void HLTBJet::analyze(
00132         const edm::Handle<edm::View<reco::Jet> >  & rawBJets,
00133         const edm::Handle<edm::View<reco::Jet> >  & correctedBJets,
00134         const edm::Handle<reco::JetTagCollection> & lifetimeBJetsL25,
00135         const edm::Handle<reco::JetTagCollection> & lifetimeBJetsL3,
00136         const edm::Handle<reco::JetTagCollection> & lifetimeBJetsL25SingleTrack,
00137         const edm::Handle<reco::JetTagCollection> & lifetimeBJetsL3SingleTrack,
00138         const edm::Handle<reco::JetTagCollection> & performanceBJetsL25,
00139         const edm::Handle<reco::JetTagCollection> & performanceBJetsL3,
00140         TTree * tree) 
00141 {
00142   // reset the tree variables
00143   clear();
00144   
00145   // if the required collections are available, fill the corresponding tree branches
00146   if (rawBJets.isValid())
00147     analyseJets(* rawBJets);
00148 
00149   if (correctedBJets.isValid())
00150     analyseCorrectedJets(* correctedBJets);
00151  
00152   if (correctedBJets.isValid() and lifetimeBJetsL25.isValid() and lifetimeBJetsL3.isValid())
00153     analyseLifetime(* correctedBJets, * lifetimeBJetsL25, * lifetimeBJetsL3);
00154 
00155   if (correctedBJets.isValid() and lifetimeBJetsL25SingleTrack.isValid() and lifetimeBJetsL3SingleTrack.isValid())
00156     analyseLifetimeSingleTrack(* correctedBJets, * lifetimeBJetsL25SingleTrack, * lifetimeBJetsL3SingleTrack);
00157 
00158   if (correctedBJets.isValid() and performanceBJetsL25.isValid() and performanceBJetsL3.isValid())
00159     analysePerformance(* correctedBJets, * performanceBJetsL25, * performanceBJetsL3);
00160 
00161 }
00162 
00163 void HLTBJet::analyseJets(const edm::View<reco::Jet> & jets)
00164 {
00165   // the jets need to be persistable, so .size() returns an 'unsigned int' to be stable across the architectures
00166   // so, for the comparison, we cast back to size_t
00167   size_t size = std::min(kMaxBJets, size_t(jets.size()) ); 
00168   NohBJetL2 = size;
00169   for (size_t i = 0; i < size; ++i) {
00170     ohBJetL2Energy[i] = jets[i].energy();
00171     ohBJetL2Et[i]     = jets[i].et();
00172     ohBJetL2Pt[i]     = jets[i].pt();
00173     ohBJetL2Eta[i]    = jets[i].eta();
00174     ohBJetL2Phi[i]    = jets[i].phi();
00175   }
00176 }
00177 
00178 void HLTBJet::analyseCorrectedJets(const edm::View<reco::Jet> & jets)
00179 {
00180   // the jets need to be persistable, so .size() returns an 'unsigned int' to be stable across the architectures
00181   // so, for the comparison, we cast back to size_t
00182   size_t size = std::min(kMaxBJets, size_t(jets.size()) );
00183   NohBJetL2Corrected = size;
00184   for (size_t i = 0; i < size; ++i) {
00185     ohBJetL2CorrectedEnergy[i] = jets[i].energy();
00186     ohBJetL2CorrectedEt[i]     = jets[i].et();
00187     ohBJetL2CorrectedPt[i]     = jets[i].pt();
00188     ohBJetL2CorrectedEta[i]    = jets[i].eta();
00189     ohBJetL2CorrectedPhi[i]    = jets[i].phi();
00190   }
00191 }
00192 
00193 void HLTBJet::analyseLifetime(
00194     const edm::View<reco::Jet>   & jets, 
00195     const reco::JetTagCollection & tagsL25, 
00196     const reco::JetTagCollection & tagsL3)
00197 {
00198   if (tagsL25.size() != jets.size()) {
00199     edm::LogWarning("OpenHLT") << kBTagLifetimeBJetsL25 << " collection has " << tagsL25.size() << " elements, but " << jets.size() << " where expected from L2" << std::endl;
00200     return;
00201   }
00202   if (tagsL3.size() != jets.size()) {
00203     edm::LogWarning("OpenHLT") << kBTagLifetimeBJetsL3 << " collection has " << tagsL3.size() << " elements, but " << jets.size() << " where expected from L2" << std::endl;
00204     return;
00205   }
00206   // the jets need to be persistable, so .size() returns an 'unsigned int' to be stable across the architectures
00207   // so, for the comparison, we cast back to size_t
00208   size_t size = std::min(kMaxBJets, size_t(jets.size()) );
00209   for (size_t i = 0; i < size; i++) {
00210     ohBJetIPL25Tag[i] = tagsL25[i].second;
00211     ohBJetIPL3Tag[i]  = tagsL3[i].second;
00212   }
00213 }
00214 
00215 void HLTBJet::analyseLifetimeSingleTrack(
00216     const edm::View<reco::Jet>   & jets, 
00217     const reco::JetTagCollection & tagsL25, 
00218     const reco::JetTagCollection & tagsL3)
00219 {
00220   if (tagsL25.size() != jets.size()) {
00221     edm::LogWarning("OpenHLT") << kBTagLifetimeBJetsL25SingleTrack << " collection has " << tagsL25.size() << " elements, but " << jets.size() << " where expected from L2" << std::endl;
00222     return;
00223   }
00224   if (tagsL3.size() != jets.size()) {
00225     edm::LogWarning("OpenHLT") << kBTagLifetimeBJetsL3SingleTrack << " collection has " << tagsL3.size() << " elements, but " << jets.size() << " where expected from L2" << std::endl;
00226     return;
00227   }
00228   // the jets need to be persistable, so .size() returns an 'unsigned int' to be stable across the architectures
00229   // so, for the comparison, we cast back to size_t
00230   size_t size = std::min(kMaxBJets, size_t(jets.size()) );
00231   for (size_t i = 0; i < size; i++) {
00232     ohBJetIPL25TagSingleTrack[i] = tagsL25[i].second;
00233     ohBJetIPL3TagSingleTrack[i]  = tagsL3[i].second;
00234   }
00235 }
00236 
00237 void HLTBJet::analysePerformance(
00238     const edm::View<reco::Jet>   & jets, 
00239     const reco::JetTagCollection & tagsL25, 
00240     const reco::JetTagCollection & tagsL3)
00241 {
00242   if (tagsL25.size() != jets.size()) {
00243     edm::LogWarning("OpenHLT") << kBTagPerformanceBJetsL25 << " collection has " << tagsL25.size() << " elements, but " << jets.size() << " where expected from L2" << std::endl;
00244     return;
00245   }
00246   if (tagsL3.size() != jets.size()) {
00247     edm::LogWarning("OpenHLT") << kBTagPerformanceBJetsL3 << " collection has " << tagsL3.size() << " elements, but " << jets.size() << " where expected from L2" << std::endl;
00248     return;
00249   }
00250   // the jets need to be persistable, so .size() returns an 'unsigned int' to be stable across the architectures
00251   // so, for the comparison, we cast back to size_t
00252   size_t size = std::min(kMaxBJets, size_t(jets.size()) );
00253   for (size_t i = 0; i < size; i++) {
00254     ohBJetPerfL25Tag[i] = (tagsL25[i].second > 0.) ? 1 : 0;
00255     ohBJetPerfL3Tag[i]  = (tagsL3[i].second  > 0.) ? 1 : 0;
00256   }
00257 }