CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TauolappInterface.cc
Go to the documentation of this file.
1 /* this is the code for the new Tauola++ */
2 
3 #include <iostream>
4 
6 
7 #include "Tauola/Tauola.h"
8 #include "Tauola/TauolaHepMCEvent.h"
9 #include "Tauola/Log.h"
10 
14 
15 #include "CLHEP/Random/RandomEngine.h"
16 
17 #include "HepMC/GenEvent.h"
18 #include "HepMC/IO_HEPEVT.h"
19 #include "HepMC/HEPEVT_Wrapper.h"
20 
22 
23 
24 using namespace gen;
25 using namespace edm;
26 using namespace std;
27 
28 CLHEP::HepRandomEngine* TauolappInterface::fRandomEngine = nullptr;
29 
30 extern "C" {
31 
32  void gen::ranmar_( float *rvec, int *lenv )
33  {
34  for(int i = 0; i < *lenv; i++)
35  *rvec++ = TauolappInterface::flat();
36  return;
37  }
38 
39  void gen::rmarin_( int*, int*, int* )
40  {
41  return;
42  }
43 
44 }
45 
47  fPolarization(false),
48  fPSet(0),
49  fIsInitialized(false),
50  fMDTAU(-1),
51  fSelectDecayByEvent(false)
52 {
53  if ( fPSet != 0 ) throw cms::Exception("TauolappInterfaceError") << "Attempt to override Tauola an existing ParameterSet\n" << std::endl;
54  fPSet = new ParameterSet(pset);
55 }
56 
58 
60  if ( fIsInitialized ) return; // do init only once
61  if ( fPSet == 0 ) throw cms::Exception("TauolappInterfaceError") << "Attempt to initialize Tauola with an empty ParameterSet\n" << std::endl;
62 
63  fIsInitialized = true;
64 
65  es.getData( fPDGTable ) ;
66 
67  Tauolapp::Tauola::setDecayingParticle(15);
68  // --> ??? Tauola::setRadiation(false);
69 
70  // polarization switch
71  //
72  // fPolarization = fPSet->getParameter<bool>("UseTauolaPolarization") ? 1 : 0 ;
73  fPolarization = fPSet->getParameter<bool>("UseTauolaPolarization");
74 
75  // read tau decay mode switches
76  //
77  ParameterSet cards = fPSet->getParameter< ParameterSet >("InputCards");
78 
79  fMDTAU = cards.getParameter< int >( "mdtau" );
80 
81  if ( fMDTAU == 0 || fMDTAU == 1 )
82  {
83  Tauolapp::Tauola::setSameParticleDecayMode( cards.getParameter< int >( "pjak1" ) ) ;
84  Tauolapp::Tauola::setOppositeParticleDecayMode( cards.getParameter< int >( "pjak2" ) ) ;
85  }
86 
87  Tauolapp::Tauola::spin_correlation.setAll(fPolarization);
88 
89  // some more options, copied over from an example
90  // Default values
91  //Tauola::setEtaK0sPi(0,0,0); // switches to decay eta K0_S and pi0 1/0 on/off.
92 
93  const HepPDT::ParticleData* PData = fPDGTable->particle(HepPDT::ParticleID( abs(Tauolapp::Tauola::getDecayingParticle()) )) ;
94  double lifetime = PData->lifetime().value();
95  Tauolapp::Tauola::setTauLifetime( lifetime );
96  //std::cout << "Setting lifetime: ctau=" << lifetime << "m or tau=" << lifetime/(299792458.0*1000.0) << "s" << std::endl;
97 
98  fPDGs.push_back( Tauolapp::Tauola::getDecayingParticle() );
99 
100  Tauolapp::Tauola::setRandomGenerator(gen::TauolappInterface::flat);
102 
103  Tauolapp::Tauola::spin_correlation.setAll(fPolarization);// Tauola switches this on during Tauola::initialise(); so we add this here to keep it on/off
104 
105  if (fPSet->exists("parameterSets")){
106  std::vector<std::string> par = fPSet->getParameter< std::vector<std::string> >("parameterSets");
107  for (unsigned int ip=0; ip<par.size(); ++ip ){
108  std::string curSet = par[ip];
109 
110  if(curSet=="setNewCurrents") Tauolapp::Tauola::setNewCurrents(fPSet->getParameter<int>(curSet));
111  if(curSet=="spinCorrelationSetAll") Tauolapp::Tauola::spin_correlation.setAll(fPSet->getParameter<bool>(curSet));
112  if(curSet=="spinCorrelationGAMMA") Tauolapp::Tauola::spin_correlation.GAMMA=fPSet->getParameter<bool>(curSet);
113  if(curSet=="spinCorrelationZ0") Tauolapp::Tauola::spin_correlation.Z0=fPSet->getParameter<bool>(curSet);
114  if(curSet=="spinCorrelationHIGGS") Tauolapp::Tauola::spin_correlation.HIGGS=fPSet->getParameter<bool>(curSet);
115  if(curSet=="spinCorrelationHIGGSH") Tauolapp::Tauola::spin_correlation.HIGGS_H=fPSet->getParameter<bool>(curSet);
116  if(curSet=="spinCorrelationHIGGSA") Tauolapp::Tauola::spin_correlation.HIGGS_A=fPSet->getParameter<bool>(curSet);
117  if(curSet=="spinCorrelationHIGGSPLUS") Tauolapp::Tauola::spin_correlation.HIGGS_PLUS=fPSet->getParameter<bool>(curSet);
118  if(curSet=="spinCorrelationHIGGSMINUS") Tauolapp::Tauola::spin_correlation.HIGGS_MINUS=fPSet->getParameter<bool>(curSet);
119  if(curSet=="spinCorrelationWPLUS") Tauolapp::Tauola::spin_correlation.W_PLUS=fPSet->getParameter<bool>(curSet);
120  if(curSet=="spinCorrelationWMINUS") Tauolapp::Tauola::spin_correlation.W_MINUS=fPSet->getParameter<bool>(curSet);
121 
122  if(curSet=="setHiggsScalarPseudoscalarPDG") Tauolapp::Tauola::setHiggsScalarPseudoscalarPDG(fPSet->getParameter<int>(curSet));
123  if(curSet=="setHiggsScalarPseudoscalarMixingAngle") Tauolapp::Tauola::setHiggsScalarPseudoscalarMixingAngle(fPSet->getParameter<double>(curSet));
124 
125  if(curSet=="setRadiation") Tauolapp::Tauola::setRadiation(fPSet->getParameter<bool>(curSet));
126  if(curSet=="setRadiationCutOff") Tauolapp::Tauola::setRadiationCutOff(fPSet->getParameter<double>(curSet));
127 
128  if(curSet=="setEtaK0sPi"){
129  std::vector<int> vpar = fPSet->getParameter<std::vector<int> >(curSet);
130  if(vpar.size()==3) Tauolapp::Tauola::setEtaK0sPi(vpar[0],vpar[1],vpar[2]);
131  else {std::cout << "WARNING invalid size for setEtaK0sPi: " << vpar.size() << " Require 3 elements " << std::endl;}
132  }
133 
134  if(curSet=="setTaukle"){
135  std::vector<double> vpar = fPSet->getParameter<std::vector<double> >(curSet);
136  if(vpar.size()==4) Tauolapp::Tauola::setTaukle(vpar[0], vpar[1], vpar[2], vpar[3]);
137  else {std::cout << "WARNING invalid size for setTaukle: " << vpar.size() << " Require 4 elements " << std::endl;}
138  }
139 
140  if(curSet=="setTauBr"){
142  std::vector<int> vJAK = cfg.getParameter<std::vector<int> >("JAK");
143  std::vector<double> vBR = cfg.getParameter<std::vector<double> >("BR");
144  if(vJAK.size() == vBR.size()){
145  for(unsigned int i=0;i<vJAK.size();i++) Tauolapp::Tauola::setTauBr(vJAK[i],vBR[i]);
146  }
147  else {std::cout << "WARNING invalid size for setTauBr - JAK: " << vJAK.size() << " BR: " << vBR.size() << std::endl;}
148  }
149  }
150  }
151 
152  // override decay modes if needs be
153  //
154  // we have to do it AFTER init because otherwises branching ratios are NOT filled in
155  //
156  if ( fMDTAU != 0 && fMDTAU != 1 )
157  {
158  decodeMDTAU( fMDTAU );
159  }
160 
161  Tauolapp::Log::LogWarning(false);
162 
163  return;
164 }
165 
167 {
168  if ( !fRandomEngine ) {
169  throw cms::Exception("LogicError")
170  << "TauolaInterface::flat: Attempt to generate random number when engine pointer is null\n"
171  << "This might mean that the code was modified to generate a random number outside the\n"
172  << "event and beginLuminosityBlock methods, which is not allowed.\n";
173  }
174  return fRandomEngine->flat();
175 }
176 
177 HepMC::GenEvent* TauolappInterface::decay( HepMC::GenEvent* evt ){
178  if ( !fIsInitialized ) return evt;
179  Tauolapp::Tauola::setRandomGenerator(gen::TauolappInterface::flat); // rest tauola++ random number incase other modules use tauola++
180  int NPartBefore = evt->particles_size();
181  int NVtxBefore = evt->vertices_size();
182 
183  // what do we do if Hep::GenEvent size is larger than 10K ???
184  // Tauola (& Photos, BTW) can only handle up to 10K via HEPEVT,
185  // and in case of CMS, it's only up to 4K !!!
186  //
187  // if ( NPartBefore > 10000 ) return evt;
188  //
189 
190  // override decay mode if needs be
191  if ( fSelectDecayByEvent )
192  {
194  }
195 
196  //construct tmp TAUOLA event
197  //
198  auto * t_event = new Tauolapp::TauolaHepMCEvent(evt);
199 
200  // another option: if one lets Pythia or another master gen to decay taus,
201  // we have to undecay them first
202  // t_event->undecayTaus();
203 
204  // run Tauola on the tmp event - HepMC::GenEvernt will be MODIFIED !!!
205  //
206  t_event->decayTaus();
207 
208  // delet tmp Tauola event
209  //
210  delete t_event;
211 
212  // do we also need to apply the lifetime and vtx position shift ???
213  // (see TauolappInterface, for example)
214  //
215  // NOTE: the procedure ASSYMES that vertex barcoding is COUNTIUOUS/SEQUENTIAL,
216  // and that the abs(barcode) corresponds to vertex "plain indexing"
217  //
218  for ( int iv=NVtxBefore+1; iv<=evt->vertices_size(); iv++ )
219  {
220  HepMC::GenVertex* GenVtx = evt->barcode_to_vertex(-iv);
221  //
222  // now find decay products with funky barcode, weed out and replace with clones of sensible barcode
223  // we can NOT change the barcode while iterating, because iterators do depend on the barcoding
224  // thus we have to take a 2-step procedure
225  //
226  std::vector<int> BCodes;
227  BCodes.clear();
228  for (HepMC::GenVertex::particle_iterator pitr= GenVtx->particles_begin(HepMC::children);
229  pitr != GenVtx->particles_end(HepMC::children); ++pitr)
230  {
231  if ( (*pitr)->barcode() > 10000 )
232  {
233  BCodes.push_back( (*pitr)->barcode() );
234  }
235  }
236  if ( BCodes.size() > 0 )
237  {
238  for ( size_t ibc=0; ibc<BCodes.size(); ibc++ )
239  {
240  HepMC::GenParticle* p1 = evt->barcode_to_particle( BCodes[ibc] );
241  int nbc = p1->barcode() - 10000 + NPartBefore;
242  p1->suggest_barcode( nbc );
243  }
244  }
245  }
246 
247  return evt;
248 
249 }
250 
252 {
253  return;
254 }
255 
257 {
258 
259  // Note-1:
260  // I have to hack the common block directly because set<...>DecayMode(...)
261  // only changes it in the Tauola++ instance but does NOT passes it over
262  // to the Fortran core - this it does only one, via initialize() stuff...
263  //
264  // So I'll do both ways of settings, just for consistency...
265  // but I probably need to communicate it to the Tauola(++) team...
266  //
267 
268  // Note-2:
269  // originally, the 1xx settings are meant for tau's from hard event,
270  // and the 2xx settings are for any tau in the event record;
271  //
272  // later one, we'll have to take this into account...
273  // but first I'll have to sort out what happens in the 1xx case
274  // to tau's coming outside of hard event (if any in the record)
275  //
276 
277  if ( mdtau == 101 || mdtau == 201 )
278  {
279  // override with electron mode for both tau's
280  //
281  Tauolapp::jaki_.jak1 = 1;
282  Tauolapp::jaki_.jak2 = 1;
283  Tauolapp::Tauola::setSameParticleDecayMode( 1 ) ;
284  Tauolapp::Tauola::setOppositeParticleDecayMode( 1 ) ;
285  return;
286  }
287 
288  if ( mdtau == 102 || mdtau == 202 )
289  {
290  // override with muon mode for both tau's
291  //
292  Tauolapp::jaki_.jak1 = 2;
293  Tauolapp::jaki_.jak2 = 2;
294  Tauolapp::Tauola::setSameParticleDecayMode( 2 ) ;
295  Tauolapp::Tauola::setOppositeParticleDecayMode( 2 ) ;
296  return;
297  }
298 
299  if ( mdtau == 111 || mdtau == 211 )
300  {
301  // override with electron mode for 1st tau
302  // and any mode for 2nd tau
303  //
304  Tauolapp::jaki_.jak1 = 1;
305  Tauolapp::jaki_.jak2 = 0;
306  Tauolapp::Tauola::setSameParticleDecayMode( 1 ) ;
307  Tauolapp::Tauola::setOppositeParticleDecayMode( 0 ) ;
308  return;
309  }
310 
311  if ( mdtau == 112 || mdtau == 212 )
312  {
313  // override with muon mode for the 1st tau
314  // and any mode for the 2nd tau
315  //
316  Tauolapp::jaki_.jak1 = 2;
317  Tauolapp::jaki_.jak2 = 0;
318  Tauolapp::Tauola::setSameParticleDecayMode( 2 ) ;
319  Tauolapp::Tauola::setOppositeParticleDecayMode( 0 ) ;
320  return;
321  }
322 
323  if ( mdtau == 121 || mdtau == 221 )
324  {
325  // override with any mode for the 1st tau
326  // and electron mode for the 2nd tau
327  //
328  Tauolapp::jaki_.jak1 = 0;
329  Tauolapp::jaki_.jak2 = 1;
330  Tauolapp::Tauola::setSameParticleDecayMode( 0 ) ;
331  Tauolapp::Tauola::setOppositeParticleDecayMode( 1 ) ;
332  return;
333  }
334 
335  if ( mdtau == 122 || mdtau == 222 )
336  {
337  // override with any mode for the 1st tau
338  // and muon mode for the 2nd tau
339  //
340  Tauolapp::jaki_.jak1 = 0;
341  Tauolapp::jaki_.jak2 = 2;
342  Tauolapp::Tauola::setSameParticleDecayMode( 0 ) ;
343  Tauolapp::Tauola::setOppositeParticleDecayMode( 2 ) ;
344  return;
345  }
346 
347  if ( mdtau == 140 || mdtau == 240 )
348  {
349  // override with pi+/- nutau mode for both tau's
350  //
351  Tauolapp::jaki_.jak1 = 3;
352  Tauolapp::jaki_.jak2 = 3;
353  Tauolapp::Tauola::setSameParticleDecayMode( 3 ) ;
354  Tauolapp::Tauola::setOppositeParticleDecayMode( 3 ) ;
355  return;
356  }
357 
358  if ( mdtau == 141 || mdtau == 241 )
359  {
360  // override with pi+/- nutau mode for the 1st tau
361  // and any mode for the 2nd tau
362  //
363  Tauolapp::jaki_.jak1 = 3;
364  Tauolapp::jaki_.jak2 = 0;
365  Tauolapp::Tauola::setSameParticleDecayMode( 3 ) ;
366  Tauolapp::Tauola::setOppositeParticleDecayMode( 0 ) ;
367  return;
368  }
369 
370  if ( mdtau == 142 || mdtau == 242 )
371  {
372  // override with any mode for the 1st tau
373  // and pi+/- nutau mode for 2nd tau
374  //
375  Tauolapp::jaki_.jak1 = 0;
376  Tauolapp::jaki_.jak2 = 3;
377  Tauolapp::Tauola::setSameParticleDecayMode( 0 ) ;
378  Tauolapp::Tauola::setOppositeParticleDecayMode( 3 ) ;
379  return;
380  }
381 
382  // OK, we come here for semi-inclusive modes
383  //
384 
385  // First of all, leptons and hadron modes sums
386  //
387  // re-scale branching ratios, just in case...
388  //
389  double sumBra = 0;
390 
391  // the number of decay modes is hardcoded at 22 because that's what it is right now in Tauola
392  // in the future, perhaps an asscess method would be useful - communicate to Tauola team...
393  //
394 
395  for ( int i=0; i<22; i++ )
396  {
397  sumBra += Tauolapp::taubra_.gamprt[i];
398  }
399  if ( sumBra == 0. ) return ; // perhaps need to throw ?
400  for ( int i=0; i<22; i++ )
401  {
402  double newBra = Tauolapp::taubra_.gamprt[i] / sumBra;
403  Tauolapp::Tauola::setTauBr( i+1, newBra );
404  }
405  sumBra = 1.0;
406 
407  double sumLeptonBra = Tauolapp::taubra_.gamprt[0] + Tauolapp::taubra_.gamprt[1];
408  double sumHadronBra = sumBra - sumLeptonBra;
409 
410  for ( int i=0; i<2; i++ )
411  {
412  fLeptonModes.push_back( i+1 );
413  fScaledLeptonBrRatios.push_back( (Tauolapp::taubra_.gamprt[i]/sumLeptonBra) );
414  }
415  for ( int i=2; i<22; i++ )
416  {
417  fHadronModes.push_back( i+1 );
418  fScaledHadronBrRatios.push_back( (Tauolapp::taubra_.gamprt[i]/sumHadronBra) );
419  }
420 
421  fSelectDecayByEvent = true;
422  return;
423 
424 }
425 
427 {
428 
429 
430  if ( fMDTAU == 100 || fMDTAU == 200 )
431  {
432  int mode = selectLeptonic();
433  Tauolapp::jaki_.jak1 = mode;
434  Tauolapp::Tauola::setSameParticleDecayMode( mode );
435  mode = selectLeptonic();
436  Tauolapp::jaki_.jak2 = mode;
437  Tauolapp::Tauola::setOppositeParticleDecayMode( mode );
438  return ;
439  }
440 
441  int modeL = selectLeptonic();
442  int modeH = selectHadronic();
443 
444  if ( fMDTAU == 110 || fMDTAU == 210 )
445  {
446  Tauolapp::jaki_.jak1 = modeL;
447  Tauolapp::jaki_.jak2 = 0;
448  Tauolapp::Tauola::setSameParticleDecayMode( modeL );
449  Tauolapp::Tauola::setOppositeParticleDecayMode( 0 );
450  return ;
451  }
452 
453  if ( fMDTAU == 120 || fMDTAU == 22 )
454  {
455  Tauolapp::jaki_.jak1 = 0;
456  Tauolapp::jaki_.jak2 = modeL;
457  Tauolapp::Tauola::setSameParticleDecayMode( 0 );
458  Tauolapp::Tauola::setOppositeParticleDecayMode( modeL );
459  return;
460  }
461 
462  if ( fMDTAU == 114 || fMDTAU == 214 )
463  {
464  Tauolapp::jaki_.jak1 = modeL;
465  Tauolapp::jaki_.jak2 = modeH;
466  Tauolapp::Tauola::setSameParticleDecayMode( modeL );
467  Tauolapp::Tauola::setOppositeParticleDecayMode( modeH );
468  return;
469  }
470 
471  if ( fMDTAU == 124 || fMDTAU == 224 )
472  {
473  Tauolapp::jaki_.jak1 = modeH;
474  Tauolapp::jaki_.jak2 = modeL;
475  Tauolapp::Tauola::setSameParticleDecayMode( modeH );
476  Tauolapp::Tauola::setOppositeParticleDecayMode( modeL );
477  return;
478  }
479 
480  if ( fMDTAU == 115 || fMDTAU == 215 )
481  {
482  Tauolapp::jaki_.jak1 = 1;
483  Tauolapp::jaki_.jak2 = modeH;
484  Tauolapp::Tauola::setSameParticleDecayMode( 1 );
485  Tauolapp::Tauola::setOppositeParticleDecayMode( modeH );
486  return;
487  }
488 
489  if ( fMDTAU == 125 || fMDTAU == 225 )
490  {
491  Tauolapp::jaki_.jak1 = modeH;
492  Tauolapp::jaki_.jak2 = 1;
493  Tauolapp::Tauola::setSameParticleDecayMode( modeH );
494  Tauolapp::Tauola::setOppositeParticleDecayMode( 1 );
495  return;
496  }
497 
498  if ( fMDTAU == 116 || fMDTAU == 216 )
499  {
500  Tauolapp::jaki_.jak1 = 2;
501  Tauolapp::jaki_.jak2 = modeH;
502  Tauolapp::Tauola::setSameParticleDecayMode( 2 );
503  Tauolapp::Tauola::setOppositeParticleDecayMode( modeH );
504  return;
505  }
506 
507  if ( fMDTAU == 126 || fMDTAU == 226 )
508  {
509  Tauolapp::jaki_.jak1 = modeH;
510  Tauolapp::jaki_.jak2 = 2;
511  Tauolapp::Tauola::setSameParticleDecayMode( modeH );
512  Tauolapp::Tauola::setOppositeParticleDecayMode( 2 );
513  return;
514  }
515 
516  if ( fMDTAU == 130 || fMDTAU == 230 )
517  {
518  Tauolapp::jaki_.jak1 = modeH;
519  Tauolapp::jaki_.jak2 = selectHadronic();
520  Tauolapp::Tauola::setSameParticleDecayMode( modeH );
521  Tauolapp::Tauola::setOppositeParticleDecayMode( Tauolapp::jaki_.jak2 );
522  return;
523  }
524 
525  if ( fMDTAU == 131 || fMDTAU == 231 )
526  {
527  Tauolapp::jaki_.jak1 = modeH;
528  Tauolapp::jaki_.jak2 = 0;
529  Tauolapp::Tauola::setSameParticleDecayMode( modeH );
530  Tauolapp::Tauola::setOppositeParticleDecayMode( 0 );
531  return;
532  }
533 
534  if ( fMDTAU == 132 || fMDTAU == 232 )
535  {
536  Tauolapp::jaki_.jak1 = 0;
537  Tauolapp::jaki_.jak2 = modeH;
538  Tauolapp::Tauola::setSameParticleDecayMode( 0 );
539  Tauolapp::Tauola::setOppositeParticleDecayMode( modeH );
540  return;
541  }
542 
543  // unlikely that we get here on unknown mdtau
544  // - there's a protection earlier
545  // but if we do, just set defaults
546  // probably need to spit a warning...
547  //
548  Tauolapp::Tauola::setSameParticleDecayMode( 0 );
549  Tauolapp::Tauola::setOppositeParticleDecayMode( 0 );
550 
551  return;
552 
553 
554 }
555 
557 {
558 
559  float prob = flat();
560 
561  if ( prob > 0. && prob <= fScaledLeptonBrRatios[0] )
562  {
563  return 1;
564  }
565  else if ( prob > fScaledLeptonBrRatios[1] && prob <=1. )
566  {
567  return 2;
568  }
569 
570  return 0;
571 }
572 
574 {
575 
576  float prob = 0.;
577  int len = 1;
578  ranmar_(&prob,&len);
579 
580  double sumBra = fScaledHadronBrRatios[0];
581  if ( prob > 0. && prob <= sumBra )
582  {
583  return fHadronModes[0];
584  }
585  else
586  {
587  int NN = fScaledHadronBrRatios.size();
588  for ( int i=1; i<NN; i++ )
589  {
590  if ( prob > sumBra && prob <= (sumBra+fScaledHadronBrRatios[i]) )
591  {
592  return fHadronModes[i];
593  }
594  sumBra += fScaledHadronBrRatios[i];
595  }
596  }
597 
598  return 0;
599 
600 }
601 
602 
T getParameter(std::string const &) const
static AlgebraicMatrix initialize()
TauolappInterface(const edm::ParameterSet &)
int i
Definition: DBlmapReader.cc:9
edm::ParameterSet * fPSet
std::vector< double > fScaledHadronBrRatios
bool exists(std::string const &parameterName) const
checks if a parameter exists
void rmarin_(int *, int *, int *)
void getData(T &iHolder) const
Definition: EventSetup.h:67
float gamprt[30]
Definition: TauolaWrapper.h:78
int jak2
Definition: TauolaWrapper.h:67
std::vector< int > fPDGs
std::vector< int > fLeptonModes
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void ranmar_(float *rvec, int *lenv)
HepPDT::ParticleData ParticleData
void init(const edm::EventSetup &)
HepMC::GenEvent * decay(HepMC::GenEvent *)
std::vector< double > fScaledLeptonBrRatios
static CLHEP::HepRandomEngine * fRandomEngine
return(e1-e2)*(e1-e2)+dp *dp
struct @413 taubra_
std::vector< int > fHadronModes
double p1[4]
Definition: TauolaWrapper.h:89
tuple cout
Definition: gather_cfg.py:121
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
volatile std::atomic< bool > shutdown_flag false
int mdtau
Definition: TauolaWrapper.h:59