CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | Static Private Attributes | Friends
gen::TauolaInterface Class Reference

#include <TauolaInterface.h>

Public Member Functions

HepMC::GenEvent * decay (HepMC::GenEvent *)
 
void disablePolarization ()
 
void enablePolarization ()
 
void init (const edm::EventSetup &)
 
const std::vector< int > & operatesOnParticles ()
 
void setPSet (const edm::ParameterSet &)
 
void setRandomEngine (CLHEP::HepRandomEngine *v)
 
void statistics ()
 
 ~TauolaInterface ()
 

Static Public Member Functions

static TauolaInterfacegetInstance ()
 

Private Member Functions

void decodeMDTAU (int)
 
float flat ()
 
void selectDecayByMDTAU ()
 
int selectHadronic ()
 
int selectLeptonic ()
 
 TauolaInterface ()
 

Private Attributes

std::vector< int > fHadronModes
 
bool fIsInitialized
 
std::vector< int > fLeptonModes
 
int fMDTAU
 
std::vector< int > fPDGs
 
edm::ESHandle
< HepPDT::ParticleDataTable
fPDGTable
 
bool fPolarization
 
edm::ParameterSetfPSet
 
CLHEP::HepRandomEngine * fRandomEngine
 
std::vector< double > fScaledHadronBrRatios
 
std::vector< double > fScaledLeptonBrRatios
 
bool fSelectDecayByEvent
 

Static Private Attributes

static TauolaInterfacefInstance = 0
 

Friends

void gen::ranmar_ (float *rvec, int *lenv)
 
double gen::TauolappInterface_RandGetter ()
 

Detailed Description

Definition at line 64 of file TauolaInterface.h.

Constructor & Destructor Documentation

TauolaInterface::~TauolaInterface ( )

Definition at line 366 of file TauolaInterface.cc.

References fInstance, and fPSet.

367 {
368 
369  if ( fPSet != 0 ) delete fPSet;
370  if ( fInstance == this ) fInstance = 0;
371 
372 }
static TauolaInterface * fInstance
edm::ParameterSet * fPSet
TauolaInterface::TauolaInterface ( )
private

Definition at line 320 of file TauolaInterface.cc.

Referenced by getInstance().

321  : fRandomEngine(nullptr), fPolarization(false), fPSet(0),
322  fIsInitialized(false), fMDTAU(-1), fSelectDecayByEvent(false)
323 {
324 }
CLHEP::HepRandomEngine * fRandomEngine
edm::ParameterSet * fPSet

Member Function Documentation

HepMC::GenEvent * TauolaInterface::decay ( HepMC::GenEvent *  evt)

Definition at line 491 of file TauolaInterface.cc.

References funct::abs(), fIsInitialized, fPDGTable, fSelectDecayByEvent, configurableAnalysis::GenParticle, create_public_lumi_plots::log, p1, RecoTau_DiTaus_pt_20-420_cfg::ParticleID, mix_2012_Summer_inTimeOnly_cff::prob, gen::ranmar_(), and selectDecayByMDTAU().

Referenced by gen::ExternalDecayDriver::decay(), and ParticleReplacerZtautau::produce().

