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 
91  firstEvent_(true),
92  abortOnUnknownPDGCode_( cfg.getUntrackedParameter<bool>( "abortOnUnknownPDGCode", true ) ),
93  saveBarCodes_( cfg.getUntrackedParameter<bool>( "saveBarCodes", false ) ),
94  chargeP_( PDGCacheMax, 0 ), chargeM_( PDGCacheMax, 0 ),
95  doSubEvent_(cfg.getUntrackedParameter<bool>( "doSubEvent", false )),
96  useCF_(cfg.getUntrackedParameter<bool>( "useCrossingFrame", false ))
97 {
98  produces<GenParticleCollection>();
99  if( saveBarCodes_ ) {
100  std::string alias( cfg.getParameter<std::string>( "@module_label" ) );
101  produces<vector<int> >().setBranchAlias( alias + "BarCodes" );
102  }
103 
104  if(useCF_) mixToken_ = mayConsume<CrossingFrame<HepMCProduct> >(InputTag(cfg.getParameter<std::string>( "mix" ),"generator"));
105  else srcToken_ = mayConsume<HepMCProduct>(cfg.getParameter<InputTag>( "src" ));
106 }
107 
109 }
110 
112  if( std::abs( id ) < PDGCacheMax )
113  return id > 0 ? chargeP_[ id ] : chargeM_[ - id ];
114  map<int, int>::const_iterator f = chargeMap_.find( id );
115  if ( f == chargeMap_.end() ) {
118  << "invalid PDG id: " << id << endl;
119  else
120  return HepPDT::ParticleID(id).threeCharge();
121  }
122  return f->second;
123 }
124 
126 
127  if (firstEvent_) {
129  es.getData( pdt );
130  for( HepPDT::ParticleDataTable::const_iterator p = pdt->begin(); p != pdt->end(); ++ p ) {
131  const HepPDT::ParticleID & id = p->first;
132  int pdgId = id.pid(), apdgId = std::abs( pdgId );
133  int q3 = id.threeCharge();
134  if ( apdgId < PDGCacheMax && pdgId > 0 ) {
135  chargeP_[ apdgId ] = q3;
136  chargeM_[ apdgId ] = -q3;
137  } else if ( apdgId < PDGCacheMax ) {
138  chargeP_[ apdgId ] = -q3;
139  chargeM_[ apdgId ] = q3;
140  } else {
141  chargeMap_[ pdgId ] = q3;
142  chargeMap_[ -pdgId ] = -q3;
143  }
144  }
145  firstEvent_ = false;
146  }
147 
148  barcodes_.clear();
149 
150  size_t totalSize = 0;
151  const GenEvent * mc = 0;
152  MixCollection<HepMCProduct>* cfhepmcprod = 0;
153  size_t npiles = 1;
154 
155  if(useCF_){
157  evt.getByToken(mixToken_,cf);
158  cfhepmcprod = new MixCollection<HepMCProduct>(cf.product());
159  npiles = cfhepmcprod->size();
160  LogDebug("GenParticleProducer")<<"npiles : "<<npiles<<endl;
161  for(unsigned int icf = 0; icf < npiles; ++icf){
162  LogDebug("GenParticleProducer")<<"subSize : "<<cfhepmcprod->getObject(icf).GetEvent()->particles_size()<<endl;
163  totalSize += cfhepmcprod->getObject(icf).GetEvent()->particles_size();
164  }
165  LogDebug("GenParticleProducer")<<"totalSize : "<<totalSize<<endl;
166  }else{
168  evt.getByToken( srcToken_, mcp );
169  mc = mcp->GetEvent();
170  if( mc == 0 )
172  << "HepMC has null pointer to GenEvent" << endl;
173  totalSize = mc->particles_size();
174  }
175 
176  // initialise containers
177  const size_t size = totalSize;
178  vector<const HepMC::GenParticle *> particles( size );
179  auto_ptr<GenParticleCollection> candsPtr( new GenParticleCollection( size ) );
180  auto_ptr<vector<int> > barCodeVector( new vector<int>( size ) );
182  GenParticleCollection & cands = * candsPtr;
183  size_t offset = 0;
184  size_t suboffset = 0;
185 
187  if(doSubEvent_ || useCF_){
188  for(size_t ipile = 0; ipile < npiles; ++ipile){
189  LogDebug("GenParticleProducer")<<"mixed object ipile : "<<ipile<<endl;
190  barcodes_.clear();
191  if(useCF_) mc = cfhepmcprod->getObject(ipile).GetEvent();
192 
193  //Look whether heavy ion/signal event
194  bool isHI = false;
195  const HepMC::HeavyIon * hi = mc->heavy_ion();
196  if(hi && hi->Ncoll_hard() > 1) isHI = true;
197  size_t num_particles = mc->particles_size();
198  LogDebug("GenParticleProducer")<<"num_particles : "<<num_particles<<endl;
199  fillIndices(mc, particles, *barCodeVector, offset);
200  // fill output collection and save association
201  for( size_t ipar = offset; ipar < offset + num_particles; ++ ipar ) {
202  const HepMC::GenParticle * part = particles[ ipar ];
203  reco::GenParticle & cand = cands[ ipar ];
204  // convert HepMC::GenParticle to new reco::GenParticle
205  convertParticle(cand, part);
206  cand.resetDaughters( ref_.id() );
207  }
208 
209  for( size_t d = offset; d < offset + num_particles; ++ d ) {
210  const HepMC::GenParticle * part = particles[ d ];
211  const GenVertex * productionVertex = part->production_vertex();
212  int sub_id = 0;
213  if ( productionVertex != 0 ) {
214  sub_id = productionVertex->id();
215  if(!isHI) sub_id = 0;
216  // search barcode map and attach daughters
217  fillDaughters(cands,part,d);
218  }else{
219  const GenVertex * endVertex = part->end_vertex();
220  if(endVertex != 0) sub_id = endVertex->id();
221  else throw cms::Exception( "SubEventID" )<<"SubEvent not determined. Particle has no production and no end vertex!"<<endl;
222  }
223  if(sub_id < 0) sub_id = 0;
224  int new_id = sub_id + suboffset;
225  GenParticleRef dref( ref_, d );
226  cands[d].setCollisionId(new_id); // For new GenParticle
227  LogDebug("VertexId")<<"SubEvent offset 3 : "<<suboffset;
228  }
229  int nsub = -2;
230  if(isHI){
231  nsub = hi->Ncoll_hard()+1;
232  suboffset += nsub;
233  }else{
234  suboffset += 1;
235  }
236  offset += num_particles;
237  }
238  }else{
239  fillIndices(mc, particles, *barCodeVector, 0);
240 
241  // fill output collection and save association
242  for( size_t i = 0; i < particles.size(); ++ i ) {
243  const HepMC::GenParticle * part = particles[ i ];
244  reco::GenParticle & cand = cands[ i ];
245  // convert HepMC::GenParticle to new reco::GenParticle
246  convertParticle(cand, part);
247  cand.resetDaughters( ref_.id() );
248  }
249 
250  // fill references to daughters
251  for( size_t d = 0; d < cands.size(); ++ d ) {
252  const HepMC::GenParticle * part = particles[ d ];
253  const GenVertex * productionVertex = part->production_vertex();
254  // search barcode map and attach daughters
255  if ( productionVertex != 0 ) fillDaughters(cands,part,d);
256  cands[d].setCollisionId(0);
257  }
258  }
259 
260  evt.put( candsPtr );
261  if(saveBarCodes_) evt.put( barCodeVector );
262  if(cfhepmcprod) delete cfhepmcprod;
263 
264 }
265 
267  Candidate::LorentzVector p4( part->momentum() );
268  int pdgId = part->pdg_id();
269  cand.setThreeCharge( chargeTimesThree( pdgId ) );
270  cand.setPdgId( pdgId );
271  cand.setStatus( part->status() );
272  cand.setP4( p4 );
273  cand.setCollisionId(0);
274  const GenVertex * v = part->production_vertex();
275  if ( v != 0 ) {
276  ThreeVector vtx = v->point3d();
277  Candidate::Point vertex( vtx.x() * mmToCm, vtx.y() * mmToCm, vtx.z() * mmToCm );
278  cand.setVertex( vertex );
279  } else {
280  cand.setVertex( Candidate::Point( 0, 0, 0 ) );
281  }
283  return true;
284 }
285 
287 
288  const GenVertex * productionVertex = part->production_vertex();
289  size_t numberOfMothers = productionVertex->particles_in_size();
290  if ( numberOfMothers > 0 ) {
291  GenVertex::particles_in_const_iterator motherIt = productionVertex->particles_in_const_begin();
292  for( ; motherIt != productionVertex->particles_in_const_end(); motherIt++) {
293  const HepMC::GenParticle * mother = * motherIt;
294  size_t m = barcodes_.find( mother->barcode() )->second;
295  cands[ m ].addDaughter( GenParticleRef( ref_, index ) );
296  cands[ index ].addMother( GenParticleRef( ref_, m ) );
297  }
298  }
299 
300  return true;
301 }
302 
303 bool GenParticleProducer::fillIndices(const HepMC::GenEvent * mc, vector<const HepMC::GenParticle*>& particles, vector<int>& barCodeVector, int offset){
304  size_t idx = offset;
305  HepMC::GenEvent::particle_const_iterator begin = mc->particles_begin(), end = mc->particles_end();
306  for( HepMC::GenEvent::particle_const_iterator p = begin; p != end; ++ p ) {
307  const HepMC::GenParticle * particle = * p;
308  size_t barCode_this_event = particle->barcode();
309  size_t barCode = barCode_this_event + offset;
310  if( barcodes_.find(barCode) != barcodes_.end() )
311  throw cms::Exception( "WrongReference" )
312  << "barcodes are duplicated! " << endl;
313  particles[idx] = particle;
314  barCodeVector[idx] = barCode;
315  barcodes_.insert( make_pair(barCode_this_event, idx ++) );
316  }
317  return true;
318 }
319 
320 
322 
324 
#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
tuple cfg
Definition: looper.py:259
bool fillDaughters(reco::GenParticleCollection &cand, const HepMC::GenParticle *part, size_t index)
void setCollisionId(int s)
Definition: GenParticle.h:38
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_
virtual void setStatus(int status)
set status word
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::Ref< GenParticleCollection > GenParticleRef
persistent reference to a GenParticle
std::vector< edm::EDGetTokenT< edm::HepMCProduct > > vectorSrcTokens_
virtual void setP4(const LorentzVector &p4)
set 4-momentum
~GenParticleProducer()
destructor
virtual void produce(edm::Event &e, const edm::EventSetup &) override
process one event
edm::EDGetTokenT< CrossingFrame< edm::HepMCProduct > > mixToken_
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:78
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
virtual void setThreeCharge(Charge qx3)
set electric charge
Definition: LeafCandidate.h:97
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:115
double p4[4]
Definition: TauolaWrapper.h:92
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:135
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)
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
tuple idx
DEBUGGING if hasattr(process,&quot;trackMonIterativeTracking2012&quot;): print &quot;trackMonIterativeTracking2012 D...
part
Definition: HCALResponse.h:20
std::map< int, int > chargeMap_
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
ProductID id() const
Accessor for product ID.
Definition: RefProd.h:137
virtual void setPdgId(int pdgId)
math::XYZPoint Point
point in the space
Definition: Candidate.h:41
bool saveBarCodes_
save bar-codes
volatile std::atomic< bool > shutdown_flag false
tuple size
Write out results.
static const int PDGCacheMax