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 
2 #include <iostream>
3 
4 // #include "GeneratorInterface/Pythia6Interface/interface/Pythia6Service.h"
5 
7 
9 // #include "FWCore/Utilities/interface/RandomNumberGenerator.h"
11 
12 #include "HepMC/GenEvent.h"
13 #include "HepMC/IO_HEPEVT.h"
14 #include "HepMC/HEPEVT_Wrapper.h"
15 
16 using namespace gen;
17 using namespace edm;
18 using namespace std;
19 
20 // const int NMXHEP = HepMC::HEPEVT_Wrapper::max_number_entries() ;
21 
22 extern "C"{
23 
24  void phoini_( void );
25  void photos_( int& );
26 
27  double phoran_(int *idummy)
28  {
29  return decayRandomEngine->flat();
30  }
31 /*
32  double phoranc_(int *idummy)
33  {
34  return decayRandomEngine->flat();
35  }
36 */
37 
38  extern struct {
39  // bool qedrad[NMXHEP];
40  bool qedrad[4000]; // hardcoded for now...
41  } phoqed_;
42 
43 }
44 
45 
47  : fOnlyPDG(-1)
48 {
49  fSpecialSettings.push_back("QED-brem-off:all");
51  fIsInitialized = false;
52 }
53 
55  : fOnlyPDG(-1)
56 {
57  fSpecialSettings.push_back("QED-brem-off:all");
58  fIsInitialized = false;
59 }
60 
61 /* -->
62 void PhotosInterface::configureOnlyFor( int ipdg )
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 */
75 
77 {
78 
79  if ( fIsInitialized ) return; // do init only once
80 
81  phoini_();
82 
83  fIsInitialized = true;
84 
85  return;
86 }
87 
88 HepMC::GenEvent* PhotosInterface::apply( HepMC::GenEvent* evt )
89 {
90 
91  // event record convertor
92  // ...well, I'm not sure we need it here,
93  // as we do it by hands, only part of the record
94  //
95  // HepMC::IO_HEPEVT conv;
96 
97  if ( !fIsInitialized ) return evt; // conv.read_next_event();
98 
99  // cross-check printout HepMC::GenEvent
100  //
101  //evt->print();
102 
103  // int numPartBefore = HepMC::HEPEVT_Wrapper::number_entries();
104  // HepMC::HEPEVT_Wrapper::print_hepevt();
105 
106  // loop over HepMC::GenEvent, find vertices
107 
108  // std::vector<int> barcodes;
109  // std::vector<Scaling> scaleFactors;
110 
111  fScaleFactors.clear();
112 
113  for ( int ip=0; ip<evt->particles_size(); ip++ )
114  {
115  fScaleFactors.push_back( Scaling(HepMC::ThreeVector(1.,1.,1.),1) );
116  phoqed_.qedrad[ip]=true;
117  }
118 
119 
120  // variables for special treatment of tau leptonic decays
121  //
122  bool tau_leptonic_decay = false;
123  int iTauDescCounter = 0;
124  int nTauDesc = 0;
125 
126  for ( int iv=1; iv<=evt->vertices_size(); iv++ )
127  {
128 
129  HepMC::GenVertex* vtx = evt->barcode_to_vertex( -iv ) ;
130  if ( vtx->particles_in_size() != 1 ) continue; // more complex than we need
131  if ( vtx->particles_out_size() <= 0 ) continue; // no outcoming particles
132  // --> if ( fOnlyPDG !=-1 && (*(vtx->particles_in_const_begin()))->pdg_id() == fOnlyPDG ) continue;
133  // now find at least one "legal" daughter
134  bool legalVtx = false;
135  for ( HepMC::GenVertex::particle_iterator pitr=vtx->particles_begin(HepMC::children);
136  pitr != vtx->particles_end(HepMC::children); ++pitr)
137  {
138  if ( (*pitr)->status() == 1 || (*pitr)->end_vertex() )
139  {
140  // OK, legal already !
141  legalVtx = true;
142  break;
143  }
144  }
145 
146  if ( !legalVtx ) continue;
147 
148  // now do all the loops again
149  //
150  // first, flush out HEPEVT & tmp barcode storage
151  //
152  HepMC::HEPEVT_Wrapper::zero_everything();
153  fBarcodes.clear();
154 
155  // add incoming particle
156  //
157  int index = 1;
158  HepMC::HEPEVT_Wrapper::set_id( index, (*(vtx->particles_in_const_begin()))->pdg_id() );
159  HepMC::FourVector vec4;
160  vec4 = (*(vtx->particles_in_const_begin()))->momentum();
161  HepMC::HEPEVT_Wrapper::set_momentum( index, vec4.x(), vec4.y(), vec4.z(), vec4.e() );
162  HepMC::HEPEVT_Wrapper::set_mass( index, (*(vtx->particles_in_const_begin()))->generated_mass() );
163  HepMC::HEPEVT_Wrapper::set_position( index, vtx->position().x(), vtx->position().y(),
164  vtx->position().z(), vtx->position().t() );
165  HepMC::HEPEVT_Wrapper::set_status( index, (*(vtx->particles_in_const_begin()))->status() );
166  HepMC::HEPEVT_Wrapper::set_parents( index, 0, 0 );
167  fBarcodes.push_back( (*(vtx->particles_in_const_begin()))->barcode() );
168 
169  // special case: avoid tau leptonic decays
170  //
171 
172  if ( fAvoidTauLeptonicDecays && !tau_leptonic_decay && abs((*(vtx->particles_in_const_begin()))->pdg_id()) == 15 )
173  {
174  for ( HepMC::GenVertex::particle_iterator pitr=vtx->particles_begin(HepMC::children);
175  pitr != vtx->particles_end(HepMC::children); ++pitr)
176  {
177  if ( abs((*pitr)->pdg_id()) == 11 || abs((*pitr)->pdg_id()) == 13 ) // leptonic decay !!!
178  // do brem off tau but NOT off decay products
179  {
180  tau_leptonic_decay = true;
181  break;
182  }
183  }
184  if ( vtx->particles_begin(HepMC::children) == vtx->particles_begin(HepMC::descendants) &&
185  vtx->particles_end(HepMC::children) == vtx->particles_end(HepMC::descendants) )
186  {
187  nTauDesc = vtx->particles_out_size();
188  }
189  else
190  {
191  for ( HepMC::GenVertex::particle_iterator pitr1=vtx->particles_begin(HepMC::children);
192  pitr1 != vtx->particles_end(HepMC::children); ++pitr1)
193  {
194  nTauDesc++;
195  }
196  }
197  // this is just the 1st tau in the branch, so it's allowed to emit
198  phoqed_.qedrad[index-1]=true;
199  iTauDescCounter = 0;
200  }
201 
202  // check if mother has ever been "altered" !
203  //
204  int mbcode = (*(vtx->particles_in_const_begin()))->barcode();
205 
206  // add outcoming particles (decay products)
207  //
208  int lastDau = 1;
209  for ( HepMC::GenVertex::particle_iterator pitr=vtx->particles_begin(HepMC::children);
210  pitr != vtx->particles_end(HepMC::children); ++pitr)
211  {
212  if ( fScaleFactors[mbcode-1].flag != 1. )
213  {
214  // yes, mother has been changed - adjust daughters
215 
216  vec4 = (*pitr)->momentum();
217  double mass2 = vec4.m2();
218  double pxn = vec4.px() * fScaleFactors[mbcode-1].weights.x();
219  double pyn = vec4.py() * fScaleFactors[mbcode-1].weights.y();
220  double pzn = vec4.pz() * fScaleFactors[mbcode-1].weights.z();
221  double en = sqrt( pxn*pxn + pyn*pyn + pzn*pzn + mass2 );
222  (*pitr)->set_momentum( HepMC::FourVector(pxn,pyn,pzn,en) );
223  int curbcode = (*pitr)->barcode();
224  double scale = fScaleFactors[curbcode-1].weights.x();
225  fScaleFactors[curbcode-1].weights.setX( scale*fScaleFactors[mbcode-1].weights.x() );
226  scale = fScaleFactors[curbcode-1].weights.y();
227  fScaleFactors[curbcode-1].weights.setY( scale*fScaleFactors[mbcode-1].weights.y() );
228  scale = fScaleFactors[curbcode-1].weights.z();
229  fScaleFactors[curbcode-1].weights.setZ( scale*fScaleFactors[mbcode-1].weights.z() );
230  fScaleFactors[curbcode-1].flag = 0;
231  }
232 
233  if ( (*pitr)->status() == 1 || (*pitr)->end_vertex() )
234  {
235  index++;
236  vec4 = (*pitr)->momentum();
237  HepMC::HEPEVT_Wrapper::set_id( index, (*pitr)->pdg_id() );
238  HepMC::HEPEVT_Wrapper::set_momentum( index, vec4.x(), vec4.y(), vec4.z(), vec4.e() );
239  HepMC::HEPEVT_Wrapper::set_mass( index, (*pitr)->generated_mass() );
240  vec4 = (*pitr)->production_vertex()->position();
241  HepMC::HEPEVT_Wrapper::set_position( index, vec4.x(), vec4.y(), vec4.z(), vec4.t() );
242  HepMC::HEPEVT_Wrapper::set_status( index, (*pitr)->status() );
243  HepMC::HEPEVT_Wrapper::set_parents( index, 1, 1 );
244  fBarcodes.push_back( (*pitr)->barcode() );
245  lastDau++;
246  if ( fAvoidTauLeptonicDecays && tau_leptonic_decay )
247  {
248  phoqed_.qedrad[index-1]=false;
249  iTauDescCounter++;
250  }
251  }
252  }
253 
254  // store, further to set NHEP in HEPEVT
255  //
256  int nentries = index;
257 
258  // reset master pointer to mother
259  index = 1;
260  HepMC::HEPEVT_Wrapper::set_children ( index, 2, lastDau ); // FIXME: need to check
261  // if last daughter>=2 !!!
262 
263  // finally, set number of entries (NHEP) in HEPEVT
264  //
265  HepMC::HEPEVT_Wrapper::set_number_entries( nentries );
266 
267  // cross-check printout HEPEVT
268  // HepMC::HEPEVT_Wrapper::print_hepevt();
269 
270  // OK, 1-level vertex is formed - now, call PHOTOS
271  //
272  photos_( index ) ;
273 
274  // another cross-check printout HEPEVT - after photos
275  // HepMC::HEPEVT_Wrapper::print_hepevt();
276 
277 
278  // now check if something has been generated
279  //
280  attachParticles( evt, vtx, nentries );
281 
282  // ugh, done with this vertex !
283 
284  // now some resets
285  //
286  if ( fAvoidTauLeptonicDecays && tau_leptonic_decay && iTauDescCounter == nTauDesc ) // treated tau leptonic decay and have come to the last descendent
287  {
288  tau_leptonic_decay = false;
289  }
290 
291  }
292 
293  // restore event number in HEPEVT (safety measure, somehow needed by Hw6)
294  HepMC::HEPEVT_Wrapper::set_event_number( evt->event_number() );
295 
296  // cross-check printout MODIFIED HepMC::GenEvent
297  // evt->print();
298 
299  // return conv.read_next_event();
300  return evt;
301 
302 }
303 
304 void PhotosInterface::attachParticles( HepMC::GenEvent* evt, HepMC::GenVertex* vtx, int nentries )
305 {
306 
307  if ( HepMC::HEPEVT_Wrapper::number_entries() > nentries )
308  {
309  // yes, need all corrections and adjustments -
310  // figure out how many photons and what particles in
311  // the decay branch have changes;
312  // also, follow up each one and correct accordingly;
313  // at the same time, add photon(s) to the GenVertex
314  //
315  int largestBarcode = -1;
316  int Nbcodes = fBarcodes.size();
317  for ( int ip=1; ip<Nbcodes; ip++ )
318  {
319  int bcode = fBarcodes[ip];
320  HepMC::GenParticle* prt = evt->barcode_to_particle( bcode );
321  if ( bcode > largestBarcode ) largestBarcode = bcode;
322  double px = HepMC::HEPEVT_Wrapper::px(ip+1);
323  double py = HepMC::HEPEVT_Wrapper::py(ip+1);
324  double pz = HepMC::HEPEVT_Wrapper::pz(ip+1);
325  double e = HepMC::HEPEVT_Wrapper::e(ip+1);
326  HepMC::FourVector mom4 = prt->momentum();
327  //double porg = sqrt( mom4.px()*mom4.px()
328  // + mom4.py()*mom4.py()
329  // + mom4.pz()*mom4.pz() ) ;
330  //double pnew = sqrt( px*px + py*py + pz*pz );
331  double scale = fScaleFactors[bcode-1].weights.x();
332  fScaleFactors[bcode-1].weights.setX( scale*(px/mom4.px()) );
333  scale = fScaleFactors[bcode-1].weights.y();
334  fScaleFactors[bcode-1].weights.setY( scale*(py/mom4.py()) );
335  scale = fScaleFactors[bcode-1].weights.z();
336  fScaleFactors[bcode-1].weights.setZ( scale*(pz/mom4.pz()) );
337  fScaleFactors[bcode-1].flag = 0;
338 
339  prt->set_momentum( HepMC::FourVector(px,py,pz,e) );
340 
341  // we do NOT adjust chaldren, etc., here - because we do it
342  // above, based on whether mother (incoming particles) has
343  // ever been modified
344 
345  }
346 
347  int newlyGen = HepMC::HEPEVT_Wrapper::number_entries() - nentries;
348 
349  if ( largestBarcode < evt->particles_size() )
350  {
351  // need to adjust barcodes down from the affected vertex/particles
352  // such that we "free up" barcodes for newly generated photons
353  // in the middle of the event record
354  //
355  for ( int ipp=evt->particles_size(); ipp>largestBarcode; ipp-- )
356  {
357  (evt->barcode_to_particle(ipp))->suggest_barcode( ipp+newlyGen );
358  }
359  }
360 
361  // now attach new generated photons to the vertex
362  //
363  for ( int ipnw=1; ipnw<=newlyGen; ipnw++ )
364  {
365  int nbcode = largestBarcode+ipnw;
366  int pdg_id = HepMC::HEPEVT_Wrapper::id( nentries+ipnw );
367  int status = HepMC::HEPEVT_Wrapper::status( nentries+ipnw );
368  double px = HepMC::HEPEVT_Wrapper::px( nentries+ipnw );
369  double py = HepMC::HEPEVT_Wrapper::py( nentries+ipnw );
370  double pz = HepMC::HEPEVT_Wrapper::pz( nentries+ipnw );
371  double e = HepMC::HEPEVT_Wrapper::e( nentries+ipnw );
372  double m = HepMC::HEPEVT_Wrapper::m( nentries+ipnw );
373 
374  HepMC::GenParticle* NewPart = new HepMC::GenParticle( HepMC::FourVector(px,py,pz,e),
375  pdg_id, status);
376  NewPart->set_generated_mass( m );
377  NewPart->suggest_barcode( nbcode );
378  vtx->add_particle_out( NewPart ) ;
379  // add/shift scale factors towards the end of the list
380  fScaleFactors.push_back( Scaling(HepMC::ThreeVector(1.,1.,1.),1) );
381  }
382  } // end of if-statement
383 
384  return;
385 }
long int flag
Definition: mlp_lapack.h:47
void phoini_(void)
void attachParticles(HepMC::GenEvent *, HepMC::GenVertex *, int)
#define abs(x)
Definition: mlp_lapack.h:159
double phoran_(int *idummy)
void photos_(int &)
std::vector< Scaling > fScaleFactors
std::vector< int > fBarcodes
T sqrt(T t)
Definition: SSEVec.h:28
struct @251 phoqed_
std::vector< std::string > fSpecialSettings
bool qedrad[4000]
CLHEP::HepRandomEngine * decayRandomEngine
tuple status
Definition: ntuplemaker.py:245
HepMC::GenEvent * apply(HepMC::GenEvent *)