492 {
493 
494  if ( !fIsInitialized ) return evt;
495 
496  int NPartBefore = evt->particles_size();
497  int NVtxBefore = evt->vertices_size();
498 
499  // what do we do if Hep::GenEvent size is larger than 10K ???
500  // Tauola (& Photos, BTW) can only handle up to 10K via HEPEVT,
501  // and in case of CMS, it's only up to 4K !!!
502  //
503  // if ( NPartBefore > 10000 ) return evt;
504  //
505 
506  // override decay mode if needs be
507  if ( fSelectDecayByEvent )
508  {
510  }
511 
512  //construct tmp TAUOLA event
513  //
514  auto * t_event = new Tauolapp::TauolaHepMCEvent(evt);
515 
516  // another option: if one lets Pythia or another master gen to decay taus,
517  // we have to undecay them first
518  // t_event->undecayTaus();
519 
520  // run Tauola on the tmp event - HepMC::GenEvernt will be MODIFIED !!!
521  //
522  t_event->decayTaus();
523 
524  // delet tmp Tauola event
525  //
526  delete t_event;
527 
528  // do we also need to apply the lifetime and vtx position shift ???
529  // (see TauolaInterface, for example)
530  //
531  // NOTE: the procedure ASSYMES that vertex barcoding is COUNTIUOUS/SEQUENTIAL,
532  // and that the abs(barcode) corresponds to vertex "plain indexing"
533  //
534  for ( int iv=NVtxBefore+1; iv<=evt->vertices_size(); iv++ )
535  {
536  HepMC::GenVertex* GenVtx = evt->barcode_to_vertex(-iv);
537  HepMC::GenParticle* GenPart = *(GenVtx->particles_in_const_begin());
538  HepMC::GenVertex* ProdVtx = GenPart->production_vertex();
539  HepMC::FourVector PMom = GenPart->momentum();
540  double mass = GenPart->generated_mass();
541  const HepPDT::ParticleData*
542  PData = fPDGTable->particle(HepPDT::ParticleID(abs(GenPart->pdg_id()))) ;
543  double lifetime = PData->lifetime().value();
544  float prob = 0.;
545  int length=1;
546  ranmar_(&prob,&length);
547  double ct = -lifetime * std::log(prob);
548  double VxDec = GenVtx->position().x();
549  VxDec += ct * (PMom.px()/mass);
550  VxDec += ProdVtx->position().x();
551  double VyDec = GenVtx->position().y();
552  VyDec += ct * (PMom.py()/mass);
553  VyDec += ProdVtx->position().y();
554  double VzDec = GenVtx->position().z();
555  VzDec += ct * (PMom.pz()/mass);
556  VzDec += ProdVtx->position().z();
557  double VtDec = GenVtx->position().t();
558  VtDec += ct * (PMom.e()/mass);
559  VtDec += ProdVtx->position().t();
560  GenVtx->set_position( HepMC::FourVector(VxDec,VyDec,VzDec,VtDec) );
561  //
562  // now find decay products with funky barcode, weed out and replace with clones of sensible barcode
563  // we can NOT change the barcode while iterating, because iterators do depend on the barcoding
564  // thus we have to take a 2-step procedure
565  //
566  std::vector<int> BCodes;
567  BCodes.clear();
568  for (HepMC::GenVertex::particle_iterator pitr= GenVtx->particles_begin(HepMC::children);
569  pitr != GenVtx->particles_end(HepMC::children); ++pitr)
570  {
571  if ( (*pitr)->barcode() > 10000 )
572  {
573  BCodes.push_back( (*pitr)->barcode() );
574  }
575  }
576  if ( BCodes.size() > 0 )
577  {
578  for ( size_t ibc=0; ibc<BCodes.size(); ibc++ )
579  {
580  HepMC::GenParticle* p1 = evt->barcode_to_particle( BCodes[ibc] );
581  int nbc = p1->barcode() - 10000 + NPartBefore;
582  p1->suggest_barcode( nbc );
583  }
584  }
585  }
586 
587  return evt;
588 
589 }
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void ranmar_(float *rvec, int *lenv)
HepPDT::ParticleData ParticleData
double p1[4]
Definition: TauolaWrapper.h:89
void TauolaInterface::decodeMDTAU ( int  mdtau)
private

Definition at line 596 of file TauolaInterface.cc.

References fHadronModes, fLeptonModes, fScaledHadronBrRatios, fScaledLeptonBrRatios, fSelectDecayByEvent, gamprt, i, reco::return(), and taubra_.

Referenced by init().

