Documentation of MARTY
A Modern ARtificial Theoretical phYsicist
Data Structures | Public Types | Public Member Functions | Protected Member Functions | Static Protected Member Functions | Protected Attributes
mty::Spectrum Class Reference

Class handling numerical diagonalization for HEP models. More...

#include <spectrum.h>

Data Structures

struct  MatrixEl
 Small structure containing a mass term (without fields) and its position in the mass matrix. More...
 

Public Types

template<class T >
using matrix = boost::numeric::ublas::matrix< T >
 Using for matrices in boost libraries. More...
 
template<class T >
using id_matrix = boost::numeric::ublas::identity_matrix< T >
 Using for identity matrices in boost libraries. More...
 
using FloatType = double
 Type of floating point variable for numerical diagonalization.
 

Public Member Functions

 Spectrum (std::vector< mty::Particle > const &t_particles, std::vector< mty::Particle > const &t_newParticles, std::vector< csl::Expr > const &terms, std::vector< std::vector< csl::Expr >> const &t_mixing)
 Constructor with 3 parameters, for standard diagonalization. More...
 
 Spectrum (std::vector< mty::Particle > const &t_newParticles, std::vector< std::vector< csl::Expr >> const &t_mass, std::vector< std::vector< csl::Expr >> const &t_mixing, std::vector< std::vector< csl::Expr >> const &t_mixing2=std::vector< std::vector< csl::Expr >>())
 Constructor with 3 parameters, for standard diagonalization. More...
 
 Spectrum (std::vector< mty::Particle > const &partA, std::vector< mty::Particle > const &partB, std::vector< mty::Particle > const &newPartA, std::vector< mty::Particle > const &newPartB, std::vector< csl::Expr > const &terms, std::vector< std::vector< csl::Expr >> const &mixingA, std::vector< std::vector< csl::Expr >> const &mixingB)
 Constructor with 6 parameters, for bi-diagonalization. More...
 
 Spectrum ()=default
 Default constructor.
 
 ~Spectrum ()=default
 Default destructor.
 
 Spectrum (Spectrum const &)=default
 Default copy constructor.
 
 Spectrum (Spectrum &&)=default
 Default move constructor.
 
Spectrumoperator= (Spectrum const &)=default
 Default copy assignement operator.
 
Spectrumoperator= (Spectrum &&)=default
 Default move assignement operator.
 
std::vector< mty::Particle > const & getParticles () const
 
matrix< csl::Expr > const & getMassTerms () const
 
matrix< csl::Expr > const & getMixings () const
 
matrix< csl::Expr > const & getMixings2 () const
 
bool getReplaceMasses () const
 
bool getReplaceMixings () const
 
void setReplaceMasses (bool t_replaceMasses)
 Sets replaceMasses. More...
 
void setReplaceMixings (bool t_replaceMixings)
 Sets replaceMixings. More...
 
bool isDiagonalized () const
 Tells if the mass matrix has been numerically diagonalized. More...
 
void setBlock (std::vector< csl::Expr > const &terms, std::vector< std::vector< csl::Expr >> const &t_mixing)
 Instantiates the Spectrum object with a mass terms and a mixing matrix for standard diagonalization. More...
 
void setBlock (std::vector< csl::Expr > const &terms, std::vector< std::vector< csl::Expr >> const &mixingA, std::vector< std::vector< csl::Expr >> const &mixingB)
 Instantiates the Spectrum object with a mass terms and two mixing matrices for bi-diagonalization. More...
 
void updateData ()
 Updates the mass matrix, if the elements have values (for example) after a lha file loading. More...
 
void applyOn (csl::Expr &expr) const
 Applies the diagonalization on an expression, replacing masses and mixings. More...
 
void applyOn (std::vector< csl::Expr > &expr) const
 Applies the diagonalization on a vector of expressions, replacing masses and mixings. More...
 

Protected Member Functions

size_t getFieldPos (QuantumField const &field) const
 Returns the position of the field field in the list. More...
 
MatrixEl getMassFromTerm (csl::Expr const &term) const
 Returns a MatrixEl object from a Lagrangian term. More...
 
