ROOTPWA
Classes | Public Types | Public Member Functions | Static Private Member Functions | Private Attributes | Friends | List of all members
rpwa::nBodyPhaseSpaceGen Class Reference

#include <nBodyPhaseSpaceGen.h>

Classes

class  rndGen

Public Types

enum  kinematicsTypeEnum { BLOCK = 1, RAUBOLD_LYNCH = 2 }
enum  weightTypeEnum {
  S_U_CHUNG = 1, NUPHAZ = 2, GENBOD = 3, FLAT = 4,
  IMPORTANCE = 5
}

Public Member Functions

 nBodyPhaseSpaceGen ()
virtual ~nBodyPhaseSpaceGen ()
bool setDecay (const std::vector< double > &daughterMasses)
 sets decay constants and prepares internal variables
bool setDecay (const unsigned int nmbOfDaughters, const double *daughterMasses)
void setSeed (const unsigned int seed)
 sets seed of random number generator
unsigned int seed ()
 returns seed of random number generator
double random ()
 returns number from internal random generator; intended for higher-level generators
void setProposalBW (double mass, double width)
 set proposal function for importance sampling in the (n-1) isobar mass as simple Breit-Wigner-Shape (to be used with IMPORTANCE weighting option)
double generateDecay (const TLorentzVector &nBody)
 generates full event with certain n-body mass and momentum and returns event weight
bool generateDecayAccepted (const TLorentzVector &nBody, const double maxWeight=0)
 generates full event with certain n-body mass and momentum only when event is accepted (return value = true) this function is more efficient, if only weighted events are needed
void pickMasses (const double nBodyMass)
 randomly choses the (n - 2) effective masses of the respective (i + 1)-body systems
double calcWeight ()
 computes event weight and breakup momenta operates on vector of intermediate two-body masses prepared by pickMasses()
void pickAngles ()
 randomly choses the (n - 1) polar and (n - 1) azimuthal angles in the respective (i + 1)-body RFs
void calcEventKinematics (const TLorentzVector &nBody)
 calculates full event kinematics from the effective masses of the (i + 1)-body systems and the Lorentz vector of the decaying system uses the break-up momenta calculated by calcWeight() and angles from pickAngles()
void setKinematicsType (const kinematicsTypeEnum kinematicsType)
 selects algorithm used to calculate event kinematics
kinematicsTypeEnum kinematicsType () const
 returns algorithm used to calculate event kinematics
void setVerbose (bool flag)
void setWeightType (const weightTypeEnum weightType)
 selects formula used for weight calculation
weightTypeEnum weightType () const
 returns formula used for weight calculation
void setMaxWeight (const double maxWeight)
 sets maximum weight used for hit-miss MC
double maxWeight () const
 returns maximum weight used for hit-miss MC
double normalization () const
 returns normalization used in weight calculation
double eventWeight () const
 returns weight of generated event
double maxWeightObserved () const
 returns maximum observed weight since instantiation
void resetMaxWeightObserved ()
 sets maximum observed weight back to zero
double estimateMaxWeight (const double nBodyMass, const unsigned int nmbOfIterations=10000)
 estimates maximum weight for given n-body mass
bool eventAccepted (const double maxWeight=0)
 applies event weight in form of hit-miss MC assumes that event weight has been already calculated by calcWeight() if maxWeight > 0 value is used as maximum weight, otherwise _maxWeight value is used
const TLorentzVector & daughter (const int index) const
 returns Lorentz vector of daughter at index
const std::vector
< TLorentzVector > & 
daughters () const
 returns Lorentz vectors of all daughters
unsigned int nmbOfDaughters () const
 returns number of daughters
double daughterMass (const int index) const
 returns invariant mass of daughter at index
double intermediateMass (const int index) const
 returns intermediate mass of (index + 1)-body system
double breakupMom (const int index) const
 returns breakup momentum in (index + 1)-body RF
double cosTheta (const int index) const
 returns polar angle in (index + 1)-body RF
