CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
GenParticles2HepMCConverter Class Reference
Inheritance diagram for GenParticles2HepMCConverter:
edm::stream::EDProducer<> edm::stream::EDProducerBase edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

 GenParticles2HepMCConverter (const edm::ParameterSet &pset)
 
void produce (edm::Event &event, const edm::EventSetup &eventSetup) override
 
 ~GenParticles2HepMCConverter ()
 
- Public Member Functions inherited from edm::stream::EDProducer<>
 EDProducer ()=default
 
- Public Member Functions inherited from edm::stream::EDProducerBase
 EDProducerBase ()
 
ModuleDescription const & moduleDescription () const
 
 ~EDProducerBase () override
 
- Public Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
std::vector< edm::ProductResolverIndex > const & indiciesForPutProducts (BranchType iBranchType) const
 
 ProducerBase ()
 
std::vector< edm::ProductResolverIndex > const & putTokenIndexToProductResolverIndex () const
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription const &)> registrationCallback () const
 used by the fwk to register list of products More...
 
void resolvePutIndicies (BranchType iBranchType, ModuleToResolverIndicies const &iIndicies, std::string const &moduleLabel)
 
virtual ~ProducerBase () noexcept(false)
 
- Public Member Functions inherited from edm::EDConsumerBase
std::vector< ConsumesInfoconsumesInfo () const
 
void convertCurrentProcessAlias (std::string const &processName)
 Convert "@currentProcess" in InputTag process names to the actual current process name. More...
 
 EDConsumerBase ()
 
 EDConsumerBase (EDConsumerBase const &)=delete
 
 EDConsumerBase (EDConsumerBase &&)=default
 
ProductResolverIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductResolverIndexAndSkipBit > &) const
 
std::vector< ProductResolverIndexAndSkipBit > const & itemsToGetFrom (BranchType iType) const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
void modulesWhoseProductsAreConsumed (std::vector< ModuleDescription const * > &modules, ProductRegistry const &preg, std::map< std::string, ModuleDescription const * > const &labelsToDesc, std::string const &processName) const
 
EDConsumerBase const & operator= (EDConsumerBase const &)=delete
 
EDConsumerBaseoperator= (EDConsumerBase &&)=default
 
bool registeredToConsume (ProductResolverIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductResolverIndexHelper const &, bool iPrefetchMayGet)
 
virtual ~EDConsumerBase () noexcept(false)
 

Private Member Functions

HepMC::FourVector FourVector (const reco::Candidate::Point &point)
 
HepMC::FourVector FourVector (const reco::Candidate::LorentzVector &lvec)
 

Private Attributes

const double cmEnergy_
 
edm::EDGetTokenT< GenEventInfoProductgenEventInfoToken_
 
edm::EDGetTokenT< reco::CandidateViewgenParticlesToken_
 
edm::ESHandle< ParticleDataTablepTable_
 
std::vector< int > signalParticlePdgIds_
 

Additional Inherited Members

- Public Types inherited from edm::stream::EDProducer<>
typedef CacheContexts< T... > CacheTypes
 
typedef CacheTypes::GlobalCache GlobalCache
 
typedef AbilityChecker< T... > HasAbility
 
typedef CacheTypes::LuminosityBlockCache LuminosityBlockCache
 
typedef LuminosityBlockContextT< LuminosityBlockCache, RunCache, GlobalCacheLuminosityBlockContext
 
typedef CacheTypes::LuminosityBlockSummaryCache LuminosityBlockSummaryCache
 
typedef CacheTypes::RunCache RunCache
 
typedef RunContextT< RunCache, GlobalCacheRunContext
 
typedef CacheTypes::RunSummaryCache RunSummaryCache
 
- Public Types inherited from edm::stream::EDProducerBase
typedef EDProducerAdaptorBase ModuleType
 
- Public Types inherited from edm::ProducerBase
using ModuleToResolverIndicies = std::unordered_multimap< std::string, std::tuple< edm::TypeID const *, const char *, edm::ProductResolverIndex >>
 
typedef ProductRegistryHelper::TypeLabelList TypeLabelList
 
- Public Types inherited from edm::EDConsumerBase
typedef ProductLabels Labels
 
- Static Public Member Functions inherited from edm::stream::EDProducerBase
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &descriptions)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Definition at line 27 of file GenParticles2HepMCConverter.cc.

Constructor & Destructor Documentation

GenParticles2HepMCConverter::GenParticles2HepMCConverter ( const edm::ParameterSet pset)
explicit

Definition at line 58 of file GenParticles2HepMCConverter.cc.

