BOSS 7.0.5
BESIII Offline Software System
Loading...
Searching...
No Matches
EvtParticle Class Referenceabstract

#include <EvtParticle.hh>

+ Inheritance diagram for EvtParticle:

Public Member Functions

 EvtParticle ()
 
virtual ~EvtParticle ()
 
virtual EvtVector4C epsParent (int i) const
 
virtual EvtVector4C eps (int i) const
 
virtual EvtVector4C epsParentPhoton (int i)
 
virtual EvtVector4C epsPhoton (int i)
 
virtual EvtDiracSpinor spParent (int) const
 
virtual EvtDiracSpinor sp (int) const
 
virtual EvtDiracSpinor spParentNeutrino () const
 
virtual EvtDiracSpinor spNeutrino () const
 
virtual EvtTensor4C epsTensorParent (int i) const
 
virtual EvtTensor4C epsTensor (int i) const
 
virtual void init (EvtId part_n, const EvtVector4R &p4)=0
 
void addDaug (EvtParticle *node)
 
void decay ()
 
void deleteTree ()
 
void deleteDaughters (bool keepChannel=false)
 
void setChannel (int i)
 
void setGeneratorFlag (int flag)
 
int getGeneratorFlag ()
 
void setIntFlag (std::vector< int > vi)
 
std::vector< int > getIntFlag ()
 
void makeDaughters (int ndaug, EvtId *id)
 
double initializePhaseSpace (int numdaughter, EvtId *daughters, double poleSize=-1., int whichTwo1=0, int whichTwo2=1)
 
EvtParticlegetDaug (int i)
 
EvtParticlenextIter (EvtParticle *rootOfTree=0)
 
void makeStdHep (EvtStdHep &stdhep, EvtSecondary &secondary, EvtId *stable_parent_ihep)
 
void makeStdHep (EvtStdHep &stdhep)
 
EvtVector4R getP4Lab ()
 
EvtVector4R getP4Restframe ()
 
EvtVector4R get4Pos ()
 
EvtParticlegetParent ()
 
void insertDaugPtr (int idaug, EvtParticle *partptr)
 
double mass () const
 
int firstornot () const
 
void setFirstOrNot ()
 
void resetFirstOrNot ()
 
EvtId getId () const
 
EvtSpinType::spintype getSpinType () const
 
int getSpinStates () const
 
const EvtVector4RgetP4 () const
 
void setP4 (const EvtVector4R &p4)
 
int getChannel () const
 
int getNDaug () const
 
void resetNDaug ()
 
void printTree () const
 
void printTreeRec (int level) const
 
std::string writeTreeRec (std::string) const
 
void dumpTree () const
 
void dumpTreeRec (int level, int dj) const
 
std::string treeStr () const
 
std::string treeStrRec (int level) const
 
void printParticle () const
 
void setLifetime (double tau)
 
void setLifetime ()
 
double getLifetime ()
 
void setDiagonalSpinDensity ()
 
void setVectorSpinDensity ()
 
void setPolarizedSpinDensity (double r00, double r11, double r22)
 
void setSpinDensityForward (const EvtSpinDensity &rho)
 
void setSpinDensityForwardHelicityBasis (const EvtSpinDensity &rho)
 
void setSpinDensityForwardHelicityBasis (const EvtSpinDensity &rho, double alpha, double beta, double gamma)
 
virtual EvtSpinDensity rotateToHelicityBasis () const =0
 
virtual EvtSpinDensity rotateToHelicityBasis (double alpha, double beta, double gamma) const =0
 
EvtSpinDensity getSpinDensityForward ()
 
void setSpinDensityBackward (const EvtSpinDensity &rho)
 
EvtSpinDensity getSpinDensityBackward ()
 
void noLifeTime ()
 
void setId (EvtId id)
 
void initDecay (bool useMinMass=false)
 
void generateMassTree ()
 
double compMassProb ()
 
void setMass (double m)
 
bool isInitialized ()
 
bool hasValidP4 ()
 
bool isDecayed ()
 
double * decayProb ()
 
void setDecayProb (double p)
 
void setInclusiveMode (int im)
 
int getInclusiveMode ()
 

Protected Member Functions

void setp (double e, double px, double py, double pz)
 
void setp (const EvtVector4R &p4)
 
void setpart_num (EvtId particle_number)
 

Protected Attributes

bool _validP4
 

Detailed Description

Definition at line 42 of file EvtParticle.hh.

Constructor & Destructor Documentation

◆ EvtParticle()

EvtParticle::EvtParticle ( )

Default constructor.

Definition at line 60 of file EvtParticle.cc.

60 {
61 _ndaug=0;
62 _parent=0;
63 _channel=-10;
64 _t=0.0;
65 _genlifetime=1;
66 _first=1;
67 _isInit=false;
68 _validP4=false;
69 _isDecayed=false;
70 _decayProb=0;
71 // _mix=false;
72}

◆ ~EvtParticle()

EvtParticle::~EvtParticle ( )
virtual

Destructor.

Definition at line 56 of file EvtParticle.cc.

56 {
57 delete _decayProb;
58}

Member Function Documentation

◆ addDaug()

void EvtParticle::addDaug ( EvtParticle node)

Add another daughter to the particle

Definition at line 104 of file EvtParticle.cc.

104 {
105 node->_daug[node->_ndaug++]=this;
106 _ndaug=0;
107 _parent=node;
108}

Referenced by EvtKstarstargamma::decay(), EvtPHOTOS::doRadCorr(), and makeDaughters().

◆ compMassProb()

double EvtParticle::compMassProb ( )

Definition at line 504 of file EvtParticle.cc.

504 {
505
506 EvtParticle *p=this;
507 //report(INFO,"EvtGen") << "compMassProb " << endl;
508 //p->printTree();
509 double mass=p->mass();
510 double parMass=0.;
511 if ( p->getParent()) {
512 parMass=p->getParent()->mass();
513 }
514
515 int nDaug=p->getNDaug();
516 double *dMasses=0;
517
518 int i;
519 if ( nDaug>0 ) {
520 dMasses=new double[nDaug];
521 for (i=0; i<nDaug; i++) dMasses[i]=p->getDaug(i)->mass();
522 }
523
524 double temp=1.0;
525 temp=EvtPDL::getMassProb(p->getId(), mass, parMass, nDaug, dMasses);
526
527 //report(INFO,"EvtGen") << temp << " " << EvtPDL::name(p->getId()) << endl;
528 //If the particle already has a mass, we dont need to include
529 //it in the probability calculation
530 if ( (!p->getParent() || _validP4 ) && temp>0.0 ) temp=1.;
531
532 delete [] dMasses;
533 // if ( temp < 0.9999999 )
534 for (i=0; i<nDaug; i++) {
535 temp*=p->getDaug(i)->compMassProb();
536 }
537 return temp;
538}
static double getMassProb(EvtId i, double mass, double massPar, int nDaug, double *massDau)
Definition: EvtPDL.hh:48
EvtId getId() const
Definition: EvtParticle.cc:113
EvtParticle * getParent()
Definition: EvtParticle.cc:87
double compMassProb()
Definition: EvtParticle.cc:504
int getNDaug() const
Definition: EvtParticle.cc:125
EvtParticle * getDaug(int i)
Definition: EvtParticle.cc:85
double mass() const
Definition: EvtParticle.cc:127

Referenced by compMassProb(), and generateMassTree().

◆ decay()

void EvtParticle::decay ( )

Decay particle

Definition at line 404 of file EvtParticle.cc.

404 {
405
406 //P is particle to decay, typically 'this' but sometime
407 //modified by mixing
408 EvtParticle* p=this;
409 //Did it mix?
410 //if ( p->getMixed() ) {
411 //should take C(p) - this should only
412 //happen the first time we call decay for this
413 //particle
414 //p->takeCConj();
415 // p->setUnMixed();
416 //}
417 EvtSpinDensity myRho; //pingrg test code
418 EvtDecayBase *decayer;
419 decayer = EvtDecayTable::getDecayFunc(p);
420 // if ( decayer ) {
421 // report(INFO,"EvtGen") << "calling decay for " << EvtPDL::name(p->getId()) << " " << p->mass() << " " << p->getP4() << " " << p->getNDaug() << " " << p << endl;
422 // report(INFO,"EvtGen") << "NDaug= " << decayer->getNDaug() << endl;
423 // int ti;
424 // for ( ti=0; ti<decayer->getNDaug(); ti++)
425 // report(INFO,"EvtGen") << "Daug " << ti << " " << EvtPDL::name(decayer->getDaug(ti)) << endl;
426 // }
427 //if (p->_ndaug>0) {
428 // report(INFO,"EvtGen") <<"Is decaying particle with daughters!!!!!"<<endl;
429 // ::abort();
430 //return;
431 //call initdecay first - April 29,2002 - Lange
432 //}
433
434 //if there are already daughters, then this step is already done!
435 // figure out the masses
436 if ( _ndaug == 0 ) generateMassTree();
437 static EvtId BS0=EvtPDL::getId("B_s0");
438 static EvtId BSB=EvtPDL::getId("anti-B_s0");
439 static EvtId BD0=EvtPDL::getId("B0");
440 static EvtId BDB=EvtPDL::getId("anti-B0");
441 if ( _ndaug==1 && (getId()==BS0||getId()==BSB||getId()==BD0||getId()==BDB) )
442 {getDaug(0)->decay();
443 std::cout<<"if ( _ndaug==1 && (getId()==BS0||getId()==BSB||getId()==BD0||getId()==BDB) )"<<endl;
444 }
445
446 else{
447 //now we have accepted a set of masses - time
448 if ( decayer ) {
449 decayer->makeDecay(p);
450 //p->printTree(); //for debugging
451 }
452 else{
453 p->_rhoBackward.SetDiag(p->getSpinStates());
454
455 }
456 }
457 _isDecayed=true;
458 return;
459}
virtual void makeDecay(EvtParticle *p)=0
static EvtDecayBase * getDecayFunc(EvtParticle *)
Definition: EvtId.hh:27
static EvtId getId(const std::string &name)
Definition: EvtPDL.cc:287
void decay()
Definition: EvtParticle.cc:404
int getSpinStates() const
Definition: EvtParticle.cc:118
void generateMassTree()
Definition: EvtParticle.cc:461
void SetDiag(int n)

Referenced by decay(), EvtJetSet::decay(), EvtLunda::decay(), EvtLundCharm::decay(), EvtPythia::decay(), EvtTauola::decay(), EvtGen::generateDecay(), EvtDecayAmp::makeDecay(), EvtDecayIncoherent::makeDecay(), and EvtDecayProb::makeDecay().

◆ decayProb()

double * EvtParticle::decayProb ( )
inline

Definition at line 390 of file EvtParticle.hh.

390{return _decayProb;}

◆ deleteDaughters()

void EvtParticle::deleteDaughters ( bool  keepChannel = false)

Definition at line 540 of file EvtParticle.cc.

540 {
541 int i;
542
543 for(i=0;i<_ndaug;i++){
544 _daug[i]->deleteTree();
545 }
546
547 _ndaug=0;
548 //if ( keepChannel ) report(INFO,"EvtGen") << "keeping \n";
549 if ( !keepChannel) _channel=-10;
550 //_t=0.0;
551 //_genlifetime=1;
552 _first=1;
553 _isInit=false;
554 //report(INFO,"EvtGen") << "calling deletedaughters " << EvtPDL::name(this->getId()) <<endl;
555}
void deleteTree()
Definition: EvtParticle.cc:557

Referenced by EvtJetSet::decay(), EvtLunda::decay(), EvtLundCharm::decay(), EvtPythia::decay(), EvtTauola::decay(), EvtVSSBMixCPT::decay(), deleteTree(), EvtConExc::gamHXSection(), and initializePhaseSpace().