double phi (const int index) const
 returns azimuth in (index + 1)-body RF
double impWeight () const
std::ostream & print (std::ostream &out=std::cout) const
 prints generator status

Static Private Member Functions

static double F (const double x, const double y)
 reduced 2-body phase-space factor = 2 * q / M needed for NUPHAS weight calculation

Private Attributes

std::vector< double > _m
 masses of daughter particles
unsigned int _n
 number of daughter particles
std::vector< double > _M
 effective masses of (i + 1)-body systems
std::vector< double > _cosTheta
 cosine of polar angle of the 2-body decay of the (i + 1)-body system
std::vector< double > _phi
 azimuthal angle of the 2-body decay of the (i + 1)-body system
std::vector< double > _mSum
 sums of daughter particle masses
std::vector< double > _breakupMom
 breakup momenta for the two-body decays: (i + 1)-body –> daughter_(i + 1) + i-body
std::vector< TLorentzVector > _daughters
 Lorentz vectors of the daughter particles.
weightTypeEnum _weightType
 switches between different weight formulas
double _norm
 normalization value
double _weight
 phase space weight of generated event
double _impweight
 importance sampling weight
double _maxWeightObserved
 maximum event weight calculated processing the input data
double _maxWeight
 maximum weight used to weight events in hit-miss MC
kinematicsTypeEnum _kinematicsType
 switches between different ways of calculating event kinematics
bool _verbose
double _isoBWMass
 Breit-Wigner mass for importance sampling proposal.
double _isoBWWidth
 Breit-Wigner width for importance sampling proposal.
rndGen _rnd
 random number generator instance

Friends

std::ostream & operator<< (std::ostream &out, const nBodyPhaseSpaceGen &gen)

Detailed Description

Definition at line 108 of file nBodyPhaseSpaceGen.h.

Member Enumeration Documentation

Enumerator:
BLOCK 
RAUBOLD_LYNCH 

Definition at line 154 of file nBodyPhaseSpaceGen.h.

Enumerator:
S_U_CHUNG 
NUPHAZ 
GENBOD 
FLAT 
IMPORTANCE 

Definition at line 164 of file nBodyPhaseSpaceGen.h.

Constructor & Destructor Documentation

nBodyPhaseSpaceGen::nBodyPhaseSpaceGen ( )

Definition at line 103 of file nBodyPhaseSpaceGen.cc.

nBodyPhaseSpaceGen::~nBodyPhaseSpaceGen ( )
virtual

Definition at line 120 of file nBodyPhaseSpaceGen.cc.

Member Function Documentation

double rpwa::nBodyPhaseSpaceGen::breakupMom ( const int  index) const
inline

returns breakup momentum in (index + 1)-body RF

Definition at line 197 of file nBodyPhaseSpaceGen.h.

References _breakupMom.

void nBodyPhaseSpaceGen::calcEventKinematics ( const TLorentzVector &  nBody)

calculates full event kinematics from the effective masses of the (i + 1)-body systems and the Lorentz vector of the decaying system uses the break-up momenta calculated by calcWeight() and angles from pickAngles()

Definition at line 414 of file nBodyPhaseSpaceGen.cc.

References _breakupMom, _cosTheta, _daughters, _kinematicsType, _m, _M, _n, _phi, BLOCK, daughter(), i, and RAUBOLD_LYNCH.

Referenced by rpwa::diffractivePhaseSpace::event(), generateDecay(), and generateDecayAccepted().

double nBodyPhaseSpaceGen::calcWeight ( )

computes event weight and breakup momenta operates on vector of intermediate two-body masses prepared by pickMasses()

!! warning: produces distorted angular distribution

Definition at line 357 of file nBodyPhaseSpaceGen.cc.

References _breakupMom, _impweight, _m, _M, _maxWeightObserved, _mSum, _n, _norm, _weight, _weightType, F(), FLAT, GENBOD, i, IMPORTANCE, M2, NUPHAZ, and S_U_CHUNG.

