CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PhotosInterface.cc
Go to the documentation of this file.
1 #include <iostream>
6 
8 
9 #include "HepMC/GenEvent.h"
10 #include "HepMC/IO_HEPEVT.h"
11 #include "HepMC/HEPEVT_Wrapper.h"
12 
13 using namespace gen;
14 using namespace edm;
15 using namespace std;
16 
17 
18 CLHEP::HepRandomEngine* PhotosInterface::fRandomEngine = nullptr;
19 
20 extern "C"{
21 
22  void phoini_( void );
23  void photos_( int& );
24 
25  double phoran_(int *idummy)
26  {
27  return PhotosInterface::flat();
28  }
29 
30  extern struct {
31  // bool qedrad[NMXHEP];
32  bool qedrad[4000]; // hardcoded for now...
33  } phoqed_;
34 
35 }
36 
37 
38 
40  : fOnlyPDG(-1)
41 {
42  fSpecialSettings.push_back("QED-brem-off:all");
44  fIsInitialized = false;
45 }
46 
48  : fOnlyPDG(-1)
49 {
50  fSpecialSettings.push_back("QED-brem-off:all");
51  fIsInitialized = false;
52 }
53 
55 {}
56 
59 }
60 
61 
63 {
64 
65  fOnlyPDG = ipdg;
66 // std::ostringstream command;
67 // command << "QED-brem-off:" << fOnlyPDG ;
68  fSpecialSettings.clear();
69 // fSpecialSettings.push_back( command.str() );
70 
71  return;
72 
73 }
74 
76 {
77 
78  if ( fIsInitialized ) return; // do init only once
79 
80  phoini_();
81 
82  fIsInitialized = true;
83 
84  return;
85 }
86 
87 HepMC::GenEvent* PhotosInterface::apply( HepMC::GenEvent* evt )
88 {
89 
90  if ( !fIsInitialized ) return evt; // conv.read_next_event();
91 
92  // loop over HepMC::GenEvent, find vertices
93 
94  // for ( int ip=0; ip<evt->particles_size(); ip++ )
95  for ( int ip=0; ip<4000; ip++ ) // 4000 is the max size of the array
96  {
97  phoqed_.qedrad[ip]=true;
98  }
99 
100  //
101  // now do actual job
102  //
103 
104  for ( int iv=1; iv<=evt->vertices_size(); iv++ )
105  {
106 
107  bool legalVtx = false;
108 
109  fSecVtxStore.clear();
110 
111  HepMC::GenVertex* vtx = evt->barcode_to_vertex( -iv ) ;
112 
113  if ( vtx->particles_in_size() != 1 ) continue; // more complex than we need
114  if ( vtx->particles_out_size() <= 1 ) continue; // no outcoming particles
115 
116  if ( (*(vtx->particles_in_const_begin()))->pdg_id() == 111 ) continue; // pi0 decay vtx - no point to try
117 
118  if ( fOnlyPDG != 1 && (*(vtx->particles_in_const_begin()))->pdg_id() != fOnlyPDG )
119  {
120  continue;
121  }
122  else
123  {
124  // requested for specific PDG ID only, typically tau (15)
125  //
126  // first check if a brem vertex, where outcoming are the same pdg id and a photon
127  //
128  bool same = false;
129  for ( HepMC::GenVertex::particle_iterator pitr=vtx->particles_begin(HepMC::children);
130  pitr != vtx->particles_end(HepMC::children); ++pitr)
131  {
132  if ( (*pitr)->pdg_id() == fOnlyPDG )
133  {
134  same = true;
135  break;
136  }
137  }
138  if ( same ) continue;
139 
140  // OK, we get here if incoming fOnlyPDG and something else outcoming
141  // call it for the whole branch starting at vtx barcode iv, and go on
142  // NOTE: theoretically, it has a danger of double counting in vertices
143  // down the decay branch originating from fOnlyPDG, but in practice
144  // it's unlikely that down the branchg there'll be more fOnlyPDG's
145 
146  // cross-check printout
147  // vtx->print();
148 
149  // overprotection...
150  //
151  if ( fOnlyPDG == 15 && fAvoidTauLeptonicDecays && isTauLeptonicDecay( vtx ) ) continue;
152 
153  applyToBranch( evt, -iv );
154  continue;
155  }
156 
157  // configured for all types of particles
158  //
159  for ( HepMC::GenVertex::particle_iterator pitr=vtx->particles_begin(HepMC::children);
160  pitr != vtx->particles_end(HepMC::children); ++pitr)
161  {
162 
163  // quark or gluon out of this vertex - no good !
164  if ( abs((*pitr)->pdg_id()) >=1 && abs((*pitr)->pdg_id()) <=8 ) break;
165  if ( abs((*pitr)->pdg_id()) == 21 ) break;
166 
167  if ( (*pitr)->status() == 1 || (*pitr)->end_vertex() )
168  {
169  // OK, legal already !
170  legalVtx = true;
171  break;
172  }
173  }
174 
175  if ( !legalVtx ) continue;
176 
177  applyToVertex( evt, -iv );
178 
179  } // end of master loop
180 
181  // restore event number in HEPEVT (safety measure, somehow needed by Hw6)
182  HepMC::HEPEVT_Wrapper::set_event_number( evt->event_number() );
183 
184  return evt;
185 
186 }
187 
188 void PhotosInterface::applyToVertex( HepMC::GenEvent* evt, int vtxbcode )
189 {
190 
191  HepMC::GenVertex* vtx = evt->barcode_to_vertex( vtxbcode );
192 
193  if ( fAvoidTauLeptonicDecays && isTauLeptonicDecay( vtx ) ) return;
194 
195  // cross-check printout
196  //
197  // vtx->print();
198 
199  // first, flush out HEPEVT & tmp barcode storage
200  //
201  HepMC::HEPEVT_Wrapper::zero_everything();
202  fBarcodes.clear();
203 
204  // add incoming particle
205  //
206  int index = 1;
207  HepMC::HEPEVT_Wrapper::set_id( index, (*(vtx->particles_in_const_begin()))->pdg_id() );
208  HepMC::FourVector vec4;
209  vec4 = (*(vtx->particles_in_const_begin()))->momentum();
210  HepMC::HEPEVT_Wrapper::set_momentum( index, vec4.x(), vec4.y(), vec4.z(), vec4.e() );
211  HepMC::HEPEVT_Wrapper::set_mass( index, (*(vtx->particles_in_const_begin()))->generated_mass() );
212  HepMC::HEPEVT_Wrapper::set_position( index, vtx->position().x(), vtx->position().y(),
213  vtx->position().z(), vtx->position().t() );
214  HepMC::HEPEVT_Wrapper::set_status( index, (*(vtx->particles_in_const_begin()))->status() );
215  HepMC::HEPEVT_Wrapper::set_parents( index, 0, 0 );
216  fBarcodes.push_back( (*(vtx->particles_in_const_begin()))->barcode() );
217 
218  // add outcoming particles (decay products)
219  //
220  int lastDau = 1;
221  for ( HepMC::GenVertex::particle_iterator pitr=vtx->particles_begin(HepMC::children);
222  pitr != vtx->particles_end(HepMC::children); ++pitr)
223  {
224 
225  if ( (*pitr)->status() == 1 || (*pitr)->end_vertex() )
226  {
227  index++;
228  vec4 = (*pitr)->momentum();
229  HepMC::HEPEVT_Wrapper::set_id( index, (*pitr)->pdg_id() );
230  HepMC::HEPEVT_Wrapper::set_momentum( index, vec4.x(), vec4.y(), vec4.z(), vec4.e() );
231  HepMC::HEPEVT_Wrapper::set_mass( index, (*pitr)->generated_mass() );
232  vec4 = (*pitr)->production_vertex()->position();
233  HepMC::HEPEVT_Wrapper::set_position( index, vec4.x(), vec4.y(), vec4.z(), vec4.t() );
234  HepMC::HEPEVT_Wrapper::set_status( index, (*pitr)->status() );
235  HepMC::HEPEVT_Wrapper::set_parents( index, 1, 1 );
236  fBarcodes.push_back( (*pitr)->barcode() );
237  lastDau++;
238  }
239  if ( (*pitr)->end_vertex() )
240  {
241  fSecVtxStore.push_back( (*pitr)->end_vertex()->barcode() );
242  }
243  }
244 
245  // store, further to set NHEP in HEPEVT
246  //
247  int nentries = index;
248 
249  // reset master pointer to mother
250  index = 1;
251  HepMC::HEPEVT_Wrapper::set_children ( index, 2, lastDau ); // FIXME: need to check
252  // if last daughter>=2 !!!
253 
254  // finally, set number of entries (NHEP) in HEPEVT
255  //
256  HepMC::HEPEVT_Wrapper::set_number_entries( nentries );
257 
258  // cross-check printout HEPEVT
259  // HepMC::HEPEVT_Wrapper::print_hepevt();
260 
261  // OK, 1-level vertex is formed - now, call PHOTOS
262  //
263  photos_( index ) ;
264 
265  // another cross-check printout HEPEVT - after photos
266  // HepMC::HEPEVT_Wrapper::print_hepevt();
267 
268 
269  // now check if something has been generated
270  // and make all adjustments to underlying vtx/parts
271  //
272  attachParticles( evt, vtx, nentries );
273 
274  // ugh, done with this vertex !
275 
276  return;
277 
278 }
279 
280 void PhotosInterface::applyToBranch( HepMC::GenEvent* evt, int vtxbcode )
281 {
282 
283 
284  fSecVtxStore.clear();
285 
286  // 1st level vertex
287  //
288  applyToVertex( evt, vtxbcode );
289 
290  // now look down the branch for more vertices, if any
291  //
292  // Note: fSecVtxStore gets filled up in applyToVertex, if necessary
293  //
294  unsigned int vcounter = 0;
295 
296  while ( vcounter < fSecVtxStore.size() )
297  {
298  applyToVertex( evt, fSecVtxStore[vcounter] );
299  vcounter++;
300  }
301 
302  fSecVtxStore.clear();
303 
304  return;
305 
306 }
307 
308 void PhotosInterface::attachParticles( HepMC::GenEvent* evt, HepMC::GenVertex* vtx, int nentries )
309 {
310 
311  if ( HepMC::HEPEVT_Wrapper::number_entries() > nentries )
312  {
313  // yes, need all corrections and adjustments -
314  // figure out how many photons and what particles in
315  // the decay branch have changes;
316  // also, follow up each one and correct accordingly;
317  // at the same time, add photon(s) to the GenVertex
318  //
319 
320  // vtx->print();
321 
322  int largestBarcode = -1;
323  int Nbcodes = fBarcodes.size();
324 
325  for ( int ip=1; ip<Nbcodes; ip++ )
326  {
327 
328  int bcode = fBarcodes[ip];
329  HepMC::GenParticle* prt = evt->barcode_to_particle( bcode );
330  if ( bcode > largestBarcode ) largestBarcode = bcode;
331  double px = HepMC::HEPEVT_Wrapper::px(ip+1);
332  double py = HepMC::HEPEVT_Wrapper::py(ip+1);
333  double pz = HepMC::HEPEVT_Wrapper::pz(ip+1);
334  double e = HepMC::HEPEVT_Wrapper::e(ip+1);
335  double m = HepMC::HEPEVT_Wrapper::m(ip+1);
336 
337  if ( prt->end_vertex() )
338  {
339 
340  HepMC::GenVertex* endVtx = prt->end_vertex();
341 
342  std::vector<int> secVtxStorage;
343  secVtxStorage.clear();
344 
345  secVtxStorage.push_back( endVtx->barcode() );
346 
347  HepMC::FourVector mom4 = prt->momentum();
348 
349  // now rescale all descendants
350  double bet1[3], bet2[3], gam1, gam2, pb;
351  double mass = mom4.m();
352  bet1[0] = -(mom4.px()/mass);
353  bet1[1] = -(mom4.py()/mass);
354  bet1[2] = -(mom4.pz()/mass);
355  bet2[0] = px/m;
356  bet2[1] = py/m;
357  bet2[2] = pz/m;
358  gam1 = mom4.e()/mass;
359  gam2 = e/m;
360 
361  unsigned int vcounter = 0;
362 
363  while ( vcounter < secVtxStorage.size() )
364  {
365 
366  HepMC::GenVertex* theVtx = evt->barcode_to_vertex( secVtxStorage[vcounter] );
367 
368  for ( HepMC::GenVertex::particle_iterator pitr=theVtx->particles_begin(HepMC::children);
369  pitr != theVtx->particles_end(HepMC::children); ++pitr)
370  {
371 
372  if ( (*pitr)->end_vertex() )
373  {
374  secVtxStorage.push_back( (*pitr)->end_vertex()->barcode() );
375  }
376 
377  if ( theVtx->particles_out_size() == 1 && (*pitr)->pdg_id() == prt->pdg_id() )
378  {
379  // carbon copy
380  (*pitr)->set_momentum( HepMC::FourVector(px,py,pz,e) );
381  continue;
382  }
383 
384  HepMC::FourVector dmom4 = (*pitr)->momentum();
385 
386  // Boost vector to parent rest frame...
387  pb = bet1[0]*dmom4.px() + bet1[1]*dmom4.py() + bet1[2]*dmom4.pz();
388  double dpx = dmom4.px() + bet1[0] * (dmom4.e() + pb/(gam1+1.) );
389  double dpy = dmom4.py() + bet1[1] * (dmom4.e() + pb/(gam1+1.) );
390  double dpz = dmom4.pz() + bet1[2] * (dmom4.e() + pb/(gam1+1.) );
391  double de = gam1*dmom4.e() + pb;
392  // ...and boost back to modified parent frame
393  pb = bet2[0]*dpx + bet2[1]*dpy + bet2[2]*dpz;
394  dpx += bet2[0] * ( de + pb/(gam2+1.) );
395  dpy += bet2[1] * ( de + pb/(gam2+1.) );
396  dpz += bet2[2] * ( de + pb/(gam2+1.) );
397  de *= gam2;
398  de += pb;
399 
400  (*pitr)->set_momentum( HepMC::FourVector(dpx,dpy,dpz,de) );
401 
402  }
403  vcounter++;
404  }
405 
406  secVtxStorage.clear();
407  }
408 
409  prt->set_momentum( HepMC::FourVector(px,py,pz,e) );
410 
411  } // ok, all affected particles update, but the photon(s) still not inserted
412 
413  int newlyGen = HepMC::HEPEVT_Wrapper::number_entries() - nentries;
414 
415  if ( largestBarcode < evt->particles_size() )
416  {
417  // need to adjust barcodes down from the affected vertex/particles
418  // such that we "free up" barcodes for newly generated photons
419  // in the middle of the event record
420  //
421  for ( int ipp=evt->particles_size(); ipp>largestBarcode; ipp-- )
422  {
423  (evt->barcode_to_particle(ipp))->suggest_barcode( ipp+newlyGen );
424  }
425  }
426 
427  // now attach new generated photons to the vertex
428  //
429  for ( int ipnw=1; ipnw<=newlyGen; ipnw++ )
430  {
431  int nbcode = largestBarcode+ipnw;
432  int pdg_id = HepMC::HEPEVT_Wrapper::id( nentries+ipnw );
433  int status = HepMC::HEPEVT_Wrapper::status( nentries+ipnw );
434  double px = HepMC::HEPEVT_Wrapper::px( nentries+ipnw );
435  double py = HepMC::HEPEVT_Wrapper::py( nentries+ipnw );
436  double pz = HepMC::HEPEVT_Wrapper::pz( nentries+ipnw );
437  double e = HepMC::HEPEVT_Wrapper::e( nentries+ipnw );
438  double m = HepMC::HEPEVT_Wrapper::m( nentries+ipnw );
439 
440  HepMC::GenParticle* NewPart = new HepMC::GenParticle( HepMC::FourVector(px,py,pz,e),
441  pdg_id, status);
442  NewPart->set_generated_mass( m );
443  NewPart->suggest_barcode( nbcode );
444  vtx->add_particle_out( NewPart ) ;
445  }
446 
447  //vtx->print();
448  //std::cout << " leaving attachParticles() " << std::endl;
449 
450  } // end of global if-statement
451 
452  return;
453 }
454 
455 bool PhotosInterface::isTauLeptonicDecay( HepMC::GenVertex* vtx )
456 {
457 
458  if ( abs((*(vtx->particles_in_const_begin()))->pdg_id()) != 15 ) return false;
459 
460  for ( HepMC::GenVertex::particle_iterator pitr=vtx->particles_begin(HepMC::children);
461  pitr != vtx->particles_end(HepMC::children); ++pitr)
462  {
463  if ( abs((*pitr)->pdg_id()) == 11 || abs((*pitr)->pdg_id()) == 13 )
464  {
465  return true;
466  }
467  }
468 
469  return false;
470 
471 }
472 
474 {
475  if ( !fRandomEngine ) {
476  throw cms::Exception("LogicError")
477  << "PhotosInterface::flat: Attempt to generate random number when engine pointer is null\n"
478  << "This might mean that the code was modified to generate a random number outside the\n"
479  << "event and beginLuminosityBlock methods, which is not allowed.\n";
480  }
481  return fRandomEngine->flat();
482 }
483 
484 
void attachParticles(HepMC::GenEvent *, HepMC::GenVertex *, int)
CLHEP::HepRandomEngine * decayRandomEngine
void photos_(int &)
std::vector< int > fBarcodes
void setRandomEngine(CLHEP::HepRandomEngine *decayRandomEngine)
void phoini_(void)
bool qedrad[4000]
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
static CLHEP::HepRandomEngine * fRandomEngine
std::vector< vec3 > vec4
Definition: HCALResponse.h:18
void applyToBranch(HepMC::GenEvent *, int)
struct @494 phoqed_
std::vector< int > fSecVtxStore
void applyToVertex(HepMC::GenEvent *, int)
#define DEFINE_EDM_PLUGIN(factory, type, name)
tuple status
Definition: ntuplemaker.py:245
bool isTauLeptonicDecay(HepMC::GenVertex *)
static double flat()
double phoran_(int *idummy)
HepMC::GenEvent * apply(HepMC::GenEvent *)
std::vector< std::string > fSpecialSettings