![]() |
![]() |
#include <TrackingPFG/PileUp/src/MCVerticesAnalyzer.cc>
Public Member Functions | |
MCVerticesAnalyzer (const edm::ParameterSet &) | |
~MCVerticesAnalyzer () | |
Private Member Functions | |
virtual void | analyze (const edm::Event &, const edm::EventSetup &) |
virtual void | beginJob () |
virtual void | beginRun (const edm::Run &, const edm::EventSetup &) |
virtual void | endJob () |
virtual void | endRun (const edm::Run &, const edm::EventSetup &) |
Private Attributes | |
TH1F * | m_hlumi |
TH1F * | m_hmainvtxx |
TH1F * | m_hmainvtxy |
TH1F * | m_hmainvtxz |
TH1F * | m_hnvtx |
TH2F * | m_hnvtxvslumi |
TH1F * | m_hnvtxweight |
TProfile * | m_hnvtxweightprof |
TH1F * | m_hpileupvtxz |
edm::InputTag | m_mctruthcollection |
edm::InputTag | m_pileupcollection |
const bool | m_useweight |
edm::InputTag | m_weight |
Description: <one line="" class="" summary>="">
Implementation: <Notes on="" implementation>="">
Definition at line 58 of file MCVerticesAnalyzer.cc.
MCVerticesAnalyzer::MCVerticesAnalyzer | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 102 of file MCVerticesAnalyzer.cc.
References m_hlumi, m_hmainvtxx, m_hmainvtxy, m_hmainvtxz, m_hnvtx, m_hnvtxvslumi, m_hnvtxweight, m_hnvtxweightprof, m_hpileupvtxz, and m_useweight.
: m_pileupcollection(iConfig.getParameter<edm::InputTag>("pileupSummaryCollection")), m_mctruthcollection(iConfig.getParameter<edm::InputTag>("mcTruthCollection")), m_useweight(iConfig.getParameter<bool>("useWeight")), m_weight(iConfig.getParameter<edm::InputTag>("weightProduct")) { //now do what ever initialization is needed edm::Service<TFileService> tfserv; m_hnvtx = tfserv->make<TH1F>("nvtx","Number of pileup vertices",60,-0.5,59.5); m_hnvtx->GetXaxis()->SetTitle("Number of Interactions"); m_hlumi = tfserv->make<TH1F>("lumi","BX luminosity*xsect",200,0.,50.); m_hlumi->GetXaxis()->SetTitle("Average Number of Interactions"); m_hnvtxvslumi = tfserv->make<TH2F>("nvtxvslumi","Npileup vs BX luminosity*xsect",200,0.,50.,60,-0.5,59.5); m_hnvtxvslumi->GetXaxis()->SetTitle("Average Number of Interactions"); m_hnvtxvslumi->GetYaxis()->SetTitle("Number of Interactions"); if(m_useweight) { m_hnvtxweight = tfserv->make<TH1F>("nvtxweight","Number of pileup vertices (1-w)",60,-0.5,59.5); m_hnvtxweight->GetXaxis()->SetTitle("Number of Interactions"); m_hnvtxweightprof = tfserv->make<TProfile>("nvtxweightprof","Mean (1-w) vs Number of pileup interactions",60,-0.5,59.5); m_hnvtxweightprof->GetXaxis()->SetTitle("Number of Interactions"); } m_hmainvtxx = tfserv->make<TH1F>("mainvtxx","Main vertex x position",200,-.5,.5); m_hmainvtxx->GetXaxis()->SetTitle("X (cm)"); m_hmainvtxy = tfserv->make<TH1F>("mainvtxy","Main vertex y position",200,-.5,.5); m_hmainvtxy->GetXaxis()->SetTitle("Y (cm)"); m_hmainvtxz = tfserv->make<TH1F>("mainvtxz","Main vertex z position",600,-30.,30.); m_hmainvtxz->GetXaxis()->SetTitle("Z (cm)"); m_hpileupvtxz = tfserv->make<TH1F>("pileupvtxz","PileUp vertices z position",600,-30.,30.); m_hpileupvtxz->GetXaxis()->SetTitle("Z (cm)"); }
MCVerticesAnalyzer::~MCVerticesAnalyzer | ( | ) |
Definition at line 144 of file MCVerticesAnalyzer.cc.
{ // do anything here that needs to be done at desctruction time // (e.g. close files, deallocate resources etc.) }
void MCVerticesAnalyzer::analyze | ( | const edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [private, virtual] |
Implements edm::EDAnalyzer.
Definition at line 159 of file MCVerticesAnalyzer.cc.
References edm::Event::getByLabel(), m_hlumi, m_hmainvtxx, m_hmainvtxy, m_hmainvtxz, m_hnvtx, m_hnvtxvslumi, m_hnvtxweight, m_hnvtxweightprof, m_hpileupvtxz, m_mctruthcollection, m_pileupcollection, m_useweight, m_weight, and histoStyle::weight.
{ using namespace edm; double weight = 1.; if(m_useweight) { Handle<double> weightprod; iEvent.getByLabel(m_weight,weightprod); weight = *weightprod; } Handle<std::vector<PileupSummaryInfo> > pileupinfos; iEvent.getByLabel(m_pileupcollection,pileupinfos); // if(pileupinfos.isValid()) { // look for the intime PileupSummaryInfo std::vector<PileupSummaryInfo>::const_iterator pileupinfo; for(pileupinfo = pileupinfos->begin(); pileupinfo != pileupinfos->end() ; ++pileupinfo) { if(pileupinfo->getBunchCrossing()==0) break; } // if(pileupinfo->getBunchCrossing()!=0) { edm::LogError("NoInTimePileUpInfo") << "Cannot find the in-time pileup info " << pileupinfo->getBunchCrossing(); } else { m_hlumi->Fill(pileupinfo->getTrueNumInteractions(),weight); m_hnvtx->Fill(pileupinfo->getPU_NumInteractions(),weight); m_hnvtxvslumi->Fill(pileupinfo->getTrueNumInteractions(),pileupinfo->getPU_NumInteractions(),weight); if(m_useweight) { m_hnvtxweight->Fill(pileupinfo->getPU_NumInteractions(),1.-weight); m_hnvtxweightprof->Fill(pileupinfo->getPU_NumInteractions(),1.-weight); } const std::vector<float>& zpositions = pileupinfo->getPU_zpositions(); for(std::vector<float>::const_iterator zpos = zpositions.begin() ; zpos != zpositions.end() ; ++zpos) { m_hpileupvtxz->Fill(*zpos,weight); } } } // main interaction part Handle< HepMCProduct > EvtHandle ; iEvent.getByLabel(m_mctruthcollection, EvtHandle ) ; if(EvtHandle.isValid()) { const HepMC::GenEvent* Evt = EvtHandle->GetEvent(); // get the first vertex if(Evt->vertices_begin() != Evt->vertices_end()) { m_hmainvtxx->Fill((*Evt->vertices_begin())->point3d().x()/10.,weight); m_hmainvtxy->Fill((*Evt->vertices_begin())->point3d().y()/10.,weight); m_hmainvtxz->Fill((*Evt->vertices_begin())->point3d().z()/10.,weight); } } }
void MCVerticesAnalyzer::beginJob | ( | void | ) | [private, virtual] |
void MCVerticesAnalyzer::beginRun | ( | const edm::Run & | iRun, |
const edm::EventSetup & | |||
) | [private, virtual] |
void MCVerticesAnalyzer::endJob | ( | void | ) | [private, virtual] |
void MCVerticesAnalyzer::endRun | ( | const edm::Run & | iRun, |
const edm::EventSetup & | |||
) | [private, virtual] |
TH1F* MCVerticesAnalyzer::m_hlumi [private] |
Definition at line 80 of file MCVerticesAnalyzer.cc.
Referenced by analyze(), and MCVerticesAnalyzer().
TH1F* MCVerticesAnalyzer::m_hmainvtxx [private] |
Definition at line 84 of file MCVerticesAnalyzer.cc.
Referenced by analyze(), and MCVerticesAnalyzer().
TH1F* MCVerticesAnalyzer::m_hmainvtxy [private] |
Definition at line 85 of file MCVerticesAnalyzer.cc.
Referenced by analyze(), and MCVerticesAnalyzer().
TH1F* MCVerticesAnalyzer::m_hmainvtxz [private] |
Definition at line 86 of file MCVerticesAnalyzer.cc.
Referenced by analyze(), and MCVerticesAnalyzer().
TH1F* MCVerticesAnalyzer::m_hnvtx [private] |
Definition at line 79 of file MCVerticesAnalyzer.cc.
Referenced by analyze(), and MCVerticesAnalyzer().
TH2F* MCVerticesAnalyzer::m_hnvtxvslumi [private] |
Definition at line 81 of file MCVerticesAnalyzer.cc.
Referenced by analyze(), and MCVerticesAnalyzer().
TH1F* MCVerticesAnalyzer::m_hnvtxweight [private] |
Definition at line 82 of file MCVerticesAnalyzer.cc.
Referenced by analyze(), and MCVerticesAnalyzer().
TProfile* MCVerticesAnalyzer::m_hnvtxweightprof [private] |
Definition at line 83 of file MCVerticesAnalyzer.cc.
Referenced by analyze(), and MCVerticesAnalyzer().
TH1F* MCVerticesAnalyzer::m_hpileupvtxz [private] |
Definition at line 87 of file MCVerticesAnalyzer.cc.
Referenced by analyze(), and MCVerticesAnalyzer().
Definition at line 75 of file MCVerticesAnalyzer.cc.
Referenced by analyze().
Definition at line 74 of file MCVerticesAnalyzer.cc.
Referenced by analyze().
const bool MCVerticesAnalyzer::m_useweight [private] |
Definition at line 76 of file MCVerticesAnalyzer.cc.
Referenced by analyze(), and MCVerticesAnalyzer().
edm::InputTag MCVerticesAnalyzer::m_weight [private] |
Definition at line 77 of file MCVerticesAnalyzer.cc.
Referenced by analyze().