Referenced by estimateMaxWeight(), rpwa::diffractivePhaseSpace::event(), generateDecay(), and generateDecayAccepted().

double rpwa::nBodyPhaseSpaceGen::cosTheta ( const int  index) const
inline

returns polar angle in (index + 1)-body RF

Definition at line 198 of file nBodyPhaseSpaceGen.h.

References _cosTheta.

const TLorentzVector& rpwa::nBodyPhaseSpaceGen::daughter ( const int  index) const
inline
double rpwa::nBodyPhaseSpaceGen::daughterMass ( const int  index) const
inline

returns invariant mass of daughter at index

Definition at line 195 of file nBodyPhaseSpaceGen.h.

References _m.

const std::vector<TLorentzVector>& rpwa::nBodyPhaseSpaceGen::daughters ( ) const
inline

returns Lorentz vectors of all daughters

Definition at line 193 of file nBodyPhaseSpaceGen.h.

References _daughters.

double nBodyPhaseSpaceGen::estimateMaxWeight ( const double  nBodyMass,
const unsigned int  nmbOfIterations = 10000 
)

estimates maximum weight for given n-body mass

Definition at line 497 of file nBodyPhaseSpaceGen.cc.

References _weight, calcWeight(), i, maxWeight(), and pickMasses().

Referenced by rpwa::diffractivePhaseSpace::BuildDaughterList(), compareAmplitudes(), and testNBodyPhaseSpaceGen().

bool rpwa::nBodyPhaseSpaceGen::eventAccepted ( const double  maxWeight = 0)
inline

applies event weight in form of hit-miss MC assumes that event weight has been already calculated by calcWeight() if maxWeight > 0 value is used as maximum weight, otherwise _maxWeight value is used

Definition at line 285 of file nBodyPhaseSpaceGen.h.

Referenced by generateDecayAccepted().

double rpwa::nBodyPhaseSpaceGen::eventWeight ( ) const
inline

returns weight of generated event

Definition at line 177 of file nBodyPhaseSpaceGen.h.

References _weight.

double rpwa::nBodyPhaseSpaceGen::F ( const double  x,
const double  y 
)
inlinestaticprivate

reduced 2-body phase-space factor = 2 * q / M needed for NUPHAS weight calculation

Definition at line 302 of file nBodyPhaseSpaceGen.h.

Referenced by calcWeight(), and pickMasses().

double nBodyPhaseSpaceGen::generateDecay ( const TLorentzVector &  nBody)

generates full event with certain n-body mass and momentum and returns event weight

Definition at line 198 of file nBodyPhaseSpaceGen.cc.

References _mSum, _n, _weight, calcEventKinematics(), calcWeight(), pickAngles(), pickMasses(), and print().

Referenced by testNBodyPhaseSpaceGen().

bool nBodyPhaseSpaceGen::generateDecayAccepted ( const TLorentzVector &  nBody,
const double  maxWeight = 0 
)

generates full event with certain n-body mass and momentum only when event is accepted (return value = true) this function is more efficient, if only weighted events are needed

Definition at line 225 of file nBodyPhaseSpaceGen.cc.

References _mSum, _n, calcEventKinematics(), calcWeight(), eventAccepted(), pickAngles(), and pickMasses().

Referenced by compareAmplitudes(), and testNBodyPhaseSpaceGen().

double rpwa::nBodyPhaseSpaceGen::impWeight ( ) const
inline

Definition at line 203 of file nBodyPhaseSpaceGen.h.

References _impweight.

Referenced by rpwa::diffractivePhaseSpace::impWeight().

double rpwa::nBodyPhaseSpaceGen::intermediateMass ( const int  index) const
inline

returns intermediate mass of (index + 1)-body system

Definition at line 196 of file nBodyPhaseSpaceGen.h.

References _M.

kinematicsTypeEnum rpwa::nBodyPhaseSpaceGen::kinematicsType ( ) const
inline

