CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/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 soft-muon-based b-tag
00047   ohBJetMuL25Tag            = new int[kMaxBJets];          // do not optimize
00048   ohBJetMuL3Tag             = new float[kMaxBJets];        // SoftMuonbyPt
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(ohBJetMuL25Tag,            '\0', kMaxBJets * sizeof(int));
00076   std::memset(ohBJetMuL3Tag,             '\0', kMaxBJets * sizeof(float));
00077   std::memset(ohBJetPerfL25Tag,          '\0', kMaxBJets * sizeof(int));
00078   std::memset(ohBJetPerfL3Tag,           '\0', kMaxBJets * sizeof(int));
00079 }
00080 
00081 HLTBJet::~HLTBJet() 
00082 { }
00083 
00084 void HLTBJet::setup(const edm::ParameterSet & config, TTree * tree)
00085 {
00086   // create the TTree branches
00087   if (tree) {
00088     tree->Branch("NohBJetL2",               & NohBJetL2,             "NohBJetL2/I");
00089     tree->Branch("ohBJetL2Energy",          ohBJetL2Energy,          "ohBJetL2Energy[NohBJetL2]/F");
00090     tree->Branch("ohBJetL2Et",              ohBJetL2Et,              "ohBJetL2Et[NohBJetL2]/F");
00091     tree->Branch("ohBJetL2Pt",              ohBJetL2Pt,              "ohBJetL2Pt[NohBJetL2]/F");
00092     tree->Branch("ohBJetL2Eta",             ohBJetL2Eta,             "ohBJetL2Eta[NohBJetL2]/F");
00093     tree->Branch("ohBJetL2Phi",             ohBJetL2Phi,             "ohBJetL2Phi[NohBJetL2]/F");
00094                 
00095     tree->Branch("NohBJetL2Corrected",      & NohBJetL2Corrected,    "NohBJetL2Corrected/I");
00096     tree->Branch("ohBJetL2CorrectedEnergy", ohBJetL2CorrectedEnergy, "ohBJetL2CorrectedEnergy[NohBJetL2Corrected]/F");
00097     tree->Branch("ohBJetL2CorrectedEt",     ohBJetL2CorrectedEt,     "ohBJetL2CorrectedEt[NohBJetL2Corrected]/F");
00098     tree->Branch("ohBJetL2CorrectedPt",     ohBJetL2CorrectedPt,     "ohBJetL2CorrectedPt[NohBJetL2Corrected]/F");
00099     tree->Branch("ohBJetL2CorrectedEta",    ohBJetL2CorrectedEta,    "ohBJetL2CorrectedEta[NohBJetL2Corrected]/F");
00100     tree->Branch("ohBJetL2CorrectedPhi",    ohBJetL2CorrectedPhi,    "ohBJetL2CorrectedPhi[NohBJetL2Corrected]/F");
00101     
00102     tree->Branch("ohBJetIPL25Tag",          ohBJetIPL25Tag,          "ohBJetIPL25Tag[NohBJetL2]/F");
00103     tree->Branch("ohBJetIPL3Tag",           ohBJetIPL3Tag,           "ohBJetIPL3Tag[NohBJetL2]/F");
00104 
00105     tree->Branch("ohBJetMuL25Tag",          ohBJetMuL25Tag,          "ohBJetMuL25Tag[NohBJetL2]/I");
00106     tree->Branch("ohBJetMuL3Tag",           ohBJetMuL3Tag,           "ohBJetMuL3Tag[NohBJetL2]/F");
00107     tree->Branch("ohBJetPerfL25Tag",        ohBJetPerfL25Tag,        "ohBJetPerfL25Tag[NohBJetL2]/I");
00108     tree->Branch("ohBJetPerfL3Tag",         ohBJetPerfL3Tag,         "ohBJetPerfL3Tag[NohBJetL2]/I");
00109   }
00110 }
00111 
00112 void HLTBJet::analyze(
00113         const edm::Handle<edm::View<reco::Jet> >  & rawBJets,
00114         const edm::Handle<edm::View<reco::Jet> >  & correctedBJets,
00115         const edm::Handle<reco::JetTagCollection> & lifetimeBJetsL25,
00116         const edm::Handle<reco::JetTagCollection> & lifetimeBJetsL3,
00117         const edm::Handle<reco::JetTagCollection> & softmuonBJetsL25,
00118         const edm::Handle<reco::JetTagCollection> & softmuonBJetsL3,
00119         const edm::Handle<reco::JetTagCollection> & performanceBJetsL25,
00120         const edm::Handle<reco::JetTagCollection> & performanceBJetsL3,
00121         TTree * tree) 
00122 {
00123   // reset the tree variables
00124   clear();
00125   
00126   // if the required collections are available, fill the corresponding tree branches
00127   if (rawBJets.isValid())
00128     analyseJets(* rawBJets);
00129 
00130   if (correctedBJets.isValid())
00131     analyseCorrectedJets(* correctedBJets);
00132  
00133   if (correctedBJets.isValid() and lifetimeBJetsL25.isValid() and lifetimeBJetsL3.isValid())
00134     analyseLifetime(* correctedBJets, * lifetimeBJetsL25, * lifetimeBJetsL3);
00135 
00136   if (correctedBJets.isValid() and softmuonBJetsL25.isValid() and softmuonBJetsL3.isValid())
00137     analyseSoftmuon(* correctedBJets, * softmuonBJetsL25, * softmuonBJetsL3);
00138   
00139   if (correctedBJets.isValid() and performanceBJetsL25.isValid() and performanceBJetsL3.isValid())
00140     analysePerformance(* correctedBJets, * performanceBJetsL25, * performanceBJetsL3);
00141 
00142 }
00143 
00144 void HLTBJet::analyseJets(const edm::View<reco::Jet> & jets)
00145 {
00146   // the jets need to be persistable, so .size() returns an 'unsigned int' to be stable across the architectures
00147   // so, for the comparison, we cast back to size_t
00148   size_t size = std::min(kMaxBJets, size_t(jets.size()) ); 
00149   NohBJetL2 = size;
00150   for (size_t i = 0; i < size; ++i) {
00151     ohBJetL2Energy[i] = jets[i].energy();
00152     ohBJetL2Et[i]     = jets[i].et();
00153     ohBJetL2Pt[i]     = jets[i].pt();
00154     ohBJetL2Eta[i]    = jets[i].eta();
00155     ohBJetL2Phi[i]    = jets[i].phi();
00156   }
00157 }
00158 
00159 void HLTBJet::analyseCorrectedJets(const edm::View<reco::Jet> & jets)
00160 {
00161   // the jets need to be persistable, so .size() returns an 'unsigned int' to be stable across the architectures
00162   // so, for the comparison, we cast back to size_t
00163   size_t size = std::min(kMaxBJets, size_t(jets.size()) );
00164   NohBJetL2Corrected = size;
00165   for (size_t i = 0; i < size; ++i) {
00166     ohBJetL2CorrectedEnergy[i] = jets[i].energy();
00167     ohBJetL2CorrectedEt[i]     = jets[i].et();
00168     ohBJetL2CorrectedPt[i]     = jets[i].pt();
00169     ohBJetL2CorrectedEta[i]    = jets[i].eta();
00170     ohBJetL2CorrectedPhi[i]    = jets[i].phi();
00171   }
00172 }
00173 
00174 void HLTBJet::analyseLifetime(
00175     const edm::View<reco::Jet>   & jets, 
00176     const reco::JetTagCollection & tagsL25, 
00177     const reco::JetTagCollection & tagsL3)
00178 {
00179   if (tagsL25.size() != jets.size()) {
00180     edm::LogWarning("OpenHLT") << kBTagLifetimeBJetsL25 << " collection has " << tagsL25.size() << " elements, but " << jets.size() << " where expected from L2" << std::endl;
00181     return;
00182   }
00183   if (tagsL3.size() != jets.size()) {
00184     edm::LogWarning("OpenHLT") << kBTagLifetimeBJetsL3 << " collection has " << tagsL3.size() << " elements, but " << jets.size() << " where expected from L2" << std::endl;
00185     return;
00186   }
00187   // the jets need to be persistable, so .size() returns an 'unsigned int' to be stable across the architectures
00188   // so, for the comparison, we cast back to size_t
00189   size_t size = std::min(kMaxBJets, size_t(jets.size()) );
00190   for (size_t i = 0; i < size; i++) {
00191     ohBJetIPL25Tag[i] = tagsL25[i].second;
00192     ohBJetIPL3Tag[i]  = tagsL3[i].second;
00193   }
00194 }
00195 
00196 void HLTBJet::analyseSoftmuon(
00197     const edm::View<reco::Jet>   & jets, 
00198     const reco::JetTagCollection & tagsL25, 
00199     const reco::JetTagCollection & tagsL3)
00200 {
00201   if (tagsL25.size() != jets.size()) {
00202     edm::LogWarning("OpenHLT") << kBTagSoftmuonBJetsL25 << " collection has " << tagsL25.size() << " elements, but " << jets.size() << " where expected from L2" << std::endl;
00203     return;
00204   }
00205   if (tagsL3.size() != jets.size()) {
00206     edm::LogWarning("OpenHLT") << kBTagSoftmuonBJetsL3 << " collection has " << tagsL3.size() << " elements, but " << jets.size() << " where expected from L2" << std::endl;
00207     return;
00208   }
00209   // the jets need to be persistable, so .size() returns an 'unsigned int' to be stable across the architectures
00210   // so, for the comparison, we cast back to size_t
00211   size_t size = std::min(kMaxBJets, size_t(jets.size()) );
00212   for (size_t i = 0; i < size; i++) {
00213     ohBJetMuL25Tag[i] = (tagsL25[i].second > 0.) ? 1 : 0;
00214     ohBJetMuL3Tag[i]  = tagsL3[i].second;
00215   }
00216 }
00217 
00218 void HLTBJet::analysePerformance(
00219     const edm::View<reco::Jet>   & jets, 
00220     const reco::JetTagCollection & tagsL25, 
00221     const reco::JetTagCollection & tagsL3)
00222 {
00223   if (tagsL25.size() != jets.size()) {
00224     edm::LogWarning("OpenHLT") << kBTagPerformanceBJetsL25 << " collection has " << tagsL25.size() << " elements, but " << jets.size() << " where expected from L2" << std::endl;
00225     return;
00226   }
00227   if (tagsL3.size() != jets.size()) {
00228     edm::LogWarning("OpenHLT") << kBTagPerformanceBJetsL3 << " collection has " << tagsL3.size() << " elements, but " << jets.size() << " where expected from L2" << std::endl;
00229     return;
00230   }
00231   // the jets need to be persistable, so .size() returns an 'unsigned int' to be stable across the architectures
00232   // so, for the comparison, we cast back to size_t
00233   size_t size = std::min(kMaxBJets, size_t(jets.size()) );
00234   for (size_t i = 0; i < size; i++) {
00235     ohBJetPerfL25Tag[i] = (tagsL25[i].second > 0.) ? 1 : 0;
00236     ohBJetPerfL3Tag[i]  = (tagsL3[i].second  > 0.) ? 1 : 0;
00237   }
00238 }