597 {
598 
599  // Note-1:
600  // I have to hack the common block directly because set<...>DecayMode(...)
601  // only changes it in the Tauola++ instance but does NOT passes it over
602  // to the Fortran core - this it does only one, via initialize() stuff...
603  //
604  // So I'll do both ways of settings, just for consistency...
605  // but I probably need to communicate it to the Tauola(++) team...
606  //
607 
608  // Note-2:
609  // originally, the 1xx settings are meant for tau's from hard event,
610  // and the 2xx settings are for any tau in the event record;
611  //
612  // later one, we'll have to take this into account...
613  // but first I'll have to sort out what happens in the 1xx case
614  // to tau's coming outside of hard event (if any in the record)
615  //
616 
617  if ( mdtau == 101 || mdtau == 201 )
618  {
619  // override with electron mode for both tau's
620  //
621  Tauolapp::jaki_.jak1 = 1;
622  Tauolapp::jaki_.jak2 = 1;
623  Tauolapp::Tauola::setSameParticleDecayMode( 1 ) ;
624  Tauolapp::Tauola::setOppositeParticleDecayMode( 1 ) ;
625  return;
626  }
627 
628  if ( mdtau == 102 || mdtau == 202 )
629  {
630  // override with muon mode for both tau's
631  //
632  Tauolapp::jaki_.jak1 = 2;
633  Tauolapp::jaki_.jak2 = 2;
634  Tauolapp::Tauola::setSameParticleDecayMode( 2 ) ;
635  Tauolapp::Tauola::setOppositeParticleDecayMode( 2 ) ;
636  return;
637  }
638 
639  if ( mdtau == 111 || mdtau == 211 )
640  {
641  // override with electron mode for 1st tau
642  // and any mode for 2nd tau
643  //
644  Tauolapp::jaki_.jak1 = 1;
645  Tauolapp::jaki_.jak2 = 0;
646  Tauolapp::Tauola::setSameParticleDecayMode( 1 ) ;
647  Tauolapp::Tauola::setOppositeParticleDecayMode( 0 ) ;
648  return;
649  }
650 
651  if ( mdtau == 112 || mdtau == 212 )
652  {
653  // override with muon mode for the 1st tau
654  // and any mode for the 2nd tau
655  //
656  Tauolapp::jaki_.jak1 = 2;
657  Tauolapp::jaki_.jak2 = 0;
658  Tauolapp::Tauola::setSameParticleDecayMode( 2 ) ;
659  Tauolapp::Tauola::setOppositeParticleDecayMode( 0 ) ;
660  return;
661  }
662 
663  if ( mdtau == 121 || mdtau == 221 )
664  {
665  // override with any mode for the 1st tau
666  // and electron mode for the 2nd tau
667  //
668  Tauolapp::jaki_.jak1 = 0;
669  Tauolapp::jaki_.jak2 = 1;
670  Tauolapp::Tauola::setSameParticleDecayMode( 0 ) ;
671  Tauolapp::Tauola::setOppositeParticleDecayMode( 1 ) ;
672  return;
673  }
674 
675  if ( mdtau == 122 || mdtau == 222 )
676  {
677  // override with any mode for the 1st tau
678  // and muon mode for the 2nd tau
679  //
680  Tauolapp::jaki_.jak1 = 0;
681  Tauolapp::jaki_.jak2 = 2;
682  Tauolapp::Tauola::setSameParticleDecayMode( 0 ) ;
683  Tauolapp::Tauola::setOppositeParticleDecayMode( 2 ) ;
684  return;
685  }
686 
687  if ( mdtau == 140 || mdtau == 240 )
688  {
689  // override with pi+/- nutau mode for both tau's
690  //
691  Tauolapp::jaki_.jak1 = 3;
692  Tauolapp::jaki_.jak2 = 3;
693  Tauolapp::Tauola::setSameParticleDecayMode( 3 ) ;
694  Tauolapp::Tauola::setOppositeParticleDecayMode( 3 ) ;
695  return;
696  }
697 
698  if ( mdtau == 141 || mdtau == 241 )
699  {
700  // override with pi+/- nutau mode for the 1st tau
701  // and any mode for the 2nd tau
702  //
703  Tauolapp::jaki_.jak1 = 3;
704  Tauolapp::jaki_.jak2 = 0;
705  Tauolapp::Tauola::setSameParticleDecayMode( 3 ) ;
706  Tauolapp::Tauola::setOppositeParticleDecayMode( 0 ) ;
707  return;
708  }
709 
710  if ( mdtau == 142 || mdtau == 242 )
711  {
712  // override with any mode for the 1st tau
713  // and pi+/- nutau mode for 2nd tau
714  //
715  Tauolapp::jaki_.jak1 = 0;
716  Tauolapp::jaki_.jak2 = 3;
717  Tauolapp::Tauola::setSameParticleDecayMode( 0 ) ;
718  Tauolapp::Tauola::setOppositeParticleDecayMode( 3 ) ;
719  return;
720  }
721 
722  // OK, we come here for semi-inclusive modes
723  //
724 
725  // First of all, leptons and hadron modes sums
726  //
727  // re-scale branching ratios, just in case...
728  //
729  double sumBra = 0;
730 
731  // the number of decay modes is hardcoded at 22 because that's what it is right now in Tauola
732  // in the future, perhaps an asscess method would be useful - communicate to Tauola team...
733  //
734 
735  for ( int i=0; i<22; i++ )
736  {
737  sumBra += Tauolapp::taubra_.gamprt[i];
738  }
739  if ( sumBra == 0. ) return ; // perhaps need to throw ?
740  for ( int i=0; i<22; i++ )
741  {
742  double newBra = Tauolapp::taubra_.gamprt[i] / sumBra;
743  Tauolapp::Tauola::setTauBr( i+1, newBra );
744  }
745  sumBra = 1.0;
746 
747  double sumLeptonBra = Tauolapp::taubra_.gamprt[0] + Tauolapp::taubra_.gamprt[1];
748  double sumHadronBra = sumBra - sumLeptonBra;
749 
750  for ( int i=0; i<2; i++ )
751  {
752  fLeptonModes.push_back( i+1 );
753  fScaledLeptonBrRatios.push_back( (Tauolapp::taubra_.gamprt[i]/sumLeptonBra) );
754  }
755  for ( int i=2; i<22; i++ )
756  {
757  fHadronModes.push_back( i+1 );
758  fScaledHadronBrRatios.push_back( (Tauolapp::taubra_.gamprt[i]/sumHadronBra) );
759  }
760 
761  fSelectDecayByEvent = true;
762  return;
763 
764 }
int i
Definition: DBlmapReader.cc:9
std::vector< double > fScaledLeptonBrRatios
float gamprt[30]
Definition: TauolaWrapper.h:78
std::vector< int > fLeptonModes
return(e1-e2)*(e1-e2)+dp *dp
std::vector< double > fScaledHadronBrRatios
std::vector< int > fHadronModes
int mdtau
Definition: TauolaWrapper.h:59
struct @357 taubra_
void gen::TauolaInterface::disablePolarization ( )
inline

