EvtRelBreitWignerBarrierFact Class Reference

#include <EvtRelBreitWignerBarrierFact.hh>

Inheritance diagram for EvtRelBreitWignerBarrierFact:

EvtAbsLineShape List of all members.

Public Member Functions

 EvtRelBreitWignerBarrierFact ()
 EvtRelBreitWignerBarrierFact (double mass, double width, double maxRange, EvtSpinType::spintype sp)
 ~EvtRelBreitWignerBarrierFact ()
EvtRelBreitWignerBarrierFactoperator= (const EvtRelBreitWignerBarrierFact &x)
 EvtRelBreitWignerBarrierFact (const EvtRelBreitWignerBarrierFact &x)
EvtAbsLineShapeclone ()
double getMassProb (double mass, double massPar, int nDaug, double *massDau)
double getRandMass (EvtId *parId, int nDaug, EvtId *dauId, EvtId *othDaugId, double maxMass, double *dauMasses)
virtual void reSetBlatt (double blatt)
double getMass ()
double getMassMin ()
double getMassMax ()
double getMaxRange ()
double getWidth ()
EvtSpinType::spintype getSpinType ()
virtual double rollMass ()
void reSetMass (double mass)
void reSetWidth (double width)
void reSetMassMin (double mass)
void reSetMassMax (double mass)
void includeBirthFactor (bool yesno)
void addFactorPn (double factor=0.)
void includeDecayFactor (bool yesno)
void setPWForDecay (int spin, EvtId d1, EvtId d2)
void setPWForBirthL (int spin, EvtId par, EvtId othD)
void fixForSP8 ()

Protected Attributes

double _blatt
bool _errorCond
bool _includeDecayFact
bool _includeBirthFact
double _addFactorPn
double _mass
double _massMin
double _massMax
double _width
double _maxRange
std::vector< EvtId_userSetPWD1
std::vector< EvtId_userSetPWD2
std::vector< int > _userSetPW
std::vector< EvtId_userSetBirthPar
std::vector< EvtId_userSetBirthOthD
std::vector< int > _userSetBirthPW
EvtSpinType::spintype _spin
bool _applyFixForSP8

Detailed Description

Definition at line 27 of file EvtRelBreitWignerBarrierFact.hh.


Constructor & Destructor Documentation

EvtRelBreitWignerBarrierFact::EvtRelBreitWignerBarrierFact (  ) 

Definition at line 36 of file EvtRelBreitWignerBarrierFact.cc.

Referenced by clone().

00036                                                            {
00037 
00038 }

EvtRelBreitWignerBarrierFact::EvtRelBreitWignerBarrierFact ( double  mass,
double  width,
double  maxRange,
EvtSpinType::spintype  sp 
)

Definition at line 43 of file EvtRelBreitWignerBarrierFact.cc.

References _blatt, _errorCond, EvtAbsLineShape::_includeBirthFact, EvtAbsLineShape::_includeDecayFact, EvtAbsLineShape::_mass, EvtAbsLineShape::_massMax, EvtAbsLineShape::_massMin, EvtAbsLineShape::_maxRange, EvtAbsLineShape::_spin, and EvtAbsLineShape::_width.

00043                                                                                                                              :
00044   EvtAbsLineShape(mass,width,maxRange,sp)
00045 { // double mDaug1, double mDaug2, int l) {
00046 
00047   _includeDecayFact=true;
00048   _includeBirthFact=true;
00049   _mass=mass;
00050   _width=width;
00051   _spin=sp;
00052   _blatt=3.0;
00053   _maxRange=maxRange;
00054   _errorCond=false;
00055 
00056   double maxdelta = 15.0*width;
00057 
00058   if ( maxRange > 0.00001 ) {
00059     _massMax=mass+maxdelta;
00060     _massMin=mass-maxRange;
00061   }
00062   else{
00063     _massMax=mass+maxdelta;
00064     _massMin=mass-15.0*width;
00065   }
00066  
00067   _massMax=mass+maxdelta;
00068   if ( _massMin< 0. ) _massMin=0.;
00069 }

EvtRelBreitWignerBarrierFact::~EvtRelBreitWignerBarrierFact (  ) 