void addMassTerm (MatrixEl &&matrixEl)
 Adds a mass term to the Spectrum from a MatrixEl object. More...
 
void diagonalize ()
 Diagonalizes the mass matrix numerically.
 
void bidiagonalize ()
 Bi-diagonalizes the mass mattrix numerically.
 
void applyDiagonalizationOn (csl::Expr &expr) const
 Applies the result of the diagonalization on an expression. More...
 
void applyBidiagonalizationOn (csl::Expr &expr) const
 Applies the result of the bi-diagonalization on an expression. More...
 

Static Protected Member Functions

static void replace (csl::Expr &expr, std::vector< csl::Expr > const &oldExpr, std::vector< csl::Expr > const &newExpr)
 Replaces corresponding mixings / masses in an expression, if the numerical diagonalization has been performed. More...
 

Protected Attributes

bool bidiagonalization
 Tells if this class is used for a standard diagonalization (false) or a bi-diagonalization (true). Default is false.
 
bool replaceMasses
 Tells if masses must be directly replaced by their numerical values when applying diagonalization (default = true).
 
bool replaceMixings
 Tells if mixings must be directly replaced by their numerical values when applying diagonalization (default = true).
 
std::vector< mty::Particlefields
 List of fields in the mass block.
 
std::vector< mty::ParticlenewFields
 List of fields introduced after the rotation.
 
matrix< csl::Exprmass
 Mass matrix. This object should contain only quantities with defined numerical when calling the diagonalization.
 
matrix< csl::Exprmixing
 Initial mixing matrix, contains the symbolic variables that enter expressions before the diagonalization is done.
 
matrix< csl::Exprmixing2
 Second mixing matrix, used for bi-diagonalization.
 
csl::Expr massData
 Numerical mass matrix, once values have been given and updateData() has been called.
 
csl::Expr diagonal = nullptr
 Diagonalized mass matrix if the diagonalization has been done.
 
csl::Expr transfer = nullptr
 Second transfer matrix used for bi-diagonalization.
 
csl::Expr transfer2 = nullptr
 Transfer matrix if the diagonalization has been done.
 

Detailed Description

Class handling numerical diagonalization for HEP models.

For some mass mixings, the standard way to proceed in MARTY (when a simple diagonalization is not possible) is to create a symbolic mixing matrix, fully general, that enters amplitudes. Calculations are done with these mixings. If later one wants to specify numerically the initial mass matrix, this class may diagonalize it and recover the mixing matrix from it. Then, the actual (numerical) mixings may be applied on an expression (for example an amplitude or Wilson coefficient), replacing old masses by numerical eigenvalues and symbolic mixings by their values.

There is three parts when using this object. First: the initialization from a MassBlock and mixing terms. It needs the initial mass matrix, and the mixing matrix that you instored to symbolically diagonalize it in full generality. Then, once all the matrix elements of the initial mass matrix have a numerical value you must call updateData() to initialize the numerical matrix. Finally, you may call applyOn() to apply this diagonalization to an expression (or several).

Two types of diagonalization are possible. For a given quantum field \( \Phi \) (bosonic or fermionic), you may have something like \( \Phi ^\dagger M \Phi \) in which case \( M \) is an hermitian matrix as the Lagrangian must be real. A standard diagonalization must then be performed, i.e.

\[ \Phi \mapsto \Phi ^\prime = U\Phi,\\ \Phi ^\dagger \mapsto \Phi ^{\prime\dagger} = \Phi ^\dagger U^\dagger,\\ UMU^\dagger \equiv D, \text{ the diagonal mass matrix. } \]

Another possibility are terms of the form \( \Phi _L^\dagger M \Phi _R + \Phi _R^\dagger M^\dagger \Phi_L \), with \( \Phi _L \) and \( \Phi _R \) different fields. In that case, a bi-diagonalization must be performed, i.e.

\[ \Phi _R \mapsto \Phi _R^\prime = U\Phi _R,\\ \Phi _L \mapsto \Phi _L^\prime = V\Phi _L,\\ V^\dagger MU = D, \text{ the diagonal mass matrix, }\\ V^\dagger MM^\dagger V = D^2,\\ U^\dagger M^\dagger M U = D^2\\ \]