Definition at line 75 of file TauolaInterface.h.

References fPolarization.

75 { fPolarization = false; return; }
void gen::TauolaInterface::enablePolarization ( )
inline

Definition at line 74 of file TauolaInterface.h.

References fPolarization.

74 { fPolarization = true; return; }
float TauolaInterface::flat ( )
private

Definition at line 463 of file TauolaInterface.cc.

References edm::hlt::Exception, fIsInitialized, fPSet, and fRandomEngine.

Referenced by gen::ranmar_(), selectLeptonic(), and gen::TauolappInterface_RandGetter().

464 {
465 
466  if ( !fPSet )
467  {
468  // throw
469  throw cms::Exception("TauolaInterfaceError")
470  << "Attempt to run random number generator of un-initialized Tauola\n"
471  << std::endl;
472  }
473 
474  if ( !fIsInitialized )
475  {
476  // throw
477  throw cms::Exception("TauolaInterfaceError")
478  << "Attempt to run random number generator of un-initialized Tauola\n"
479  << std::endl;
480  }
481 
482  if ( !fRandomEngine ) {
483  throw cms::Exception("LogicError")
484  << "TauolaInterface::flat: Attempt to generate random number when engine pointer is null\n"
485  << "This might mean that the code was modified to generate a random number outside the\n"
486  << "event and beginLuminosityBlock methods, which is not allowed.\n";
487  }
488  return fRandomEngine->flat();
489 }
CLHEP::HepRandomEngine * fRandomEngine
edm::ParameterSet * fPSet
TauolaInterface * TauolaInterface::getInstance ( )
static