Definition at line 40 of file EvtRelBreitWignerBarrierFact.cc.

00040                                                             {
00041 }

EvtRelBreitWignerBarrierFact::EvtRelBreitWignerBarrierFact ( const EvtRelBreitWignerBarrierFact x  ) 

Definition at line 71 of file EvtRelBreitWignerBarrierFact.cc.

References _blatt, _errorCond, EvtAbsLineShape::_includeBirthFact, EvtAbsLineShape::_includeDecayFact, EvtAbsLineShape::_massMax, EvtAbsLineShape::_massMin, EvtAbsLineShape::_maxRange, and x.

00071                                                                                                 :
00072   EvtAbsLineShape(x)
00073 {
00074   _massMax=x._massMax;
00075   _massMin=x._massMin;
00076   _blatt=x._blatt;
00077   _maxRange=x._maxRange;
00078   _includeDecayFact=x._includeDecayFact;
00079   _includeBirthFact=x._includeBirthFact;
00080   _errorCond=x._errorCond;
00081 
00082 }


Member Function Documentation

void EvtAbsLineShape::addFactorPn ( double  factor = 0.  )  [inline, inherited]

Definition at line 71 of file EvtAbsLineShape.hh.

References EvtAbsLineShape::_addFactorPn.

Referenced by EvtPartProp::addFactorPn(), and getRandMass().

00071 { _addFactorPn = factor;}

EvtAbsLineShape * EvtRelBreitWignerBarrierFact::clone (  )  [virtual]

Reimplemented from EvtAbsLineShape.

Definition at line 99 of file EvtRelBreitWignerBarrierFact.cc.

References EvtRelBreitWignerBarrierFact().

00099                                                      {
00100 
00101   return new EvtRelBreitWignerBarrierFact(*this);
00102 }

void EvtAbsLineShape::fixForSP8 (  )  [inline, inherited]

Definition at line 87 of file EvtAbsLineShape.hh.

References EvtAbsLineShape::_applyFixForSP8.

Referenced by EvtPartProp::fixLSForSP8().

00087 { _applyFixForSP8=true;}

double EvtAbsLineShape::getMass (  )  [inline, inherited]

Definition at line 41 of file EvtAbsLineShape.hh.

References EvtAbsLineShape::_mass.

Referenced by EvtPartProp::getMass(), and EvtPartProp::newLineShape().

00041 {return _mass;}

double EvtAbsLineShape::getMassMax (  )  [inline, inherited]

Definition at line 43 of file EvtAbsLineShape.hh.

References EvtAbsLineShape::_massMax.

Referenced by EvtPartProp::getMassMax().

00043 {return _massMax;}

double EvtAbsLineShape::getMassMin (  )  [inline, inherited]

Definition at line 42 of file EvtAbsLineShape.hh.

References EvtAbsLineShape::_massMin.

Referenced by EvtPartProp::getMassMin().

00042 {return _massMin;} 

double EvtRelBreitWignerBarrierFact::getMassProb ( double  mass,
double  massPar,
int  nDaug,
double *  massDau 
) [virtual]

Reimplemented from EvtAbsLineShape.

Definition at line 105 of file EvtRelBreitWignerBarrierFact.cc.

References _errorCond, EvtAbsLineShape::_width, EvtAbsLineShape::getMassProb(), and genRecEmupikp::i.

00105                                                                                                        {
00106 
00107   _errorCond=false;
00108   //return EvtAbsLineShape::getMassProb(mass,massPar,nDaug,massDau);
00109   if (nDaug!=2) return EvtAbsLineShape::getMassProb(mass,massPar,nDaug,massDau);
00110 
00111   double dTotMass=0.;
00112 
00113   int i;
00114   for (i=0; i<nDaug; i++) {
00115     dTotMass+=massDau[i];
00116   }
00117   //report(INFO,"EvtGen") << mass << " " << massPar << " " << dTotMass << " "<< endl;
00118   //    if ( (mass-dTotMass)<0.0001 ) return 0.;
00119   //report(INFO,"EvtGen") << mass << " " << dTotMass << endl;
00120   if ( (mass<dTotMass) ) return 0.;
00121 
00122   if ( _width< 0.0001) return 1.;
00123 
00124   if ( massPar>0.0000000001 ) {
00125     if ( mass > massPar) return 0.;
00126   }
00127 
00128   if ( _errorCond ) return 0.;
00129 
00130   // we did all the work in getRandMass
00131   return 1.;
00132 }

