CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HepMCSplitter.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: HepMCSplitter
4 // Class: HepMCSplitter
5 //
13 //
14 // Original Author: Tomasz Maciej Frueboes
15 //
16 //
17 
18 /*
19 Usage.
20 
21 1) Add to prodFilterSequence two modules:
22 
23 process.genUE = cms.EDProducer("HepMCSplitter",
24  input = cms.InputTag("generator"),
25  todo = cms.vstring("UE")
26 )
27 
28 process.genZtautau = cms.EDProducer("HepMCSplitter",
29  input = cms.InputTag("generator"),
30  todo = cms.vstring("Ztautau")
31 )
32 
33 process.ProductionFilterSequence = cms.Sequence(process.generator+process.tfFilter + process.genUE + process.genZtautau)
34 
35 2) Replace "generator" with desired part, eg:
36 
37 process.mix.mixObjects.mixHepMC.input = cms.VInputTag(cms.InputTag("genZtautau"))
38 process.g4SimHits.Generator.HepMCProductLabel = cms.string("genZtautau")
39 process.mergedtruth.HepMCDataLabels.remove("generator")
40 process.mergedtruth.HepMCDataLabels.extend(["genZtautau"])
41 # Note - above was done for python wo vtx smearing
42 
43 */
44 
45 
46 // system include files
47 #include <memory>
48 
49 // user include files
52 
55 
60 
65 
67 #include <HepMC/GenEvent.h>
68 #include <HepMC/GenVertex.h>
69 #include <HepMC/GenParticle.h>
70 
71 // class decleration
72 //
73 
75  public:
76  explicit HepMCSplitter(const edm::ParameterSet&);
78 
79  private:
80  virtual void beginJob() ;
81  virtual void produce(edm::Event&, const edm::EventSetup&);
82  virtual void endJob() ;
83 
84 
85 
86  int cntStableParticles(const HepMC::GenEvent * evtUE);
87  // ----------member data ---------------------------
89  std::vector<std::string> _todo;
90  bool _doUe;
91  bool _doZtautau;
92 };
93 
94 //
95 // constants, enums and typedefs
96 //
97 
98 
99 //
100 // static data member definitions
101 //
102 
103 //
104 // constructors and destructor
105 //
107  _input(iConfig.getParameter<edm::InputTag>("input") ),
108  _todo(iConfig.getParameter< std::vector<std::string> >("todo"))
109 {
110 
111  // Oscar producer uses string, not input tag :(
112 
113 
114  _doUe = std::find(_todo.begin(), _todo.end(), "UE")!=_todo.end();
115  _doZtautau = std::find(_todo.begin(), _todo.end(), "Ztautau")!=_todo.end();
116 
117 
118  if (_todo.size()==1 ) {
119  produces<edm::HepMCProduct>();
120 
121  if (!_doUe && !_doZtautau )
122  throw cms::Exception("") << "todo wrong - select from Ztautau and UE \n";
123  } else if (_todo.size()==2 ) {
124  produces<edm::HepMCProduct>("Ztautau");
125  produces<edm::HepMCProduct>("UE");
126  if (!_doUe || !_doZtautau )
127  throw cms::Exception("") << "todo wrong - select from Ztautau and UE \n";
128  } else {
129  throw cms::Exception("") << "todo wrong\n";
130  }
131 }
132 
134 {
135 
136  // do anything here that needs to be done at desctruction time
137  // (e.g. close files, deallocate resources etc.)
138 
139 }
140 
141 
142 //
143 // member functions
144 //
145 
146 // ------------ method called to produce the data ------------
147 void
149 {
150  using namespace edm;
151  using namespace reco;
152 
153  // Get the product
154 
156  iEvent.getByLabel(_input, prodIn);
157 
158 
159 
160  // find Z decay and prod vertex
161  int ZprodVtx, ZdecayVtx = 0;
162  bool vtxFound = false;
163  HepMC::GenEvent::particle_const_iterator parIt, parItE;
164  parIt = prodIn->GetEvent()->particles_begin();
165  parItE = prodIn->GetEvent()->particles_end();
166  for ( ;parIt != parItE; ++parIt){
167  int pdg = abs( (*parIt)->pdg_id() ) ;
168  if (pdg == 23){
169  ZdecayVtx = (*parIt)->end_vertex()->barcode();
170  ZprodVtx = (*parIt)->production_vertex()->barcode();
171  vtxFound=true;
172  break;
173  }
174  }
175 
176  if ( !vtxFound) throw cms::Exception("Z vtx not found");
177 // std::cout << " Found Z - prod vtx " << ZprodVtx
178 // << " decay vtx " << ZdecayVtx
179 // << std::endl;
180 
181 
182 
183  // iterate over all vertices, mark those comming from Z
184  std::set<int> barcodesZ;
185  barcodesZ.insert(ZdecayVtx);
186 
187  HepMC::GenEvent::vertex_const_iterator itVtx, itVtxE;
188  itVtx = prodIn->GetEvent()->vertices_begin();
189  itVtxE = prodIn->GetEvent()->vertices_end();
190 
191  for(;itVtx!=itVtxE;++itVtx){
192  if ( barcodesZ.find( (*itVtx)->barcode() )!=barcodesZ.end()){
193 
194  // iterate over out particles. Copy their decay vertices.to Z vertices
195  HepMC::GenVertex::particles_out_const_iterator parIt, parItE;
196  parIt=(*itVtx)->particles_out_const_begin();
197  parItE=(*itVtx)->particles_out_const_end();
198  for ( ; parIt!=parItE; ++parIt){
199  if ((*parIt)->end_vertex()){
200  barcodesZ.insert((*parIt)->end_vertex()->barcode());
201  //std::cout << " Inserted " << (*parIt)->end_vertex()->barcode() << std::endl;
202  }
203  }
204 
205  }
206  }
207 
208  barcodesZ.insert(ZprodVtx);
209  barcodesZ.insert(ZdecayVtx);
210 
211  //std::cout << " Z vertices identified " << barcodesZ.size() << std::endl;
212 
213  // prepare final products
214  // std::auto_ptr<edm::HepMCProduct> prodTauTau(new edm::HepMCProduct());
215 
216  std::auto_ptr<HepMC::GenEvent> evtUE(new HepMC::GenEvent());
217  std::auto_ptr<HepMC::GenEvent> evtTauTau(new HepMC::GenEvent());
218 
219  // Copy the vertices
220  itVtx = prodIn->GetEvent()->vertices_begin();
221  itVtxE = prodIn->GetEvent()->vertices_end();
222  int newBarcode = 1;
223 
224  std::map<int, HepMC::GenVertex*> barcodeMapUE;
225  std::map<int, HepMC::GenVertex*> barcodeMapTauTau;
226  for(;itVtx!=itVtxE;++itVtx){
227  bool isZ = barcodesZ.find( (*itVtx)->barcode() )!=barcodesZ.end();
228 
229  bool leadToZ = false;
230  if ( (*parIt)->production_vertex() ) {
231  leadToZ = std::abs( (*itVtx)->barcode()) < std::abs(ZprodVtx);
232  }
233 
234 
235  // special case - Z production vertex. Put in both
236  if (isZ || (*itVtx)->barcode()==ZprodVtx || leadToZ ){
237  HepMC::GenVertex* newvertex = new HepMC::GenVertex((*itVtx)->position(),
238  (*itVtx)->id(),
239  (*itVtx)->weights() );
240 
241  newvertex->suggest_barcode( -(++newBarcode) );
242  //std::cout << "Adding vtx " << (*itVtx)->barcode() << std::endl;
243  barcodeMapTauTau[(*itVtx)->barcode()] = newvertex;
244  evtTauTau->add_vertex(newvertex);
245  }
246 
247  if (!isZ || (*itVtx)->barcode()==ZprodVtx){
248  HepMC::GenVertex* newvertex = new HepMC::GenVertex((*itVtx)->position(),
249  (*itVtx)->id(),
250  (*itVtx)->weights() );
251 
252  newvertex->suggest_barcode( -(++newBarcode) );
253  barcodeMapUE[ (*itVtx)->barcode()] = newvertex;
254  evtUE->add_vertex(newvertex);
255  }
256  }
257 
258 
259 
260  //std::cout << " Vertices copied " << std::endl;
261 
262 
263 
264  // iterate over particles again. Copy and assign
265  parIt = prodIn->GetEvent()->particles_begin();
266  parItE = prodIn->GetEvent()->particles_end();
267  for ( ;parIt != parItE; ++parIt){
268 
269  // check if prod vertex is Z vertex.
270  bool fromZprodVtx = false;
271  if ( (*parIt)->production_vertex() ){
272  fromZprodVtx = ( (*parIt)->production_vertex()->barcode() == ZprodVtx);
273  }
274 
275  // check if from Z
276  bool fromZ = false;
277  if ( (*parIt)->production_vertex() ){
278  fromZ = barcodesZ.find( (*parIt)->production_vertex()->barcode() )!=barcodesZ.end();
279  }
280 
281 
282  // incoming particle (proton) or interaction that lead to Z boson. Copy to both events
283  bool leadToZ = false;
284  //if ( (*parIt)->production_vertex() && (*parIt)->end_vertex()) {
285  if ( (*parIt)->end_vertex()) {
286  //leadToZ = std::abs((*parIt)->production_vertex()->barcode()) < std::abs(ZprodVtx);
287  leadToZ = std::abs((*parIt)->end_vertex()->barcode()) < std::abs(ZdecayVtx);
288  }
289 
290  //if ( !(*parIt)->production_vertex() || leadToZ ) {
291  if ( leadToZ ) {
292 
293  HepMC::GenParticle* newparticle = new HepMC::GenParticle(*(*parIt));
294  newparticle->suggest_barcode(newBarcode++);
295  if ( (*parIt)->end_vertex()){
296  //XX Add xcheck if this barcode is in map
297  barcodeMapUE[ (*parIt)->end_vertex()->barcode()]->add_particle_in(newparticle);
298  } else {
299  //std::cout << "XXX particle without vertices!" << std::endl;
300  }
301 
302  newparticle = new HepMC::GenParticle(**parIt);
303  newparticle->suggest_barcode(newBarcode++);
304  if ( (*parIt)->end_vertex()){
305  //XX Add xcheck if this barcode is in map
306  //std::cout << "Trying: " << (*parIt)->end_vertex()->barcode() << std::endl;
307  barcodeMapTauTau[ (*parIt)->end_vertex()->barcode()]->add_particle_in(newparticle);
308  } else {
309  //std::cout << "XXX particle without vertices!" << std::endl;
310  }
311 
312  } else if (fromZprodVtx) { // Special case - particle comes from Z prod vertex
313  bool isZ = (*parIt)->pdg_id();
314  HepMC::GenParticle* newparticle = new HepMC::GenParticle(**parIt);
315  newparticle->suggest_barcode(newBarcode++);
316  if (isZ){
317  //XX Add xcheck if this barcode is in map
318  if ( (*parIt)->end_vertex())
319  barcodeMapTauTau[(*parIt)->end_vertex()->barcode()]->add_particle_in(newparticle);
320  if ((*parIt)->production_vertex())
321  barcodeMapTauTau[(*parIt)->production_vertex()->barcode()]->add_particle_out(newparticle);
322  } else {
323  if ((*parIt)->end_vertex())
324  barcodeMapUE[(*parIt)->end_vertex()->barcode()]->add_particle_in(newparticle);
325  if ((*parIt)->production_vertex())
326  barcodeMapUE[(*parIt)->production_vertex()->barcode()]->add_particle_out(newparticle);
327  }
328  } else {
329  // Generic case
330  HepMC::GenParticle* newparticle = new HepMC::GenParticle(**parIt);
331  newparticle->suggest_barcode(newBarcode++);
332  if (fromZ){
333  //XX Add xcheck if this barcode is in map
334  if ((*parIt)->end_vertex())
335  barcodeMapTauTau[(*parIt)->end_vertex()->barcode()]->add_particle_in(newparticle);
336  if ((*parIt)->production_vertex())
337  barcodeMapTauTau[(*parIt)->production_vertex()->barcode()]->add_particle_out(newparticle);
338  } else {
339  if ((*parIt)->end_vertex())
340  // std::cout << "trying " << (*parIt)->end_vertex()->barcode() << std::endl;
341  barcodeMapUE[(*parIt)->end_vertex()->barcode()]->add_particle_in(newparticle);
342  if ((*parIt)->production_vertex())
343  barcodeMapUE[(*parIt)->production_vertex()->barcode()]->add_particle_out(newparticle);
344  }
345  }
346  } // end particles iteration
347 
348 // std::cout << "XXXXXXXXXXXXXXXXXX" << std::endl;
349 
350  int orgCnt = cntStableParticles(prodIn->GetEvent());
351  int tautauCnt = cntStableParticles(evtTauTau.get());
352  int UEcnt = cntStableParticles(evtUE.get());
353  if (orgCnt != tautauCnt + UEcnt){
354  std::cout << " XXX Warning - wrong num particles - ORG: "
355  << orgCnt
356  << " tautau " << tautauCnt
357  << " UE " << UEcnt
358  << std::endl;
359  }
360 
361 // std::cout << "=====================================\n";
362 // std::cout << "EV Tau \n";
363 // //evtTauTau->print();
364 // std::cout << "=====================================\n";
365 // std::cout << "EV org \n";
366 // //prodIn->GetEvent()->print();
367 //
368 // std::cout << "=====================================\n";
369 // std::cout << "EV EU \n";
370  // evtUE->print();
371 
372  if (_doUe){
373  std::auto_ptr<edm::HepMCProduct> prodUE(new edm::HepMCProduct());
374  prodUE->addHepMCData( evtUE.release() );
375  if (_doZtautau) iEvent.put( prodUE, "UE");
376  else iEvent.put( prodUE);
377  }
378 
379  if (_doZtautau) {
380  std::auto_ptr<edm::HepMCProduct> prodTauTau(new edm::HepMCProduct());
381 // std::cout << "===================================== " << iEvent.id() << std::endl;
382 // evtTauTau->print();
383  prodTauTau->addHepMCData( evtTauTau.release() );
384  if (_doUe ) iEvent.put( prodTauTau, "Ztautau");
385  else iEvent.put( prodTauTau);
386  }
387 
388 
389 }
390 
391 // ------------ method called once each job just before starting event loop ------------
392 void
394 {
395 }
396 
397 // ------------ method called once each job just after ending the event loop ------------
398 void
400 }
401 
402 int HepMCSplitter::cntStableParticles(const HepMC::GenEvent * evt){
403 
404  int ret = 0;
405  HepMC::GenEvent::particle_const_iterator parIt, parItE;
406  parIt = evt->particles_begin();
407  parItE = evt->particles_end();
408  for ( ;parIt != parItE; ++parIt){
409  if ( (*parIt)->status()==1) ++ ret;
410 
411  }
412 
413  return ret;
414 }
415 
416 //define this as a plug-in
std::vector< std::string > _todo
HepMCSplitter(const edm::ParameterSet &)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
virtual void endJob()
int cntStableParticles(const HepMC::GenEvent *evtUE)
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
virtual void produce(edm::Event &, const edm::EventSetup &)
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:405
edm::InputTag _input
tuple cout
Definition: gather_cfg.py:121
virtual void beginJob()