Definition at line 357 of file TauolaInterface.cc.

References fInstance, and TauolaInterface().

Referenced by gen::ExternalDecayDriver::ExternalDecayDriver(), and gen::TauolappInterface_RandGetter().

358 {
359 
360  if ( fInstance == 0 ) fInstance = new TauolaInterface() ;
361  return fInstance;
362 
363 }
static TauolaInterface * fInstance
void TauolaInterface::init ( const edm::EventSetup es)

Definition at line 390 of file TauolaInterface.cc.

References decodeMDTAU(), edm::hlt::Exception, fIsInitialized, fMDTAU, fPDGs, fPDGTable, fPolarization, fPSet, edm::EventSetup::getData(), edm::ParameterSet::getParameter(), and gen::TauolappInterface_RandGetter().

Referenced by ParticleReplacerZtautau::beginRun(), and gen::ExternalDecayDriver::init().

391 {
392 
393  if ( fIsInitialized ) return; // do init only once
394 
395  if ( fPSet == 0 )
396  {
397 
398  throw cms::Exception("TauolaInterfaceError")
399  << "Attempt to initialize Tauola with an empty ParameterSet\n"
400  << std::endl;
401  }
402 
403  fIsInitialized = true;
404 
405  es.getData( fPDGTable ) ;
406 
407  Tauolapp::Tauola::setDecayingParticle(15);
408  // --> ??? Tauola::setRadiation(false);
409 
410  // polarization switch
411  //
412  // fPolarization = fPSet->getParameter<bool>("UseTauolaPolarization") ? 1 : 0 ;
413  fPolarization = fPSet->getParameter<bool>("UseTauolaPolarization");
414 
415  // read tau decay mode switches
416  //
417  ParameterSet cards = fPSet->getParameter< ParameterSet >("InputCards");
418 
419  fMDTAU = cards.getParameter< int >( "mdtau" );
420 
421  if ( fMDTAU == 0 || fMDTAU == 1 )
422  {
423  Tauolapp::Tauola::setSameParticleDecayMode( cards.getParameter< int >( "pjak1" ) ) ;
424  Tauolapp::Tauola::setOppositeParticleDecayMode( cards.getParameter< int >( "pjak2" ) ) ;
425  }
426 
427  Tauolapp::Tauola::setTauLifetime(0.0);
428  Tauolapp::Tauola::spin_correlation.setAll(fPolarization);
429 
430  // some more options, copied over from an example
431  // - maybe will use later...
432  //
433  //Tauola::setEtaK0sPi(0,0,0); // switches to decay eta K0_S and pi0 1/0 on/off.
434  //
435 
436 //
437 // const HepPDT::ParticleData*
438 // PData = fPDGTable->particle(HepPDT::ParticleID( abs(Tauola::getDecayingParticle()) )) ;
439 // double lifetime = PData->lifetime().value();
440 // Tauola::setTauLifetime( lifetime );
441 
442  fPDGs.push_back( Tauolapp::Tauola::getDecayingParticle() );
443 
444  Tauolapp::Tauola::setRandomGenerator(&gen::TauolappInterface_RandGetter);
445  Tauolapp::Tauola::initialize();
446 
447  Tauolapp::Tauola::spin_correlation.setAll(fPolarization);// Tauola switches this on during Tauola::initialise(); so we add this here to keep it on/off
448 
449  // override decay modes if needs be
450  //
451  // we have to do it AFTER init because otherwises branching ratios are NOT filled in
452  //
453  if ( fMDTAU != 0 && fMDTAU != 1 )
454  {
455  decodeMDTAU( fMDTAU );
456  }
457 
458  Tauolapp::Log::LogWarning(false);
459 
460  return;
461 }
T getParameter(std::string const &) const
void getData(T &iHolder) const
Definition: EventSetup.h:67
double TauolappInterface_RandGetter()
edm::ESHandle< HepPDT::ParticleDataTable > fPDGTable
edm::ParameterSet * fPSet
std::vector< int > fPDGs
const std::vector<int>& gen::TauolaInterface::operatesOnParticles ( )
inline