◆ deleteTree()

void EvtParticle::deleteTree ( )

Delete a decay chain

Definition at line 557 of file EvtParticle.cc.

557 {
558
559 this->deleteDaughters();
560
561 delete this;
562
563}
void deleteDaughters(bool keepChannel=false)
Definition: EvtParticle.cc:540

Referenced by EvtbTosllAmp::CalcMaxProb(), EvtSemiLeptonicAmp::CalcMaxProb(), deleteDaughters(), EvtGen::generateDecay(), EvtGen::generateEvent(), EvtPsi3Sdecay::PHSPDecay(), and EvtConExc::~EvtConExc().

◆ dumpTree()

void EvtParticle::dumpTree ( ) const

Definition at line 978 of file EvtParticle.cc.

978 { //pingrg, Mar. 25,2008
979
980 report(INFO,"EvtGen") << "This is the current decay chain to dump"<<endl;
981 report(INFO,"") << "This top particle is "<<
982 EvtPDL::name(_id).c_str()<<endl;
983
984 this->dumpTreeRec(0,0);
985 report(INFO,"EvtGen") << "End of decay chain."<<endl;
986
987}
ostream & report(Severity severity, const char *facility)
Definition: EvtReport.cc:36
@ INFO
Definition: EvtReport.hh:52
static std::string name(EvtId i)
Definition: EvtPDL.hh:64
void dumpTreeRec(int level, int dj) const
Definition: EvtParticle.cc:948

◆ dumpTreeRec()

void EvtParticle::dumpTreeRec ( int  level,
int  dj 
) const

Definition at line 948 of file EvtParticle.cc.

948 { //pingrg, Mar. 25,2008
949
950 int newlevel,i;
951 newlevel = level +1;
952
953
954 if (_ndaug!=0) {
955
956 int Ids= EvtPDL::getStdHep(_id);
957 std::string c1,cid;
958 if(Ids>0) c1="p";
959 if(Ids<0) {
960 c1="m";Ids=-1*Ids;
961 }
962
963 cid=c1+IntToStr(Ids);
964
965 report(INFO,"") <<newlevel<<" "<< cid<<" "<<dj;
966 report(INFO,"") << " = ";
967
968 int Nchannel=this->getChannel()+1;
969 report(INFO,"") <<Nchannel<<endl;
970
971 for(i=0;i<_ndaug;i++){
972 _daug[i]->dumpTreeRec(newlevel,i);
973 }
974 }
975}
std::string IntToStr(int a)
static int getStdHep(EvtId id)
Definition: EvtPDL.hh:56
int getChannel() const
Definition: EvtParticle.cc:123

Referenced by dumpTree(), and dumpTreeRec().

◆ eps()

EvtVector4C EvtParticle::eps ( int  i) const
virtual

Returns polarization vector in the particles own restframe.

Reimplemented in EvtVectorParticle.

Definition at line 576 of file EvtParticle.cc.

576 {
577 EvtVector4C temp;
579 report(ERROR,"EvtGen") << "and you have asked for the:"<<i
580 <<"th polarization vector."
581 <<" I.e. you thought it was a"
582 <<" vector particle!" << endl;
583 ::abort();
584 return temp;
585}
@ ERROR
Definition: EvtReport.hh:49
void printParticle() const

Referenced by EvtEDM::decay(), EvtJpipi::decay(), EvtOmegaDalitz::decay(), EvtSVPCP::decay(), EvtSVPHelAmp::decay(), EvtVll::decay(), EvtVPHOtoVISR::decay(), EvtVPHOtoVISRHi::decay(), EvtVSPPwave::decay(), EvtVSS::decay(), EvtVSSBMixCPT::decay(), EvtVSSMix::decay(), EvtVVP::decay(), EvtVVpipi::decay(), EvtVVPIPI_WEIGHTED::decay(), EvtVVSPwave::decay(), and EvtSVVHelAmp::SVVHel().

◆ epsParent()

EvtVector4C EvtParticle::epsParent ( int  i) const
virtual

Returns polarization vector in the parents restframe.

Reimplemented in EvtVectorParticle.

Definition at line 565 of file EvtParticle.cc.

565 {
566 EvtVector4C temp;
568 report(ERROR,"EvtGen") << "and you have asked for the:"<<i
569 <<"th polarization vector."
570 <<" I.e. you thought it was a"
571 <<" vector particle!" << endl;
572 ::abort();
573 return temp;
574}

Referenced by EvtbTosllVectorAmp::CalcAmp(), EvtSemiLeptonicVectorAmp::CalcAmp(), EvtBHadronic::decay(), EvtKstarnunu::decay(), EvtSSDCP::decay(), EvtVectorIsr::decay(), EvtVPHOtoVISR::decay(), EvtVPHOtoVISRHi::decay(), and EvtVVPIPI_WEIGHTED::decay().

◆ epsParentPhoton()

EvtVector4C EvtParticle::epsParentPhoton ( int  i)
virtual

Returns polarization vector in the parents restframe for a photon.

Reimplemented in EvtPhotonParticle.

Definition at line 587 of file EvtParticle.cc.

587 {
588 EvtVector4C temp;
590 report(ERROR,"EvtGen") << "and you have asked for the:"<<i
591 <<"th polarization vector of photon."
592 <<" I.e. you thought it was a"
593 <<" photon particle!" << endl;
594 ::abort();
595 return temp;
596}

Referenced by EvtLNuGamma::decay(), EvtSVPCP::decay(), EvtSVPHelAmp::decay(), EvtVectorIsr::decay(), EvtVSPPwave::decay(), and EvtVVP::decay().

◆ epsPhoton()

EvtVector4C EvtParticle::epsPhoton ( int  i)
virtual

Returns polarization vector in the particles own restframe for a photon.

Reimplemented in EvtPhotonParticle.

Definition at line 598 of file EvtParticle.cc.

598 {
599 EvtVector4C temp;
601 report(ERROR,"EvtGen") << "and you have asked for the:"<<i
602 <<"th polarization vector of a photon."
603 <<" I.e. you thought it was a"
604 <<" photon particle!" << endl;
605 ::abort();
606 return temp;
607}

◆ epsTensor()

EvtTensor4C EvtParticle::epsTensor ( int  i) const
virtual

Returns tensor in the particles own restframe for a spin 2 particle.

Reimplemented in EvtTensorParticle.

Definition at line 670 of file EvtParticle.cc.

670 {
671 int temp;
672 temp = i;
673 EvtTensor4C tempC;
675 report(ERROR,"EvtGen") << "and you have asked for the:"<<i
676 <<"th tensor."
677 <<" I.e. you thought it was a"
678 <<" Tensor particle!" << endl;
679 ::abort();
680 return tempC;
681}

Referenced by EvtTSS::decay(), and EvtTVSPwave::decay().

◆ epsTensorParent()

EvtTensor4C EvtParticle::epsTensorParent ( int  i) const
virtual

Returns tensor in the parents restframe for a spin 2 particle.

Reimplemented in EvtTensorParticle.

Definition at line 657 of file EvtParticle.cc.

657 {
658 int temp;
659 temp = i;
660 EvtTensor4C tempC;
662 report(ERROR,"EvtGen") << "and you have asked for the:"<<i
663 <<"th tensor."
664 <<" I.e. you thought it was a"
665 <<" Tensor particle!" << endl;
666 ::abort();
667 return tempC;
668}

Referenced by EvtSemiLeptonicTensorAmp::CalcAmp(), EvtBHadronic::decay(), EvtSSDCP::decay(), EvtSTS::decay(), and EvtSTSCP::decay().

◆ firstornot()

int EvtParticle::firstornot ( ) const

Used internally to decide if first time particle is decayed.

Definition at line 111 of file EvtParticle.cc.

111{ return _first;}

Referenced by EvtDecayBase::findMasses().

◆ generateMassTree()

void EvtParticle::generateMassTree ( )

Definition at line 461 of file EvtParticle.cc.

461 {
462 double massProb=1.;
463 double ranNum=2.;
464 int counter=0;
465 EvtParticle *p=this;
466 while (massProb<ranNum) {
467 //check it out the first time.
468 p->initDecay();
469 //report(INFO,"EvtGen") << "calling massProb \n";
470 massProb=p->compMassProb();
471 ranNum=EvtRandom::Flat();
472 // report(INFO,"EvtGen") << "end of iter " << massProb <<endl;
473 counter++;
474
475 if ( counter > 10000 ) {
476 if ( counter == 10001 ) {
477 report(INFO,"EvtGen") << "Too many iterations to determine the mass tree. Parent mass= "<< p->mass() << _p<<" " << massProb <<endl;
478 p->printTree();
479 report(INFO,"EvtGen") << "will take next combo with non-zero likelihood\n";
480 }
481 if ( massProb>0. ) massProb=2.0;
482 if ( counter > 20000 ) {
483 // one last try - take the minimum masses
484 p->initDecay(true);
485 p->printTree();
486 massProb=p->compMassProb();
487 if ( massProb>0. ) {
488 massProb=2.0;
489 report(INFO,"EvtGen") << "Taking the minimum mass of all particles in the chain\n";
490 }
491 else {
492 report(INFO,"EvtGen") << "Sorry, no luck finding a valid set of masses. This may be a pathological combo\n";
493 std::cout<<EvtPDL::name(p->getId())<<": Parent mass "<<p->getP4().mass()<<" with momentum "<<p->getP4()<<std::endl;
494 if(EvtGlobalSet::ConExcPythia){EvtGlobalSet::ConExcPythia=0;return;}else{abort();}
495 //assert(0);
496 }
497 }
498 }
499 }
500 // report(INFO,"EvtGen") << counter << endl;
501 // p->printTree();
502}
static bool ConExcPythia
Definition: EvtGlobalSet.hh:20
void initDecay(bool useMinMass=false)
Definition: EvtParticle.cc:237
const EvtVector4R & getP4() const
Definition: EvtParticle.cc:121
void printTree() const
Definition: EvtParticle.cc:897
static double Flat()
Definition: EvtRandom.cc:74
double mass() const
Definition: EvtVector4R.cc:39

Referenced by decay(), EvtBtoKD3P::decay(), and initializePhaseSpace().

◆ get4Pos()

EvtVector4R EvtParticle::get4Pos ( )

Returns the 4position of the particle in the lab frame.

Definition at line 706 of file EvtParticle.cc.

706 {
707
708 EvtVector4R temp,mom;
709 EvtParticle *ptemp;
710
711 temp.set(0.0,0.0,0.0,0.0);
712 ptemp=getParent();
713
714 if (ptemp==0) return temp;
715
716 temp=(ptemp->_t/ptemp->mass())*(ptemp->getP4());
717
718 while (ptemp->getParent()!=0) {
719 ptemp=ptemp->getParent();
720 mom=ptemp->getP4();
721 temp=boostTo(temp,mom);
722 temp=temp+(ptemp->_t/ptemp->mass())*(ptemp->getP4());
723 }
724
725 return temp;
726}
EvtDiracSpinor boostTo(const EvtDiracSpinor &sp, const EvtVector4R p4)
void set(int i, double d)
Definition: EvtVector4R.hh:183

Referenced by makeStdHep().

◆ getChannel()

int EvtParticle::getChannel ( ) const

◆ getDaug()

EvtParticle * EvtParticle::getDaug ( int  i)

Get pointer the the i:th daugther.

Definition at line 85 of file EvtParticle.cc.

85{ return _daug[i]; }