double EvtAbsLineShape::getMaxRange (  )  [inline, inherited]

Definition at line 44 of file EvtAbsLineShape.hh.

References EvtAbsLineShape::_maxRange.

Referenced by EvtPartProp::getMaxRange(), and EvtPartProp::newLineShape().

00044 {return _maxRange;}

double EvtRelBreitWignerBarrierFact::getRandMass ( EvtId parId,
int  nDaug,
EvtId dauId,
EvtId othDaugId,
double  maxMass,
double *  dauMasses 
) [virtual]

Reimplemented from EvtAbsLineShape.

Definition at line 134 of file EvtRelBreitWignerBarrierFact.cc.

References EvtAbsLineShape::_addFactorPn, EvtAbsLineShape::_applyFixForSP8, _blatt, _errorCond, EvtAbsLineShape::_includeBirthFact, EvtAbsLineShape::_includeDecayFact, EvtAbsLineShape::_mass, EvtAbsLineShape::_massMax, EvtAbsLineShape::_massMin, EvtAbsLineShape::_spin, EvtAbsLineShape::_userSetBirthOthD, EvtAbsLineShape::_userSetBirthPar, EvtAbsLineShape::_userSetBirthPW, EvtAbsLineShape::_userSetPW, EvtAbsLineShape::_userSetPWD1, EvtAbsLineShape::_userSetPWD2, EvtAbsLineShape::_width, EvtMassAmp::addBirthFact(), EvtMassAmp::addBirthFactFF(), EvtMassAmp::addDeathFact(), EvtMassAmp::addDeathFactFF(), EvtMassAmp::addFactorPn(), EvtAbsLineShape::addFactorPn(), EvtPdf< T >::evaluate(), EvtPDL::getMeanMass(), EvtAbsLineShape::getRandMass(), EvtSpinType::getSpin2(), EvtPDL::getSpinType(), max, EvtTwoBodyVertex::set_f(), EvtMassAmp::setBirthVtx(), and EvtPoint1D::value().

