CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/PhysicsTools/TagAndProbe/src/BaseTreeFiller.cc

Go to the documentation of this file.
00001 #include "PhysicsTools/TagAndProbe/interface/BaseTreeFiller.h"
00002 #include "FWCore/ServiceRegistry/interface/Service.h"
00003 #include "CommonTools/UtilAlgos/interface/TFileService.h"
00004 
00005 #include "DataFormats/METReco/interface/MET.h"
00006 #include "DataFormats/METReco/interface/METCollection.h"
00007 #include "DataFormats/METReco/interface/CaloMET.h"
00008 #include "DataFormats/METReco/interface/CaloMETCollection.h"
00009 #include "DataFormats/METReco/interface/PFMET.h"
00010 #include "DataFormats/METReco/interface/PFMETCollection.h"
00011 #include "DataFormats/VertexReco/interface/Vertex.h"
00012 #include "DataFormats/VertexReco/interface/VertexFwd.h"
00013 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
00014 #include "PhysicsTools/TagAndProbe/interface/ColinsSoperVariables.h"
00015 
00016 #include <TList.h>
00017 #include <TObjString.h>
00018 
00019 tnp::ProbeVariable::~ProbeVariable() {}
00020 
00021 tnp::ProbeFlag::~ProbeFlag() {}
00022 
00023 void tnp::ProbeFlag::init(const edm::Event &iEvent) const {
00024     if (external_) {
00025         edm::Handle<edm::View<reco::Candidate> > view;
00026         iEvent.getByLabel(src_, view);
00027         passingProbes_.clear();
00028         for (size_t i = 0, n = view->size(); i < n; ++i) passingProbes_.push_back(view->refAt(i));
00029     }
00030 
00031 }
00032 
00033 void tnp::ProbeFlag::fill(const reco::CandidateBaseRef &probe) const {
00034     if (external_) {
00035         value_ = (std::find(passingProbes_.begin(), passingProbes_.end(), probe) != passingProbes_.end());
00036     } else {
00037         value_ = bool(cut_(*probe));
00038     }
00039 }
00040 
00041 tnp::BaseTreeFiller::BaseTreeFiller(const char *name, const edm::ParameterSet iConfig) {
00042     // make trees as requested
00043     edm::Service<TFileService> fs;
00044     tree_ = fs->make<TTree>(name,name);
00045 
00046     // add the branches
00047     addBranches_(tree_, iConfig, "");
00048 
00049     // set up weights, if needed
00050     if (iConfig.existsAs<double>("eventWeight")) { 
00051         weightMode_ = Fixed;
00052         weight_ = iConfig.getParameter<double>("eventWeight");
00053     } else if (iConfig.existsAs<edm::InputTag>("eventWeight")) { 
00054         weightMode_ = External;
00055         weightSrc_ = iConfig.getParameter<edm::InputTag>("eventWeight");
00056     } else {
00057         weightMode_ = None;
00058     }
00059     if (weightMode_ != None) {
00060         tree_->Branch("weight", &weight_, "weight/F");
00061     }
00062 
00063     addRunLumiInfo_ = iConfig.existsAs<bool>("addRunLumiInfo") ? iConfig.getParameter<bool>("addRunLumiInfo") : false;
00064     if (addRunLumiInfo_) {
00065          tree_->Branch("run",  &run_,  "run/i");
00066          tree_->Branch("lumi", &lumi_, "lumi/i");
00067          tree_->Branch("event", &event_, "event/i");
00068     }
00069     addEventVariablesInfo_ = iConfig.existsAs<bool>("addEventVariablesInfo") ? iConfig.getParameter<bool>("addEventVariablesInfo") : false;
00070     if (addEventVariablesInfo_) {      
00071       tree_->Branch("event_nPV"        ,&mNPV_                 ,"mNPV/I");
00072       tree_->Branch("event_met_calomet"    ,&mMET_                ,"mMET/F");
00073       tree_->Branch("event_met_calosumet"  ,&mSumET_              ,"mSumET/F");
00074       tree_->Branch("event_met_calometsignificance",&mMETSign_    ,"mMETSign/F");
00075       tree_->Branch("event_met_tcmet"    ,&mtcMET_                ,"mtcMET/F");
00076       tree_->Branch("event_met_tcsumet"  ,&mtcSumET_              ,"mtcSumET/F");
00077       tree_->Branch("event_met_tcmetsignificance",&mtcMETSign_    ,"mtcMETSign/F");
00078       tree_->Branch("event_met_pfmet"    ,&mpfMET_                ,"mpfMET/F");
00079       tree_->Branch("event_met_pfsumet"  ,&mpfSumET_              ,"mpfSumET/F");
00080       tree_->Branch("event_met_pfmetsignificance",&mpfMETSign_    ,"mpfMETSign/F");
00081       tree_->Branch("event_PrimaryVertex_x"  ,&mPVx_              ,"mPVx/F");
00082       tree_->Branch("event_PrimaryVertex_y"  ,&mPVy_              ,"mPVy/F");
00083       tree_->Branch("event_PrimaryVertex_z"  ,&mPVz_              ,"mPVz/F");
00084       tree_->Branch("event_BeamSpot_x"       ,&mBSx_              ,"mBSx/F");
00085       tree_->Branch("event_BeamSpot_y"       ,&mBSy_              ,"mBSy/F");
00086       tree_->Branch("event_BeamSpot_z"       ,&mBSz_              ,"mBSz/F");
00087     }
00088 
00089     ignoreExceptions_ = iConfig.existsAs<bool>("ignoreExceptions") ? iConfig.getParameter<bool>("ignoreExceptions") : false;
00090 }
00091 
00092 tnp::BaseTreeFiller::BaseTreeFiller(BaseTreeFiller &main, const edm::ParameterSet &iConfig, const std::string &branchNamePrefix) :
00093     addEventVariablesInfo_(false),
00094     tree_(0)
00095 {
00096     addBranches_(main.tree_, iConfig, branchNamePrefix);
00097 }
00098 
00099 void
00100 tnp::BaseTreeFiller::addBranches_(TTree *tree, const edm::ParameterSet &iConfig, const std::string &branchNamePrefix) {
00101     // set up variables
00102     edm::ParameterSet variables = iConfig.getParameter<edm::ParameterSet>("variables");
00103     //.. the ones that are strings
00104     std::vector<std::string> stringVars = variables.getParameterNamesForType<std::string>();
00105     for (std::vector<std::string>::const_iterator it = stringVars.begin(), ed = stringVars.end(); it != ed; ++it) {
00106         vars_.push_back(tnp::ProbeVariable(branchNamePrefix + *it, variables.getParameter<std::string>(*it)));
00107     }
00108     //.. the ones that are InputTags
00109     std::vector<std::string> inputTagVars = variables.getParameterNamesForType<edm::InputTag>();
00110     for (std::vector<std::string>::const_iterator it = inputTagVars.begin(), ed = inputTagVars.end(); it != ed; ++it) {
00111         vars_.push_back(tnp::ProbeVariable(branchNamePrefix + *it, variables.getParameter<edm::InputTag>(*it)));
00112     }
00113     // set up flags
00114     edm::ParameterSet flags = iConfig.getParameter<edm::ParameterSet>("flags");
00115     //.. the ones that are strings
00116     std::vector<std::string> stringFlags = flags.getParameterNamesForType<std::string>();
00117     for (std::vector<std::string>::const_iterator it = stringFlags.begin(), ed = stringFlags.end(); it != ed; ++it) {
00118         flags_.push_back(tnp::ProbeFlag(branchNamePrefix + *it, flags.getParameter<std::string>(*it)));
00119     }
00120     //.. the ones that are InputTags
00121     std::vector<std::string> inputTagFlags = flags.getParameterNamesForType<edm::InputTag>();
00122     for (std::vector<std::string>::const_iterator it = inputTagFlags.begin(), ed = inputTagFlags.end(); it != ed; ++it) {
00123         flags_.push_back(tnp::ProbeFlag(branchNamePrefix + *it, flags.getParameter<edm::InputTag>(*it)));
00124     }
00125 
00126     // then make all the variables in the trees
00127     for (std::vector<tnp::ProbeVariable>::iterator it = vars_.begin(), ed = vars_.end(); it != ed; ++it) {
00128         tree->Branch(it->name().c_str(), it->address(), (it->name()+"/F").c_str());
00129     }
00130     
00131     for (std::vector<tnp::ProbeFlag>::iterator it = flags_.begin(), ed = flags_.end(); it != ed; ++it) {
00132         tree->Branch(it->name().c_str(), it->address(), (it->name()+"/I").c_str());
00133     }
00134     
00135 }
00136 
00137 tnp::BaseTreeFiller::~BaseTreeFiller() { }
00138 
00139 void tnp::BaseTreeFiller::init(const edm::Event &iEvent) const {
00140     run_  = iEvent.id().run();
00141     lumi_ = iEvent.id().luminosityBlock();
00142     event_ = iEvent.id().event(); 
00143 
00144     for (std::vector<tnp::ProbeVariable>::const_iterator it = vars_.begin(), ed = vars_.end(); it != ed; ++it) {
00145         it->init(iEvent);
00146     }
00147     for (std::vector<tnp::ProbeFlag>::const_iterator it = flags_.begin(), ed = flags_.end(); it != ed; ++it) {
00148         it->init(iEvent);
00149     }
00150     if (weightMode_ == External) {
00151         edm::Handle<double> weight;
00152         iEvent.getByLabel(weightSrc_, weight);
00153         weight_ = *weight;
00154     }
00155 
00156     if (addEventVariablesInfo_) {
00159         edm::Handle<reco::VertexCollection> recVtxs;
00160         iEvent.getByLabel("offlinePrimaryVertices",recVtxs);
00161         mNPV_ = 0;
00162         mPVx_ =  100.0;
00163         mPVy_ =  100.0;
00164         mPVz_ =  100.0;
00165 
00166         for(unsigned int ind=0;ind<recVtxs->size();ind++) {
00167           if (!((*recVtxs)[ind].isFake()) && ((*recVtxs)[ind].ndof()>4) 
00168               && (fabs((*recVtxs)[ind].z())<=24.0) &&  
00169               ((*recVtxs)[ind].position().Rho()<=2.0) ) {
00170             mNPV_++;
00171             if(mNPV_==1) { // store the first good primary vertex
00172               mPVx_ = (*recVtxs)[ind].x();
00173               mPVy_ = (*recVtxs)[ind].y();
00174               mPVz_ = (*recVtxs)[ind].z();
00175             }
00176           }
00177         }
00178 
00179 
00181         edm::Handle<reco::BeamSpot> beamSpot;
00182         iEvent.getByLabel("offlineBeamSpot", beamSpot);
00183         mBSx_ = beamSpot->position().X();
00184         mBSy_ = beamSpot->position().Y();
00185         mBSz_ = beamSpot->position().Z();
00186 
00187 
00189         edm::Handle<reco::CaloMETCollection> met;
00190         iEvent.getByLabel("met",met);
00191         if (met->size() == 0) {
00192           mMET_   = -1;
00193           mSumET_ = -1;
00194           mMETSign_ = -1;
00195         }
00196         else {
00197           mMET_   = (*met)[0].et();
00198           mSumET_ = (*met)[0].sumEt();
00199           mMETSign_ = (*met)[0].significance();
00200         }
00201 
00203         edm::Handle<reco::METCollection> tcmet;
00204         iEvent.getByLabel("tcMet", tcmet);
00205         if (tcmet->size() == 0) {
00206           mtcMET_   = -1;
00207           mtcSumET_ = -1;
00208           mtcMETSign_ = -1;
00209         }
00210         else {
00211           mtcMET_   = (*tcmet)[0].et();
00212           mtcSumET_ = (*tcmet)[0].sumEt();
00213           mtcMETSign_ = (*tcmet)[0].significance();
00214         }
00215 
00217         edm::Handle<reco::PFMETCollection> pfmet;
00218         iEvent.getByLabel("pfMet", pfmet);
00219         if (pfmet->size() == 0) {
00220           mpfMET_   = -1;
00221           mpfSumET_ = -1;
00222           mpfMETSign_ = -1;
00223         }
00224         else {
00225           mpfMET_   = (*pfmet)[0].et();
00226           mpfSumET_ = (*pfmet)[0].sumEt();
00227           mpfMETSign_ = (*pfmet)[0].significance();
00228         }
00229     }
00230 
00231 }
00232 
00233 void tnp::BaseTreeFiller::fill(const reco::CandidateBaseRef &probe) const {
00234     for (std::vector<tnp::ProbeVariable>::const_iterator it = vars_.begin(), ed = vars_.end(); it != ed; ++it) {
00235         if (ignoreExceptions_)  {
00236             try{ it->fill(probe); } catch(cms::Exception &ex ){}
00237         } else {
00238             it->fill(probe);
00239         }
00240     }
00241 
00242     for (std::vector<tnp::ProbeFlag>::const_iterator it = flags_.begin(), ed = flags_.end(); it != ed; ++it) {
00243         if (ignoreExceptions_)  {
00244             try{ it->fill(probe); } catch(cms::Exception &ex ){}
00245         } else {
00246             it->fill(probe);
00247         }
00248     }
00249     if (tree_) tree_->Fill();
00250 }
00251 void tnp::BaseTreeFiller::writeProvenance(const edm::ParameterSet &pset) const {
00252     TList *list = tree_->GetUserInfo();
00253     list->Add(new TObjString(pset.dump().c_str()));
00254 }