returns algorithm used to calculate event kinematics

Definition at line 157 of file nBodyPhaseSpaceGen.h.

References _kinematicsType.

Referenced by setKinematicsType().

double rpwa::nBodyPhaseSpaceGen::maxWeight ( ) const
inline

returns maximum weight used for hit-miss MC

Definition at line 175 of file nBodyPhaseSpaceGen.h.

References _maxWeight.

Referenced by rpwa::diffractivePhaseSpace::BuildDaughterList(), estimateMaxWeight(), rpwa::diffractivePhaseSpace::event(), and setMaxWeight().

double rpwa::nBodyPhaseSpaceGen::maxWeightObserved ( ) const
inline

returns maximum observed weight since instantiation

Definition at line 178 of file nBodyPhaseSpaceGen.h.

References _maxWeightObserved.

unsigned int rpwa::nBodyPhaseSpaceGen::nmbOfDaughters ( ) const
inline

returns number of daughters

Definition at line 194 of file nBodyPhaseSpaceGen.h.

References _n.

Referenced by setDecay().

double rpwa::nBodyPhaseSpaceGen::normalization ( ) const
inline

returns normalization used in weight calculation

Definition at line 176 of file nBodyPhaseSpaceGen.h.

References _norm.

double rpwa::nBodyPhaseSpaceGen::phi ( const int  index) const
inline

returns azimuth in (index + 1)-body RF

Definition at line 199 of file nBodyPhaseSpaceGen.h.

References _phi.

void rpwa::nBodyPhaseSpaceGen::pickAngles ( )
inline

randomly choses the (n - 1) polar and (n - 1) azimuthal angles in the respective (i + 1)-body RFs

Definition at line 274 of file nBodyPhaseSpaceGen.h.

References _cosTheta, _n, _phi, i, and random().

Referenced by rpwa::diffractivePhaseSpace::event(), generateDecay(), and generateDecayAccepted().

void nBodyPhaseSpaceGen::pickMasses ( const double  nBodyMass)

randomly choses the (n - 2) effective masses of the respective (i + 1)-body systems

Definition at line 249 of file nBodyPhaseSpaceGen.cc.

References _impweight, _isoBWMass, _isoBWWidth, _m, _M, _mSum, _n, _weight, _weightType, F(), i, IMPORTANCE, NUPHAZ, and random().

Referenced by estimateMaxWeight(), rpwa::diffractivePhaseSpace::event(), generateDecay(), and generateDecayAccepted().

ostream & nBodyPhaseSpaceGen::print ( std::ostream &  out = std::cout) const

prints generator status

Definition at line 513 of file nBodyPhaseSpaceGen.cc.

References _daughters, _kinematicsType, _maxWeight, _maxWeightObserved, _n, _norm, _weight, _weightType, and i.

Referenced by generateDecay().

double rpwa::nBodyPhaseSpaceGen::random ( )
inline

returns number from internal random generator; intended for higher-level generators

Definition at line 124 of file nBodyPhaseSpaceGen.h.

References _rnd, and rpwa::nBodyPhaseSpaceGen::rndGen::pick().

Referenced by rpwa::diffractivePhaseSpace::event(), pickAngles(), and pickMasses().

void rpwa::nBodyPhaseSpaceGen::resetMaxWeightObserved ( )
inline

sets maximum observed weight back to zero

Definition at line 179 of file nBodyPhaseSpaceGen.h.

References _maxWeightObserved.

Referenced by setDecay().

unsigned int rpwa::nBodyPhaseSpaceGen::seed ( )
inline

returns seed of random number generator

Definition at line 123 of file nBodyPhaseSpaceGen.h.

References _rnd, and rpwa::nBodyPhaseSpaceGen::rndGen::seed().

bool nBodyPhaseSpaceGen::setDecay ( const std::vector< double > &  daughterMasses)
bool nBodyPhaseSpaceGen::setDecay ( const unsigned int  nmbOfDaughters,
const double *  daughterMasses 
)

