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::setTauLifetime(0.0);
88  Tauolapp::Tauola::spin_correlation.setAll(fPolarization);
89 
90  // some more options, copied over from an example
91  // - maybe will use later...
92  //
93  //Tauola::setEtaK0sPi(0,0,0); // switches to decay eta K0_S and pi0 1/0 on/off.
94  //
95 
96 //
97 // const HepPDT::ParticleData*
98 // PData = fPDGTable->particle(HepPDT::ParticleID( abs(Tauola::getDecayingParticle()) )) ;
99 // double lifetime = PData->lifetime().value();
100 // Tauola::setTauLifetime( lifetime );
101 
102  fPDGs.push_back( Tauolapp::Tauola::getDecayingParticle() );
103 
104  Tauolapp::Tauola::setRandomGenerator(gen::TauolappInterface::flat);
106 
107  Tauolapp::Tauola::spin_correlation.setAll(fPolarization);// Tauola switches this on during Tauola::initialise(); so we add this here to keep it on/off
108 
109  // override decay modes if needs be
110  //
111  // we have to do it AFTER init because otherwises branching ratios are NOT filled in
112  //
113  if ( fMDTAU != 0 && fMDTAU != 1 )
114  {
115  decodeMDTAU( fMDTAU );
116  }
117 
118  Tauolapp::Log::LogWarning(false);
119 
120  return;
121 }
122 
124 {
125  if ( !fRandomEngine ) {
126  throw cms::Exception("LogicError")
127  << "TauolaInterface::flat: Attempt to generate random number when engine pointer is null\n"
128  << "This might mean that the code was modified to generate a random number outside the\n"
129  << "event and beginLuminosityBlock methods, which is not allowed.\n";
130  }
131  return fRandomEngine->flat();
132 }
133 
134 HepMC::GenEvent* TauolappInterface::decay( HepMC::GenEvent* evt ){
135  if ( !fIsInitialized ) return evt;
136  Tauolapp::Tauola::setRandomGenerator(gen::TauolappInterface::flat); // rest tauola++ random number incase other modules use tauola++
137  int NPartBefore = evt->particles_size();
138  int NVtxBefore = evt->vertices_size();
139 
140  // what do we do if Hep::GenEvent size is larger than 10K ???
141  // Tauola (& Photos, BTW) can only handle up to 10K via HEPEVT,
142  // and in case of CMS, it's only up to 4K !!!
143  //
144  // if ( NPartBefore > 10000 ) return evt;
145  //
146 
147  // override decay mode if needs be
148  if ( fSelectDecayByEvent )
149  {
151  }
152 
153  //construct tmp TAUOLA event
154  //
155  auto * t_event = new Tauolapp::TauolaHepMCEvent(evt);
156 
157  // another option: if one lets Pythia or another master gen to decay taus,
158  // we have to undecay them first
159  // t_event->undecayTaus();
160 
161  // run Tauola on the tmp event - HepMC::GenEvernt will be MODIFIED !!!
162  //
163  t_event->decayTaus();
164 
165  // delet tmp Tauola event
166  //
167  delete t_event;
168 
169  // do we also need to apply the lifetime and vtx position shift ???
170  // (see TauolappInterface, for example)
171  //
172  // NOTE: the procedure ASSYMES that vertex barcoding is COUNTIUOUS/SEQUENTIAL,
173  // and that the abs(barcode) corresponds to vertex "plain indexing"
174  //
175  for ( int iv=NVtxBefore+1; iv<=evt->vertices_size(); iv++ )
176  {
177  HepMC::GenVertex* GenVtx = evt->barcode_to_vertex(-iv);
178  HepMC::GenParticle* GenPart = *(GenVtx->particles_in_const_begin());
179  HepMC::GenVertex* ProdVtx = GenPart->production_vertex();
180  HepMC::FourVector PMom = GenPart->momentum();
181  double mass = GenPart->generated_mass();
182  const HepPDT::ParticleData*
183  PData = fPDGTable->particle(HepPDT::ParticleID(abs(GenPart->pdg_id()))) ;
184  double lifetime = PData->lifetime().value();
185  float prob = 0.;
186  int length=1;
187  ranmar_(&prob,&length);
188  double ct = -lifetime * std::log(prob);
189  double VxDec = GenVtx->position().x();
190  VxDec += ct * (PMom.px()/mass);
191  VxDec += ProdVtx->position().x();
192  double VyDec = GenVtx->position().y();
193  VyDec += ct * (PMom.py()/mass);
194  VyDec += ProdVtx->position().y();
195  double VzDec = GenVtx->position().z();
196  VzDec += ct * (PMom.pz()/mass);
197  VzDec += ProdVtx->position().z();
198  double VtDec = GenVtx->position().t();
199  VtDec += ct * (PMom.e()/mass);
200  VtDec += ProdVtx->position().t();
201  GenVtx->set_position( HepMC::FourVector(VxDec,VyDec,VzDec,VtDec) );
202  //
203  // now find decay products with funky barcode, weed out and replace with clones of sensible barcode
204  // we can NOT change the barcode while iterating, because iterators do depend on the barcoding
205  // thus we have to take a 2-step procedure
206  //
207  std::vector<int> BCodes;
208  BCodes.clear();
209  for (HepMC::GenVertex::particle_iterator pitr= GenVtx->particles_begin(HepMC::children);
210  pitr != GenVtx->particles_end(HepMC::children); ++pitr)
211  {
212  if ( (*pitr)->barcode() > 10000 )
213  {
214  BCodes.push_back( (*pitr)->barcode() );
215  }
216  }
217  if ( BCodes.size() > 0 )
218  {
219  for ( size_t ibc=0; ibc<BCodes.size(); ibc++ )
220  {
221  HepMC::GenParticle* p1 = evt->barcode_to_particle( BCodes[ibc] );
222  int nbc = p1->barcode() - 10000 + NPartBefore;
223  p1->suggest_barcode( nbc );
224  }
225  }
226  }
227 
228  return evt;
229 
230 }
231 
233 {
234  return;
235 }
236 
238 {
239 
240  // Note-1:
241  // I have to hack the common block directly because set<...>DecayMode(...)
242  // only changes it in the Tauola++ instance but does NOT passes it over
243  // to the Fortran core - this it does only one, via initialize() stuff...
244  //
245  // So I'll do both ways of settings, just for consistency...
246  // but I probably need to communicate it to the Tauola(++) team...
247  //
248 
249  // Note-2:
250  // originally, the 1xx settings are meant for tau's from hard event,
251  // and the 2xx settings are for any tau in the event record;
252  //
253  // later one, we'll have to take this into account...
254  // but first I'll have to sort out what happens in the 1xx case
255  // to tau's coming outside of hard event (if any in the record)
256  //
257 
258  if ( mdtau == 101 || mdtau == 201 )
259  {
260  // override with electron mode for both tau's
261  //
262  Tauolapp::jaki_.jak1 = 1;
263  Tauolapp::jaki_.jak2 = 1;
264  Tauolapp::Tauola::setSameParticleDecayMode( 1 ) ;
265  Tauolapp::Tauola::setOppositeParticleDecayMode( 1 ) ;
266  return;
267  }
268 
269  if ( mdtau == 102 || mdtau == 202 )
270  {
271  // override with muon mode for both tau's
272  //
273  Tauolapp::jaki_.jak1 = 2;
274  Tauolapp::jaki_.jak2 = 2;
275  Tauolapp::Tauola::setSameParticleDecayMode( 2 ) ;
276  Tauolapp::Tauola::setOppositeParticleDecayMode( 2 ) ;
277  return;
278  }
279 
280  if ( mdtau == 111 || mdtau == 211 )
281  {
282  // override with electron mode for 1st tau
283  // and any mode for 2nd tau
284  //
285  Tauolapp::jaki_.jak1 = 1;
286  Tauolapp::jaki_.jak2 = 0;
287  Tauolapp::Tauola::setSameParticleDecayMode( 1 ) ;
288  Tauolapp::Tauola::setOppositeParticleDecayMode( 0 ) ;
289  return;
290  }
291 
292  if ( mdtau == 112 || mdtau == 212 )
293  {
294  // override with muon mode for the 1st tau
295  // and any mode for the 2nd tau
296  //
297  Tauolapp::jaki_.jak1 = 2;
298  Tauolapp::jaki_.jak2 = 0;
299  Tauolapp::Tauola::setSameParticleDecayMode( 2 ) ;
300  Tauolapp::Tauola::setOppositeParticleDecayMode( 0 ) ;
301  return;
302  }
303 
304  if ( mdtau == 121 || mdtau == 221 )
305  {
306  // override with any mode for the 1st tau
307  // and electron mode for the 2nd tau
308  //
309  Tauolapp::jaki_.jak1 = 0;
310  Tauolapp::jaki_.jak2 = 1;
311  Tauolapp::Tauola::setSameParticleDecayMode( 0 ) ;
312  Tauolapp::Tauola::setOppositeParticleDecayMode( 1 ) ;
313  return;
314  }
315 
316  if ( mdtau == 122 || mdtau == 222 )
317  {
318  // override with any mode for the 1st tau
319  // and muon mode for the 2nd tau
320  //
321  Tauolapp::jaki_.jak1 = 0;
322  Tauolapp::jaki_.jak2 = 2;
323  Tauolapp::Tauola::setSameParticleDecayMode( 0 ) ;
324  Tauolapp::Tauola::setOppositeParticleDecayMode( 2 ) ;
325  return;
326  }
327 
328  if ( mdtau == 140 || mdtau == 240 )
329  {
330  // override with pi+/- nutau mode for both tau's
331  //
332  Tauolapp::jaki_.jak1 = 3;
333  Tauolapp::jaki_.jak2 = 3;
334  Tauolapp::Tauola::setSameParticleDecayMode( 3 ) ;
335  Tauolapp::Tauola::setOppositeParticleDecayMode( 3 ) ;
336  return;
337  }
338 
339  if ( mdtau == 141 || mdtau == 241 )
340  {
341  // override with pi+/- nutau mode for the 1st tau
342  // and any mode for the 2nd tau
343  //
344  Tauolapp::jaki_.jak1 = 3;
345  Tauolapp::jaki_.jak2 = 0;
346  Tauolapp::Tauola::setSameParticleDecayMode( 3 ) ;
347  Tauolapp::Tauola::setOppositeParticleDecayMode( 0 ) ;
348  return;
349  }
350 
351  if ( mdtau == 142 || mdtau == 242 )
352  {
353  // override with any mode for the 1st tau
354  // and pi+/- nutau mode for 2nd tau
355  //
356  Tauolapp::jaki_.jak1 = 0;
357  Tauolapp::jaki_.jak2 = 3;
358  Tauolapp::Tauola::setSameParticleDecayMode( 0 ) ;
359  Tauolapp::Tauola::setOppositeParticleDecayMode( 3 ) ;
360  return;
361  }
362 
363  // OK, we come here for semi-inclusive modes
364  //
365 
366  // First of all, leptons and hadron modes sums
367  //
368  // re-scale branching ratios, just in case...
369  //
370  double sumBra = 0;
371 
372  // the number of decay modes is hardcoded at 22 because that's what it is right now in Tauola
373  // in the future, perhaps an asscess method would be useful - communicate to Tauola team...
374  //
375 
376  for ( int i=0; i<22; i++ )
377  {
378  sumBra += Tauolapp::taubra_.gamprt[i];
379  }
380  if ( sumBra == 0. ) return ; // perhaps need to throw ?
381  for ( int i=0; i<22; i++ )
382  {
383  double newBra = Tauolapp::taubra_.gamprt[i] / sumBra;
384  Tauolapp::Tauola::setTauBr( i+1, newBra );
385  }
386  sumBra = 1.0;
387 
388  double sumLeptonBra = Tauolapp::taubra_.gamprt[0] + Tauolapp::taubra_.gamprt[1];
389  double sumHadronBra = sumBra - sumLeptonBra;
390 
391  for ( int i=0; i<2; i++ )
392  {
393  fLeptonModes.push_back( i+1 );
394  fScaledLeptonBrRatios.push_back( (Tauolapp::taubra_.gamprt[i]/sumLeptonBra) );
395  }
396  for ( int i=2; i<22; i++ )
397  {
398  fHadronModes.push_back( i+1 );
399  fScaledHadronBrRatios.push_back( (Tauolapp::taubra_.gamprt[i]/sumHadronBra) );
400  }
401 
402  fSelectDecayByEvent = true;
403  return;
404 
405 }
406 
408 {
409 
410 
411  if ( fMDTAU == 100 || fMDTAU == 200 )
412  {
413  int mode = selectLeptonic();
414  Tauolapp::jaki_.jak1 = mode;
415  Tauolapp::Tauola::setSameParticleDecayMode( mode );
416  mode = selectLeptonic();
417  Tauolapp::jaki_.jak2 = mode;
418  Tauolapp::Tauola::setOppositeParticleDecayMode( mode );
419  return ;
420  }
421 
422  int modeL = selectLeptonic();
423  int modeH = selectHadronic();
424 
425  if ( fMDTAU == 110 || fMDTAU == 210 )
426  {
427  Tauolapp::jaki_.jak1 = modeL;
428  Tauolapp::jaki_.jak2 = 0;
429  Tauolapp::Tauola::setSameParticleDecayMode( modeL );
430  Tauolapp::Tauola::setOppositeParticleDecayMode( 0 );
431  return ;
432  }
433 
434  if ( fMDTAU == 120 || fMDTAU == 22 )
435  {
436  Tauolapp::jaki_.jak1 = 0;
437  Tauolapp::jaki_.jak2 = modeL;
438  Tauolapp::Tauola::setSameParticleDecayMode( 0 );
439  Tauolapp::Tauola::setOppositeParticleDecayMode( modeL );
440  return;
441  }
442 
443  if ( fMDTAU == 114 || fMDTAU == 214 )
444  {
445  Tauolapp::jaki_.jak1 = modeL;
446  Tauolapp::jaki_.jak2 = modeH;
447  Tauolapp::Tauola::setSameParticleDecayMode( modeL );
448  Tauolapp::Tauola::setOppositeParticleDecayMode( modeH );
449  return;
450  }
451 
452  if ( fMDTAU == 124 || fMDTAU == 224 )
453  {
454  Tauolapp::jaki_.jak1 = modeH;
455  Tauolapp::jaki_.jak2 = modeL;
456  Tauolapp::Tauola::setSameParticleDecayMode( modeH );
457  Tauolapp::Tauola::setOppositeParticleDecayMode( modeL );
458  return;
459  }
460 
461  if ( fMDTAU == 115 || fMDTAU == 215 )
462  {
463  Tauolapp::jaki_.jak1 = 1;
464  Tauolapp::jaki_.jak2 = modeH;
465  Tauolapp::Tauola::setSameParticleDecayMode( 1 );
466  Tauolapp::Tauola::setOppositeParticleDecayMode( modeH );
467  return;
468  }
469 
470  if ( fMDTAU == 125 || fMDTAU == 225 )
471  {
472  Tauolapp::jaki_.jak1 = modeH;
473  Tauolapp::jaki_.jak2 = 1;
474  Tauolapp::Tauola::setSameParticleDecayMode( modeH );
475  Tauolapp::Tauola::setOppositeParticleDecayMode( 1 );
476  return;
477  }
478 
479  if ( fMDTAU == 116 || fMDTAU == 216 )
480  {
481  Tauolapp::jaki_.jak1 = 2;
482  Tauolapp::jaki_.jak2 = modeH;
483  Tauolapp::Tauola::setSameParticleDecayMode( 2 );
484  Tauolapp::Tauola::setOppositeParticleDecayMode( modeH );
485  return;
486  }
487 
488  if ( fMDTAU == 126 || fMDTAU == 226 )
489  {
490  Tauolapp::jaki_.jak1 = modeH;
491  Tauolapp::jaki_.jak2 = 2;
492  Tauolapp::Tauola::setSameParticleDecayMode( modeH );
493  Tauolapp::Tauola::setOppositeParticleDecayMode( 2 );
494  return;
495  }
496 
497  if ( fMDTAU == 130 || fMDTAU == 230 )
498  {
499  Tauolapp::jaki_.jak1 = modeH;
500  Tauolapp::jaki_.jak2 = selectHadronic();
501  Tauolapp::Tauola::setSameParticleDecayMode( modeH );
502  Tauolapp::Tauola::setOppositeParticleDecayMode( Tauolapp::jaki_.jak2 );
503  return;
504  }
505 
506  if ( fMDTAU == 131 || fMDTAU == 231 )
507  {
508  Tauolapp::jaki_.jak1 = modeH;
509  Tauolapp::jaki_.jak2 = 0;
510  Tauolapp::Tauola::setSameParticleDecayMode( modeH );
511  Tauolapp::Tauola::setOppositeParticleDecayMode( 0 );
512  return;
513  }
514 
515  if ( fMDTAU == 132 || fMDTAU == 232 )
516  {
517  Tauolapp::jaki_.jak1 = 0;
518  Tauolapp::jaki_.jak2 = modeH;
519  Tauolapp::Tauola::setSameParticleDecayMode( 0 );
520  Tauolapp::Tauola::setOppositeParticleDecayMode( modeH );
521  return;
522  }
523 
524  // unlikely that we get here on unknown mdtau
525  // - there's a protection earlier
526  // but if we do, just set defaults
527  // probably need to spit a warning...
528  //
529  Tauolapp::Tauola::setSameParticleDecayMode( 0 );
530  Tauolapp::Tauola::setOppositeParticleDecayMode( 0 );
531 
532  return;
533 
534 
535 }
536 
538 {
539 
540  float prob = flat();
541 
542  if ( prob > 0. && prob <= fScaledLeptonBrRatios[0] )
543  {
544  return 1;
545  }
546  else if ( prob > fScaledLeptonBrRatios[1] && prob <=1. )
547  {
548  return 2;
549  }
550 
551  return 0;
552 }
553 
555 {
556 
557  float prob = 0.;
558  int len = 1;
559  ranmar_(&prob,&len);
560 
561  double sumBra = fScaledHadronBrRatios[0];
562  if ( prob > 0. && prob <= sumBra )
563  {
564  return fHadronModes[0];
565  }
566  else
567  {
568  int NN = fScaledHadronBrRatios.size();
569  for ( int i=1; i<NN; i++ )
570  {
571  if ( prob > sumBra && prob <= (sumBra+fScaledHadronBrRatios[i]) )
572  {
573  return fHadronModes[i];
574  }
575  sumBra += fScaledHadronBrRatios[i];
576  }
577  }
578 
579  return 0;
580 
581 }
582 
583 
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
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
struct @407 taubra_
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
std::vector< int > fHadronModes
double p1[4]
Definition: TauolaWrapper.h:89
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
volatile std::atomic< bool > shutdown_flag false
int mdtau
Definition: TauolaWrapper.h:59