CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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  */
17 
18 #include <vector>
19 #include <string>
20 #include <map>
21 #include <set>
22 
23 namespace edm { class ParameterSet; }
24 namespace HepMC { class GenParticle; class GenEvent; }
25 
27  public:
32 
34  virtual void produce( edm::Event& e, const edm::EventSetup&) override;
35  int chargeTimesThree( int ) const;
37  bool fillDaughters(reco::GenParticleCollection& cand, const HepMC::GenParticle * part, size_t index);
38  bool fillIndices(const HepMC::GenEvent * mc, std::vector<const HepMC::GenParticle*>& particles, std::vector<int>& barCodeVector, int offset);
39  std::map<int, size_t> barcodes_;
41 
42  private:
45  std::vector<edm::EDGetTokenT<edm::HepMCProduct> > vectorSrcTokens_;
47 
55  std::vector<int> chargeP_, chargeM_;
56  std::map<int, int> chargeMap_;
57 
60  bool useCF_;
61 
64 
65 };
66 
69 
70 //#include "SimDataFormats/HiGenData/interface/SubEventMap.h"
72 
80 #include <fstream>
81 #include <algorithm>
82 using namespace edm;
83 using namespace reco;
84 using namespace std;
85 using namespace HepMC;
86 
87 static const int PDGCacheMax = 32768;
88 static const double mmToCm = 0.1;
89 static const double mmToNs = 1.0/299792458e-6;
90 
92  firstEvent_(true),
93  abortOnUnknownPDGCode_( cfg.getUntrackedParameter<bool>( "abortOnUnknownPDGCode", true ) ),
94  saveBarCodes_( cfg.getUntrackedParameter<bool>( "saveBarCodes", false ) ),
95  chargeP_( PDGCacheMax, 0 ), chargeM_( PDGCacheMax, 0 ),
96  doSubEvent_(cfg.getUntrackedParameter<bool>( "doSubEvent", false )),
97  useCF_(cfg.getUntrackedParameter<bool>( "useCrossingFrame", false ))
98 {
99  produces<GenParticleCollection>();
100  produces<math::XYZPointF>("xyz0");
101  produces<float>("t0");
102  if( saveBarCodes_ ) {
103  std::string alias( cfg.getParameter<std::string>( "@module_label" ) );
104  produces<vector<int> >().setBranchAlias( alias + "BarCodes" );
105  }
106 
107  if(useCF_) mixToken_ = mayConsume<CrossingFrame<HepMCProduct> >(InputTag(cfg.getParameter<std::string>( "mix" ),"generatorSmeared"));
108  else srcToken_ = mayConsume<HepMCProduct>(cfg.getParameter<InputTag>( "src" ));
109 }
110 
112 }
113 
115  if( std::abs( id ) < PDGCacheMax )
116  return id > 0 ? chargeP_[ id ] : chargeM_[ - id ];
117  map<int, int>::const_iterator f = chargeMap_.find( id );
118  if ( f == chargeMap_.end() ) {
121  << "invalid PDG id: " << id << endl;
122  else
123  return HepPDT::ParticleID(id).threeCharge();
124  }
125  return f->second;
126 }
127 
129 
130  if (firstEvent_) {
132  es.getData( pdt );
133  for( HepPDT::ParticleDataTable::const_iterator p = pdt->begin(); p != pdt->end(); ++ p ) {
134  const HepPDT::ParticleID & id = p->first;
135  int pdgId = id.pid(), apdgId = std::abs( pdgId );
136  int q3 = id.threeCharge();
137  if ( apdgId < PDGCacheMax && pdgId > 0 ) {
138  chargeP_[ apdgId ] = q3;
139  chargeM_[ apdgId ] = -q3;
140  } else if ( apdgId < PDGCacheMax ) {
141  chargeP_[ apdgId ] = -q3;
142  chargeM_[ apdgId ] = q3;
143  } else {
144  chargeMap_[ pdgId ] = q3;
145  chargeMap_[ -pdgId ] = -q3;
146  }
147  }
148  firstEvent_ = false;
149  }
150 
151  barcodes_.clear();
152 
153  size_t totalSize = 0;
154  const GenEvent * mc = 0;
155  MixCollection<HepMCProduct>* cfhepmcprod = 0;
156  size_t npiles = 1;
157 
158  if(useCF_){
160  evt.getByToken(mixToken_,cf);
161  cfhepmcprod = new MixCollection<HepMCProduct>(cf.product());
162  npiles = cfhepmcprod->size();
163  LogDebug("GenParticleProducer")<<"npiles : "<<npiles<<endl;
164  for(unsigned int icf = 0; icf < npiles; ++icf){
165  LogDebug("GenParticleProducer")<<"subSize : "<<cfhepmcprod->getObject(icf).GetEvent()->particles_size()<<endl;
166  totalSize += cfhepmcprod->getObject(icf).GetEvent()->particles_size();
167  }
168  LogDebug("GenParticleProducer")<<"totalSize : "<<totalSize<<endl;
169  }else{
171  evt.getByToken( srcToken_, mcp );
172  mc = mcp->GetEvent();
173  if( mc == 0 )
175  << "HepMC has null pointer to GenEvent" << endl;
176  totalSize = mc->particles_size();
177  }
178 
179  // initialise containers
180  const size_t size = totalSize;
181  vector<const HepMC::GenParticle *> particles( size );
182  auto candsPtr = std::make_unique<GenParticleCollection>(size);
183  auto barCodeVector = std::make_unique<vector<int>>(size);
184  std::unique_ptr<math::XYZPointF> xyz0Ptr(new math::XYZPointF(0.,0.,0.));
185  std::unique_ptr<float> t0Ptr(new float(0.f));
187  GenParticleCollection & cands = * candsPtr;
188  size_t offset = 0;
189  size_t suboffset = 0;
190 
192  if(doSubEvent_ || useCF_){
193  for(size_t ipile = 0; ipile < npiles; ++ipile){
194  LogDebug("GenParticleProducer")<<"mixed object ipile : "<<ipile<<endl;
195  barcodes_.clear();
196  if(useCF_) mc = cfhepmcprod->getObject(ipile).GetEvent();
197 
198  //Look whether heavy ion/signal event
199  bool isHI = false;
200  const HepMC::HeavyIon * hi = mc->heavy_ion();
201  if(hi && hi->Ncoll_hard() > 1) isHI = true;
202  size_t num_particles = mc->particles_size();
203  LogDebug("GenParticleProducer")<<"num_particles : "<<num_particles<<endl;
204  if (ipile == 0) {
205  auto origin = (*mc->vertices_begin())->position();
206  xyz0Ptr->SetXYZ(origin.x() * mmToCm, origin.y() * mmToCm, origin.z() * mmToCm);
207  *t0Ptr = origin.t() * mmToNs;
208  }
209  fillIndices(mc, particles, *barCodeVector, offset);
210  // fill output collection and save association
211  for( size_t ipar = offset; ipar < offset + num_particles; ++ ipar ) {
212  const HepMC::GenParticle * part = particles[ ipar ];
213  reco::GenParticle & cand = cands[ ipar ];
214  // convert HepMC::GenParticle to new reco::GenParticle
215  convertParticle(cand, part);
216  cand.resetDaughters( ref_.id() );
217  }
218 
219  for( size_t d = offset; d < offset + num_particles; ++ d ) {
220  const HepMC::GenParticle * part = particles[ d ];
221  const GenVertex * productionVertex = part->production_vertex();
222  int sub_id = 0;
223  if ( productionVertex != 0 ) {
224  sub_id = productionVertex->id();
225  if(!isHI) sub_id = 0;
226  // search barcode map and attach daughters
227  fillDaughters(cands,part,d);
228  }else{
229  const GenVertex * endVertex = part->end_vertex();
230  if(endVertex != 0) sub_id = endVertex->id();
231  else throw cms::Exception( "SubEventID" )<<"SubEvent not determined. Particle has no production and no end vertex!"<<endl;
232  }
233  if(sub_id < 0) sub_id = 0;
234  int new_id = sub_id + suboffset;
235  GenParticleRef dref( ref_, d );
236  cands[d].setCollisionId(new_id); // For new GenParticle
237  LogDebug("VertexId")<<"SubEvent offset 3 : "<<suboffset;
238  }
239  int nsub = -2;
240  if(isHI){
241  nsub = hi->Ncoll_hard()+1;
242  suboffset += nsub;
243  }else{
244  suboffset += 1;
245  }
246  offset += num_particles;
247  }
248  }else{
249  auto origin = (*mc->vertices_begin())->position();
250  xyz0Ptr->SetXYZ(origin.x() * mmToCm, origin.y() * mmToCm, origin.z() * mmToCm);
251  *t0Ptr = origin.t() * mmToNs;
252  fillIndices(mc, particles, *barCodeVector, 0);
253 
254  // fill output collection and save association
255  for( size_t i = 0; i < particles.size(); ++ i ) {
256  const HepMC::GenParticle * part = particles[ i ];
257  reco::GenParticle & cand = cands[ i ];
258  // convert HepMC::GenParticle to new reco::GenParticle
259  convertParticle(cand, part);
260  cand.resetDaughters( ref_.id() );
261  }
262 
263  // fill references to daughters
264  for( size_t d = 0; d < cands.size(); ++ d ) {
265  const HepMC::GenParticle * part = particles[ d ];
266  const GenVertex * productionVertex = part->production_vertex();
267  // search barcode map and attach daughters
268  if ( productionVertex != 0 ) fillDaughters(cands,part,d);
269  cands[d].setCollisionId(0);
270  }
271  }
272 
273  evt.put(std::move(candsPtr));
274  if(saveBarCodes_) evt.put(std::move(barCodeVector));
275  if(cfhepmcprod) delete cfhepmcprod;
276  evt.put(std::move(xyz0Ptr),"xyz0");
277  evt.put(std::move(t0Ptr),"t0");
278 }
279 
281  Candidate::LorentzVector p4( part->momentum() );
282  int pdgId = part->pdg_id();
283  cand.setThreeCharge( chargeTimesThree( pdgId ) );
284  cand.setPdgId( pdgId );
285  cand.setStatus( part->status() );
286  cand.setP4( p4 );
287  cand.setCollisionId(0);
288  const GenVertex * v = part->production_vertex();
289  if ( v != 0 ) {
290  ThreeVector vtx = v->point3d();
291  Candidate::Point vertex( vtx.x() * mmToCm, vtx.y() * mmToCm, vtx.z() * mmToCm );
292  cand.setVertex( vertex );
293  } else {
294  cand.setVertex( Candidate::Point( 0, 0, 0 ) );
295  }
297  return true;
298 }
299 
301 
302  const GenVertex * productionVertex = part->production_vertex();
303  size_t numberOfMothers = productionVertex->particles_in_size();
304  if ( numberOfMothers > 0 ) {
305  GenVertex::particles_in_const_iterator motherIt = productionVertex->particles_in_const_begin();
306  for( ; motherIt != productionVertex->particles_in_const_end(); motherIt++) {
307  const HepMC::GenParticle * mother = * motherIt;
308  size_t m = barcodes_.find( mother->barcode() )->second;
309  cands[ m ].addDaughter( GenParticleRef( ref_, index ) );
310  cands[ index ].addMother( GenParticleRef( ref_, m ) );
311  }
312  }
313 
314  return true;
315 }
316 
317 bool GenParticleProducer::fillIndices(const HepMC::GenEvent * mc, vector<const HepMC::GenParticle*>& particles, vector<int>& barCodeVector, int offset){
318  size_t idx = offset;
319  HepMC::GenEvent::particle_const_iterator begin = mc->particles_begin(), end = mc->particles_end();
320  for( HepMC::GenEvent::particle_const_iterator p = begin; p != end; ++ p ) {
321  const HepMC::GenParticle * particle = * p;
322  size_t barCode_this_event = particle->barcode();
323  size_t barCode = barCode_this_event + offset;
324  if( barcodes_.find(barCode) != barcodes_.end() )
325  throw cms::Exception( "WrongReference" )
326  << "barcodes are duplicated! " << endl;
327  particles[idx] = particle;
328  barCodeVector[idx] = barCode;
329  barcodes_.insert( make_pair(barCode_this_event, idx ++) );
330  }
331  return true;
332 }
333 
334 
336 
338 
#define LogDebug(id)
GenParticleProducer(const edm::ParameterSet &)
constructor
std::vector< GenParticle > GenParticleCollection
collection of GenParticles
T getParameter(std::string const &) const
bool abortOnUnknownPDGCode_
unknown code treatment flag
MCTruthHelper< HepMC::GenParticle > mcTruthHelper_
void fillGenStatusFlags(const P &p, reco::GenStatusFlags &statusFlags)
int i
Definition: DBlmapReader.cc:9
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
tuple cfg
Definition: looper.py:293
bool fillDaughters(reco::GenParticleCollection &cand, const HepMC::GenParticle *part, size_t index)
void setCollisionId(int s)
Definition: GenParticle.h:38
static const double mmToNs
std::vector< int > chargeP_
charge indices
MCTruthHelper< reco::GenParticle > mcTruthHelperGenParts_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:457
std::map< int, size_t > barcodes_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void setPdgId(int pdgId) final
edm::Ref< GenParticleCollection > GenParticleRef
persistent reference to a GenParticle
std::vector< edm::EDGetTokenT< edm::HepMCProduct > > vectorSrcTokens_
~GenParticleProducer()
destructor
virtual void produce(edm::Event &e, const edm::EventSetup &) override
process one event
edm::EDGetTokenT< CrossingFrame< edm::HepMCProduct > > mixToken_
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float > > XYZPointF
point in space with cartesian internal representation
Definition: Point3D.h:10
int chargeTimesThree(int) const
edm::EDGetTokenT< edm::HepMCProduct > srcToken_
source collection name
int size() const
Definition: MixCollection.h:24
void getData(T &iHolder) const
Definition: EventSetup.h:79
U second(std::pair< T, U > const &p)
reco::GenParticleRefProd ref_
tuple d
Definition: ztail.py:151
void resetDaughters(const edm::ProductID &id)
set daughters product ID
bool firstEvent_
whether the first event was looked at
double p4[4]
Definition: TauolaWrapper.h:92
def move
Definition: eostools.py:510
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
double f[11][100]
#define end
Definition: vmac.h:37
virtual void setVertex(const Point &vertex)
set vertex
RefProd< PROD > getRefBeforePut()
Definition: Event.h:134
const GenStatusFlags & statusFlags() const
Definition: GenParticle.h:41
bool fillIndices(const HepMC::GenEvent *mc, std::vector< const HepMC::GenParticle * > &particles, std::vector< int > &barCodeVector, int offset)
bool convertParticle(reco::GenParticle &cand, const HepMC::GenParticle *part)
virtual void setStatus(int status) final
set status word
const T & getObject(unsigned int ip) const
Definition: MixCollection.h:31
const HepMC::GenEvent * GetEvent() const
Definition: HepMCProduct.h:35
T const * product() const
Definition: Handle.h:81
static const double mmToCm
part
Definition: HCALResponse.h:20
std::map< int, int > chargeMap_
virtual void setP4(const LorentzVector &p4) final
set 4-momentum
math::XYZTLorentzVector LorentzVector
Lorentz vector.
Definition: Candidate.h:37
bool doSubEvent_
input &amp; output modes
std::vector< int > chargeM_
#define begin
Definition: vmac.h:30
static int position[264][3]
Definition: ReadPGInfo.cc:509
ProductID id() const
Accessor for product ID.
Definition: RefProd.h:137
math::XYZPoint Point
point in the space
Definition: Candidate.h:41
bool saveBarCodes_
save bar-codes
volatile std::atomic< bool > shutdown_flag false
virtual void setThreeCharge(Charge qx3) final
set electric charge
Definition: LeafCandidate.h:97
tuple size
Write out results.
static const int PDGCacheMax