00134                                                                                                                                           {
00135   if ( nDaug!=2) {
00136     if(fabs(_addFactorPn) >0.00000001) EvtAbsLineShape::addFactorPn(_addFactorPn);
00137     double themass = EvtAbsLineShape::getRandMass(parId,nDaug,dauId,othDaugId,maxMass,dauMasses);
00138     return themass;
00139   }
00140 
00141   if ( _width< 0.00001) return _mass;
00142 
00143   //first figure out L - take the lowest allowed.
00144 
00145   EvtSpinType::spintype spinD1=EvtPDL::getSpinType(dauId[0]);
00146   EvtSpinType::spintype spinD2=EvtPDL::getSpinType(dauId[1]);
00147 
00148   int t1=EvtSpinType::getSpin2(spinD1);
00149   int t2=EvtSpinType::getSpin2(spinD2);
00150   int t3=EvtSpinType::getSpin2(_spin);
00151 
00152   int Lmin=-10;
00153 
00154 
00155   // the user has overridden the partial wave to use.
00156   for ( int vC=0; vC<_userSetPW.size(); vC++) {
00157     if ( dauId[0]==_userSetPWD1[vC] &&  dauId[1]==_userSetPWD2[vC] ) Lmin=2*_userSetPW[vC]; 
00158     if ( dauId[0]==_userSetPWD2[vC] &&  dauId[1]==_userSetPWD1[vC] ) Lmin=2*_userSetPW[vC]; 
00159   }
00160   
00161   // allow for special cases.
00162   if (Lmin<-1 ) {
00163     
00164     //There are some things I don't know how to deal with
00165     if ( t3>4) return EvtAbsLineShape::getRandMass(parId,nDaug,dauId,othDaugId,maxMass,dauMasses);
00166     if ( t1>4) return EvtAbsLineShape::getRandMass(parId,nDaug,dauId,othDaugId,maxMass,dauMasses);
00167     if ( t2>4) return EvtAbsLineShape::getRandMass(parId,nDaug,dauId,othDaugId,maxMass,dauMasses);
00168     
00169     //figure the min and max allowwed "spins" for the daughters state
00170     Lmin=std::max(t3-t2-t1,std::max(t2-t3-t1,t1-t3-t2));
00171     if (Lmin<0) Lmin=0;
00172     assert(Lmin==0||Lmin==2||Lmin==4);
00173   }
00174 
00175   //double massD1=EvtPDL::getMeanMass(dauId[0]);
00176   //double massD2=EvtPDL::getMeanMass(dauId[1]);
00177   double massD1=dauMasses[0];
00178   double massD2=dauMasses[1];
00179 
00180   // I'm not sure how to define the vertex factor here - so retreat to nonRel code.
00181   if ( (massD1+massD2)> _mass ) return  EvtAbsLineShape::getRandMass(parId,nDaug,dauId,othDaugId,maxMass,dauMasses);
00182 
00183   //parent vertex factor not yet implemented
00184   double massOthD=-10.;
00185   double massParent=-10.;
00186   int birthl=-10;
00187   if ( othDaugId) {
00188     EvtSpinType::spintype spinOth=EvtPDL::getSpinType(*othDaugId);
00189     EvtSpinType::spintype spinPar=EvtPDL::getSpinType(*parId);
00190     
00191     int tt1=EvtSpinType::getSpin2(spinOth);
00192     int tt2=EvtSpinType::getSpin2(spinPar);
00193     int tt3=EvtSpinType::getSpin2(_spin);
00194     
00195     
00196     //figure the min and max allowwed "spins" for the daughters state
00197     if ( (tt1<=4) && ( tt2<=4) ) {
00198       birthl=std::max(tt3-tt2-tt1,std::max(tt2-tt3-tt1,tt1-tt3-tt2));
00199       if (birthl<0) birthl=0;
00200     
00201       massOthD=EvtPDL::getMeanMass(*othDaugId);
00202       massParent=EvtPDL::getMeanMass(*parId);
00203     
00204     }
00205 
00206 
00207      // allow user to override
00208      for (size_t vC=0; vC<_userSetBirthPW.size(); vC++) {
00209        if ( *othDaugId==_userSetBirthOthD[vC] && *parId==_userSetBirthPar[vC]){
00210         birthl=2*_userSetBirthPW[vC];
00211        } 
00212      }
00213 
00214   }
00215   double massM=_massMax;
00216   if ( (maxMass > -0.5) && (maxMass < massM) ) massM=maxMass;
00217 
00218   //special case... if the parent mass is _fixed_ we can do a little better
00219   //and only for a two body decay as that seems to be where we have problems
00220 
00221   // Define relativistic propagator amplitude
00222 
00223   EvtTwoBodyVertex vd(massD1,massD2,_mass,Lmin/2);
00224   vd.set_f(_blatt);
00225   EvtPropBreitWignerRel bw(_mass,_width);
00226   EvtMassAmp amp(bw,vd);
00227   if ( _includeDecayFact) {
00228     amp.addDeathFact();
00229     amp.addDeathFactFF();
00230   }
00231 
00232   if(fabs(_addFactorPn) >0.00000001){
00233     //    std::cout<<"EvtRelBreitWignerBarrierFact "<< _addFactorPn<<std::endl;
00234     amp.addFactorPn( _addFactorPn);
00235   }
00236   if ( massParent>-1.) {
00237     if ( _includeBirthFact ) {
00238 
00239       EvtTwoBodyVertex vb(_mass,massOthD,massParent,birthl/2);
00240       if ( _applyFixForSP8 ) vb.set_f(_blatt);
00241       amp.setBirthVtx(vb);
00242       amp.addBirthFact();
00243       amp.addBirthFactFF();
00244     }
00245   }
00246 
00247 
00248 
00249   EvtAmpPdf<EvtPoint1D> pdf(amp);
00250 
00251   // Estimate maximum and create predicate for accept reject
00252 
00253 
00254   double tempMaxLoc=_mass;
00255   if ( maxMass>-0.5 && maxMass<_mass) tempMaxLoc=maxMass;
00256   double tempMax=_massMax;
00257   if ( maxMass>-0.5 && maxMass<_massMax) tempMax=maxMass;
00258   double tempMinMass=_massMin;
00259   if ( massD1+massD2 > _massMin) tempMinMass=massD1+massD2;
00260 
00261   //redo sanity check - is there a solution to our problem.
00262   //if not return an error condition that is caught by the
00263   //mass prob calculation above.
00264   if ( tempMinMass > tempMax ) {
00265     _errorCond=true;
00266     return tempMinMass;
00267   }
00268   
00269   if ( tempMaxLoc < tempMinMass) tempMaxLoc=tempMinMass;
00270 
00271   double safetyFactor=1.2;
00272   if ( _applyFixForSP8 ) safetyFactor=1.4;
00273 
00274   EvtPdfMax<EvtPoint1D> max(safetyFactor*pdf.evaluate(EvtPoint1D(tempMinMass,tempMax,tempMaxLoc)));
00275 
00276   EvtPdfPred<EvtPoint1D> pred(pdf);
00277   pred.setMax(max);
00278 
00279   EvtIntervalFlatPdf flat(tempMinMass,tempMax);
00280   EvtPdfGen<EvtPoint1D> gen(flat);
00281   EvtPredGen<EvtPdfGen<EvtPoint1D>,EvtPdfPred<EvtPoint1D> > predgen(gen,pred);
00282 
00283   EvtPoint1D point = predgen();
00284   return point.value();
00285 
00286 }