In this case there is actually two diagonalizations to perform, for \( MM^\dagger \) and \( M^\dagger M \).

See also
csl::Diagonalizer

Member Typedef Documentation

◆ id_matrix

Using for identity matrices in boost libraries.

Template Parameters
TType of the value type in the matrix.

◆ matrix

Using for matrices in boost libraries.

Template Parameters
TType of the value type in the matrix.

Constructor & Destructor Documentation

◆ Spectrum() [1/3]

mty::Spectrum::Spectrum ( std::vector< mty::Particle > const &  t_particles,
std::vector< mty::Particle > const &  t_newParticles,
std::vector< csl::Expr > const &  terms,
std::vector< std::vector< csl::Expr >> const &  t_mixing 
)

Constructor with 3 parameters, for standard diagonalization.

Parameters
t_particlesParticles to diagonalize.
t_newParticlesNew particles introduces after rotation.
termsVector of mass terms for particles.
t_mixingMixing introduced to symbolically diagonalize the mass matrix (2D matrix as vector).

◆ Spectrum() [2/3]

mty::Spectrum::Spectrum ( std::vector< mty::Particle > const &  t_newParticles,
std::vector< std::vector< csl::Expr >> const &  t_mass,
std::vector< std::vector< csl::Expr >> const &  t_mixing,
std::vector< std::vector< csl::Expr >> const &  t_mixing2 = std::vector<std::vector<csl::Expr>>() 
)

Constructor with 3 parameters, for standard diagonalization.

This constructor will not do anything apart from reading parameters. These parameters are considered to have already been treated by a Spectrum object. It us used in particular to re-build Spectrum object in a different program (Lagrangian code generation).

Parameters
t_newParticlesNew particles introduces after rotation.
t_massInitial mass matrix for the Spectrum object.
t_mixingMixing introduced to symbolically diagonalize the mass matrix (2D matrix as vector).
t_mixingMixing introduced to symbolically diagonalize the mass matrix (2D matrix as vector) for right particles, if there is.

◆ Spectrum() [3/3]

mty::Spectrum::Spectrum ( std::vector< mty::Particle > const &  partA,
std::vector< mty::Particle > const &  partB,
std::vector< mty::Particle > const &  newPartA,
std::vector< mty::Particle > const &  newPartB,
std::vector< csl::Expr > const &  terms,
std::vector< std::vector< csl::Expr >> const &  mixingA,
std::vector< std::vector< csl::Expr >> const &  mixingB 
)

Constructor with 6 parameters, for bi-diagonalization.

Parameters
partALeft block of particles to diagonalize
partBRight block of particles to diagonalize
newPartANew particles introduces after rotation, left block.
newPartBNew particles introduces after rotation, right block.
termsVector of mass terms for the particles.
mixingAMixing introduced to symbolically diagonalize the mass matrix (2D matrix as vector) for left particles.
mixingBMixing introduced to symbolically diagonalize the mass matrix (2D matrix as vector) for right particles.

Member Function Documentation

◆ addMassTerm()

void mty::Spectrum::addMassTerm ( MatrixEl &&  matrixEl)
protected

Adds a mass term to the Spectrum from a MatrixEl object.

Parameters
matrixElMatrix element containing the term and its position in the mass matrix.

◆ applyBidiagonalizationOn()

void mty::Spectrum::applyBidiagonalizationOn ( csl::Expr expr) const
protected

Applies the result of the bi-diagonalization on an expression.

Parameters
exprcsl::Expr to modify.

◆ applyDiagonalizationOn()

void mty::Spectrum::applyDiagonalizationOn ( csl::Expr expr) const
protected

Applies the result of the diagonalization on an expression.

Parameters
exprcsl::Expr to modify.

◆ applyOn() [1/2]

void mty::Spectrum::applyOn ( csl::Expr expr) const

Applies the diagonalization on an expression, replacing masses and mixings.

This function should be called if updateData() has already been called and if all elements of the initial mass matrix have defined values.

Parameters
exprExpression on which the diagonalization is applied.

◆ applyOn() [2/2]

void mty::Spectrum::applyOn ( std::vector< csl::Expr > &  expr) const

Applies the diagonalization on a vector of expressions, replacing masses and mixings.

