CMS 3D CMS Logo

GenParticleProducer.cc
Go to the documentation of this file.
1 /* \class GenParticleProducer
2  *
3  * \author Luca Lista, INFN
4  *
5  * Convert HepMC GenEvent format into a collection of type
6  * CandidateCollection containing objects of type GenParticle
7  *
8  *
9  */
18 
19 #include <vector>
20 #include <string>
21 #include <unordered_map>
22 
23 namespace edm { class ParameterSet; }
24 namespace HepMC { class GenParticle; class GenEvent; }
25 
26 static constexpr int PDGCacheMax = 32768;
27 
28 namespace {
29  struct IDto3Charge {
30  IDto3Charge(HepPDT::ParticleDataTable const&,
31  bool abortOnUnknownPDGCode);
32 
33  int chargeTimesThree(int) const;
34 
35  private:
36  std::vector<int> chargeP_, chargeM_;
37  std::unordered_map<int, int> chargeMap_;
38  bool abortOnUnknownPDGCode_;
39 
40  };
41 
42  IDto3Charge::IDto3Charge(HepPDT::ParticleDataTable const& iTable,
43  bool iAbortOnUnknownPDGCode):
44  chargeP_( PDGCacheMax, 0 ), chargeM_( PDGCacheMax, 0 ),
45  abortOnUnknownPDGCode_(iAbortOnUnknownPDGCode) {
46  for( auto const& p: iTable ) {
47  const HepPDT::ParticleID & id = p.first;
48  int pdgId = id.pid(), apdgId = std::abs( pdgId );
49  int q3 = id.threeCharge();
50  if ( apdgId < PDGCacheMax && pdgId > 0 ) {
51  chargeP_[ apdgId ] = q3;
52  chargeM_[ apdgId ] = -q3;
53  } else if ( apdgId < PDGCacheMax ) {
54  chargeP_[ apdgId ] = -q3;
55  chargeM_[ apdgId ] = q3;
56  } else {
57  chargeMap_.emplace( pdgId, q3);
58  chargeMap_.emplace( -pdgId, -q3);
59  }
60  }
61  }
62 
63  int IDto3Charge::chargeTimesThree( int id ) const {
64  if( std::abs( id ) < PDGCacheMax )
65  return id > 0 ? chargeP_[ id ] : chargeM_[ - id ];
66  auto f = chargeMap_.find( id );
67  if ( f == chargeMap_.end() ) {
68  if ( abortOnUnknownPDGCode_ )
70  << "invalid PDG id: " << id;
71  else
72  return HepPDT::ParticleID(id).threeCharge();
73  }
74  return f->second;
75 }
76 
77 }
78 
79 class GenParticleProducer : public edm::global::EDProducer<edm::RunCache<IDto3Charge>> {
80  public:
84  ~GenParticleProducer() override;
85 
87  void produce( edm::StreamID, edm::Event& e, const edm::EventSetup&) const override;
88  std::shared_ptr<IDto3Charge> globalBeginRun(const edm::Run&, const edm::EventSetup&) const override;
89  void globalEndRun(edm::Run const&, edm::EventSetup const&) const override {};
90 
91  bool convertParticle(reco::GenParticle& cand, const HepMC::GenParticle * part, const IDto3Charge& id2Charge) const;
92  bool fillDaughters(reco::GenParticleCollection& cand, const HepMC::GenParticle * part, reco::GenParticleRefProd const& ref, size_t index, std::unordered_map<int, size_t>& barcodes) const;
93  bool fillIndices(const HepMC::GenEvent * mc, std::vector<const HepMC::GenParticle*>& particles, std::vector<int>& barCodeVector, int offset,std::unordered_map<int, size_t>& barcodes) const;
94 
95  private:
98  std::vector<edm::EDGetTokenT<edm::HepMCProduct> > vectorSrcTokens_;
100 
105 
108  bool useCF_;
109 
112 
113 };
114 
116 
117 //#include "SimDataFormats/HiGenData/interface/SubEventMap.h"
119 
128 #include <fstream>
129 #include <algorithm>
130 using namespace edm;
131 using namespace reco;
132 using namespace std;
133 using namespace HepMC;
134 
135 static const double mmToCm = 0.1;
136 static const double mmToNs = 1.0/299792458e-6;
137 
139  abortOnUnknownPDGCode_( cfg.getUntrackedParameter<bool>( "abortOnUnknownPDGCode", true ) ),
140  saveBarCodes_( cfg.getUntrackedParameter<bool>( "saveBarCodes", false ) ),
141  doSubEvent_(cfg.getUntrackedParameter<bool>( "doSubEvent", false )),
142  useCF_(cfg.getUntrackedParameter<bool>( "useCrossingFrame", false ))
143 {
144  produces<GenParticleCollection>();
145  produces<math::XYZPointF>("xyz0");
146  produces<float>("t0");
147  if( saveBarCodes_ ) {
148  std::string alias( cfg.getParameter<std::string>( "@module_label" ) );
149  produces<vector<int> >().setBranchAlias( alias + "BarCodes" );
150  }
151 
152  if(useCF_) mixToken_ = mayConsume<CrossingFrame<HepMCProduct> >(InputTag(cfg.getParameter<std::string>( "mix" ),"generatorSmeared"));
153  else srcToken_ = mayConsume<HepMCProduct>(cfg.getParameter<InputTag>( "src" ));
154 }
155 
157 }
158 
159 std::shared_ptr<IDto3Charge> GenParticleProducer::globalBeginRun(const Run&, const EventSetup&es) const {
161  es.getData( pdt );
162  return std::make_shared<IDto3Charge>(*pdt, abortOnUnknownPDGCode_);
163 }
164 
165 void GenParticleProducer::produce( StreamID, Event& evt, const EventSetup& es ) const {
166 
167  std::unordered_map<int, size_t> barcodes;
168 
169  size_t totalSize = 0;
170  const GenEvent * mc = nullptr;
171  MixCollection<HepMCProduct>* cfhepmcprod = nullptr;
172  size_t npiles = 1;
173 
174  if(useCF_){
176  evt.getByToken(mixToken_,cf);
177  cfhepmcprod = new MixCollection<HepMCProduct>(cf.product());
178  npiles = cfhepmcprod->size();
179  LogDebug("GenParticleProducer")<<"npiles : "<<npiles<<endl;
180  for(unsigned int icf = 0; icf < npiles; ++icf){
181  LogDebug("GenParticleProducer")<<"subSize : "<<cfhepmcprod->getObject(icf).GetEvent()->particles_size()<<endl;
182  totalSize += cfhepmcprod->getObject(icf).GetEvent()->particles_size();
183  }
184  LogDebug("GenParticleProducer")<<"totalSize : "<<totalSize<<endl;
185  }else{
187  evt.getByToken( srcToken_, mcp );
188  mc = mcp->GetEvent();
189  if( mc == nullptr )
191  << "HepMC has null pointer to GenEvent" << endl;
192  totalSize = mc->particles_size();
193  }
194 
195  // initialise containers
196  const size_t size = totalSize;
197  vector<const HepMC::GenParticle *> particles( size );
198  auto candsPtr = std::make_unique<GenParticleCollection>(size);
199  auto barCodeVector = std::make_unique<vector<int>>(size);
200  std::unique_ptr<math::XYZPointF> xyz0Ptr(new math::XYZPointF(0.,0.,0.));
201  std::unique_ptr<float> t0Ptr(new float(0.f));
203  GenParticleCollection & cands = * candsPtr;
204  size_t offset = 0;
205  size_t suboffset = 0;
206 
207  IDto3Charge const& id2Charge = *runCache(evt.getRun().index());
209  if(doSubEvent_ || useCF_){
210  for(size_t ipile = 0; ipile < npiles; ++ipile){
211  LogDebug("GenParticleProducer")<<"mixed object ipile : "<<ipile<<endl;
212  barcodes.clear();
213  if(useCF_) mc = cfhepmcprod->getObject(ipile).GetEvent();
214 
215  //Look whether heavy ion/signal event
216  bool isHI = false;
217  const HepMC::HeavyIon * hi = mc->heavy_ion();
218  if(hi && hi->Ncoll_hard() > 1) isHI = true;
219  size_t num_particles = mc->particles_size();
220  LogDebug("GenParticleProducer")<<"num_particles : "<<num_particles<<endl;
221  if (ipile == 0) {
222  auto origin = (*mc->vertices_begin())->position();
223  xyz0Ptr->SetXYZ(origin.x() * mmToCm, origin.y() * mmToCm, origin.z() * mmToCm);
224  *t0Ptr = origin.t() * mmToNs;
225  }
226  fillIndices(mc, particles, *barCodeVector, offset,barcodes);
227  // fill output collection and save association
228  for( size_t ipar = offset; ipar < offset + num_particles; ++ ipar ) {
229  const HepMC::GenParticle * part = particles[ ipar ];
230  reco::GenParticle & cand = cands[ ipar ];
231  // convert HepMC::GenParticle to new reco::GenParticle
232  convertParticle(cand, part,id2Charge);
233  cand.resetDaughters( ref.id() );
234  }
235 
236  for( size_t d = offset; d < offset + num_particles; ++ d ) {
237  const HepMC::GenParticle * part = particles[ d ];
238  const GenVertex * productionVertex = part->production_vertex();
239  int sub_id = 0;
240  if ( productionVertex != nullptr ) {
241  sub_id = productionVertex->id();
242  if(!isHI) sub_id = 0;
243  // search barcode map and attach daughters
244  fillDaughters(cands,part,ref,d,barcodes);
245  }else{
246  const GenVertex * endVertex = part->end_vertex();
247  if(endVertex != nullptr) sub_id = endVertex->id();
248  else throw cms::Exception( "SubEventID" )<<"SubEvent not determined. Particle has no production and no end vertex!"<<endl;
249  }
250  if(sub_id < 0) sub_id = 0;
251  int new_id = sub_id + suboffset;
252  GenParticleRef dref( ref, d );
253  cands[d].setCollisionId(new_id); // For new GenParticle
254  LogDebug("VertexId")<<"SubEvent offset 3 : "<<suboffset;
255  }
256  int nsub = -2;
257  if(isHI){
258  nsub = hi->Ncoll_hard()+1;
259  suboffset += nsub;
260  }else{
261  suboffset += 1;
262  }
263  offset += num_particles;
264  }
265  }else{
266  auto origin = (*mc->vertices_begin())->position();
267  xyz0Ptr->SetXYZ(origin.x() * mmToCm, origin.y() * mmToCm, origin.z() * mmToCm);
268  *t0Ptr = origin.t() * mmToNs;
269  fillIndices(mc, particles, *barCodeVector, 0,barcodes);
270 
271  // fill output collection and save association
272  for( size_t i = 0; i < particles.size(); ++ i ) {
273  const HepMC::GenParticle * part = particles[ i ];
274  reco::GenParticle & cand = cands[ i ];
275  // convert HepMC::GenParticle to new reco::GenParticle
276  convertParticle(cand, part,id2Charge);
277  cand.resetDaughters( ref.id() );
278  }
279 
280  // fill references to daughters
281  for( size_t d = 0; d < cands.size(); ++ d ) {
282  const HepMC::GenParticle * part = particles[ d ];
283  const GenVertex * productionVertex = part->production_vertex();
284  // search barcode map and attach daughters
285  if ( productionVertex != nullptr ) fillDaughters(cands,part,ref,d,barcodes);
286  cands[d].setCollisionId(0);
287  }
288  }
289 
290  evt.put(std::move(candsPtr));
291  if(saveBarCodes_) evt.put(std::move(barCodeVector));
292  if(cfhepmcprod) delete cfhepmcprod;
293  evt.put(std::move(xyz0Ptr),"xyz0");
294  evt.put(std::move(t0Ptr),"t0");
295 }
296 
297 bool GenParticleProducer::convertParticle(reco::GenParticle& cand, const HepMC::GenParticle * part, IDto3Charge const& id2Charge) const {
298  Candidate::LorentzVector p4( part->momentum() );
299  int pdgId = part->pdg_id();
300  cand.setThreeCharge( id2Charge.chargeTimesThree( pdgId ) );
301  cand.setPdgId( pdgId );
302  cand.setStatus( part->status() );
303  cand.setP4( p4 );
304  cand.setCollisionId(0);
305  const GenVertex * v = part->production_vertex();
306  if ( v != nullptr ) {
307  ThreeVector vtx = v->point3d();
308  Candidate::Point vertex( vtx.x() * mmToCm, vtx.y() * mmToCm, vtx.z() * mmToCm );
309  cand.setVertex( vertex );
310  } else {
311  cand.setVertex( Candidate::Point( 0, 0, 0 ) );
312  }
314  return true;
315 }
316 
317 bool GenParticleProducer::fillDaughters(reco::GenParticleCollection& cands, const HepMC::GenParticle * part, reco::GenParticleRefProd const& ref, size_t index, std::unordered_map<int, size_t>& barcodes) const {
318 
319  const GenVertex * productionVertex = part->production_vertex();
320  size_t numberOfMothers = productionVertex->particles_in_size();
321  if ( numberOfMothers > 0 ) {
322  GenVertex::particles_in_const_iterator motherIt = productionVertex->particles_in_const_begin();
323  for( ; motherIt != productionVertex->particles_in_const_end(); motherIt++) {
324  const HepMC::GenParticle * mother = * motherIt;
325  size_t m = barcodes.find( mother->barcode() )->second;
326  cands[ m ].addDaughter( GenParticleRef( ref, index ) );
327  cands[ index ].addMother( GenParticleRef( ref, m ) );
328  }
329  }
330 
331  return true;
332 }
333 
334 bool GenParticleProducer::fillIndices(const HepMC::GenEvent * mc, vector<const HepMC::GenParticle*>& particles, vector<int>& barCodeVector, int offset, std::unordered_map<int, size_t>& barcodes) const {
335  size_t idx = offset;
336  HepMC::GenEvent::particle_const_iterator begin = mc->particles_begin(), end = mc->particles_end();
337  for( HepMC::GenEvent::particle_const_iterator p = begin; p != end; ++ p ) {
338  const HepMC::GenParticle * particle = * p;
339  size_t barCode_this_event = particle->barcode();
340  size_t barCode = barCode_this_event + offset;
341  if( barcodes.find(barCode) != barcodes.end() )
342  throw cms::Exception( "WrongReference" )
343  << "barcodes are duplicated! " << endl;
344  particles[idx] = particle;
345  barCodeVector[idx] = barCode;
346  barcodes.insert( make_pair(barCode_this_event, idx ++) );
347  }
348  return true;
349 }
350 
351 
353 
355 
#define LogDebug(id)
size
Write out results.
GenParticleProducer(const edm::ParameterSet &)
constructor
void produce(edm::StreamID, edm::Event &e, const edm::EventSetup &) const override
process one event
void globalEndRun(edm::Run const &, edm::EventSetup const &) const override
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
T getParameter(std::string const &) const
bool abortOnUnknownPDGCode_
unknown code treatment flag
MCTruthHelper< HepMC::GenParticle > mcTruthHelper_
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:125
void setCollisionId(int s)
Definition: GenParticle.h:38
static const double mmToNs
bool fillDaughters(reco::GenParticleCollection &cand, const HepMC::GenParticle *part, reco::GenParticleRefProd const &ref, size_t index, std::unordered_map< int, size_t > &barcodes) const
MCTruthHelper< reco::GenParticle > mcTruthHelperGenParts_
HepPDT::ParticleDataTable ParticleDataTable
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
std::vector< edm::EDGetTokenT< edm::HepMCProduct > > vectorSrcTokens_
Run const & getRun() const
Definition: Event.cc:99
edm::EDGetTokenT< CrossingFrame< edm::HepMCProduct > > mixToken_
bool fillIndices(const HepMC::GenEvent *mc, std::vector< const HepMC::GenParticle * > &particles, std::vector< int > &barCodeVector, int offset, std::unordered_map< int, size_t > &barcodes) const
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float > > XYZPointF
point in space with cartesian internal representation
Definition: Point3D.h:10
edm::EDGetTokenT< edm::HepMCProduct > srcToken_
source collection name
bool convertParticle(reco::GenParticle &cand, const HepMC::GenParticle *part, const IDto3Charge &id2Charge) const
int size() const
Definition: MixCollection.h:20
void fillGenStatusFlags(const P &p, reco::GenStatusFlags &statusFlags) const
void setVertex(const Point &vertex) override
set vertex
bool getData(T &iHolder) const
Definition: EventSetup.h:111
U second(std::pair< T, U > const &p)
void resetDaughters(const edm::ProductID &id)
set daughters product ID
edm::Ref< edm::HepMCProduct, HepMC::GenParticle > GenParticleRef
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
double p4[4]
Definition: TauolaWrapper.h:92
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
RunIndex index() const
Definition: Run.cc:21
double f[11][100]
static int PDGCacheMax
#define end
Definition: vmac.h:39
RefProd< PROD > getRefBeforePut()
Definition: Event.h:150
const GenStatusFlags & statusFlags() const
Definition: GenParticle.h:41
ProductID id() const
Accessor for product ID.
Definition: RefProd.h:137
const T & getObject(unsigned int ip) const
Definition: MixCollection.h:27
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:38
T const * product() const
Definition: Handle.h:74
void setThreeCharge(Charge qx3) final
set electric charge
Definition: LeafCandidate.h:97
static const double mmToCm
part
Definition: HCALResponse.h:20
~GenParticleProducer() override
destructor
std::shared_ptr< IDto3Charge > globalBeginRun(const edm::Run &, const edm::EventSetup &) const override
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
bool doSubEvent_
input & output modes
fixed size matrix
#define begin
Definition: vmac.h:32
HLT enums.
static int position[264][3]
Definition: ReadPGInfo.cc:509
math::XYZPoint Point
point in the space
Definition: Candidate.h:41
bool saveBarCodes_
save bar-codes
void setStatus(int status) final
set status word
void setPdgId(int pdgId) final
void setP4(const LorentzVector &p4) final
set 4-momentum
def move(src, dest)
Definition: eostools.py:511
#define constexpr
Definition: Run.h:45