EvtSpinType::spintype EvtAbsLineShape::getSpinType (  )  [inline, inherited]

Definition at line 46 of file EvtAbsLineShape.hh.

References EvtAbsLineShape::_spin.

Referenced by EvtPartProp::newLineShape().

00046 {return _spin;}

double EvtAbsLineShape::getWidth (  )  [inline, inherited]

Definition at line 45 of file EvtAbsLineShape.hh.

References EvtAbsLineShape::_width.

Referenced by EvtPartProp::getWidth(), and EvtPartProp::newLineShape().

00045 {return _width;} 

void EvtAbsLineShape::includeBirthFactor ( bool  yesno  )  [inline, inherited]

Definition at line 70 of file EvtAbsLineShape.hh.

References EvtAbsLineShape::_includeBirthFact.

Referenced by EvtPartProp::includeBirthFactor().

00070 { _includeBirthFact = yesno; }

void EvtAbsLineShape::includeDecayFactor ( bool  yesno  )  [inline, inherited]

Definition at line 72 of file EvtAbsLineShape.hh.

References EvtAbsLineShape::_includeDecayFact.

Referenced by EvtPartProp::includeDecayFactor().

00072 { _includeDecayFact = yesno; }

EvtRelBreitWignerBarrierFact & EvtRelBreitWignerBarrierFact::operator= ( const EvtRelBreitWignerBarrierFact x  ) 

Definition at line 84 of file EvtRelBreitWignerBarrierFact.cc.

References _blatt, _errorCond, EvtAbsLineShape::_includeBirthFact, EvtAbsLineShape::_includeDecayFact, EvtAbsLineShape::_mass, EvtAbsLineShape::_massMax, EvtAbsLineShape::_massMin, EvtAbsLineShape::_maxRange, EvtAbsLineShape::_spin, EvtAbsLineShape::_width, and x.

00084                                                                                                            {
00085   _mass=x._mass;
00086   _width=x._width;
00087   _spin=x._spin;
00088   _massMax=x._massMax;
00089   _massMin=x._massMin;
00090   _blatt=x._blatt;
00091   _maxRange=x._maxRange;
00092   _includeDecayFact=x._includeDecayFact;
00093   _includeBirthFact=x._includeBirthFact;
00094   _errorCond=x._errorCond;
00095 
00096   return *this;
00097 }

virtual void EvtRelBreitWignerBarrierFact::reSetBlatt ( double  blatt  )  [inline, virtual]