Referenced by EvtRexc::angularSampling(), EvtConExc::angularSampling(), EvtHypNonLepton::calcAmp(), EvtbTosllScalarAmp::CalcAmp(), EvtbTosllVectorAmp::CalcAmp(), EvtSemiLeptonicBaryonAmp::CalcAmp(), EvtSemiLeptonicTensorAmp::CalcAmp(), EvtSemiLeptonicVectorAmp::CalcAmp(), EvtbTosllAmp::CalcMaxProb(), EvtSemiLeptonicAmp::CalcMaxProb(), EvtConExc::checkdecay(), EvtPsi3Sdecay::choseDecay(), compMassProb(), EvtSecondary::createSecondary(), decay(), EvtAngH2::decay(), EvtAngSam::decay(), EvtAngSam3::decay(), EvtAngSamLab::decay(), EvtAngSamX::decay(), EvtAV2GV::decay(), EvtBHadronic::decay(), EvtBody3::decay(), EvtBsquark::decay(), EvtBTo3piCP::decay(), EvtBTo4piCP::decay(), EvtBtoKD3P::decay(), EvtBToKpipiCP::decay(), EvtBtoXsEtap::decay(), EvtBtoXsgamma::decay(), EvtBtoXsll::decay(), EvtCalHelAmp::decay(), EvtCBTo3piMPP::decay(), EvtCBTo3piP00::decay(), EvtChi0BB1::decay(), EvtChi0BB2::decay(), EvtChi1BB1::decay(), EvtChi1BB2::decay(), EvtChi2BB1::decay(), EvtChi2BB2::decay(), EvtConExc::decay(), EvtD0mixDalitz::decay(), EvtD0ToKpipi0pi0::decay(), EvtD0ToKpipipi::decay(), EvtDDalitz::decay(), EvtDeBD::decay(), EvtDIY::decay(), EvtDToKpienu::decay(), EvtDToKSpipipi::decay(), EvtDTopipienu::decay(), EvtEDM::decay(), EvtEtaDalitz::decay(), EvtEtap2gpipi::decay(), EvtEtap2gpipiB::decay(), EvtFlatQ2::decay(), EvtHAngSam3::decay(), EvtHelPPJ::decay(), EvtHypWK::decay(), EvtIntervalDecayAmp< T >::decay(), EvtJ2BB1::decay(), EvtJ2BB2::decay(), EvtJ2BB3::decay(), EvtJetSet::decay(), EvtJPE::decay(), EvtJpipi::decay(), EvtJscont::decay(), EvtJTO3P::decay(), EvtKstarnunu::decay(), EvtKstarstargamma::decay(), EvtLambdac2pKpi::decay(), EvtLambdaP_BarGamma::decay(), EvtLNuGamma::decay(), EvtLunda::decay(), EvtLundCharm::decay(), EvtMassH1::decay(), EvtMassH2::decay(), EvtMBody3::decay(), EvtmH2::decay(), EvtmPhsp::decay(), EvtMultibody::decay(), EvtOmegaDalitz::decay(), EvtOpenCharm::decay(), EvtP2GC0::decay(), EvtP2GC1::decay(), EvtP2GC2::decay(), EvtPhiDalitz::decay(), EvtPhokhara::decay(), EvtPhokhara_4pi::decay(), EvtPhokhara_K0K0::decay(), EvtPhokhara_KK::decay(), EvtPhokhara_LLB::decay(), EvtPhokhara_nnbar::decay(), EvtPhokhara_pi0pi0pipi::decay(), EvtPhokhara_pipi::decay(), EvtPhokhara_pipieta::decay(), EvtPhokhara_pipipi0::decay(), EvtPhokhara_ppbar::decay(), EvtPhsp::decay(), EvtPi0Dalitz::decay(), EvtPycont::decay(), EvtPyGaGa::decay(), EvtPythia::decay(), EvtRexc::decay(), EvtRhoPi::decay(), EvtSingleParticle::decay(), EvtSingleParticle2::decay(), EvtSinglePoint::decay(), EvtSll::decay(), EvtSLN::decay(), EvtSPL::decay(), EvtSSDCP::decay(), EvtSTS::decay(), EvtSTSCP::decay(), EvtSVPCP::decay(), EvtSVPHelAmp::decay(), EvtSVS::decay(), EvtSVSCP::decay(), EvtSVSCPiso::decay(), EvtSVSCPLH::decay(), EvtSVSNONCPEIGEN::decay(), EvtT2GV::decay(), EvtTauGamMu::decay(), EvtTauHadnu::decay(), EvtTaulnunu::decay(), EvtTauola::decay(), EvtTauScalarnu::decay(), EvtTauVectornu::decay(), EvtTrackGen::decay(), EvtTSS::decay(), EvtTVSPwave::decay(), EvtVectorIsr::decay(), EvtVll::decay(), EvtVPHOtoVISR::decay(), EvtVPHOtoVISRHi::decay(), EvtVSPPwave::decay(), EvtVSS::decay(), EvtVSSBMixCPT::decay(), EvtVSSMix::decay(), EvtVub::decay(), EvtVubHybrid::decay(), EvtVubNLO::decay(), EvtVVP::decay(), EvtVVpipi::decay(), EvtVVPIPI_WEIGHTED::decay(), EvtVVSPwave::decay(), EvtConExc::difgamXs(), EvtPHOTOS::doRadCorr(), EvtEvalHelAmp::evalAmp(), EvtLunda::ExclusiveDecay(), EvtDecayBase::findMass(), EvtDecayBase::findMasses(), EvtDecayBase::findMaxMass(), EvtConExc::findMaxXS(), EvtPsi3Sdecay::findMode(), EvtFSPick::FSPick(), EvtConExc::gamHXSection(), EvtGen::generateDecay(), initDecay(), initializePhaseSpace(), EvtDecayAmp::makeDecay(), EvtDecayIncoherent::makeDecay(), EvtDecayProb::makeDecay(), EvtDecayTag::makeTag(), EvtCPUtil::OtherB(), EvtConExc::photonSampling(), EvtPsi3Sdecay::PHSPDecay(), EvtConExc::Rad1difXs(), EvtConExc::Rad2difXs(), EvtConExc::SetP4(), EvtConExc::SetP4Rvalue(), and EvtSVVHelAmp::SVVHel().

◆ getGeneratorFlag()

int EvtParticle::getGeneratorFlag ( )
inline

get generator information; pingrg-2011-1-6

Definition at line 146 of file EvtParticle.hh.

146{ return _generatorFlag ;}

◆ getId()

EvtId EvtParticle::getId ( ) const

Returns Id of particle.

Definition at line 113 of file EvtParticle.cc.

113{ return _id;}

Referenced by EvtRexc::angularSampling(), EvtConExc::angularSampling(), EvtbTosllScalarAmp::CalcAmp(), EvtbTosllVectorAmp::CalcAmp(), EvtSemiLeptonicBaryonAmp::CalcAmp(), EvtSemiLeptonicTensorAmp::CalcAmp(), EvtSemiLeptonicVectorAmp::CalcAmp(), EvtConExc::checkdecay(), EvtPsi3Sdecay::choseDecay(), compMassProb(), EvtSecondary::createSecondary(), decay(), EvtAngH2::decay(), EvtBsquark::decay(), EvtBTo3piCP::decay(), EvtBtoXsll::decay(), EvtCalHelAmp::decay(), EvtCBTo3piMPP::decay(), EvtCBTo3piP00::decay(), EvtConExc::decay(), EvtDDalitz::decay(), EvtDeBD::decay(), EvtDMix::decay(), EvtDToKpienu::decay(), EvtDTopipienu::decay(), EvtEtap2gpipi::decay(), EvtEtap2gpipiB::decay(), EvtHypWK::decay(), EvtJetSet::decay(), EvtJscont::decay(), EvtKstarstargamma::decay(), EvtLambdac2pKpi::decay(), EvtLNuGamma::decay(), EvtLunda::decay(), EvtLundCharm::decay(), EvtOpenCharm::decay(), EvtPhokhara::decay(), EvtPhokhara_4pi::decay(), EvtPhokhara_K0K0::decay(), EvtPhokhara_KK::decay(), EvtPhokhara_LLB::decay(), EvtPhokhara_nnbar::decay(), EvtPhokhara_pi0pi0pipi::decay(), EvtPhokhara_pipi::decay(), EvtPhokhara_pipieta::decay(), EvtPhokhara_pipipi0::decay(), EvtPhokhara_ppbar::decay(), EvtPycont::decay(), EvtPyGaGa::decay(), EvtPythia::decay(), EvtRexc::decay(), EvtSVSCPiso::decay(), EvtSVSNONCPEIGEN::decay(), EvtSVVNONCPEIGEN::decay(), EvtTauHadnu::decay(), EvtTaulnunu::decay(), EvtTauola::decay(), EvtTauScalarnu::decay(), EvtTauVectornu::decay(), EvtTrackGen::decay(), EvtVPHOtoVISR::decay(), EvtVPHOtoVISRHi::decay(), EvtVSSBMixCPT::decay(), EvtVSSMix::decay(), EvtPHOTOS::doRadCorr(), EvtLunda::ExclusiveDecay(), EvtDecayBase::findMass(), EvtDecayBase::findMasses(), EvtDecayBase::findMaxMass(), EvtConExc::findMaxXS(), EvtPsi3Sdecay::findMode(), EvtFSPick::FSPick(), EvtConExc::gamHXSection(), generateMassTree(), EvtDecayTable::getDecayFunc(), EvtParticleDecayList::getDecayModel(), initDecay(), initializePhaseSpace(), EvtDecayAmp::makeDecay(), EvtDecayProb::makeDecay(), makeStdHep(), EvtDecayTag::makeTag(), EvtCPUtil::OtherB(), EvtConExc::photonSampling(), EvtPsi3Sdecay::PHSPDecay(), printTreeRec(), EvtDiracParticle::rotateToHelicityBasis(), EvtHighSpinParticle::rotateToHelicityBasis(), EvtRaritaSchwingerParticle::rotateToHelicityBasis(), EvtConExc::selectMode(), setLifetime(), EvtConExc::SetP4(), EvtConExc::SetP4Rvalue(), treeStrRec(), and writeTreeRec().

◆ getInclusiveMode()

int EvtParticle::getInclusiveMode ( )
inline

Definition at line 394 of file EvtParticle.hh.

394{return _inclusiveMode;}

◆ getIntFlag()

std::vector< int > EvtParticle::getIntFlag ( )
inline

get int flag for ConExc: pingrg-2015-2-7

Definition at line 158 of file EvtParticle.hh.

158 {
159 return _intFlag;
160 }

Referenced by EvtDecayTag::getModeTag().

◆ getLifetime()

double EvtParticle::getLifetime ( )

Returns the lifetime.

Definition at line 99 of file EvtParticle.cc.

99 {
100
101 return _t;
102}

Referenced by EvtVSSBMixCPT::decay(), EvtVSSMix::decay(), and EvtCPUtil::OtherB().

◆ getNDaug()

◆ getP4()

const EvtVector4R & EvtParticle::getP4 ( ) const

Returns 4momentum in parents restframe.

Definition at line 121 of file EvtParticle.cc.

121{ return _p;}

