#include <CFWriter.h>
Public Member Functions | |
void | beginJob () |
void | beginRun (const edm::Run &run, const edm::EventSetup &es) |
CFWriter (const edm::ParameterSet &conf) | |
virtual void | produce (edm::Event &e, const edm::EventSetup &c) |
virtual void | put (edm::Event &e) |
virtual | ~CFWriter () |
Private Types | |
typedef std::vector < edm::HepMCProduct > | HepMCProductContainer |
Private Member Functions | |
virtual void | branchesActivate (const std::string &friendlyName, std::string subdet, InputTag &tag, std::string &label) |
PCrossingFrame< SimTrack > | fctTest (PCrossingFrame< SimTrack > p) |
Private Attributes | |
bool | flagHepMCProduct_ |
bool | flagPCaloHit_ |
bool | flagPSimHit_ |
bool | flagSimTrack_ |
bool | flagSimVertex_ |
std::vector< std::string > | labCaloHit |
std::vector< std::string > | labSimHit |
bool | useCurrentProcessOnly_ |
std::vector< std::string > | wantedBranches_ |
Definition at line 31 of file CFWriter.h.
typedef std::vector<edm::HepMCProduct> edm::CFWriter::HepMCProductContainer [private] |
Definition at line 56 of file CFWriter.h.
CFWriter::CFWriter | ( | const edm::ParameterSet & | conf | ) | [explicit] |
Definition at line 38 of file CFWriter.cc.
References branchesActivate(), edm::InputTag::encode(), edm::ParameterSet::exists(), flagHepMCProduct_, flagPCaloHit_, flagPSimHit_, flagSimTrack_, flagSimVertex_, edm::ParameterSet::getParameter(), edm::ParameterSet::getParameterNames(), labCaloHit, label, labSimHit, h::names, subdets, GlobalPosition_Frontier_DevDB_cff::tag, and o2o::tags.
: flagSimTrack_(false), flagSimVertex_(false), flagHepMCProduct_(false), flagPCaloHit_(false), flagPSimHit_(false) { //register your products ParameterSet ps=iConfig.getParameter<ParameterSet>("mixObjects"); std::vector<std::string> names = ps.getParameterNames(); for (std::vector<std::string>::iterator it=names.begin();it!= names.end();++it) { ParameterSet pset=ps.getParameter<ParameterSet>((*it)); if (!pset.exists("type")) continue; //to allow replacement by empty pset std::string object = pset.getParameter<std::string>("type"); std::vector<InputTag> tags=pset.getParameter<std::vector<InputTag> >("input"); //SimTracks if (object=="SimTrack") { flagSimTrack_ = true; InputTag tag; if (tags.size()>0) tag=tags[0]; std::string label; branchesActivate(TypeID(typeid(std::vector<SimTrack>)).friendlyClassName(),std::string(""),tag,label); produces<PCrossingFrame<SimTrack> >(label); LogInfo("MixingModule") <<"Add PCrossingFrame<SimTrack> "<<object<<"s with InputTag= "<<tag.encode()<<", label will be "<<label; } //SimVertices else if (object=="SimVertex") { flagSimVertex_ = true; InputTag tag; if (tags.size()>0) tag=tags[0]; std::string label; branchesActivate(TypeID(typeid(std::vector<SimVertex>)).friendlyClassName(),std::string(""),tag,label); produces<PCrossingFrame<SimVertex> >(label); LogInfo("MixingModule") <<"Add SimVertexContainer "<<object<<"s with InputTag= "<<tag.encode()<<", label will be "<<label; } // PCaloHit else if (object=="PCaloHit"){ flagPCaloHit_ = true; std::vector<std::string> subdets=pset.getParameter<std::vector<std::string> >("subdets"); for (unsigned int ii=0;ii<subdets.size();ii++) { InputTag tag; if (tags.size()==1) tag=tags[0]; else if(tags.size()>1) tag=tags[ii]; std::string label; branchesActivate(TypeID(typeid(std::vector<PCaloHit>)).friendlyClassName(),subdets[ii],tag,label); produces<PCrossingFrame<PCaloHit> >(label); LogInfo("MixingModule") <<"Add PCrossingFrame<PCaloHit> "<<object<<"s with InputTag= "<<tag.encode()<<", label will be "<<label; // fill table with labels labCaloHit.push_back(label); } } // PSimHit else if (object=="PSimHit"){ flagPSimHit_ = true; std::vector<std::string> subdets=pset.getParameter<std::vector<std::string> >("subdets"); for (unsigned int ii=0;ii<subdets.size();ii++) { InputTag tag; if (tags.size()==1) tag=tags[0]; else if(tags.size()>1) tag=tags[ii]; std::string label; branchesActivate(TypeID(typeid(std::vector<PSimHit>)).friendlyClassName(),subdets[ii],tag,label); produces<PCrossingFrame<PSimHit> >(label); LogInfo("MixingModule") <<"Add PSimHitContainer "<<object<<"s with InputTag= "<<tag.encode()<<", label will be "<<label; // fill table with labels labSimHit.push_back(label); } //end for } // HepMCProduct else if (object=="HepMCProduct"){ flagHepMCProduct_ = true; InputTag tag; if (tags.size()>0) tag=tags[0]; std::string label; branchesActivate(TypeID(typeid(HepMCProduct)).friendlyClassName(),std::string(""),tag,label); produces<PCrossingFrame<edm::HepMCProduct> >(label); LogInfo("MixingModule") <<"Add HepMCProduct "<<object<<"s with InputTag= "<<tag.encode()<<", label will be "<<label; } else LogWarning("MixingModule") <<"You did not mix a type of object("<<object<<")."; }//end for }
CFWriter::~CFWriter | ( | ) | [virtual] |
Definition at line 145 of file CFWriter.cc.
{}
void edm::CFWriter::beginJob | ( | void | ) | [inline, virtual] |
void CFWriter::beginRun | ( | const edm::Run & | run, |
const edm::EventSetup & | es | ||
) |
Definition at line 140 of file CFWriter.cc.
{ }
void CFWriter::branchesActivate | ( | const std::string & | friendlyName, |
std::string | subdet, | ||
InputTag & | tag, | ||
std::string & | label | ||
) | [private, virtual] |
Definition at line 241 of file CFWriter.cc.
References EgammaHLTValidationUtils::getProcessName(), edm::InputTag::instance(), edm::InputTag::label(), useCurrentProcessOnly_, and wantedBranches_.
Referenced by CFWriter().
{ label=tag.label()+tag.instance(); wantedBranches_.push_back(friendlyName + '_' + tag.label() + '_' + tag.instance()); //if useCurrentProcessOnly, we have to change the input tag if (useCurrentProcessOnly_) { const std::string processName = edm::Service<edm::service::TriggerNamesService>()->getProcessName(); tag = InputTag(tag.label(),tag.instance(),processName); } }
PCrossingFrame<SimTrack> edm::CFWriter::fctTest | ( | PCrossingFrame< SimTrack > | p | ) | [inline, private] |
Definition at line 45 of file CFWriter.h.
References gather_cfg::cout, and AlCaHLTBitMon_ParallelJobs::p.
void CFWriter::produce | ( | edm::Event & | e, |
const edm::EventSetup & | c | ||
) | [virtual] |
Implements edm::EDProducer.
Definition at line 148 of file CFWriter.cc.
References flagHepMCProduct_, flagPCaloHit_, flagPSimHit_, flagSimTrack_, flagSimVertex_, edm::Event::getByLabel(), labCaloHit, labSimHit, edm::Handle< T >::product(), and edm::Event::put().
{ if (flagSimTrack_){ bool gotTracks; edm::Handle<CrossingFrame<SimTrack> > cf_simtrack; gotTracks = iEvent.getByLabel("mix","g4SimHits",cf_simtrack); if (gotTracks){ PCrossingFrame<SimTrack> * PCFbis = new PCrossingFrame<SimTrack>(*cf_simtrack.product()); std::auto_ptr<PCrossingFrame<SimTrack> > pOutTrack(PCFbis); iEvent.put(pOutTrack,"g4SimHits"); } else{ LogInfo("MixingModule") << " Please, check if the object <SimTrack> has been mixed by the MixingModule!"; } }//end if flagSimTrack_ //SimVertex if (flagSimVertex_){ bool gotSimVertex; edm::Handle<CrossingFrame<SimVertex> > cf_simvtx; gotSimVertex=iEvent.getByLabel("mix","g4SimHits",cf_simvtx); if (gotSimVertex){ PCrossingFrame<SimVertex> * PCFvtx = new PCrossingFrame<SimVertex>(*cf_simvtx.product()); std::auto_ptr<PCrossingFrame<SimVertex> > pOutVertex(PCFvtx); iEvent.put(pOutVertex,"g4SimHits"); } else{ LogInfo("MixingModule") << " Please, check if the object <SimVertex> has been mixed by the MixingModule!"; } } // PCaloHit if (flagPCaloHit_){ for (unsigned int ii=0;ii<labCaloHit.size();ii++){ bool gotPCaloHit; edm::Handle<CrossingFrame<PCaloHit> > cf_calohit; gotPCaloHit=iEvent.getByLabel("mix",labCaloHit[ii],cf_calohit); if (gotPCaloHit){ PCrossingFrame<PCaloHit> * PCFPhCaloHit = new PCrossingFrame<PCaloHit>(*cf_calohit.product()); std::auto_ptr<PCrossingFrame<PCaloHit> > pOutHCalo(PCFPhCaloHit); iEvent.put(pOutHCalo,labCaloHit[ii]); } else{ LogInfo("MixingModule") << " Please, check if the object <PCaloHit> " << labCaloHit[ii] << " has been mixed by the MixingModule!"; } } } if (flagPSimHit_){ for (unsigned int ii=0;ii<labSimHit.size();ii++) { bool gotPSimHit; edm::Handle<CrossingFrame<PSimHit> > cf_simhit; gotPSimHit=iEvent.getByLabel("mix",labSimHit[ii],cf_simhit); if (gotPSimHit){ PCrossingFrame<PSimHit> * PCFSimHit = new PCrossingFrame<PSimHit>(*cf_simhit.product()); std::auto_ptr<PCrossingFrame<PSimHit> > pOutSimHit(PCFSimHit); iEvent.put(pOutSimHit,labSimHit[ii]); } else{ LogInfo("MixingModule") << " Please, check if the object <PSimHit> " << labSimHit[ii] << " has been mixed by the MixingModule!"; } } } //HepMCProduct if (flagHepMCProduct_){ bool gotHepMCProduct; edm::Handle<CrossingFrame<edm::HepMCProduct> > cf_hepmc; gotHepMCProduct=iEvent.getByLabel("mix","generator",cf_hepmc); if (gotHepMCProduct){ PCrossingFrame<edm::HepMCProduct> * PCFHepMC = new PCrossingFrame<edm::HepMCProduct>(*cf_hepmc.product()); std::auto_ptr<PCrossingFrame<edm::HepMCProduct> > pOuthepmcpr(PCFHepMC); iEvent.put(pOuthepmcpr,"generator"); } else{ LogInfo("MixingModule") << " Please, check if the object <HepMCProduct> has been mixed by the MixingModule!"; } }// end if flagHepMCProduct_ }
virtual void edm::CFWriter::put | ( | edm::Event & | e | ) | [inline, virtual] |
Definition at line 42 of file CFWriter.h.
{;}
bool edm::CFWriter::flagHepMCProduct_ [private] |
Definition at line 52 of file CFWriter.h.
Referenced by CFWriter(), and produce().
bool edm::CFWriter::flagPCaloHit_ [private] |
Definition at line 53 of file CFWriter.h.
Referenced by CFWriter(), and produce().
bool edm::CFWriter::flagPSimHit_ [private] |
Definition at line 54 of file CFWriter.h.
Referenced by CFWriter(), and produce().
bool edm::CFWriter::flagSimTrack_ [private] |
Definition at line 50 of file CFWriter.h.
Referenced by CFWriter(), and produce().
bool edm::CFWriter::flagSimVertex_ [private] |
Definition at line 51 of file CFWriter.h.
Referenced by CFWriter(), and produce().
std::vector<std::string> edm::CFWriter::labCaloHit [private] |
Definition at line 58 of file CFWriter.h.
Referenced by CFWriter(), and produce().
std::vector<std::string> edm::CFWriter::labSimHit [private] |
Definition at line 57 of file CFWriter.h.
Referenced by CFWriter(), and produce().
bool edm::CFWriter::useCurrentProcessOnly_ [private] |
Definition at line 49 of file CFWriter.h.
Referenced by branchesActivate().
std::vector<std::string> edm::CFWriter::wantedBranches_ [private] |
Definition at line 48 of file CFWriter.h.
Referenced by branchesActivate().