References genEventInfoToken_, genParticlesToken_, edm::ParameterSet::getParameter(), and signalParticlePdgIds_.

58  :
59  cmEnergy_(pset.getUntrackedParameter<double>("cmEnergy", 13000)) // dummy value to set incident proton pz for particle gun samples
60 {
61  genParticlesToken_ = consumes<reco::CandidateView>(pset.getParameter<edm::InputTag>("genParticles"));
62  genEventInfoToken_ = consumes<GenEventInfoProduct>(pset.getParameter<edm::InputTag>("genEventInfo"));
63  signalParticlePdgIds_ = pset.getParameter<std::vector<int>>("signalParticlePdgIds");
64 
65  produces<edm::HepMCProduct>("unsmeared");
66 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< GenEventInfoProduct > genEventInfoToken_
edm::EDGetTokenT< reco::CandidateView > genParticlesToken_
GenParticles2HepMCConverter::~GenParticles2HepMCConverter ( )
inline

Definition at line 31 of file GenParticles2HepMCConverter.cc.

31 {};

Member Function Documentation

HepMC::FourVector GenParticles2HepMCConverter::FourVector ( const reco::Candidate::Point point)
inlineprivate

Definition at line 44 of file GenParticles2HepMCConverter.cc.

Referenced by produce().

45  {
46  return HepMC::FourVector(10*point.x(), 10*point.y(), 10*point.z(), 0);
47  };
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
HepMC::FourVector GenParticles2HepMCConverter::FourVector ( const reco::Candidate::LorentzVector lvec)
inlineprivate

Definition at line 49 of file GenParticles2HepMCConverter.cc.

References SiStripPI::max.

50  {
51  // Avoid negative mass, set minimum m^2 = 0
52  return HepMC::FourVector(lvec.px(), lvec.py(), lvec.pz(), std::hypot(lvec.P(), std::max(0., lvec.mass())));
53  };
void GenParticles2HepMCConverter::produce ( edm::Event event,
const edm::EventSetup eventSetup 
)
override

Definition at line 68 of file GenParticles2HepMCConverter.cc.

References funct::abs(), GenEventInfoProduct::alphaQCD(), GenEventInfoProduct::alphaQED(), edm::View< T >::at(), class-composition::children, cmEnergy_, reco::Candidate::daughter(), DEFINE_FWK_MODULE, reco::Candidate::eta(), edm::EventID::event(), spr::find(), FourVector(), genEventInfoToken_, GenParticle::GenParticle, genParticlesToken_, edm::EventSetup::getData(), mps_fire::i, gen::PdfInfo::id, edm::EventBase::id(), ResonanceBuilder::mass, reco::Candidate::mass(), reco::Candidate::mother(), eostools::move(), gen::n, reco::Candidate::numberOfMothers(), or, AlCaHLTBitMon_ParallelJobs::p, reco::Candidate::p4(), GenEventInfoProduct::pdf(), common_cff::pdgId, reco::Candidate::pdgId(), pTable_, reco::Candidate::pz(), GenEventInfoProduct::qScale(), gen::PdfInfo::scalePDF, signalParticlePdgIds_, GenEventInfoProduct::signalProcessID(), edm::View< T >::size(), reco::Candidate::status(), findQualityFiles::v, reco::Candidate::vertex(), GenEventInfoProduct::weights(), gen::PdfInfo::x, and gen::PdfInfo::xPDF.

69 {
70  edm::Handle<reco::CandidateView> genParticlesHandle;
71  event.getByToken(genParticlesToken_, genParticlesHandle);
72 
73  edm::Handle<GenEventInfoProduct> genEventInfoHandle;
74  event.getByToken(genEventInfoToken_, genEventInfoHandle);
75 
76  eventSetup.getData(pTable_);
77 
78  HepMC::GenEvent* hepmc_event = new HepMC::GenEvent();
79  hepmc_event->set_event_number(event.id().event());
80  hepmc_event->set_signal_process_id(genEventInfoHandle->signalProcessID());
81  hepmc_event->set_event_scale(genEventInfoHandle->qScale());
82  hepmc_event->set_alphaQED(genEventInfoHandle->alphaQED());
83  hepmc_event->set_alphaQCD(genEventInfoHandle->alphaQCD());
84 
85  hepmc_event->weights() = genEventInfoHandle->weights();
86 
87  // Set PDF
88  const gen::PdfInfo* pdf = genEventInfoHandle->pdf();
89  if ( pdf != nullptr ) {
90  const int pdf_id1 = pdf->id.first, pdf_id2 = pdf->id.second;
91  const double pdf_x1 = pdf->x.first, pdf_x2 = pdf->x.second;
92  const double pdf_scalePDF = pdf->scalePDF;
93  const double pdf_xPDF1 = pdf->xPDF.first, pdf_xPDF2 = pdf->xPDF.second;
94  HepMC::PdfInfo hepmc_pdfInfo(pdf_id1, pdf_id2, pdf_x1, pdf_x2, pdf_scalePDF, pdf_xPDF1, pdf_xPDF2);
95  hepmc_event->set_pdf_info(hepmc_pdfInfo);
96  }
97 
98  // Prepare list of HepMC::GenParticles
99  std::map<const reco::Candidate*, HepMC::GenParticle*> genCandToHepMCMap;
100  HepMC::GenParticle* hepmc_parton1 = nullptr, * hepmc_parton2 = nullptr;
101  std::vector<HepMC::GenParticle*> hepmc_particles;
102  const reco::Candidate* parton1 = nullptr, * parton2 = nullptr;
103  for ( unsigned int i=0, n=genParticlesHandle->size(); i<n; ++i )
104  {
105  const reco::Candidate* p = &genParticlesHandle->at(i);
106  HepMC::GenParticle* hepmc_particle = new HepMC::GenParticle(FourVector(p->p4()), p->pdgId(), p->status());
107  hepmc_particle->suggest_barcode(i+1);
108 
109  // Assign particle's generated mass from the standard particle data table
110  double particleMass;
111  if ( pTable_->particle(p->pdgId()) ) particleMass = pTable_->particle(p->pdgId())->mass();
112  else particleMass = p->mass();
113 
114  hepmc_particle->set_generated_mass(particleMass);
115 
116  hepmc_particles.push_back(hepmc_particle);
117  genCandToHepMCMap[p] = hepmc_particle;
118 
119  // Find incident proton pair
120  if ( p->mother() == nullptr and std::abs(p->eta()) > 5 and std::abs(p->pz()) > 1000 ) {
121  if ( !parton1 and p->pz() > 0 ) {
122  parton1 = p;
123  hepmc_parton1 = hepmc_particle;
124  }
125  else if ( !parton2 and p->pz() < 0 ) {
126  parton2 = p;
127  hepmc_parton2 = hepmc_particle;
128  }
129  }
130  }
131 
132  HepMC::GenVertex* vertex1 = nullptr;
133  HepMC::GenVertex* vertex2 = nullptr;
134  if ( parton1 == nullptr and parton2 == nullptr ) {
135  // Particle gun samples do not have incident partons. Put dummy incident particle and prod vertex
136  // Note: leave parton1 and parton2 as nullptr since it is not used anymore after creating hepmc_parton1 and 2
137  const reco::Candidate::LorentzVector nullP4(0,0,0,0);
138  const reco::Candidate::LorentzVector beamP4(0,0,cmEnergy_/2, cmEnergy_/2);
139  vertex1 = new HepMC::GenVertex(FourVector(nullP4));
140  vertex2 = new HepMC::GenVertex(FourVector(nullP4));
141  hepmc_parton1 = new HepMC::GenParticle(FourVector(+beamP4), 2212, 4);
142  hepmc_parton2 = new HepMC::GenParticle(FourVector(-beamP4), 2212, 4);
143  }
144  else {
145  // Put incident beam particles : proton -> parton vertex
146  vertex1 = new HepMC::GenVertex(FourVector(parton1->vertex()));
147  vertex2 = new HepMC::GenVertex(FourVector(parton2->vertex()));
148  }
149  hepmc_event->add_vertex(vertex1);
150  hepmc_event->add_vertex(vertex2);
151  vertex1->add_particle_in(hepmc_parton1);
152  vertex2->add_particle_in(hepmc_parton2);
153  //hepmc_event->set_beam_particles(hepmc_parton1, hepmc_parton2);
154 
155  // Prepare vertex list
156  typedef std::map<const reco::Candidate*, HepMC::GenVertex*> ParticleToVertexMap;
157  ParticleToVertexMap particleToVertexMap;
158  particleToVertexMap[parton1] = vertex1;
159  particleToVertexMap[parton2] = vertex2;
160  for ( unsigned int i=0, n=genParticlesHandle->size(); i<n; ++i )
161  {
162  const reco::Candidate* p = &genParticlesHandle->at(i);
163  if ( p == parton1 or p == parton2 ) continue;
164 
165  // Connect mother-daughters for the other cases
166  for ( unsigned int j=0, nMothers=p->numberOfMothers(); j<nMothers; ++j )
167  {
168  // Mother-daughter hierarchy defines vertex
169  const reco::Candidate* elder = p->mother(j)->daughter(0);
170  HepMC::GenVertex* vertex;
171  if ( particleToVertexMap.find(elder) == particleToVertexMap.end() )
172  {
173  vertex = new HepMC::GenVertex(FourVector(elder->vertex()));
174  hepmc_event->add_vertex(vertex);
175  particleToVertexMap[elder] = vertex;
176  }
177  else
178  {
179  vertex = particleToVertexMap[elder];
180  }
181 
182  // Vertex is found. Now connect each other
183  const reco::Candidate* mother = p->mother(j);
184  vertex->add_particle_in(genCandToHepMCMap[mother]);
185  vertex->add_particle_out(hepmc_particles[i]);
186  }
187  }
188 
189  // Finalize HepMC event record
190  bool hasSignalVertex = false;
191  if ( !signalParticlePdgIds_.empty() ) {
192  // Loop over all vertices to assign the signal vertex, decaying to a signal particle
193  for ( auto v = hepmc_event->vertices_begin(); v != hepmc_event->vertices_end(); ++v ) {
194  for ( auto p = (*v)->particles_begin(HepMC::children);
195  p != (*v)->particles_end(HepMC::children); ++p ) {
196  const int pdgId = (*p)->pdg_id();
198  hepmc_event->set_signal_process_vertex(*v);
199  hasSignalVertex = true;
200  break;
201  }
202  }
203  if ( hasSignalVertex ) break;
204  }
205  }
206  // Set the default signal vertex if still not set
207  if ( !hasSignalVertex ) hepmc_event->set_signal_process_vertex(*(vertex1->vertices_begin()));
208 
209  std::unique_ptr<edm::HepMCProduct> hepmc_product(new edm::HepMCProduct());
210  hepmc_product->addHepMCData(hepmc_event);
211  event.put(std::move(hepmc_product), "unsmeared");
212 
213 }
EventNumber_t event() const
Definition: EventID.h:41
virtual double pz() const =0
z coordinate of momentum vector
edm::EDGetTokenT< GenEventInfoProduct > genEventInfoToken_
const PDF * pdf() const
double alphaQCD() const
virtual const Candidate * daughter(size_type i) const =0
return daughter at a given position, i = 0, ... numberOfDaughters() - 1 (read only mode) ...
std::pair< double, double > x
Definition: PdfInfo.h:13
edm::ESHandle< ParticleDataTable > pTable_
size_type size() const
virtual const Candidate * mother(size_type i=0) const =0
return pointer to mother
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
HepMC::FourVector FourVector(const reco::Candidate::Point &point)
edm::EDGetTokenT< reco::CandidateView > genParticlesToken_
void getData(T &iHolder) const
Definition: EventSetup.h:78
std::pair< double, double > xPDF
Definition: PdfInfo.h:14
virtual int status() const =0
status word
double qScale() const
virtual const LorentzVector & p4() const =0
four-momentum Lorentz vector
virtual int pdgId() const =0
PDG identifier.
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
virtual size_type numberOfMothers() const =0
number of mothers (zero or one in most of but not all the cases)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double alphaQED() const
std::pair< int, int > id
Definition: PdfInfo.h:12
virtual double eta() const =0
momentum pseudorapidity
virtual double mass() const =0
mass
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
edm::EventID id() const
Definition: EventBase.h:60
std::vector< double > & weights()
const_reference at(size_type pos) const
virtual const Point & vertex() const =0
vertex position
unsigned int signalProcessID() const
double scalePDF
Definition: PdfInfo.h:15
def move(src, dest)
Definition: eostools.py:510

Member Data Documentation

const double GenParticles2HepMCConverter::cmEnergy_
private

Definition at line 41 of file GenParticles2HepMCConverter.cc.

Referenced by produce().

edm::EDGetTokenT<GenEventInfoProduct> GenParticles2HepMCConverter::genEventInfoToken_
private

Definition at line 37 of file GenParticles2HepMCConverter.cc.

Referenced by GenParticles2HepMCConverter(), and produce().

edm::EDGetTokenT<reco::CandidateView> GenParticles2HepMCConverter::genParticlesToken_
private

Definition at line 36 of file GenParticles2HepMCConverter.cc.

Referenced by GenParticles2HepMCConverter(), and produce().

edm::ESHandle<ParticleDataTable> GenParticles2HepMCConverter::pTable_
private

Definition at line 38 of file GenParticles2HepMCConverter.cc.

Referenced by produce().

std::vector<int> GenParticles2HepMCConverter::signalParticlePdgIds_
private

Definition at line 40 of file GenParticles2HepMCConverter.cc.

Referenced by GenParticles2HepMCConverter(), and produce().