Reimplemented from EvtAbsLineShape.

Definition at line 46 of file EvtRelBreitWignerBarrierFact.hh.

References _blatt.

00046 { _blatt = blatt; }

void EvtAbsLineShape::reSetMass ( double  mass  )  [inline, inherited]

Definition at line 50 of file EvtAbsLineShape.hh.

References EvtAbsLineShape::_mass.

Referenced by EvtPartProp::reSetMass().

00050 { _mass=mass;}

void EvtAbsLineShape::reSetMassMax ( double  mass  )  [inline, inherited]

Definition at line 68 of file EvtAbsLineShape.hh.

References EvtAbsLineShape::_massMax.

Referenced by EvtPartProp::reSetMassMax().

00068 { _massMax=mass;}

void EvtAbsLineShape::reSetMassMin ( double  mass  )  [inline, inherited]

Definition at line 67 of file EvtAbsLineShape.hh.

References EvtAbsLineShape::_massMin.

Referenced by EvtPartProp::reSetMassMin().

00067 { _massMin=mass;}

void EvtAbsLineShape::reSetWidth ( double  width  )  [inline, inherited]

Definition at line 51 of file EvtAbsLineShape.hh.

References EvtAbsLineShape::_mass, EvtAbsLineShape::_massMax, EvtAbsLineShape::_massMin, EvtAbsLineShape::_maxRange, and EvtAbsLineShape::_width.

Referenced by EvtPartProp::reSetWidth().

00051                                 { 
00052           _width=width;
00053           // <--- added by L.L. Wang to fix a bug
00054           double maxdelta=15.0*width;
00055           if ( _maxRange > 0.00001 ) {
00056                   _massMax=_mass+maxdelta;
00057                   _massMin=_mass-_maxRange;
00058           }
00059           else{
00060                   _massMax=_mass+maxdelta;
00061                   _massMin=_mass-15.0*width;
00062           }
00063           if ( _massMin< 0. ) _massMin=0.;
00064           _massMax=_mass+maxdelta;
00065           // ---> //
00066   }

double EvtAbsLineShape::rollMass (  )  [virtual, inherited]

Definition at line 102 of file EvtAbsLineShape.cc.

References EvtAbsLineShape::_mass, EvtAbsLineShape::_massMax, EvtAbsLineShape::_massMin, EvtAbsLineShape::_width, EvtRandom::Flat(), tan(), and subSeperate::temp.

Referenced by EvtPartProp::rollMass().

00102                                  {
00103 
00104   double ymin, ymax;
00105   double temp;
00106 
00107   if ( _width < 0.0001 ) {
00108     return _mass;
00109   }
00110   else{
00111     ymin = atan( 2.0*(_massMin-_mass)/_width);
00112     ymax = atan( 2.0*(_massMax-_mass)/_width);
00113 
00114     temp= ( _mass + ((_width/2.0)*tan(EvtRandom::Flat(ymin,ymax))));
00115 
00116     return temp;
00117   }
00118 }

void EvtAbsLineShape::setPWForBirthL ( int  spin,
EvtId  par,
EvtId  othD 
) [inline, inherited]

Definition at line 78 of file EvtAbsLineShape.hh.

References EvtAbsLineShape::_userSetBirthOthD, EvtAbsLineShape::_userSetBirthPar, and EvtAbsLineShape::_userSetBirthPW.

Referenced by EvtPartProp::setPWForBirthL().

00078                                                         { 
00079           _userSetBirthPW.push_back(spin);
00080           _userSetBirthOthD.push_back(othD);
00081           _userSetBirthPar.push_back(par);
00082   }

void EvtAbsLineShape::setPWForDecay ( int  spin,
EvtId  d1,
EvtId  d2 
) [inline, inherited]

Definition at line 73 of file EvtAbsLineShape.hh.

References EvtAbsLineShape::_userSetPW, EvtAbsLineShape::_userSetPWD1, and EvtAbsLineShape::_userSetPWD2.

Referenced by EvtPartProp::setPWForDecay().

00073                                                     { 
00074           _userSetPW.push_back(spin);
00075           _userSetPWD1.push_back(d1);
00076           _userSetPWD2.push_back(d2);
00077   }