Referenced by EvtbTosllScalarAmp::CalcAmp(), EvtbTosllVectorAmp::CalcAmp(), EvtSemiLeptonicBaryonAmp::CalcAmp(), EvtSemiLeptonicTensorAmp::CalcAmp(), EvtSemiLeptonicVectorAmp::CalcAmp(), EvtPsi3Sdecay::choseDecay(), EvtAngSam::decay(), EvtAngSam3::decay(), EvtAngSamX::decay(), EvtAV2GV::decay(), EvtBHadronic::decay(), EvtBsquark::decay(), EvtBTo4piCP::decay(), EvtBtoKD3P::decay(), EvtCalHelAmp::decay(), EvtChi0BB1::decay(), EvtChi0BB2::decay(), EvtChi1BB1::decay(), EvtChi1BB2::decay(), EvtChi2BB1::decay(), EvtChi2BB2::decay(), EvtConExc::decay(), EvtD0mixDalitz::decay(), EvtD0ToKpipi0pi0::decay(), EvtD0ToKpipipi::decay(), EvtDDalitz::decay(), EvtDeBD::decay(), EvtDIY::decay(), EvtDToKpienu::decay(), EvtDToKSpipipi::decay(), EvtDTopipienu::decay(), EvtEDM::decay(), EvtEtaDalitz::decay(), EvtEtap2gpipi::decay(), EvtEtap2gpipiB::decay(), EvtFlatQ2::decay(), EvtHAngSam3::decay(), EvtHelPPJ::decay(), EvtHypWK::decay(), EvtJ2BB1::decay(), EvtJ2BB2::decay(), EvtJ2BB3::decay(), EvtJPE::decay(), EvtJpipi::decay(), EvtJTO3P::decay(), EvtKstarnunu::decay(), EvtKstarstargamma::decay(), EvtLNuGamma::decay(), EvtMBody3::decay(), EvtOmegaDalitz::decay(), EvtOpenCharm::decay(), EvtP2GC0::decay(), EvtP2GC1::decay(), EvtP2GC2::decay(), EvtPhiDalitz::decay(), EvtPhsp::decay(), EvtPi0Dalitz::decay(), EvtRexc::decay(), EvtRhoPi::decay(), EvtSPL::decay(), EvtSSDCP::decay(), EvtSTS::decay(), EvtSTSCP::decay(), EvtSVPCP::decay(), EvtSVPHelAmp::decay(), EvtSVS::decay(), EvtSVSCP::decay(), EvtSVSCPLH::decay(), EvtSVSNONCPEIGEN::decay(), EvtT2GV::decay(), EvtTauGamMu::decay(), EvtTauHadnu::decay(), EvtTauScalarnu::decay(), EvtTSS::decay(), EvtVSPPwave::decay(), EvtVSS::decay(), EvtVSSBMixCPT::decay(), EvtVSSMix::decay(), EvtVVpipi::decay(), EvtVVPIPI_WEIGHTED::decay(), EvtVVSPwave::decay(), EvtConExc::difgamXs(), EvtPHOTOS::doRadCorr(), EvtVectorParticle::epsParent(), EvtPhotonParticle::epsParentPhoton(), EvtTensorParticle::epsTensorParent(), EvtEvalHelAmp::evalAmp(), EvtConExc::findMaxXS(), EvtConExc::gamHXSection(), generateMassTree(), get4Pos(), getP4Lab(), EvtDecayAmp::makeDecay(), EvtPsi3Sdecay::PHSPDecay(), EvtConExc::Rad1difXs(), EvtConExc::Rad2difXs(), EvtDiracParticle::rotateToHelicityBasis(), EvtRaritaSchwingerParticle::rotateToHelicityBasis(), and EvtSVVHelAmp::SVVHel().

◆ getP4Lab()

EvtVector4R EvtParticle::getP4Lab ( )

Gets 4vector in the labframe, i.e., the frame in which the root particles momentum is measured.

Definition at line 685 of file EvtParticle.cc.

685 {
686 EvtVector4R temp,mom;
687 EvtParticle *ptemp;
688
689 temp=this->getP4();
690 ptemp=this;
691
692 while (ptemp->getParent()!=0) {
693 ptemp=ptemp->getParent();
694 mom=ptemp->getP4();
695 temp=boostTo(temp,mom);
696 }
697 return temp;
698}

Referenced by EvtRexc::angularSampling(), EvtConExc::angularSampling(), EvtAngH2::decay(), EvtBody3::decay(), EvtDeBD::decay(), EvtDIY::decay(), EvtEtap2gpipi::decay(), EvtEtap2gpipiB::decay(), EvtLambdaP_BarGamma::decay(), EvtMassH1::decay(), EvtMassH2::decay(), EvtMBody3::decay(), EvtmH2::decay(), EvtmPhsp::decay(), EvtMultibody::decay(), EvtRhoPi::decay(), EvtSingleParticle::decay(), EvtSingleParticle2::decay(), EvtSPL::decay(), makeStdHep(), and EvtConExc::photonSampling().

◆ getP4Restframe()

EvtVector4R EvtParticle::getP4Restframe ( )

Gets 4vector in the particles restframe, i.e. this functiont will return (m,0,0,0)

Definition at line 700 of file EvtParticle.cc.

700 {
701
702 return EvtVector4R(mass(),0.0,0.0,0.0);
703
704}

Referenced by EvtKstarstargamma::decay(), EvtSSDCP::decay(), EvtSVSCPLH::decay(), EvtVPHOtoVISR::decay(), and EvtVPHOtoVISRHi::decay().

◆ getParent()

EvtParticle * EvtParticle::getParent ( )

Returns pointer to parent particle.

Definition at line 87 of file EvtParticle.cc.

87{ return _parent;}

Referenced by compMassProb(), EvtDDalitz::decay(), EvtDecayBase::findMass(), EvtDecayBase::findMaxMass(), get4Pos(), getP4Lab(), initDecay(), EvtDecayAmp::makeDecay(), and EvtCPUtil::OtherB().

◆ getSpinDensityBackward()

EvtSpinDensity EvtParticle::getSpinDensityBackward ( )
inline

Get backward spin density matrix.

Definition at line 357 of file EvtParticle.hh.

357{return _rhoBackward;}

Referenced by EvtDecayAmp::makeDecay().

◆ getSpinDensityForward()

EvtSpinDensity EvtParticle::getSpinDensityForward ( )
inline

Get forward spin density matrix.

Definition at line 347 of file EvtParticle.hh.

347{return _rhoForward;}

Referenced by EvtDecayAmp::makeDecay(), and EvtDecayIncoherent::makeDecay().

◆ getSpinStates()

int EvtParticle::getSpinStates ( ) const

Returns number of spin states of the particle.

Definition at line 118 of file EvtParticle.cc.

static EvtSpinType::spintype getSpinType(EvtId i)
Definition: EvtPDL.hh:61
static int getSpinStates(spintype stype)
Definition: EvtSpinType.hh:64

Referenced by EvtbTosllAmp::CalcMaxProb(), EvtSemiLeptonicAmp::CalcMaxProb(), decay(), EvtDecayIncoherent::makeDecay(), EvtDecayProb::makeDecay(), setDiagonalSpinDensity(), setPolarizedSpinDensity(), and setVectorSpinDensity().

◆ getSpinType()

EvtSpinType::spintype EvtParticle::getSpinType ( ) const

Returns particle type.

Definition at line 115 of file EvtParticle.cc.

116 { return EvtPDL::getSpinType(_id);}

◆ hasValidP4()

bool EvtParticle::hasValidP4 ( )
inline

Definition at line 383 of file EvtParticle.hh.

383{return _validP4;}

Referenced by EvtDecayBase::findMaxMass(), EvtParticleDecayList::getDecayModel(), and initDecay().

◆ init()

◆ initDecay()

void EvtParticle::initDecay ( bool  useMinMass = false)

Definition at line 237 of file EvtParticle.cc.

237 {
238
239 EvtParticle* p=this;
240 // carefull - the parent mass might be fixed in stone..
241 EvtParticle* par=p->getParent();
242 double parMass=-1.;
243 if ( par != 0 ) {
244 if ( par->hasValidP4() ) parMass=par->mass();
245 int i;
246 for ( i=0;i<par->getNDaug();i++) {
247 EvtParticle *tDaug=par->getDaug(i);
248 if ( p != tDaug )
249 parMass-=EvtPDL::getMinMass(tDaug->getId());
250 }
251 }
252
253 if ( _isInit ) {
254 //we have already been here - just reroll the masses!
255 if ( _ndaug>0) {
256 int ii;
257 for(ii=0;ii<_ndaug;ii++){
258 if ( _ndaug==1 || EvtPDL::getWidth(p->getDaug(ii)->getId()) > 0.0000001)
259 p->getDaug(ii)->initDecay(useMinMass);
260 else p->getDaug(ii)->setMass(EvtPDL::getMeanMass(p->getDaug(ii)->getId()));
261 }
262 }
263
264 int j;
265 EvtId *dauId=0;
266 double *dauMasses=0;
267 if ( _ndaug > 0) {
268 dauId=new EvtId[_ndaug];
269 dauMasses=new double[_ndaug];
270 for (j=0;j<_ndaug;j++) {
271 dauId[j]=p->getDaug(j)->getId();
272 dauMasses[j]=p->getDaug(j)->mass();
273 }
274 }
275 EvtId *parId=0;
276 EvtId *othDauId=0;
277 EvtParticle *tempPar=p->getParent();
278 if (tempPar) {
279 parId=new EvtId(tempPar->getId());
280 if ( tempPar->getNDaug()==2 ) {
281 if ( tempPar->getDaug(0) == this ) othDauId=new EvtId(tempPar->getDaug(1)->getId());
282 else othDauId=new EvtId(tempPar->getDaug(0)->getId());
283 }
284 }
285 if ( p->getParent() && _validP4==false ) {
286 if ( !useMinMass ) p->setMass(EvtPDL::getRandMass(p->getId(),parId,_ndaug,dauId,othDauId,parMass,dauMasses));
287 else p->setMass(EvtPDL::getMinMass(p->getId()));
288 }
289 if ( parId) delete parId;
290 if ( othDauId) delete othDauId;
291 if ( dauId) delete [] dauId;
292 if ( dauMasses) delete [] dauMasses;
293 return;
294 }
295
296
297 //Will include effects of mixing here
298 //added by Lange Jan4,2000
299 static EvtId BS0=EvtPDL::getId("B_s0");
300 static EvtId BSB=EvtPDL::getId("anti-B_s0");
301 static EvtId BD0=EvtPDL::getId("B0");
302 static EvtId BDB=EvtPDL::getId("anti-B0");
303
304 //only makes sense if there is no parent particle
305 if ( (getNDaug()==0 && getParent()==0) && (getId()==BS0||getId()==BSB||getId()==BD0||getId()==BDB)){
306 double t;
307 int mix;
309 setLifetime(t);
310
311 if (mix) {
312
313 EvtScalarParticle* scalar_part;
314
315 scalar_part=new EvtScalarParticle;
316 if (getId()==BS0) {
317 EvtVector4R p_init(EvtPDL::getMass(BSB),0.0,0.0,0.0);
318 scalar_part->init(BSB,p_init);
319 }
320 else if (getId()==BSB) {
321 EvtVector4R p_init(EvtPDL::getMass(BS0),0.0,0.0,0.0);
322 scalar_part->init(BS0,p_init);
323 }
324 else if (getId()==BD0) {
325 EvtVector4R p_init(EvtPDL::getMass(BDB),0.0,0.0,0.0);
326 scalar_part->init(BDB,p_init);
327 }
328 else if (getId()==BDB) {
329 EvtVector4R p_init(EvtPDL::getMass(BD0),0.0,0.0,0.0);
330 scalar_part->init(BD0,p_init);
331 }
332
333 scalar_part->setLifetime(0);
334
335 scalar_part->setDiagonalSpinDensity();
336
337 insertDaugPtr(0,scalar_part);
338
339 _ndaug=1;
340 _isInit=true;
341 p=scalar_part;
342 p->initDecay(useMinMass);
343 return;
344
345
346 }
347 }
348 if ( _ndaug==1 ) std::cout << "hi " << EvtPDL::name(this->getId()) << std::endl;
349
350 EvtDecayBase *decayer;
351 decayer = EvtDecayTable::getDecayFunc(p);
352
353 if ( decayer ) {
354 p->makeDaughters(decayer->nRealDaughters(),decayer->getDaugs());
355 // report(INFO,"EvtGen") << "has found decay " << decayer->nRealDaughters() << endl;
356 //then loop over the daughters and init their decay
357 int i;
358 for(i=0;i<p->getNDaug();i++){
359 if ( EvtPDL::getWidth(p->getDaug(i)->getId()) > 0.0000001)
360 p->getDaug(i)->initDecay(useMinMass);
361 else p->getDaug(i)->setMass(EvtPDL::getMeanMass(p->getDaug(i)->getId()));
362 // report(INFO,"EvtGen") << "has inited " << EvtPDL::name(p->getDaug(i)->getId()) <<
363 // " " << EvtPDL::name(p->getId()) << endl;
364 }
365 }
366 //figure masses in separate step...
367 // if ( p->getParent() && _validP4==false ) EvtDecayBase::findMass(p);
368
369 int j;
370 EvtId *dauId=0;
371 double *dauMasses=0;
372 int nDaugT=p->getNDaug();
373 if ( nDaugT > 0) {
374 dauId=new EvtId[nDaugT];
375 dauMasses=new double[nDaugT];
376 for (j=0;j<nDaugT;j++) {
377 dauId[j]=p->getDaug(j)->getId();
378 dauMasses[j]=p->getDaug(j)->mass();
379 }
380 }
381
382 EvtId *parId=0;
383 EvtId *othDauId=0;
384 EvtParticle *tempPar=p->getParent();
385 if (tempPar) {
386 parId=new EvtId(tempPar->getId());
387 if ( tempPar->getNDaug()==2 ) {
388 if ( tempPar->getDaug(0) == this ) othDauId=new EvtId(tempPar->getDaug(1)->getId());
389 else othDauId=new EvtId(tempPar->getDaug(0)->getId());
390 }
391 }
392 if ( p->getParent() && p->hasValidP4()==false ) {
393 if ( !useMinMass ) p->setMass(EvtPDL::getRandMass(p->getId(),parId,p->getNDaug(),dauId,othDauId,parMass,dauMasses));
394 else p->setMass(EvtPDL::getMinMass(p->getId()));
395 }
396 if ( parId) delete parId;
397 if ( othDauId) delete othDauId;
398 if ( dauId) delete [] dauId;
399 if ( dauMasses) delete [] dauMasses;
400 _isInit=true;
401}
TTree * t
Definition: binning.cxx:23
static void incoherentMix(const EvtId id, double &t, int &mix)
Definition: EvtCPUtil.cc:294
virtual int nRealDaughters()
EvtId * getDaugs()
Definition: EvtDecayBase.hh:65
static double getWidth(EvtId i)
Definition: EvtPDL.hh:54
static double getRandMass(EvtId i, EvtId *parId, int nDaug, EvtId *dauId, EvtId *othDaugId, double maxMass, double *dauMasses)
Definition: EvtPDL.hh:47
static double getMeanMass(EvtId i)
Definition: EvtPDL.hh:45
static double getMinMass(EvtId i)
Definition: EvtPDL.hh:51
static double getMass(EvtId i)
Definition: EvtPDL.hh:46
void setMass(double m)
Definition: EvtParticle.hh:372
void makeDaughters(int ndaug, EvtId *id)
void insertDaugPtr(int idaug, EvtParticle *partptr)
Definition: EvtParticle.hh:220
bool hasValidP4()
Definition: EvtParticle.hh:383
void setLifetime()
Definition: EvtParticle.cc:93
void setDiagonalSpinDensity()
Definition: EvtParticle.cc:133
void setLifetime(double tau)
Definition: EvtParticle.cc:89
void init(EvtId part_n, double e, double px, double py, double pz)