Definition at line 77 of file TauolaInterface.h.

References fPDGs.

Referenced by gen::ExternalDecayDriver::init().

77 { return fPDGs; }
std::vector< int > fPDGs
void TauolaInterface::selectDecayByMDTAU ( )
private

Definition at line 766 of file TauolaInterface.cc.

References fMDTAU, jak2, alignBH_cfg::mode, reco::return(), selectHadronic(), and selectLeptonic().

Referenced by decay().

767 {
768 
769 
770  if ( fMDTAU == 100 || fMDTAU == 200 )
771  {
772  int mode = selectLeptonic();
773  Tauolapp::jaki_.jak1 = mode;
774  Tauolapp::Tauola::setSameParticleDecayMode( mode );
775  mode = selectLeptonic();
776  Tauolapp::jaki_.jak2 = mode;
777  Tauolapp::Tauola::setOppositeParticleDecayMode( mode );
778  return ;
779  }
780 
781  int modeL = selectLeptonic();
782  int modeH = selectHadronic();
783 
784  if ( fMDTAU == 110 || fMDTAU == 210 )
785  {
786  Tauolapp::jaki_.jak1 = modeL;
787  Tauolapp::jaki_.jak2 = 0;
788  Tauolapp::Tauola::setSameParticleDecayMode( modeL );
789  Tauolapp::Tauola::setOppositeParticleDecayMode( 0 );
790  return ;
791  }
792 
793  if ( fMDTAU == 120 || fMDTAU == 22 )
794  {
795  Tauolapp::jaki_.jak1 = 0;
796  Tauolapp::jaki_.jak2 = modeL;
797  Tauolapp::Tauola::setSameParticleDecayMode( 0 );
798  Tauolapp::Tauola::setOppositeParticleDecayMode( modeL );
799  return;
800  }
801 
802  if ( fMDTAU == 114 || fMDTAU == 214 )
803  {
804  Tauolapp::jaki_.jak1 = modeL;
805  Tauolapp::jaki_.jak2 = modeH;
806  Tauolapp::Tauola::setSameParticleDecayMode( modeL );
807  Tauolapp::Tauola::setOppositeParticleDecayMode( modeH );
808  return;
809  }
810 
811  if ( fMDTAU == 124 || fMDTAU == 224 )
812  {
813  Tauolapp::jaki_.jak1 = modeH;
814  Tauolapp::jaki_.jak2 = modeL;
815  Tauolapp::Tauola::setSameParticleDecayMode( modeH );
816  Tauolapp::Tauola::setOppositeParticleDecayMode( modeL );
817  return;
818  }
819 
820  if ( fMDTAU == 115 || fMDTAU == 215 )
821  {
822  Tauolapp::jaki_.jak1 = 1;
823  Tauolapp::jaki_.jak2 = modeH;
824  Tauolapp::Tauola::setSameParticleDecayMode( 1 );
825  Tauolapp::Tauola::setOppositeParticleDecayMode( modeH );
826  return;
827  }
828 
829  if ( fMDTAU == 125 || fMDTAU == 225 )
830  {
831  Tauolapp::jaki_.jak1 = modeH;
832  Tauolapp::jaki_.jak2 = 1;
833  Tauolapp::Tauola::setSameParticleDecayMode( modeH );
834  Tauolapp::Tauola::setOppositeParticleDecayMode( 1 );
835  return;
836  }
837 
838  if ( fMDTAU == 116 || fMDTAU == 216 )
839  {
840  Tauolapp::jaki_.jak1 = 2;
841  Tauolapp::jaki_.jak2 = modeH;
842  Tauolapp::Tauola::setSameParticleDecayMode( 2 );
843  Tauolapp::Tauola::setOppositeParticleDecayMode( modeH );
844  return;
845  }
846 
847  if ( fMDTAU == 126 || fMDTAU == 226 )
848  {
849  Tauolapp::jaki_.jak1 = modeH;
850  Tauolapp::jaki_.jak2 = 2;
851  Tauolapp::Tauola::setSameParticleDecayMode( modeH );
852  Tauolapp::Tauola::setOppositeParticleDecayMode( 2 );
853  return;
854  }
855 
856  if ( fMDTAU == 130 || fMDTAU == 230 )
857  {
858  Tauolapp::jaki_.jak1 = modeH;
859  Tauolapp::jaki_.jak2 = selectHadronic();
860  Tauolapp::Tauola::setSameParticleDecayMode( modeH );
861  Tauolapp::Tauola::setOppositeParticleDecayMode( Tauolapp::jaki_.jak2 );
862  return;
863  }
864 
865  if ( fMDTAU == 131 || fMDTAU == 231 )
866  {
867  Tauolapp::jaki_.jak1 = modeH;
868  Tauolapp::jaki_.jak2 = 0;
869  Tauolapp::Tauola::setSameParticleDecayMode( modeH );
870  Tauolapp::Tauola::setOppositeParticleDecayMode( 0 );
871  return;
872  }
873 
874  if ( fMDTAU == 132 || fMDTAU == 232 )
875  {
876  Tauolapp::jaki_.jak1 = 0;
877  Tauolapp::jaki_.jak2 = modeH;
878  Tauolapp::Tauola::setSameParticleDecayMode( 0 );
879  Tauolapp::Tauola::setOppositeParticleDecayMode( modeH );
880  return;
881  }
882 
883  // unlikely that we get here on unknown mdtau
884  // - there's a protection earlier
885  // but if we do, just set defaults
886  // probably need to spit a warning...
887  //
888  Tauolapp::Tauola::setSameParticleDecayMode( 0 );
889  Tauolapp::Tauola::setOppositeParticleDecayMode( 0 );
890 
891  return;
892 
893 
894 }
int jak2
Definition: TauolaWrapper.h:67
return(e1-e2)*(e1-e2)+dp *dp
int TauolaInterface::selectHadronic ( )
private