Member Data Documentation

double EvtAbsLineShape::_addFactorPn [protected, inherited]

Definition at line 92 of file EvtAbsLineShape.hh.

Referenced by EvtAbsLineShape::addFactorPn(), getRandMass(), and EvtAbsLineShape::getRandMass().

bool EvtAbsLineShape::_applyFixForSP8 [protected, inherited]

Definition at line 112 of file EvtAbsLineShape.hh.

Referenced by EvtAbsLineShape::EvtAbsLineShape(), EvtAbsLineShape::fixForSP8(), getRandMass(), and EvtAbsLineShape::operator=().

double EvtRelBreitWignerBarrierFact::_blatt [protected]

Definition at line 50 of file EvtRelBreitWignerBarrierFact.hh.

Referenced by EvtRelBreitWignerBarrierFact(), getRandMass(), operator=(), and reSetBlatt().

bool EvtRelBreitWignerBarrierFact::_errorCond [protected]

Definition at line 51 of file EvtRelBreitWignerBarrierFact.hh.

Referenced by EvtRelBreitWignerBarrierFact(), getMassProb(), getRandMass(), and operator=().

bool EvtAbsLineShape::_includeBirthFact [protected, inherited]

Definition at line 91 of file EvtAbsLineShape.hh.

Referenced by EvtAbsLineShape::EvtAbsLineShape(), EvtRelBreitWignerBarrierFact(), getRandMass(), EvtAbsLineShape::includeBirthFactor(), operator=(), and EvtAbsLineShape::operator=().

bool EvtAbsLineShape::_includeDecayFact [protected, inherited]

Definition at line 90 of file EvtAbsLineShape.hh.

Referenced by EvtAbsLineShape::EvtAbsLineShape(), EvtRelBreitWignerBarrierFact(), getRandMass(), EvtAbsLineShape::includeDecayFactor(), operator=(), and EvtAbsLineShape::operator=().

double EvtAbsLineShape::_mass [protected, inherited]

Definition at line 93 of file EvtAbsLineShape.hh.

Referenced by EvtAbsLineShape::EvtAbsLineShape(), EvtFlatLineShape::EvtFlatLineShape(), EvtManyDeltaFuncLineShape::EvtManyDeltaFuncLineShape(), EvtRelBreitWignerBarrierFact(), EvtAbsLineShape::getMass(), getRandMass(), EvtAbsLineShape::getRandMass(), operator=(), EvtManyDeltaFuncLineShape::operator=(), EvtFlatLineShape::operator=(), EvtAbsLineShape::operator=(), EvtAbsLineShape::reSetMass(), EvtAbsLineShape::reSetWidth(), and EvtAbsLineShape::rollMass().

double EvtAbsLineShape::_massMax [protected, inherited]

Definition at line 95 of file EvtAbsLineShape.hh.

Referenced by EvtAbsLineShape::EvtAbsLineShape(), EvtFlatLineShape::EvtFlatLineShape(), EvtManyDeltaFuncLineShape::EvtManyDeltaFuncLineShape(), EvtRelBreitWignerBarrierFact(), EvtAbsLineShape::getMassMax(), getRandMass(), EvtManyDeltaFuncLineShape::getRandMass(), EvtFlatLineShape::getRandMass(), EvtAbsLineShape::getRandMass(), operator=(), EvtManyDeltaFuncLineShape::operator=(), EvtFlatLineShape::operator=(), EvtAbsLineShape::operator=(), EvtAbsLineShape::reSetMassMax(), EvtAbsLineShape::reSetWidth(), and EvtAbsLineShape::rollMass().

double EvtAbsLineShape::_massMin [protected, inherited]

Definition at line 94 of file EvtAbsLineShape.hh.