Referenced by generateMassTree(), and initDecay().

◆ initializePhaseSpace()

double EvtParticle::initializePhaseSpace ( int  numdaughter,
EvtId daughters,
double  poleSize = -1.,
int  whichTwo1 = 0,
int  whichTwo2 = 1 
)

Similar to the routine above except that here momentum is generated according to phase space daughters are filled with this momentum.

Definition at line 1071 of file EvtParticle.cc.

1073 {
1074
1075 double m_b;
1076 int i;
1077 //lange
1078 // this->makeDaughters(numdaughter,daughters);
1079
1080 static EvtVector4R p4[100];
1081 static double mass[100];
1082
1083 m_b = this->mass();
1084
1085 //lange - Jan2,2002 - Need to check to see if the daughters of the parent
1086 // have changed. If so, delete them and start over.
1087 //report(INFO,"EvtGen") << "the parent is\n";
1088 //if ( this->getParent() ) {
1089 // if ( this->getParent()->getParent() ) this->getParent()->getParent()->printTree();
1090 // this->getParent()->printTree();
1091 //}
1092 //report(INFO,"EvtGen") << "and this is\n";
1093 //if ( this) this->printTree();
1094 bool resetDaughters=false;
1095 if ( numdaughter != this->getNDaug() && this->getNDaug() > 0 ) resetDaughters=true;
1096 if ( numdaughter == this->getNDaug() )
1097 for (i=0; i<numdaughter;i++) {
1098 if ( this->getDaug(i)->getId() != daughters[i] ) resetDaughters=true;
1099 //report(INFO,"EvtGen") << this->getDaug(i)->getId() << " " << daughters[i] << endl;
1100 }
1101
1102 if ( resetDaughters ) {
1103 // report(INFO,"EvtGen") << "reseting daughters\n";
1104 //for (i=0; i<numdaughter;i++) {
1105 // report(INFO,"EvtGen") << "reset " <<i<< " "<< EvtPDL::name(this->getDaug(i)->getId()) << " " << EvtPDL::name(daughters[i]) << endl;
1106 //}
1107 bool t1=true;
1108 //but keep the decay channel of the parent.
1109 this->deleteDaughters(t1);
1110 this->makeDaughters(numdaughter,daughters);
1111 this->generateMassTree();
1112 }
1113
1114 double weight=0.;
1115 // EvtDecayBase::findMasses( this, numdaughter, daughters, mass );
1116 //get the list of masses
1117 //report(INFO,"EvtGen") << "mpar= " << m_b << " " << this <<endl;
1118 for (i=0; i<numdaughter;i++) {
1119 mass[i]=this->getDaug(i)->mass();
1120 // report(INFO,"EvtGen") << "mass " << i << " " << mass[i] << " " << this->getDaug(i) << endl;
1121 }
1122
1123 if ( poleSize<-0.1) {
1124 EvtGenKine::PhaseSpace( numdaughter, mass, p4, m_b );
1125 for(i=0;i<numdaughter;i++){
1126 this->getDaug(i)->init(daughters[i],p4[i]);
1127 }
1128
1129 }
1130 else {
1131 if ( numdaughter != 3 ) {
1132 report(ERROR,"EvtGen") << "Only can generate pole phase space "
1133 << "distributions for 3 body final states"
1134 << endl<<"Will terminate."<<endl;
1135 ::abort();
1136 }
1137 bool ok=false;
1138 if ( (whichTwo1 == 1 && whichTwo2 == 0 ) ||
1139 (whichTwo1 == 0 && whichTwo2 == 1 ) ) {
1141 poleSize, p4);
1142 //report(INFO,"EvtGen") << "here " << weight << " " << poleSize << endl;
1143 this->getDaug(0)->init(daughters[0],p4[0]);
1144 this->getDaug(1)->init(daughters[1],p4[1]);
1145 this->getDaug(2)->init(daughters[2],p4[2]);
1146 ok=true;
1147 }
1148 if ( (whichTwo1 == 1 && whichTwo2 == 2 ) ||
1149 (whichTwo1 == 2 && whichTwo2 == 1 ) ) {
1151 poleSize, p4);
1152 this->getDaug(0)->init(daughters[0],p4[2]);
1153 this->getDaug(1)->init(daughters[1],p4[1]);
1154 this->getDaug(2)->init(daughters[2],p4[0]);
1155 ok=true;
1156 }
1157 if ( (whichTwo1 == 0 && whichTwo2 == 2 ) ||
1158 (whichTwo1 == 2 && whichTwo2 == 0 ) ) {
1160 poleSize, p4);
1161 this->getDaug(0)->init(daughters[0],p4[1]);
1162 this->getDaug(1)->init(daughters[1],p4[0]);
1163 this->getDaug(2)->init(daughters[2],p4[2]);
1164 ok=true;
1165 }
1166 if ( !ok) {
1167 report(ERROR,"EvtGen") << "Invalid pair of particle to generate a pole dist"
1168 << whichTwo1 << " " << whichTwo2
1169 << endl<<"Will terminate."<<endl;
1170 ::abort();
1171 }
1172 }
1173
1174 return weight;
1175}
*********DOUBLE PRECISION m_pi INTEGER m_lenwt !max no of aux weights INTEGER m_phmax !maximum photon multiplicity ISR FSR *DOUBLE COMPLEX m_Pauli4 DOUBLE COMPLEX m_AmpBorn DOUBLE COMPLEX m_AmpBoxy DOUBLE COMPLEX m_AmpBorn1 DOUBLE COMPLEX m_AmpBorn2 DOUBLE COMPLEX m_AmpExpo2p DOUBLE COMPLEX m_Rmat DOUBLE COMPLEX m_BoxGZut !DOUBLE COMPLEX m_F1finPair2 !DOUBLE PRECISION m_Vcut DOUBLE PRECISION m_Alfinv DOUBLE PRECISION m_Lorin1 DOUBLE PRECISION m_Lorin2 DOUBLE PRECISION m_b
Definition: GPS.h:30
*********Class see also m_nmax DOUBLE PRECISION m_MasPhot DOUBLE PRECISION m_phsu DOUBLE PRECISION m_Xenph DOUBLE PRECISION m_r2 DOUBLE PRECISION m_WtMass INTEGER m_nmax INTEGER m_Nevgen INTEGER m_IsFSR INTEGER m_MarTot *COMMON c_KarFin $ !Output file $ !Event serial number $ !alpha QED at Thomson limit $ !minimum energy at CMS for remooval $ !infrared dimensionless $ !dummy photon IR regulator $ !crude photon multiplicity enhancement factor *EVENT $ !MC crude volume of PhhSpace *Sfactors $ !YFS formfactor IR part only $ !YFS formfactor non IR finite part $ !mass weight
Definition: KarFin.h:34
static double PhaseSpacePole(double M, double m1, double m2, double m3, double a, EvtVector4R p4[10])
Definition: EvtGenKine.cc:272
static double PhaseSpace(int ndaug, double mass[30], EvtVector4R p4[30], double mp)
Definition: EvtGenKine.cc:50
virtual void init(EvtId part_n, const EvtVector4R &p4)=0