Definition at line 913 of file TauolaInterface.cc.

References fHadronModes, fScaledHadronBrRatios, i, mix_2012_Summer_inTimeOnly_cff::prob, and gen::ranmar_().

Referenced by selectDecayByMDTAU().

914 {
915 
916  float prob = 0.;
917  int len = 1;
918  ranmar_(&prob,&len);
919 
920  double sumBra = fScaledHadronBrRatios[0];
921  if ( prob > 0. && prob <= sumBra )
922  {
923  return fHadronModes[0];
924  }
925  else
926  {
927  int NN = fScaledHadronBrRatios.size();
928  for ( int i=1; i<NN; i++ )
929  {
930  if ( prob > sumBra && prob <= (sumBra+fScaledHadronBrRatios[i]) )
931  {
932  return fHadronModes[i];
933  }
934  sumBra += fScaledHadronBrRatios[i];
935  }
936  }
937 
938  return 0;
939 
940 }
int i
Definition: DBlmapReader.cc:9
void ranmar_(float *rvec, int *lenv)
std::vector< double > fScaledHadronBrRatios
std::vector< int > fHadronModes
int TauolaInterface::selectLeptonic ( )
private

Definition at line 896 of file TauolaInterface.cc.

References flat(), fScaledLeptonBrRatios, and mix_2012_Summer_inTimeOnly_cff::prob.

Referenced by selectDecayByMDTAU().

897 {
898 
899  float prob = flat();
900 
901  if ( prob > 0. && prob <= fScaledLeptonBrRatios[0] )
902  {
903  return 1;
904  }
905  else if ( prob > fScaledLeptonBrRatios[1] && prob <=1. )
906  {
907  return 2;
908  }
909 
910  return 0;
911 }
std::vector< double > fScaledLeptonBrRatios
void TauolaInterface::setPSet ( const edm::ParameterSet pset)