Definition at line 184 of file nBodyPhaseSpaceGen.cc.

References i, nmbOfDaughters(), and setDecay().

void rpwa::nBodyPhaseSpaceGen::setKinematicsType ( const kinematicsTypeEnum  kinematicsType)
inline

selects algorithm used to calculate event kinematics

Definition at line 156 of file nBodyPhaseSpaceGen.h.

References _kinematicsType, and kinematicsType().

Referenced by compareAmplitudes(), rpwa::diffractivePhaseSpace::diffractivePhaseSpace(), and testNBodyPhaseSpaceGen().

void rpwa::nBodyPhaseSpaceGen::setMaxWeight ( const double  maxWeight)
inline

sets maximum weight used for hit-miss MC

Definition at line 174 of file nBodyPhaseSpaceGen.h.

References _maxWeight, and maxWeight().

Referenced by rpwa::diffractivePhaseSpace::BuildDaughterList(), compareAmplitudes(), and testNBodyPhaseSpaceGen().

void rpwa::nBodyPhaseSpaceGen::setProposalBW ( double  mass,
double  width 
)
inline

set proposal function for importance sampling in the (n-1) isobar mass as simple Breit-Wigner-Shape (to be used with IMPORTANCE weighting option)

Definition at line 127 of file nBodyPhaseSpaceGen.h.

References _isoBWMass, _isoBWWidth, and mass.

Referenced by rpwa::diffractivePhaseSpace::SetImportanceBW().

void rpwa::nBodyPhaseSpaceGen::setSeed ( const unsigned int  seed)
inline

sets seed of random number generator

Definition at line 122 of file nBodyPhaseSpaceGen.h.

References _rnd, and rpwa::nBodyPhaseSpaceGen::rndGen::setSeed().

Referenced by compareAmplitudes(), rpwa::diffractivePhaseSpace::SetSeed(), and testNBodyPhaseSpaceGen().

void rpwa::nBodyPhaseSpaceGen::setVerbose ( bool  flag)
inline

Definition at line 160 of file nBodyPhaseSpaceGen.h.

References _verbose.

Referenced by rpwa::diffractivePhaseSpace::setVerbose().

void rpwa::nBodyPhaseSpaceGen::setWeightType ( const weightTypeEnum  weightType)
inline
weightTypeEnum rpwa::nBodyPhaseSpaceGen::weightType ( ) const
inline

returns formula used for weight calculation

Definition at line 171 of file nBodyPhaseSpaceGen.h.

References _weightType.

Referenced by setWeightType().

Friends And Related Function Documentation

std::ostream& operator<< ( std::ostream &  out,
const nBodyPhaseSpaceGen gen 
)
friend

Definition at line 205 of file nBodyPhaseSpaceGen.h.

Member Data Documentation

std::vector<double> rpwa::nBodyPhaseSpaceGen::_breakupMom
private

breakup momenta for the two-body decays: (i + 1)-body –> daughter_(i + 1) + i-body

Definition at line 225 of file nBodyPhaseSpaceGen.h.

Referenced by breakupMom(), calcEventKinematics(), calcWeight(), and setDecay().

std::vector<double> rpwa::nBodyPhaseSpaceGen::_cosTheta
private

cosine of polar angle of the 2-body decay of the (i + 1)-body system

Definition at line 222 of file nBodyPhaseSpaceGen.h.

Referenced by calcEventKinematics(), cosTheta(), pickAngles(), and setDecay().

std::vector<TLorentzVector> rpwa::nBodyPhaseSpaceGen::_daughters
private

Lorentz vectors of the daughter particles.

Definition at line 226 of file nBodyPhaseSpaceGen.h.

Referenced by calcEventKinematics(), daughter(), daughters(), print(), and setDecay().

double rpwa::nBodyPhaseSpaceGen::_impweight
private

importance sampling weight