Referenced by EvtPsi3Sdecay::choseDecay(), EvtAngH2::decay(), EvtAngSam::decay(), EvtAngSam3::decay(), EvtAngSamLab::decay(), EvtAngSamX::decay(), EvtAV2GV::decay(), EvtBHadronic::decay(), EvtBody3::decay(), EvtBsquark::decay(), EvtBto2piCPiso::decay(), EvtBTo4piCP::decay(), EvtBtoKD3P::decay(), EvtBtoKpiCPiso::decay(), EvtbTosllAli::decay(), EvtbTosllBall::decay(), EvtCalHelAmp::decay(), EvtChi0BB1::decay(), EvtChi0BB2::decay(), EvtChi1BB1::decay(), EvtChi1BB2::decay(), EvtChi2BB1::decay(), EvtChi2BB2::decay(), EvtConExc::decay(), EvtD0mixDalitz::decay(), EvtD0ToKpipi0pi0::decay(), EvtD0ToKpipipi::decay(), EvtDDalitz::decay(), EvtDeBD::decay(), EvtDIY::decay(), EvtDMix::decay(), EvtDToKpienu::decay(), EvtDToKSpipipi::decay(), EvtDTopipienu::decay(), EvtEDM::decay(), EvtEtaDalitz::decay(), EvtEtap2gpipi::decay(), EvtEtap2gpipiB::decay(), EvtFlatQ2::decay(), EvtHAngSam3::decay(), EvtHelAmp::decay(), EvtHelPPJ::decay(), EvtHQET::decay(), EvtHQET2::decay(), EvtHypNonLepton::decay(), EvtHypWK::decay(), EvtISGW::decay(), EvtISGW2::decay(), EvtJ2BB1::decay(), EvtJ2BB2::decay(), EvtJ2BB3::decay(), EvtJPE::decay(), EvtJpipi::decay(), EvtJTO3P::decay(), EvtKKLambdaC::decay(), EvtKstarnunu::decay(), EvtKstarstargamma::decay(), EvtLambdac2pKpi::decay(), EvtLambdaP_BarGamma::decay(), EvtLNuGamma::decay(), EvtMassH1::decay(), EvtMassH2::decay(), EvtMBody3::decay(), EvtMelikhov::decay(), EvtmH2::decay(), EvtmPhsp::decay(), EvtMultibody::decay(), EvtOmegaDalitz::decay(), EvtP2GC0::decay(), EvtP2GC1::decay(), EvtP2GC2::decay(), EvtPartWave::decay(), EvtPBB1::decay(), EvtPBB2::decay(), EvtPhiDalitz::decay(), EvtPhsp::decay(), EvtPi0Dalitz::decay(), EvtRexc::decay(), EvtRhoPi::decay(), EvtS2GV::decay(), EvtSLBKPole::decay(), EvtSll::decay(), EvtSLN::decay(), EvtSLPole::decay(), EvtSPL::decay(), EvtSSDCP::decay(), EvtSSSCP::decay(), EvtSSSCPpng::decay(), EvtSSSCPT::decay(), EvtSTS::decay(), EvtSTSCP::decay(), EvtSVPCP::decay(), EvtSVPHelAmp::decay(), EvtSVS::decay(), EvtSVSCP::decay(), EvtSVSCPiso::decay(), EvtSVSCPLH::decay(), EvtSVSNONCPEIGEN::decay(), EvtSVVNONCPEIGEN::decay(), EvtT2GV::decay(), EvtTauGamMu::decay(), EvtTauHadnu::decay(), EvtTaulnunu::decay(), EvtTauScalarnu::decay(), EvtTauVectornu::decay(), EvtTrackGen::decay(), EvtTSS::decay(), EvtTVSPwave::decay(), EvtVectorIsr::decay(), EvtVll::decay(), EvtVSPPwave::decay(), EvtVSS::decay(), EvtVSSBMixCPT::decay(), EvtVSSMix::decay(), EvtVub::decay(), EvtVubHybrid::decay(), EvtVubNLO::decay(), EvtVVP::decay(), EvtVVpipi::decay(), EvtVVPIPI_WEIGHTED::decay(), EvtVVSPwave::decay(), EvtLunda::ExclusiveDecay(), EvtConExc::findMaxXS(), EvtConExc::gamHXSection(), EvtPsi3Sdecay::PHSPDecay(), and EvtSVVHelAmp::SVVHel().

◆ insertDaugPtr()

void EvtParticle::insertDaugPtr ( int  idaug,
EvtParticle partptr 
)
inline

Makes partptr the idaug:th daugther.

Definition at line 220 of file EvtParticle.hh.

220 { _daug[idaug]=partptr;
221 partptr->_parent=this; }

Referenced by initDecay(), and EvtCPUtil::OtherB().

◆ isDecayed()

bool EvtParticle::isDecayed ( )
inline

Definition at line 384 of file EvtParticle.hh.

384{return _isDecayed;}

◆ isInitialized()

bool EvtParticle::isInitialized ( )
inline

Definition at line 382 of file EvtParticle.hh.

382{return _isInit;}

Referenced by EvtDecayBase::findMaxMass().

◆ makeDaughters()

void EvtParticle::makeDaughters ( int  ndaug,
EvtId id 
)

Creates the daughters in the list of ids and adds them to the parent. Note that momentum is left uninitialized, this is only creation.

Definition at line 1177 of file EvtParticle.cc.

1177 {
1178
1179 int i;
1180 if ( _channel < 0 ) {
1181 //report(INFO,"EvtGen") << "setting channel " << EvtPDL::name(this->getId()) << endl;
1182 setChannel(0);
1183 }
1184 EvtParticle* pdaug;
1185 if (_ndaug!=0 ){
1186 if (_ndaug!=ndaugstore){
1187 report(ERROR,"EvtGen") << "Asking to make a different number of "
1188 << "daughters than what was previously created."
1189 << endl<<"Will terminate."<<endl;
1190 ::abort();
1191 }
1192 }
1193 else{
1194 for(i=0;i<ndaugstore;i++){
1196 pdaug->setId(id[i]);
1197 pdaug->addDaug(this);
1198 }
1199
1200 } //else
1201} //makeDaughters
static EvtParticle * particleFactory(EvtSpinType::spintype spinType)
void setChannel(int i)
Definition: EvtParticle.cc:81
void addDaug(EvtParticle *node)
Definition: EvtParticle.cc:104
void setId(EvtId id)
Definition: EvtParticle.hh:365

Referenced by EvtbTosllAmp::CalcMaxProb(), EvtSemiLeptonicAmp::CalcMaxProb(), EvtPsi3Sdecay::choseDecay(), EvtBTo3piCP::decay(), EvtBToKpipiCP::decay(), EvtBtoXsEtap::decay(), EvtBtoXsgamma::decay(), EvtBtoXsll::decay(), EvtCBTo3piMPP::decay(), EvtCBTo3piP00::decay(), EvtConExc::decay(), EvtJetSet::decay(), EvtJscont::decay(), EvtLunda::decay(), EvtLundCharm::decay(), EvtOpenCharm::decay(), EvtPhokhara::decay(), EvtPhokhara_4pi::decay(), EvtPhokhara_K0K0::decay(), EvtPhokhara_KK::decay(), EvtPhokhara_LLB::decay(), EvtPhokhara_nnbar::decay(), EvtPhokhara_pi0pi0pipi::decay(), EvtPhokhara_pipi::decay(), EvtPhokhara_pipieta::decay(), EvtPhokhara_pipipi0::decay(), EvtPhokhara_ppbar::decay(), EvtPycont::decay(), EvtPyGaGa::decay(), EvtPythia::decay(), EvtRexc::decay(), EvtSingleParticle::decay(), EvtSingleParticle2::decay(), EvtSinglePoint::decay(), EvtTauola::decay(), EvtLunda::ExclusiveDecay(), EvtConExc::findMaxXS(), EvtConExc::gamHXSection(), initDecay(), initializePhaseSpace(), and EvtPsi3Sdecay::PHSPDecay().

◆ makeStdHep() [1/2]

void EvtParticle::makeStdHep ( EvtStdHep stdhep)

Definition at line 795 of file EvtParticle.cc.