Definition at line 374 of file TauolaInterface.cc.

References edm::hlt::Exception, and fPSet.

Referenced by gen::ExternalDecayDriver::ExternalDecayDriver(), ParticleReplacerParticleGun::ParticleReplacerParticleGun(), and ParticleReplacerZtautau::ParticleReplacerZtautau().

375 {
376 
377  if ( fPSet != 0 )
378  {
379  throw cms::Exception("TauolaInterfaceError")
380  << "Attempt to override Tauola an existing ParameterSet\n"
381  << std::endl;
382  }
383 
384  fPSet = new ParameterSet(pset);
385 
386  return;
387 
388 }
edm::ParameterSet * fPSet
void gen::TauolaInterface::setRandomEngine ( CLHEP::HepRandomEngine *  v)
inline

Definition at line 81 of file TauolaInterface.h.

References fRandomEngine, and gen::v.

Referenced by gen::ExternalDecayDriver::setRandomEngine().

81 { fRandomEngine = v; }
double v[5][pyjets_maxn]
CLHEP::HepRandomEngine * fRandomEngine
void TauolaInterface::statistics ( )

Definition at line 591 of file TauolaInterface.cc.

Referenced by ParticleReplacerZtautau::endJob(), and gen::ExternalDecayDriver::statistics().

592 {
593  return;
594 }

Friends And Related Function Documentation

void gen::ranmar_ ( float *  rvec,
int *  lenv 
)
friend

Member Data Documentation

std::vector<int> gen::TauolaInterface::fHadronModes
private

Definition at line 109 of file TauolaInterface.h.

Referenced by decodeMDTAU(), and selectHadronic().

TauolaInterface * TauolaInterface::fInstance = 0
staticprivate

Definition at line 113 of file TauolaInterface.h.

Referenced by getInstance(), and ~TauolaInterface().

bool gen::TauolaInterface::fIsInitialized
private

Definition at line 104 of file TauolaInterface.h.

Referenced by decay(), flat(), and init().

std::vector<int> gen::TauolaInterface::fLeptonModes
private

Definition at line 108 of file TauolaInterface.h.

Referenced by decodeMDTAU().

int gen::TauolaInterface::fMDTAU
private

Definition at line 106 of file TauolaInterface.h.

Referenced by init(), and selectDecayByMDTAU().

std::vector<int> gen::TauolaInterface::fPDGs
private

Definition at line 100 of file TauolaInterface.h.

Referenced by init(), and operatesOnParticles().

edm::ESHandle<HepPDT::ParticleDataTable> gen::TauolaInterface::fPDGTable
private

Definition at line 102 of file TauolaInterface.h.

Referenced by decay(), and init().

bool gen::TauolaInterface::fPolarization
private

Definition at line 101 of file TauolaInterface.h.

Referenced by disablePolarization(), enablePolarization(), and init().

edm::ParameterSet* gen::TauolaInterface::fPSet
private

Definition at line 103 of file TauolaInterface.h.

Referenced by flat(), init(), setPSet(), and ~TauolaInterface().

CLHEP::HepRandomEngine* gen::TauolaInterface::fRandomEngine
private

Definition at line 99 of file TauolaInterface.h.

Referenced by flat(), and setRandomEngine().

std::vector<double> gen::TauolaInterface::fScaledHadronBrRatios
private

Definition at line 111 of file TauolaInterface.h.

Referenced by decodeMDTAU(), and selectHadronic().

std::vector<double> gen::TauolaInterface::fScaledLeptonBrRatios
private

Definition at line 110 of file TauolaInterface.h.

Referenced by decodeMDTAU(), and selectLeptonic().

bool gen::TauolaInterface::fSelectDecayByEvent
private

Definition at line 107 of file TauolaInterface.h.

Referenced by decay(), and decodeMDTAU().