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