795 {
796
797 //first add particle to the stdhep list;
798 stdhep.createParticle(getP4Lab(),get4Pos(),-1,-1,
800
801 int i;
802 for(i=0;i<_ndaug;i++){
803 stdhep.createParticle(_daug[i]->getP4Lab(),_daug[i]->get4Pos(),0,0,
804 EvtPDL::getStdHep(_daug[i]->getId()));
805 }
806
807 for(i=0;i<_ndaug;i++){
808 _daug[i]->makeStdHepRec(1+i,1+i,stdhep);
809 }
810 return;
811
812}
EvtVector4R getP4Lab()
Definition: EvtParticle.cc:685
EvtVector4R get4Pos()
Definition: EvtParticle.cc:706
void createParticle(EvtVector4R p4, EvtVector4R x, int prntfirst, int prntlast, int id)
Definition: EvtStdHep.cc:41

◆ makeStdHep() [2/2]

void EvtParticle::makeStdHep ( EvtStdHep stdhep,
EvtSecondary secondary,
EvtId stable_parent_ihep 
)

Makes stdhep list

Definition at line 759 of file EvtParticle.cc.

760 {
761
762 //first add particle to the stdhep list;
763 stdhep.createParticle(getP4Lab(),get4Pos(),-1,-1,
765
766 int ii=0;
767
768 //lets see if this is a longlived particle and terminate the
769 //list building!
770
771 while (list_of_stable[ii]!=EvtId(-1,-1)) {
772 if (getId()==list_of_stable[ii]){
773 secondary.createSecondary(0,this);
774 return;
775 }
776 ii++;
777 }
778
779
780
781
782 int i;
783 for(i=0;i<_ndaug;i++){
784 stdhep.createParticle(_daug[i]->getP4Lab(),_daug[i]->get4Pos(),0,0,
785 EvtPDL::getStdHep(_daug[i]->getId()));
786 }
787
788 for(i=0;i<_ndaug;i++){
789 _daug[i]->makeStdHepRec(1+i,1+i,stdhep,secondary,list_of_stable);
790 }
791 return;
792
793}
void createSecondary(int stdhepindex, EvtParticle *prnt)
Definition: EvtSecondary.cc:40

Referenced by EvtGen::generateDecay(), and EvtGen::generateEvent().

◆ mass()

double EvtParticle::mass ( ) const

Returns mass of particle.

Definition at line 127 of file EvtParticle.cc.

127 {
128
129 return _p.mass();
130}

Referenced by EvtbTosllScalarAmp::CalcAmp(), EvtbTosllVectorAmp::CalcAmp(), EvtSemiLeptonicBaryonAmp::CalcAmp(), EvtSemiLeptonicTensorAmp::CalcAmp(), EvtSemiLeptonicVectorAmp::CalcAmp(), EvtbTosllAmp::CalcMaxProb(), EvtSemiLeptonicAmp::CalcMaxProb(), EvtConExc::checkdecay(), EvtPsi3Sdecay::choseDecay(), compMassProb(), EvtBHadronic::decay(), EvtBsquark::decay(), EvtBtoXsEtap::decay(), EvtBtoXsgamma::decay(), EvtBtoXsll::decay(), EvtConExc::decay(), EvtDDalitz::decay(), EvtEDM::decay(), EvtEtaDalitz::decay(), EvtHypWK::decay(), EvtJetSet::decay(), EvtJpipi::decay(), EvtJscont::decay(), EvtKstarnunu::decay(), EvtLunda::decay(), EvtLundCharm::decay(), EvtOmegaDalitz::decay(), EvtOpenCharm::decay(), EvtPhokhara::decay(), EvtPhokhara_4pi::decay(), EvtPhokhara_K0K0::decay(), EvtPhokhara_KK::decay(), EvtPhokhara_LLB::decay(), EvtPhokhara_nnbar::decay(), EvtPhokhara_pi0pi0pipi::decay(), EvtPhokhara_pipi::decay(), EvtPhokhara_pipieta::decay(), EvtPhokhara_pipipi0::decay(), EvtPhokhara_ppbar::decay(), EvtPycont::decay(), EvtPyGaGa::decay(), EvtPythia::decay(), EvtRexc::decay(), EvtSll::decay(), EvtSLN::decay(), EvtSSDCP::decay(), EvtSTS::decay(), EvtSTSCP::decay(), EvtSVS::decay(), EvtSVSCPLH::decay(), EvtSVSNONCPEIGEN::decay(), EvtTauola::decay(), EvtTVSPwave::decay(), EvtVectorIsr::decay(), EvtVll::decay(), EvtVPHOtoVISR::decay(), EvtVPHOtoVISRHi::decay(), EvtVSPPwave::decay(), EvtVub::decay(), EvtVubHybrid::decay(), EvtVubNLO::decay(), EvtVVpipi::decay(), EvtVVPIPI_WEIGHTED::decay(), EvtPHOTOS::doRadCorr(), EvtLunda::ExclusiveDecay(), EvtDecayBase::findMass(), EvtDecayBase::findMasses(), EvtDecayBase::findMaxMass(), generateMassTree(), get4Pos(), EvtParticleDecayList::getDecayModel(), getP4Restframe(), EvtRaritaSchwingerParticle::init(), EvtDiracParticle::init(), EvtPhokhara::init_evt(), EvtPhokhara::init_mode(), initDecay(), initializePhaseSpace(), EvtDecayAmp::makeDecay(), EvtDecayProb::makeDecay(), EvtPsi3Sdecay::PHSPDecay(), printTreeRec(), EvtDiracParticle::rotateToHelicityBasis(), and EvtRaritaSchwingerParticle::rotateToHelicityBasis().

◆ nextIter()

EvtParticle * EvtParticle::nextIter ( EvtParticle rootOfTree = 0)

Iterates over the particles in a decay chain.

Definition at line 729 of file EvtParticle.cc.

729 {
730
731 EvtParticle *bpart;
732 EvtParticle *current;
733
734 current=this;
735 int i;
736
737 if (_ndaug!=0) return _daug[0];
738
739 do{
740 bpart=current->_parent;
741 if (bpart==0) return 0;
742 i=0;
743 while (bpart->_daug[i]!=current) {i++;}
744
745 if ( bpart==rootOfTree ) {
746 if ( i== bpart->_ndaug-1 ) return 0;
747 }
748
749 i++;
750 current=bpart;
751
752 }while(i>=bpart->_ndaug);
753
754 return bpart->_daug[i];
755
756}

◆ noLifeTime()

void EvtParticle::noLifeTime ( )
inline

Definition at line 362 of file EvtParticle.hh.

362{ _genlifetime=0; }

Referenced by EvtbTosllAmp::CalcMaxProb(), and EvtSemiLeptonicAmp::CalcMaxProb().

◆ printParticle()

void EvtParticle::printParticle ( ) const

Prints information for the particle.

Definition at line 1000 of file EvtParticle.cc.

1000 {
1001
1002 switch (EvtPDL::getSpinType(_id)){
1004 report(INFO,"EvtGen") << "This is a scalar particle:"<<EvtPDL::name(_id).c_str()<<"\n";
1005 break;
1007 report(INFO,"EvtGen") << "This is a vector particle:"<<EvtPDL::name(_id).c_str()<<"\n";
1008 break;
1010 report(INFO,"EvtGen") << "This is a tensor particle:"<<EvtPDL::name(_id).c_str()<<"\n";
1011 break;
1012 case EvtSpinType::DIRAC:
1013 report(INFO,"EvtGen") << "This is a dirac particle:"<<EvtPDL::name(_id).c_str()<<"\n";
1014 break;
1016 report(INFO,"EvtGen") << "This is a photon:"<<EvtPDL::name(_id).c_str()<<"\n";
1017 break;
1019 report(INFO,"EvtGen") << "This is a neutrino:"<<EvtPDL::name(_id).c_str()<<"\n";
1020 break;
1022 report(INFO,"EvtGen") << "This is a raritaschwinger:"<<EvtPDL::name(_id).c_str()<<"\n";
1023 break;
1025 report(INFO,"EvtGen") << "This is a string:"<<EvtPDL::name(_id).c_str()<<"\n";
1026 break;
1027 default:
1028 report(INFO,"EvtGen") <<"Unknown particle type in EvtParticle::printParticle()"<<endl;
1029 break;
1030 }
1031 report(INFO,"EvtGen") << "Number of daughters:"<<_ndaug<<"\n";
1032
1033
1034}

Referenced by eps(), epsParent(), epsParentPhoton(), epsPhoton(), epsTensor(), epsTensorParent(), sp(), spNeutrino(), spParent(), and spParentNeutrino().

◆ printTree()

void EvtParticle::printTree ( ) const

Prints out the particle "tree" of a given particle. The tree consists of all daughters (and their daughters, etc) and their properties.

Definition at line 897 of file EvtParticle.cc.

897 {
898
899 report(INFO,"EvtGen") << "This is the current decay chain"<<endl;
900 report(INFO,"") << "This top particle is "<<
901 EvtPDL::name(_id).c_str()<<endl;
902
903 this->printTreeRec(0);
904 report(INFO,"EvtGen") << "End of decay chain."<<endl;
905
906}
void printTreeRec(int level) const
Definition: EvtParticle.cc:870

Referenced by EvtConExc::decay(), EvtDecayBase::findMass(), and generateMassTree().

◆ printTreeRec()

void EvtParticle::printTreeRec ( int  level) const

Definition at line 870 of file EvtParticle.cc.

870 {
871
872 int newlevel,i;
873 newlevel = level +1;
874
875
876 if (_ndaug!=0) {
877 if ( level > 0 ) {
878 for (i=0;i<(5*level);i++) {
879 report(INFO,"") <<" ";
880 }
881 }
882 report(INFO,"") << EvtPDL::name(_id).c_str();
883 report(INFO,"") << " -> ";
884 for(i=0;i<_ndaug;i++){
885 report(INFO,"") << EvtPDL::name(_daug[i]->getId()).c_str()<<" ";
886 }
887 for(i=0;i<_ndaug;i++){
888 report(INFO,"") << _daug[i]->mass()<<" ";
889 }
890 report(INFO,"")<<endl;
891 for(i=0;i<_ndaug;i++){
892 _daug[i]->printTreeRec(newlevel);
893 }
894 }
895}

Referenced by printTree(), and printTreeRec().

◆ resetFirstOrNot()

void EvtParticle::resetFirstOrNot ( )

Definition at line 77 of file EvtParticle.cc.

77 {
78 _first=1;
79}

Referenced by EvtGen::generateDecay().

◆ resetNDaug()

void EvtParticle::resetNDaug ( )
inline

Definition at line 269 of file EvtParticle.hh.

269{_ndaug=0; return;}

Referenced by EvtOpenCharm::decay(), and EvtGen::generateDecay().

◆ rotateToHelicityBasis() [1/2]

virtual EvtSpinDensity EvtParticle::rotateToHelicityBasis ( ) const
pure virtual

Returns a rotation matrix need to rotate the basis state to the helicity basis. The EvtSpinDensity matrix is just use as a matrix here. This function is to be implemented in each derived class.

Implemented in EvtDiracParticle, EvtHighSpinParticle, EvtNeutrinoParticle, EvtPhotonParticle, EvtRaritaSchwingerParticle, EvtScalarParticle, EvtStringParticle, EvtTensorParticle, and EvtVectorParticle.

Referenced by EvtMultibody::decay(), and setSpinDensityForwardHelicityBasis().

◆ rotateToHelicityBasis() [2/2]

virtual EvtSpinDensity EvtParticle::rotateToHelicityBasis ( double  alpha,
double  beta,
double  gamma 
) const
pure virtual

◆ setChannel()

void EvtParticle::setChannel ( int  i)

Should only be used internally.

Definition at line 81 of file EvtParticle.cc.

81 {
82 _channel=i;
83}

Referenced by EvtParticleDecayList::getDecayModel(), and makeDaughters().

◆ setDecayProb()

void EvtParticle::setDecayProb ( double  p)

Definition at line 1204 of file EvtParticle.cc.

1204 {
1205
1206 if ( _decayProb == 0 ) _decayProb=new double;
1207 *_decayProb=prob;
1208}

Referenced by EvtDecayAmp::makeDecay(), EvtDecayIncoherent::makeDecay(), and EvtDecayProb::makeDecay().

◆ setDiagonalSpinDensity()

void EvtParticle::setDiagonalSpinDensity ( )

Set diagonal spindensity matrix.

Definition at line 133 of file EvtParticle.cc.

133 {
134
135 _rhoForward.SetDiag(getSpinStates());
136}

Referenced by EvtbTosllAmp::CalcMaxProb(), EvtSemiLeptonicAmp::CalcMaxProb(), EvtJetSet::decay(), EvtLunda::decay(), EvtLundCharm::decay(), EvtPythia::decay(), EvtTauola::decay(), and initDecay().

◆ setFirstOrNot()

void EvtParticle::setFirstOrNot ( )

Definition at line 74 of file EvtParticle.cc.

74 {
75 _first=0;
76}

Referenced by EvtDecayBase::findMasses().

◆ setGeneratorFlag()

void EvtParticle::setGeneratorFlag ( int  flag)
inline

set generator information; pingrg-2011-1-6

Definition at line 141 of file EvtParticle.hh.

141{_generatorFlag = flag;}

Referenced by EvtLundCharm::decay(), and EvtOpenCharm::decay().

◆ setId()

void EvtParticle::setId ( EvtId  id)
inline

Definition at line 365 of file EvtParticle.hh.

365{ _id=id;}

Referenced by makeDaughters().

◆ setInclusiveMode()

void EvtParticle::setInclusiveMode ( int  im)
inline

Definition at line 393 of file EvtParticle.hh.

393{_inclusiveMode=im ;}

◆ setIntFlag()

void EvtParticle::setIntFlag ( std::vector< int >  vi)
inline

set int flag for ConExc: pingrg-2015-2-7

Definition at line 151 of file EvtParticle.hh.

151 {
152 _intFlag=vi;
153 }

Referenced by EvtConExc::decay().

◆ setLifetime() [1/2]

void EvtParticle::setLifetime ( )

Generate lifetime according to pure exponential.

Definition at line 93 of file EvtParticle.cc.

93 {
94 if (_genlifetime){
96 }
97}
static double getctau(EvtId i)
Definition: EvtPDL.hh:55

Referenced by EvtHighSpinParticle::init(), EvtRaritaSchwingerParticle::init(), EvtDiracParticle::init(), EvtNeutrinoParticle::init(), EvtTensorParticle::init(), EvtPhotonParticle::init(), EvtScalarParticle::init(), EvtVectorParticle::init(), and initDecay().

◆ setLifetime() [2/2]

void EvtParticle::setLifetime ( double  tau)

Set lifetime of the particle in parents restframe.

Definition at line 89 of file EvtParticle.cc.

89 {
90 _t=tau;
91}

Referenced by EvtD0mixDalitz::decay(), EvtDMix::decay(), EvtVSSBMixCPT::decay(), EvtVSSMix::decay(), EvtVub::decay(), EvtVubHybrid::decay(), initDecay(), and EvtCPUtil::OtherB().

◆ setMass()

◆ setp() [1/2]

void EvtParticle::setp ( const EvtVector4R p4)
inlineprotected

Definition at line 399 of file EvtParticle.hh.

399{ _p =p4; }

◆ setp() [2/2]

void EvtParticle::setp ( double  e,
double  px,
double  py,
double  pz 
)
inlineprotected

◆ setP4()

void EvtParticle::setP4 ( const EvtVector4R p4)
inline

Sets the 4momentum in the parents restframe.

Definition at line 258 of file EvtParticle.hh.

258{_p=p4;}

Referenced by EvtPHOTOS::doRadCorr().

◆ setpart_num()

void EvtParticle::setpart_num ( EvtId  particle_number)
inlineprotected

Definition at line 400 of file EvtParticle.hh.

401 {
402 assert(_channel==-10||
403 _id.getId()==particle_number.getId()||
404 _id.getId()==-1);
405 _id = particle_number;
406 }
int getId() const
Definition: EvtId.hh:41

Referenced by EvtHighSpinParticle::init(), EvtRaritaSchwingerParticle::init(), EvtStringParticle::init(), EvtDiracParticle::init(), EvtNeutrinoParticle::init(), EvtPhotonParticle::init(), EvtScalarParticle::init(), EvtTensorParticle::init(), and EvtVectorParticle::init().

◆ setPolarizedSpinDensity()

void EvtParticle::setPolarizedSpinDensity ( double  r00,
double  r11,
double  r22 
)

Definition at line 157 of file EvtParticle.cc.

157 {
158
159 if (getSpinStates()!=3) {
160 report(ERROR,"EvtGen")<<"Error in EvtParticle::setVectorSpinDensity"<<endl;
161 report(ERROR,"EvtGen")<<"spin_states:"<<getSpinStates()<<endl;
162 report(ERROR,"EvtGen")<<"particle:"<<EvtPDL::name(_id).c_str()<<endl;
163 ::abort();
164 }
165
166 EvtSpinDensity rho;
167
168 //Set helicity +1 and -1 to 1.
169 rho.SetDiag(getSpinStates());
170 rho.Set(0,0,EvtComplex(r00,0.));
171 rho.Set(1,1,EvtComplex(r11,0.));
172 rho.Set(2,2,EvtComplex(r22,0.));
174
175}
void setSpinDensityForwardHelicityBasis(const EvtSpinDensity &rho)
Definition: EvtParticle.cc:178
void Set(int i, int j, const EvtComplex &rhoij)

◆ setSpinDensityBackward()

void EvtParticle::setSpinDensityBackward ( const EvtSpinDensity rho)
inline

Set backward spin density matrix.

Definition at line 352 of file EvtParticle.hh.

352{_rhoBackward=rho;}

Referenced by EvtDecayAmp::makeDecay(), EvtDecayIncoherent::makeDecay(), and EvtDecayProb::makeDecay().

◆ setSpinDensityForward()

void EvtParticle::setSpinDensityForward ( const EvtSpinDensity rho)
inline

Set forward spin density matrix.

Definition at line 321 of file EvtParticle.hh.

321{_rhoForward=rho;}

Referenced by EvtVectorIsr::decay(), EvtDecayAmp::makeDecay(), EvtDecayIncoherent::makeDecay(), EvtDecayProb::makeDecay(), and EvtParticleFactory::particleFactory().

◆ setSpinDensityForwardHelicityBasis() [1/2]

void EvtParticle::setSpinDensityForwardHelicityBasis ( const EvtSpinDensity rho)

Set forward spin density matrix according to the density matrix rho in the helicity amplitude basis.

Definition at line 178 of file EvtParticle.cc.

178 {
179
181
182 assert(R.GetDim()==rho.GetDim());
183
184 int n=rho.GetDim();
185
186 _rhoForward.SetDim(n);
187
188 int i,j,k,l;
189
190 for(i=0;i<n;i++){
191 for(j=0;j<n;j++){
192 EvtComplex tmp=0.0;
193 for(k=0;k<n;k++){
194 for(l=0;l<n;l++){
195 tmp+=R.Get(l,i)*rho.Get(l,k)*conj(R.Get(k,j));
196 }
197 }
198 _rhoForward.Set(i,j,tmp);
199 }
200 }
201
202 //report(INFO,"EvtGen") << "_rhoForward:"<<_rhoForward<<endl;
203
204}
const Int_t n
Evt3Rank3C conj(const Evt3Rank3C &t2)
Definition: Evt3Rank3C.cc:175
virtual EvtSpinDensity rotateToHelicityBasis() const =0
int GetDim() const
const EvtComplex & Get(int i, int j) const
void SetDim(int n)
complex_t R(double Q2, double M2, double G, double Mp2, double Mm2)
Definition: TUtil.h:27

Referenced by setPolarizedSpinDensity(), and setVectorSpinDensity().

◆ setSpinDensityForwardHelicityBasis() [2/2]

void EvtParticle::setSpinDensityForwardHelicityBasis ( const EvtSpinDensity rho,
double  alpha,
double  beta,
double  gamma 
)

Definition at line 206 of file EvtParticle.cc.

209 {
210
212
213 assert(R.GetDim()==rho.GetDim());
214
215 int n=rho.GetDim();
216
217 _rhoForward.SetDim(n);
218
219 int i,j,k,l;
220
221 for(i=0;i<n;i++){
222 for(j=0;j<n;j++){
223 EvtComplex tmp=0.0;
224 for(k=0;k<n;k++){
225 for(l=0;l<n;l++){
226 tmp+=R.Get(l,i)*rho.Get(l,k)*conj(R.Get(k,j));
227 }
228 }
229 _rhoForward.Set(i,j,tmp);
230 }
231 }
232
233 //report(INFO,"EvtGen") << "_rhoForward:"<<_rhoForward<<endl;
234
235}
const double alpha

◆ setVectorSpinDensity()

void EvtParticle::setVectorSpinDensity ( )

Set spindensity matrix for e+e- -> V

Definition at line 138 of file EvtParticle.cc.

138 {
139
140 if (getSpinStates()!=3) {
141 report(ERROR,"EvtGen")<<"Error in EvtParticle::setVectorSpinDensity"<<endl;
142 report(ERROR,"EvtGen")<<"spin_states:"<<getSpinStates()<<endl;
143 report(ERROR,"EvtGen")<<"particle:"<<EvtPDL::name(_id).c_str()<<endl;
144 ::abort();
145 }
146
147 EvtSpinDensity rho;
148
149 //Set helicity +1 and -1 to 1.
150 rho.SetDiag(getSpinStates());
151 rho.Set(1,1,EvtComplex(0.0,0.0));
153
154}

Referenced by EvtGen::generateEvent().

◆ sp()

EvtDiracSpinor EvtParticle::sp ( int  i) const
virtual

Returns Dirac spinor in the particles own restframe for a Dirac particle.

Reimplemented in EvtDiracParticle.

Definition at line 622 of file EvtParticle.cc.

622 {
623 EvtDiracSpinor tempD;
624 int temp;
625 temp = i;
627 report(ERROR,"EvtGen") << "and you have asked for the:"<<i
628 <<"th dirac spinor."
629 <<" I.e. you thought it was a"
630 <<" Dirac particle!" << endl;
631 ::abort();
632 return tempD;
633}

Referenced by EvtHypNonLepton::calcAmp(), EvtSemiLeptonicBaryonAmp::CalcAmp(), EvtTauHadnu::decay(), EvtTaulnunu::decay(), EvtTauScalarnu::decay(), and EvtTauVectornu::decay().

◆ spNeutrino()

EvtDiracSpinor EvtParticle::spNeutrino ( ) const
virtual

Returns Dirac spinor in the particles own restframe for a Neutrino particle.

Reimplemented in EvtNeutrinoParticle.

Definition at line 646 of file EvtParticle.cc.

646 {
647 EvtDiracSpinor tempD;
649 report(ERROR,"EvtGen") << "and you have asked for the "
650 <<"dirac spinor."
651 <<" I.e. you thought it was a"
652 <<" neutrino particle!" << endl;
653 ::abort();
654 return tempD;
655}

◆ spParent()

EvtDiracSpinor EvtParticle::spParent ( int  i) const
virtual

Returns Dirac spinor in the parents restframe for a Dirac particle.

Reimplemented in EvtDiracParticle.

Definition at line 609 of file EvtParticle.cc.

609 {
610 EvtDiracSpinor tempD;
611 int temp;
612 temp = i;
614 report(ERROR,"EvtGen") << "and you have asked for the:"<<i
615 <<"th dirac spinor."
616 <<" I.e. you thought it was a"
617 <<" Dirac particle!" << endl;
618 ::abort();
619 return tempD;
620}

Referenced by EvtHypNonLepton::calcAmp(), EvtbTosllScalarAmp::CalcAmp(), EvtbTosllVectorAmp::CalcAmp(), EvtSemiLeptonicBaryonAmp::CalcAmp(), EvtSemiLeptonicTensorAmp::CalcAmp(), EvtSemiLeptonicVectorAmp::CalcAmp(), EvtBsquark::decay(), EvtEDM::decay(), EvtLNuGamma::decay(), EvtSll::decay(), EvtSLN::decay(), EvtTaulnunu::decay(), and EvtVll::decay().

◆ spParentNeutrino()

EvtDiracSpinor EvtParticle::spParentNeutrino ( ) const
virtual

Returns Dirac spinor in the parents restframe for a Neutrino particle.

Reimplemented in EvtNeutrinoParticle.

Definition at line 635 of file EvtParticle.cc.

635 {
636 EvtDiracSpinor tempD;
638 report(ERROR,"EvtGen") << "and you have asked for the "
639 <<"dirac spinor."
640 <<" I.e. you thought it was a"
641 <<" neutrino particle!" << endl;
642 ::abort();
643 return tempD;
644}

Referenced by EvtSemiLeptonicBaryonAmp::CalcAmp(), EvtSemiLeptonicTensorAmp::CalcAmp(), EvtSemiLeptonicVectorAmp::CalcAmp(), EvtKstarnunu::decay(), EvtLNuGamma::decay(), EvtSLN::decay(), EvtTauHadnu::decay(), EvtTaulnunu::decay(), EvtTauScalarnu::decay(), and EvtTauVectornu::decay().

◆ treeStr()

std::string EvtParticle::treeStr ( ) const

Definition at line 990 of file EvtParticle.cc.

990 {
991
992 std::string retval=EvtPDL::name(_id);
993 retval+=" -> ";
994
995 retval+=treeStrRec(0);
996
997 return retval;
998}
std::string treeStrRec(int level) const
Definition: EvtParticle.cc:908

◆ treeStrRec()

std::string EvtParticle::treeStrRec ( int  level) const

Definition at line 908 of file EvtParticle.cc.

908 {
909
910 int newlevel,i;
911 newlevel = level +1;
912
913 std::string retval="";
914
915 for(i=0;i<_ndaug;i++){
916 retval+=EvtPDL::name(_daug[i]->getId());
917 if ( _daug[i]->getNDaug() > 0 ) {
918 retval+= " (";
919 retval+= _daug[i]->treeStrRec(newlevel);
920 retval+= ") ";
921 }
922 else{
923 if ( i!=_ndaug-1) retval+=" ";
924 }
925 }
926
927 return retval;
928}

Referenced by treeStr(), and treeStrRec().

◆ writeTreeRec()

std::string EvtParticle::writeTreeRec ( std::string  resonance) const

Definition at line 930 of file EvtParticle.cc.

930 { //pingrg, Jun. 6, 2008
931 std::string retval="";
932
933 if (resonance == EvtPDL::name(_id).c_str() && _ndaug!=0) {
934 retval=resonance+": "+IntToStr(_ndaug)+" ";
935 for(int i=0;i<_ndaug;i++){
936 retval += EvtPDL::name(_daug[i]->getId()).c_str();
937 retval += " ";
938 }
939 }
940
941 for(int i=0;i<_ndaug;i++){
942 _daug[i]->writeTreeRec(resonance);
943 }
944 std::cout<<retval;
945 return retval;
946}
std::string writeTreeRec(std::string) const
Definition: EvtParticle.cc:930
char * c_str(Index i)
Definition: EvtCyclic3.cc:252

Referenced by writeTreeRec().

Member Data Documentation

◆ _validP4


The documentation for this class was generated from the following files: