1#ifndef AMPGEN_PARTICLE_H
2#define AMPGEN_PARTICLE_H
9#if __cplusplus >= 201703L
12 using namespace ::std;
14#elif __cplusplus >= 201402L
15 #include <experimental/optional>
17 using namespace ::std;
18 using namespace ::std::experimental;
21# error "Require c++ std >=14"
118 Particle(
const std::string& decayString,
const std::vector<std::string>& finalStates = {},
const bool& orderDaughters = true );
121 Particle(
const std::string&
name,
const std::vector<Particle>& particles);
169 std::pair<size_t,size_t>
orbitalRange(
const bool& converseParity =
true )
const;
177 stdx::optional<std::string>
attribute(
const std::string& key)
const;
220 std::string
texLabel(
const bool& printHead =
false,
const bool& recurse=
true )
const;
234 std::shared_ptr<Particle>
daughter(
const std::string&
name,
const int& maxDepth=-1)
const;
237 std::vector<std::shared_ptr<Particle>>
daughters()
const;
272 for(
auto& st : container )
if ( m_lineshape.find(st) != std::string::npos )
return true;
297 std::string m_name = {
""};
299 std::string m_lineshape = {
"BW"};
300 std::string m_vertexName = {
""};
301 std::string m_uniqueString = {
""};
303 int m_polState = {0};
304 unsigned m_index = {999};
305 unsigned m_originalIndex = {999};
306 unsigned m_orbital = {0};
307 unsigned m_spinConfigurationNumber = {0};
308 unsigned m_minL = {0};
309 bool m_usesDefaultLineshape = {
false};
310 bool m_isStateGood = {
true};
311 std::vector<std::shared_ptr<Particle>> m_daughters;
312 std::vector<std::string> m_modifiers;
313 const Particle* m_parent = {
nullptr};
314 void pdgLookup(
bool quiet=
false);
315 bool hasModifier(
const std::string& modifier )
const;
316 std::string modifierString()
const;
317 void sortDaughters();
321 std::make_pair(
"Covariant",
"[default] Covariant Tensor, based on Rarita-Schwinger constraints on the allowed covariant wavefunctions.")
322 , std::make_pair(
"Canonical",
"Canonical formulation, based on rotational properties of wavefunctions, i.e. Wigner D-matrices and Clebsch-Gordan for (L,S) expansion.") ) };
324 optionalHelpString(
"Basis to use for calculating external polarisation tensors / spinors.",
325 std::make_pair(
"Dirac",
"[default] Quantises along the z-axis")
326 , std::make_pair(
"Weyl",
"Quantises along the direction of motion") )};
327 NamedParameter<std::string> m_defaultModifier = {
"Particle::DefaultModifier",
"",
"Default modifier to use for lineshapes, for example to use normalised vs unnormalised Blatt-Weisskopf factors."};
Deals with final state configuration of events, specifically dealing with the ordering of particles i...
Wrapper class for shared_ptrs to virtual expressions for use in conjunction with operators to build e...
A parameter with value specified by the user at runtime, either in an options file or via the command...
Describes a particle, its decay process and subsequent decay products, which are also Particles.
stdx::optional< std::string > attribute(const std::string &key) const
Return the additional optional attribute keyed by variable key.
unsigned int matches(const Particle &other) const
matches Check the matching between two decay chains, according to the MatchState enum.
std::string texLabel(const bool &printHead=false, const bool &recurse=true) const
Decay descriptor formatted as LaTeX for this decay.
void setName(const std::string &name)
Set the particle name.
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 setParent(const Particle *particle)
Set the parent particle for this particle.
void setOrbital(const unsigned int &orbital)
Set the orbital quantum number 'L' for this decay.
void setDaughter(const unsigned int &index, const Particle &p)
Particle(const std::string &name, const unsigned int &index)
Constructor by name and with an index to match to the event type.
int finalStateParity() const
Returns the parity of the final state of this particle.
void conjThis()
CP conjugate 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 exc...
std::vector< std::shared_ptr< Particle > > daughters() const
Vector of decay products of this particle.
std::pair< size_t, size_t > orbitalRange(const bool &converseParity=true) const
Returns the range of orbital angular momentum between the decay products.
bool expand(const Particle &particle)
Tensor Q() const
Calculates the momentum difference between the decay products (only well defined for quasi two-body p...
std::vector< std::shared_ptr< Particle > > getFinalStateParticles(const bool &sort=true) const
Returns the final state particles for this decay process.
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 a...
int C() const
Returns the C quantum number for this decay.
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.
std::string lineshape() const
Name of the propagator to use for the decay of this particle.
std::shared_ptr< Particle > daughter(const std::string &name, const int &maxDepth=-1) const
Returns in indexth decay product of this particle (as constant)
unsigned L() const
Returns the orbital angular.
Particle quasiStableTree() const
Calculate the particle tree only including quasi-stable processes, i.e.
const Particle * parent() const
Returns the parent of the particle.
int CP() const
Returns the CP of this decay.
int parity() const
Returns the parity of this particle.
Tensor spinTensor(DebugSymbols *db=nullptr) const
Calculates the spin tensor or generalised current for this particle.
unsigned index() const
Returns the current index of the particle in event data structure. Can differ from the original index...
void clearDecayProducts()
Remove all of the decay products of this particle.
void setPolarisationState(const int &state)
Set the polarisation state of this particle, which is twice the projection of the spin along the quan...
std::shared_ptr< Particle > daughter(const size_t &index)
Returns the indexth decay product of 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 p...
bool conservesParity(unsigned int L=0) const
Check whether the decay of this particle with angular momentum L conserves parity or not.
Particle(const std::string &name, const std::vector< Particle > &particles)
Constructor that takes a set of particles.
unsigned originalIndex() const
Returns the original index of the particle.
void addModifier(const std::string &mod)
Add some modifier to the particle, such as a lineshape or a different spin state.
std::string orbitalString() const
The string that describes the spin/orbit configuration of this decay.
bool operator<(const Particle &other)
void setPolarisationState(const std::vector< int > &state)
std::string topologicalString() const
The string that describes the spin/orbital topology of this decay, i.e.
std::string makeUniqueString()
Generate the decay descriptor for this decay.
const ParticleProperties * props() const
Return the particleProperties object for this particle.
Particle()
Default Constructor.
QuarkContent daughterQuarks() const
Returns the quark content of the sum of the decay products of this particle.
Particle conj(bool invertHead=true, bool reorder=true)
(Quasi) Constructor that returns the (quasi)CP conjugated amplitude. The full behaviour of the amplit...
bool isQuasiStable() const
Check whether the particle is quasi-stable, i.e. may have some appreciable flight distance.
void setDaughters(const std::vector< Particle > &particles)
Expression propagator(DebugSymbols *db=nullptr) const
Calculates the lineshape / propagator for this particle.
bool lineshapeContains(const std::vector< std::string > &container) const
Check if lineshape contains a substring.
Expression massSq() const
Calculates the invariant mass-squared of the mass of this particle.
void setLineshape(const std::string &lineshape)
Set the lineshape for the decay of this particle.
bool operator>(const Particle &other)
Tensor P() const
Calculates the momentum sum of the decay products.
QuarkContent quarks() const
Return the quarks of this particle.
std::shared_ptr< Particle > daughter(const size_t &index) const
Returns in indexth decay product of this particle (as constant)
void setOrdering(const std::vector< size_t > &ordering)
Set some particle ordering of the decay products of this particle, mostly used internally by the symm...
Expression getExpression(DebugSymbols *db=nullptr, const std::vector< int > &={})
Calculates the total expression for this particle, including symmetrisation and the current polarisat...
std::string name() const
Name of the decaying particle.
bool isStable() const
Check whether this particle is stable, has any decay products.
bool isStateGood() const
Returns whether this particle, and its decays have been configured correctly.
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,...
int polState() const
Returns the polarisation state, i.e. twice the projection of the spin along the quantisation axis,...
void parseModifier(const std::string &mod)
Parse some set of modifiers, delimited with semicolons.
std::string decayDescriptor() const
Returns the unique string (i.e. decay descriptor) that identifies this decay, which can be parsed to ...
double mass() const
Returns the (PDG) mass of the particle.
EventType eventType() const
Return the eventType for this decay (i.e. the initial and final state particles)
double S() const
Returns the spin configuration of the decay products of the particle.
double spin() const
Returns the spin of the particle.
void setDaughter(const Particle &particle, const unsigned int &index)
Set the index'th daughter of this to particle.
Tensor transitionMatrix(DebugSymbols *db=nullptr)
Calculate the transition matrix for this decay.
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 a...
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...
std::string uniqueString() const
Returns the unique string (i.e. decay descriptor) that identifies this decay, which can be parsed to ...
void addDaughter(const std::shared_ptr< Particle > &particle)
Add a decay product.
static bool isValidDecayDescriptor(const std::string &decayDescriptor)
Class that contains the PDG properties (mass, width, charges, etc.) for a single particle species,...
#define declare_enum(name,...)
std::ostream & operator<<(std::ostream &os, const CompiledExpressionBase &expression)
std::string optionalHelpString(const std::string &header, const T &... args)
std::vector< DebugSymbol > DebugSymbols
declare_enum(coordinateType, cartesian, polar) declare_enum(angType