Definition at line 230 of file nBodyPhaseSpaceGen.h.

Referenced by calcWeight(), impWeight(), and pickMasses().

double rpwa::nBodyPhaseSpaceGen::_isoBWMass
private

Breit-Wigner mass for importance sampling proposal.

Definition at line 239 of file nBodyPhaseSpaceGen.h.

Referenced by pickMasses(), and setProposalBW().

double rpwa::nBodyPhaseSpaceGen::_isoBWWidth
private

Breit-Wigner width for importance sampling proposal.

Definition at line 240 of file nBodyPhaseSpaceGen.h.

Referenced by pickMasses(), and setProposalBW().

kinematicsTypeEnum rpwa::nBodyPhaseSpaceGen::_kinematicsType
private

switches between different ways of calculating event kinematics

Definition at line 233 of file nBodyPhaseSpaceGen.h.

Referenced by calcEventKinematics(), kinematicsType(), print(), and setKinematicsType().

std::vector<double> rpwa::nBodyPhaseSpaceGen::_m
private

masses of daughter particles

Definition at line 216 of file nBodyPhaseSpaceGen.h.

Referenced by calcEventKinematics(), calcWeight(), daughterMass(), pickMasses(), and setDecay().

std::vector<double> rpwa::nBodyPhaseSpaceGen::_M
private

effective masses of (i + 1)-body systems

Definition at line 221 of file nBodyPhaseSpaceGen.h.

Referenced by calcEventKinematics(), calcWeight(), intermediateMass(), pickMasses(), and setDecay().

double rpwa::nBodyPhaseSpaceGen::_maxWeight
private

maximum weight used to weight events in hit-miss MC

Definition at line 232 of file nBodyPhaseSpaceGen.h.

Referenced by maxWeight(), print(), and setMaxWeight().

double rpwa::nBodyPhaseSpaceGen::_maxWeightObserved
private

maximum event weight calculated processing the input data

Definition at line 231 of file nBodyPhaseSpaceGen.h.

Referenced by calcWeight(), maxWeightObserved(), print(), and resetMaxWeightObserved().

std::vector<double> rpwa::nBodyPhaseSpaceGen::_mSum
private

sums of daughter particle masses

Definition at line 224 of file nBodyPhaseSpaceGen.h.

Referenced by calcWeight(), generateDecay(), generateDecayAccepted(), pickMasses(), and setDecay().

unsigned int rpwa::nBodyPhaseSpaceGen::_n
private
double rpwa::nBodyPhaseSpaceGen::_norm
private

normalization value

Definition at line 228 of file nBodyPhaseSpaceGen.h.

Referenced by calcWeight(), normalization(), print(), and setDecay().

std::vector<double> rpwa::nBodyPhaseSpaceGen::_phi
private

azimuthal angle of the 2-body decay of the (i + 1)-body system

Definition at line 223 of file nBodyPhaseSpaceGen.h.

Referenced by calcEventKinematics(), phi(), pickAngles(), and setDecay().

rndGen rpwa::nBodyPhaseSpaceGen::_rnd
private

random number generator instance

Definition at line 263 of file nBodyPhaseSpaceGen.h.

Referenced by random(), seed(), and setSeed().

bool rpwa::nBodyPhaseSpaceGen::_verbose
private

Definition at line 237 of file nBodyPhaseSpaceGen.h.

Referenced by setVerbose().

double rpwa::nBodyPhaseSpaceGen::_weight
private

phase space weight of generated event

Definition at line 229 of file nBodyPhaseSpaceGen.h.

Referenced by calcWeight(), estimateMaxWeight(), eventWeight(), generateDecay(), pickMasses(), and print().

weightTypeEnum rpwa::nBodyPhaseSpaceGen::_weightType
private

switches between different weight formulas

Definition at line 227 of file nBodyPhaseSpaceGen.h.

Referenced by calcWeight(), pickMasses(), print(), setDecay(), setWeightType(), and weightType().


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