Encodes a multi-body decay tree structure, is largely limited to describe a sequence of quasi two-body processes (i.e. the isobar model) by the implemented lineshapes (propagators) and spin structure (vertices) Decay chains are usually constucted by parsing strings, for example:
D0{rho(770)0{pi+,pi-},K0S0}
is the tree structure for the decay of a neutral D-meson into a \(\rho\) meson and the short-lived neutral kaon, \( K^0_S \). States are lists of particles, separated by commas, and included within {} braces. As a trivial example, this decay descriptor can be used to generate some of its own documentation, in this case a formatted latex string:
std::cout << example.texLabel(true) << std::endl;
Describes a particle, its decay process and subsequent decay products, which are also Particles.
produces the latex source for the output
\[D^{0}\rightarrow \rho(770)^{0}\left[\pi^{+} \pi^{-} \right] K_{S}^{0}
\]
There are also modifiers to the amplitude that alter the details of the decay process. For example, there may multiple orbital angular momentum substates. These will be denoted by the use of [] braces, so for example
a(1)(1260)+{rho(770)0,pi+}
a(1)(1260)+[D]{rho(770)0,pi+}
correspond to the S-wave and D-wave decays of the \(a(1)(1260)^{+}\) meson. As a rule, the lowest orbital state permitted by the relevant conservation laws of the decay is used if the orbital state is not specified, so the conservation of angular momentum, and the conservation of parity if the decay proceeds via the strong or electromagnetic force.
The modifier syntax is also used to specify a different choice of lineshape for the resonance. For example, a common parameterisation for the \(\rho(770)\) meson is the Gounaris-Sakurai propagator, which accounts for dispersive corrections to the \(I=1\) dipion system. In this example
rho(770)0[GounarisSakurai]{pi+,pi-}
Multiple modifiers can be applied to the same particle, by including them in a semi-colon separated list. For example, for three-body decays of broad resonances, such as the \(a_1(1260)\), the propagator must take the evolution of the width of resonance from a numerical calculation, which can be supplied via a cubic spline (hence using the GSpline propagator). Additionally, for the decay chain \( a_1(1260)^{+} \to \rho^{0} \pi^{+} \), the decay products can either be in a relative S or D wave, and hence the two particle descriptors
a(1)(1260)+[GSpline]{rho(770)0,pi+}
a(1)(1260)+[D;GSpline]{rho(770)0,pi+}
are relevant for the decay.
Similar to other components of AmpGen, Particles will rarely be constructed in the C++ context, and will instead be instantiated dynamically at runtime from a user supplied options file.
Definition at line 102 of file Particle.h.
|
| Particle () |
| Default Constructor.
|
|
| Particle (const std::string &name, const Particle &p1, const Particle &p2) |
| Constructor that takes a pair of other particles (i.e. this particle's decay products) as arguments and looks up the properties of this particle using the particle name.
|
|
| Particle (const int &pdg_id, const Particle &p1, const Particle &p2) |
| Constructor that takes a pair of other particles (i.e. this particle's decay products) as arguments and looks up the properties of this particle using the PDG MC ID.
|
|
| Particle (const std::string &name, const unsigned int &index) |
| Constructor by name and with an index to match to the event type.
|
|
| Particle (const std::string &decayString, const std::vector< std::string > &finalStates={}, const bool &orderDaughters=true) |
| Constructor that takes a decay descriptor as an argument and a list of final state particles to match to the event type. Constructs the entire decay tree.
|
|
| Particle (const std::string &name, const std::vector< Particle > &particles) |
| Constructor that takes a set of particles.
|
|
Particle | conj (bool invertHead=true, bool reorder=true) |
| (Quasi) Constructor that returns the (quasi)CP conjugated amplitude. The full behaviour of the amplitude is made more complicated by the ordering convention.
|
|
void | conjThis () |
| CP conjugate this particle //.
|
|
void | setOrbital (const unsigned int &orbital) |
| Set the orbital quantum number 'L' for this decay.
|
|
void | setLineshape (const std::string &lineshape) |
| Set the lineshape for the decay of this particle.
|
|
void | setDaughter (const Particle &particle, const unsigned int &index) |
| Set the index'th daughter of this to particle.
|
|
void | setParent (const Particle *particle) |
| Set the parent particle for this particle.
|
|
void | setIndex (const unsigned int &index, const bool &setOri=false) |
| Set the index of this particle, i.e. where it is positioned in the event data structure.
|
|
void | clearDecayProducts () |
| Remove all of the decay products of this particle.
|
|
void | addModifier (const std::string &mod) |
| Add some modifier to the particle, such as a lineshape or a different spin state.
|
|
void | parseModifier (const std::string &mod) |
| Parse some set of modifiers, delimited with semicolons.
|
|
void | setOrdering (const std::vector< size_t > &ordering) |
| Set some particle ordering of the decay products of this particle, mostly used internally by the symmetrisation.
|
|
void | setName (const std::string &name) |
| Set the particle name.
|
|
void | addDaughter (const std::shared_ptr< Particle > &particle) |
| Add a decay product.
|
|
void | setPolarisationState (const int &state) |
| Set the polarisation state of this particle, which is twice the projection of the spin along the quantisation axis.
|
|
void | setPolarisationState (const std::vector< int > &state) |
|
std::pair< size_t, size_t > | orbitalRange (const bool &converseParity=true) const |
| Returns the range of orbital angular momentum between the decay products.
|
|
std::vector< std::pair< double, double > > | spinOrbitCouplings (const bool &conserveParity=true) const |
| Returns the set of possible spin-orbit couplings allowed by conservation of angular momentum, and if specified parity.
|
|
void | setDaughters (const std::vector< Particle > &particles) |
|
stdx::optional< std::string > | attribute (const std::string &key) const |
| Return the additional optional attribute keyed by variable key.
|
|
const ParticleProperties * | props () const |
| Return the particleProperties object for this particle.
|
|
QuarkContent | quarks () const |
| Return the quarks of this particle.
|
|
QuarkContent | daughterQuarks () const |
| Returns the quark content of the sum of the decay products of this particle.
|
|
int | parity () const |
| Returns the parity of this particle.
|
|
int | finalStateParity () const |
| Returns the parity of the final state of this particle.
|
|
int | polState () const |
| Returns the polarisation state, i.e. twice the projection of the spin along the quantisation axis, of this particle.
|
|
int | CP () const |
| Returns the CP of this decay.
|
|
int | C () const |
| Returns the C quantum number for this decay.
|
|
double | mass () const |
| Returns the (PDG) mass of the particle.
|
|
double | spin () const |
| Returns the spin of the particle.
|
|
double | S () const |
| Returns the spin configuration of the decay products of the particle.
|
|
unsigned | L () const |
| Returns the orbital angular.
|
|
bool | isHead () const |
| Returns whether if this particle is the head of the decay, i.e. has no parent.
|
|
bool | isWeakDecay () const |
| Returns whether is this particle decays weakly or not.
|
|
bool | isStateGood () const |
| Returns whether this particle, and its decays have been configured correctly.
|
|
bool | isStable () const |
| Check whether this particle is stable, has any decay products.
|
|
bool | isQuasiStable () const |
| Check whether the particle is quasi-stable, i.e. may have some appreciable flight distance.
|
|
bool | conservesParity (unsigned int L=0) const |
| Check whether the decay of this particle with angular momentum L conserves parity or not.
|
|
unsigned | index () const |
| Returns the current index of the particle in event data structure. Can differ from the original index due to symmetrisation.
|
|
unsigned | originalIndex () const |
| Returns the original index of the particle.
|
|
std::string | name () const |
| Name of the decaying particle.
|
|
std::string | lineshape () const |
| Name of the propagator to use for the decay of this particle.
|
|
std::string | uniqueString () const |
| Returns the unique string (i.e. decay descriptor) that identifies this decay, which can be parsed to generate the decay tree.
|
|
std::string | decayDescriptor () const |
| Returns the unique string (i.e. decay descriptor) that identifies this decay, which can be parsed to generate the decay tree.
|
|
std::string | topologicalString () const |
| The string that describes the spin/orbital topology of this decay, i.e.
|
|
std::string | orbitalString () const |
| The string that describes the spin/orbit configuration of this decay.
|
|
std::string | texLabel (const bool &printHead=false, const bool &recurse=true) const |
| Decay descriptor formatted as LaTeX for this decay.
|
|
EventType | eventType () const |
| Return the eventType for this decay (i.e. the initial and final state particles)
|
|
const Particle * | parent () const |
| Returns the parent of the particle.
|
|
std::shared_ptr< Particle > | daughter (const size_t &index) |
| Returns the indexth decay product of this particle.
|
|
std::shared_ptr< Particle > | daughter (const size_t &index) const |
| Returns in indexth decay product of this particle (as constant)
|
|
std::shared_ptr< Particle > | daughter (const std::string &name, const int &maxDepth=-1) const |
| Returns in indexth decay product of this particle (as constant)
|
|
std::vector< std::shared_ptr< Particle > > | daughters () const |
| Vector of decay products of this particle.
|
|
std::vector< std::vector< size_t > > | identicalDaughterOrderings () const |
| Get orderings of the final state that are identical to each other, i.e. those that only differ by exchanging identical particles.
|
|
std::vector< std::shared_ptr< Particle > > | getFinalStateParticles (const bool &sort=true) const |
| Returns the final state particles for this decay process.
|
|
Particle | quasiStableTree () const |
| Calculate the particle tree only including quasi-stable processes, i.e.
|
|
Tensor | P () const |
| Calculates the momentum sum of the decay products.
|
|
Tensor | Q () const |
| Calculates the momentum difference between the decay products (only well defined for quasi two-body processes )
|
|
Tensor | spinTensor (DebugSymbols *db=nullptr) const |
| Calculates the spin tensor or generalised current for this particle.
|
|
Tensor | externalSpinTensor (const int &polState, DebugSymbols *db=nullptr) const |
| Calculates the polarisation vector / spinor etc. of this particle, used for the initial/final state particles.
|
|
Expression | massSq () const |
| Calculates the invariant mass-squared of the mass of this particle.
|
|
Expression | propagator (DebugSymbols *db=nullptr) const |
| Calculates the lineshape / propagator for this particle.
|
|
Expression | getExpression (DebugSymbols *db=nullptr, const std::vector< int > &={}) |
| Calculates the total expression for this particle, including symmetrisation and the current polarisation state.
|
|
bool | lineshapeContains (const std::vector< std::string > &container) const |
| Check if lineshape contains a substring.
|
|
Tensor | transitionMatrix (DebugSymbols *db=nullptr) |
| Calculate the transition matrix for this decay.
|
|
bool | operator< (const Particle &other) |
|
bool | operator> (const Particle &other) |
|
unsigned int | matches (const Particle &other) const |
| matches Check the matching between two decay chains, according to the MatchState enum.
|
|
std::string | makeUniqueString () |
| Generate the decay descriptor for this decay.
|
|
Particle | clone () const |
|
void | setDaughter (const unsigned int &index, const Particle &p) |
|
bool | expand (const Particle &particle) |
|