Referenced by EvtAbsLineShape::EvtAbsLineShape(), EvtFlatLineShape::EvtFlatLineShape(), EvtManyDeltaFuncLineShape::EvtManyDeltaFuncLineShape(), EvtRelBreitWignerBarrierFact(), EvtAbsLineShape::getMassMin(), getRandMass(), EvtManyDeltaFuncLineShape::getRandMass(), EvtFlatLineShape::getRandMass(), EvtAbsLineShape::getRandMass(), operator=(), EvtManyDeltaFuncLineShape::operator=(), EvtFlatLineShape::operator=(), EvtAbsLineShape::operator=(), EvtAbsLineShape::reSetMassMin(), EvtAbsLineShape::reSetWidth(), and EvtAbsLineShape::rollMass().

double EvtAbsLineShape::_maxRange [protected, inherited]

Definition at line 97 of file EvtAbsLineShape.hh.

Referenced by EvtAbsLineShape::EvtAbsLineShape(), EvtFlatLineShape::EvtFlatLineShape(), EvtManyDeltaFuncLineShape::EvtManyDeltaFuncLineShape(), EvtRelBreitWignerBarrierFact(), EvtAbsLineShape::getMaxRange(), operator=(), EvtManyDeltaFuncLineShape::operator=(), EvtFlatLineShape::operator=(), EvtAbsLineShape::operator=(), and EvtAbsLineShape::reSetWidth().

EvtSpinType::spintype EvtAbsLineShape::_spin [protected, inherited]

Definition at line 110 of file EvtAbsLineShape.hh.

Referenced by EvtAbsLineShape::EvtAbsLineShape(), EvtFlatLineShape::EvtFlatLineShape(), EvtManyDeltaFuncLineShape::EvtManyDeltaFuncLineShape(), EvtRelBreitWignerBarrierFact(), getRandMass(), EvtAbsLineShape::getSpinType(), operator=(), EvtManyDeltaFuncLineShape::operator=(), EvtFlatLineShape::operator=(), and EvtAbsLineShape::operator=().

std::vector<EvtId> EvtAbsLineShape::_userSetBirthOthD [protected, inherited]

Definition at line 107 of file EvtAbsLineShape.hh.

Referenced by getRandMass(), and EvtAbsLineShape::setPWForBirthL().

std::vector<EvtId> EvtAbsLineShape::_userSetBirthPar [protected, inherited]

Definition at line 107 of file EvtAbsLineShape.hh.

Referenced by getRandMass(), and EvtAbsLineShape::setPWForBirthL().

std::vector<int> EvtAbsLineShape::_userSetBirthPW [protected, inherited]

Definition at line 108 of file EvtAbsLineShape.hh.

Referenced by getRandMass(), and EvtAbsLineShape::setPWForBirthL().

std::vector<int> EvtAbsLineShape::_userSetPW [protected, inherited]

Definition at line 104 of file EvtAbsLineShape.hh.

Referenced by getRandMass(), and EvtAbsLineShape::setPWForDecay().

std::vector<EvtId> EvtAbsLineShape::_userSetPWD1 [protected, inherited]

Definition at line 103 of file EvtAbsLineShape.hh.

Referenced by getRandMass(), and EvtAbsLineShape::setPWForDecay().

std::vector<EvtId> EvtAbsLineShape::_userSetPWD2 [protected, inherited]

Definition at line 103 of file EvtAbsLineShape.hh.

Referenced by getRandMass(), and EvtAbsLineShape::setPWForDecay().

double EvtAbsLineShape::_width [protected, inherited]

Definition at line 96 of file EvtAbsLineShape.hh.

Referenced by EvtAbsLineShape::EvtAbsLineShape(), EvtFlatLineShape::EvtFlatLineShape(), EvtManyDeltaFuncLineShape::EvtManyDeltaFuncLineShape(), EvtRelBreitWignerBarrierFact(), getMassProb(), EvtAbsLineShape::getMassProb(), getRandMass(), EvtManyDeltaFuncLineShape::getRandMass(), EvtAbsLineShape::getRandMass(), EvtAbsLineShape::getWidth(), operator=(), EvtManyDeltaFuncLineShape::operator=(), EvtFlatLineShape::operator=(), EvtAbsLineShape::operator=(), EvtAbsLineShape::reSetWidth(), and EvtAbsLineShape::rollMass().


Generated on Tue Nov 29 23:19:21 2016 for BOSS_7.0.2 by  doxygen 1.4.7