This function should be called if updateData() has already been called and if all elements of the initial mass matrix have defined values.

Parameters
exprExpressions on which the diagonalization is applied.

◆ getFieldPos()

size_t mty::Spectrum::getFieldPos ( QuantumField const &  field) const
protected

Returns the position of the field field in the list.

Parameters
fieldQuantumField.
Returns
The position of field in the attribute fields.

◆ getMassFromTerm()

Spectrum::MatrixEl mty::Spectrum::getMassFromTerm ( csl::Expr const &  term) const
protected

Returns a MatrixEl object from a Lagrangian term.

Parameters
termLagrangian term.
Returns
The MatrixEl containing the position of the term in the matrix, and the term itself.

◆ getReplaceMasses()

bool mty::Spectrum::getReplaceMasses ( ) const
inline
Returns
replaceMasses.

◆ getReplaceMixings()

bool mty::Spectrum::getReplaceMixings ( ) const
inline
Returns
replaceMixings.

◆ isDiagonalized()

bool mty::Spectrum::isDiagonalized ( ) const

Tells if the mass matrix has been numerically diagonalized.

If this function returns true, it means that replacement may be applied on expressions, using applyOn(). If initial mass parameters have a numerical value, consider calling updateData(), then checking that everything is going well checking that this function returns true, and finally apply the diagonalization on any expression you want using applyOn().

Returns
True if all masses and mixings have been computed numerically.
False else.
See also
applyOn(), updateData()

◆ replace()

void mty::Spectrum::replace ( csl::Expr expr,
std::vector< csl::Expr > const &  oldExpr,
std::vector< csl::Expr > const &  newExpr 
)
staticprotected

Replaces corresponding mixings / masses in an expression, if the numerical diagonalization has been performed.

Is abbreviations are encountered, the content of the abbreviation os also looked up for the replacement.

Parameters
exprExpression in which we apply the diagonalization.
oldExprOld mixings / masses.
newExprNew mixings / masses.

◆ setBlock() [1/2]

void mty::Spectrum::setBlock ( std::vector< csl::Expr > const &  terms,
std::vector< std::vector< csl::Expr >> const &  t_mixing 
)

Instantiates the Spectrum object with a mass terms and a mixing matrix for standard diagonalization.

This function simply gathers all terms in terms to create the mass matrix. The mixing matrix used to symbolically diagonalize the fields must also be given.

Parameters
termsMass terms for the fields.
t_mixingMixing matrix.

◆ setBlock() [2/2]

void mty::Spectrum::setBlock ( std::vector< csl::Expr > const &  terms,
std::vector< std::vector< csl::Expr >> const &  mixingA,
std::vector< std::vector< csl::Expr >> const &  mixingB 
)

Instantiates the Spectrum object with a mass terms and two mixing matrices for bi-diagonalization.

This function simply gathers all terms in terms to create the mass matrix. The mixing matrices used to symbolically diagonalize the fields must also be given.

Parameters
termsMass terms for the fields.
mixingAMixing matrix of left fields.
mixingBMixing matrix of right fields.

◆ setReplaceMasses()

void mty::Spectrum::setReplaceMasses ( bool  t_replaceMasses)
inline

Sets replaceMasses.

Parameters
t_replaceMassesNew value for replaceMasses.

◆ setReplaceMixings()

void mty::Spectrum::setReplaceMixings ( bool  t_replaceMixings)
inline

Sets replaceMixings.

Parameters
t_replaceMixingsNew value for replaceMixings.

◆ updateData()

void mty::Spectrum::updateData ( )

Updates the mass matrix, if the elements have values (for example) after a lha file loading.

This function should be called once all the initial elements of the mass matrix have a defined value.

To properly define a value, you must use the function csl::Abstract::setValue() of the symbolic constants. For example with a mass matrix \( \left(\begin{array}{cc}a && b \\ c && d \end{array}\right),\) one must write

M[0][0]->setValue(523.5); // good

and not

M[0][0] = 523.5; // very bad

In the latter case, the matrix element is replaced by the value 523.5, but the variable a that actually appears in the Spectrum class still has no value. The diagonalization cannot be performed then.


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