Belle II Software  light-2205-abys
BeamParameters Class Reference

This class contains the nominal beam parameters and the parameters used for smearing of the primary vertex and the beam energy and momentum. More...

#include <BeamParameters.h>

Inheritance diagram for BeamParameters:
Collaboration diagram for BeamParameters:

Public Types

enum  EGenerationFlags {
  c_generateCMS = 1 << 0 ,
  c_smearBeamEnergy = 1 << 1 ,
  c_smearBeamDirection = 1 << 2 ,
  c_smearBeam = c_smearBeamEnergy | c_smearBeamDirection ,
  c_smearVertex = 1 << 3 ,
  c_smearALL = c_smearVertex | c_smearBeam
}
 Possible Flags for initial event generation. More...
 

Public Member Functions

 BeamParameters ()
 Default constructor.
 
 BeamParameters (const BeamParameters &b)
 Copy constructor.
 
BeamParametersoperator= (const BeamParameters &b)
 Assignment operator.
 
bool operator== (const BeamParameters &b) const
 Equality operator.
 
void setGenerationFlags (int flags) override
 Set the generation flags to be used for event generation (ORed combination of EGenerationFlags). More...
 
void setCovHER (const TMatrixDSym &cov)
 Set the covariance matrix for HER (E, theta_x, theta_y) where E is the energy, theta_x is the horizontal angle between nominal direction and actual direction in spread and theta_y is the vertical angle The upper triangle will be saved.
 
void setCovLER (const TMatrixDSym &cov)
 Set the covariance matrix for LER (E, theta_x, theta_y) where E is the energy, theta_x is the horizontal angle between nominal direction and actual direction in spread and theta_y is the vertical angle. More...
 
void setCovVertex (const TMatrixDSym &cov)
 Set the covariance matrix of the vertex position. More...
 
void setHER (double energy, double angleX, double angleY, const std::vector< double > &cov)
 Set the HER FourVector and error matrix from beam energy, angle and covariance entries. More...
 
void setLER (double energy, double angleX, double angleY, const std::vector< double > &cov)
 Set the LER FourVector and error matrix from beam energy, angle and covariance entries. More...
 
void setVertex (const TVector3 &vertex, const std::vector< double > &cov)
 Set the vertex position and error matrix. More...
 
TMatrixDSym getCovHER () const
 Get the covariance matrix of HER (E, theta_x, theta_y) where E is the energy, theta_x is the horizontal angle between nominal direction and actual direction in spread and theta_y is the vertical angle.
 
TMatrixDSym getCovLER () const
 Get the covariance matrix of LER (E, theta_x, theta_y) where E is the energy, theta_x is the horizontal angle between nominal direction and actual direction in spread and theta_y is the vertical angle.
 
TMatrixDSym getCovVertex () const
 Get the covariance matrix of the vertex position.
 
void setLER (const ROOT::Math::PxPyPzEVector &ler)
 Set the Low Energy Beam 4-momentum.
 
void setHER (const ROOT::Math::PxPyPzEVector &her)
 Set the High Energy Beam 4-momentum.
 
void setVertex (const B2Vector3D &vertex)
 Set the vertex position.
 
bool operator== (const MCInitialParticles &b) const
 Equality operator.
 
void set (const ROOT::Math::PxPyPzEVector &her, const ROOT::Math::PxPyPzEVector &ler, const B2Vector3D &vertex)
 Set the initial event values, i.e. More...
 
void setHER (const ROOT::Math::PxPyPzEVector &her)
 Set the High Energy Beam 4-momentum.
 
void setLER (const ROOT::Math::PxPyPzEVector &ler)
 Set the Low Energy Beam 4-momentum.
 
void setVertex (const B2Vector3D &vertex)
 Set the vertex position.
 
void setTime (double time)
 Set collison time.
 
const ROOT::Math::PxPyPzEVector & getHER () const
 Get 4vector of the high energy beam.
 
const ROOT::Math::PxPyPzEVector & getLER () const
 Get 4vector of the low energy beam.
 
const B2Vector3DgetVertex () const
 Get the position of the collision.
 
double getTime () const
 Get collison time.
 
double getEnergy () const
 Get the the actual collision energy (in lab system)
 
double getMass () const
 Get the invariant mass of the collision (= energy in CMS)
 
const ROOT::Math::LorentzRotation & getLabToCMS () const
 Return the LorentzRotation to convert from lab to CMS frame.
 
const ROOT::Math::LorentzRotation & getCMSToLab () const
 Return the LorentzRotation to convert from CMS to lab frame.
 
bool getValidFlag () const
 Get the flag to check if a valid MCInitialParticles object was already generated and filled in an event.
 
int getGenerationFlags () const
 Get the generation flags to be used for event generation (ORed combination of EGenerationFlags)
 
bool hasGenerationFlags (int flags) const
 Check if a certain set of EGenerationFlags is set.
 
std::string getGenerationFlagString (const std::string &separator=" ") const
 Return string representation of all active flags for printing. More...
 

Protected Attributes

int m_generationFlags {0}
 Flags to be used when generating events.
 

Private Member Functions

 ClassDefOverride (BeamParameters, 3)
 ClassDef.
 
void calculateBoost () const
 Calculate the boost if necessary.
 
void resetBoost ()
 Reset cached transformations after changing parameters.
 
 ClassDef (MCInitialParticles, 4)
 ROOT Dictionary.
 

Static Private Member Functions

static ROOT::Math::PxPyPzEVector getFourVector (double energy, double angleX, double angleY)
 Calculate FourVector of a beam from energy and angle wrt the z-axis. More...
 
static void setCovMatrix (Double32_t *member, const std::vector< double > &cov, bool common)
 Set covariance matrix from vector of entries. More...
 
static void setCovMatrix (Double32_t *member, const TMatrixDSym &cov)
 Set covariance matrix from ROOT Matrix object.
 
static TMatrixDSym getCovMatrix (const Double32_t *member)
 Obtain covariance matrix from a given float array.
 

Private Attributes

Double32_t m_covHER [6]
 Covariance matrix of the high energy beam at the IP.
 
Double32_t m_covLER [6]
 Covariance matrix of the low energy beam at the IP.
 
Double32_t m_covVertex [6]
 Covariance matrix of the vertex position.
 
ROOT::Math::PxPyPzEVector m_her
 HER 4vector.
 
ROOT::Math::PxPyPzEVector m_ler
 LER 4vector.
 
B2Vector3D m_vertex
 collision position
 
double m_time = 0
 collision time
 
ROOT::Math::LorentzRotation * m_labToCMS {nullptr}
 Boost from Lab into CMS. More...
 
ROOT::Math::LorentzRotation * m_CMSToLab {nullptr}
 transient More...
 
double m_invariantMass {0.0}
 transient More...
 
bool m_validFlag = false
 transient More...
 

Detailed Description

This class contains the nominal beam parameters and the parameters used for smearing of the primary vertex and the beam energy and momentum.

It is event independent (but might be run or even sub run dependent)

Definition at line 25 of file BeamParameters.h.

Member Enumeration Documentation

◆ EGenerationFlags

enum EGenerationFlags
inherited

Possible Flags for initial event generation.

Enumerator
c_generateCMS 

generate initial event in CMS instead of lab

c_smearBeamEnergy 

smear energy of HER and LER (but not direction)

c_smearBeamDirection 

smear direction of HER and LER (but not energy)

c_smearBeam 

smear the full beam momentum (energy and direction)

c_smearVertex 

smear vertex

c_smearALL 

smear all

Definition at line 36 of file MCInitialParticles.h.

Member Function Documentation

◆ getFourVector()

ROOT::Math::PxPyPzEVector getFourVector ( double  energy,
double  angleX,
double  angleY 
)
staticprivate

Calculate FourVector of a beam from energy and angle wrt the z-axis.

Negative angles will be treated as angle = M_PI - fabs(angle)

Parameters
energybeam energy
angleXhorizontal angle wrt z-axis
angleYvertical angle wrt z-axis

Definition at line 49 of file BeamParameters.cc.

50 {
51  double p = sqrt(pow(energy, 2) - pow(Const::electronMass, 2));
52  if (angleX < 0) {
53  angleX = M_PI + angleX;
54  angleY = M_PI + angleY;
55  }
56 
57  double Sign = cos(angleX) < 0 ? -1 : 1;
58 
59  double pzTerm = 1 - pow(sin(angleX), 2) - pow(sin(angleY), 2);
60  B2ASSERT("Angular term should be positive", pzTerm >= 0);
61 
62  double pz = Sign * p * sqrt(pzTerm);
63 
64  return ROOT::Math::PxPyPzEVector(p * sin(angleX), p * sin(angleY), pz, energy);
65 
66 }
static const double electronMass
electron mass
Definition: Const.h:565

◆ getGenerationFlagString()

std::string getGenerationFlagString ( const std::string &  separator = " ") const
inherited

Return string representation of all active flags for printing.

Parameters
separatorseparation string to be put between flags

Definition at line 13 of file MCInitialParticles.cc.

◆ set()

void set ( const ROOT::Math::PxPyPzEVector &  her,
const ROOT::Math::PxPyPzEVector &  ler,
const B2Vector3D vertex 
)
inlineinherited

Set the initial event values, i.e.

the four momenta of both beams and the vertex

Parameters
her4vector of the high energy beam
ler4vector of the low energy beam
vertexposition of the actual collision vertex

Definition at line 97 of file MCInitialParticles.h.

◆ setCovLER()

void setCovLER ( const TMatrixDSym &  cov)
inline

Set the covariance matrix for LER (E, theta_x, theta_y) where E is the energy, theta_x is the horizontal angle between nominal direction and actual direction in spread and theta_y is the vertical angle.

The upper triangle will be saved.

Definition at line 97 of file BeamParameters.h.

97 { setCovMatrix(m_covLER, cov); }
Double32_t m_covLER[6]
Covariance matrix of the low energy beam at the IP.
static void setCovMatrix(Double32_t *member, const std::vector< double > &cov, bool common)
Set covariance matrix from vector of entries.

◆ setCovMatrix()

void setCovMatrix ( Double32_t *  member,
const std::vector< double > &  cov,
bool  common 
)
staticprivate

Set covariance matrix from vector of entries.

The vector for the covariance matrix can have either 0, 1, 3, 6 or 9 entries:

  • 0 entries means no error.
  • 1 entry will be treated as the common variance for all variables if common is true. Otherwise it will be taken as the variance of just the first element.
  • 3 entries will be treated as uncorrelated variances for all variables (diagonal of the matrix)
  • 6 entries will be interpreted as the upper triangle of the covariance matrix with the element order being (0, 0), (0, 1), (0, 2), (1, 1), (1, 2), (2, 2)
  • 9 entries for the full covariance matrix in row-major order. The symmetry of the matrix will not be checked but just the lower triangle will be used.
Parameters
memberto the member which contains the matrix
coventries for the covariance matrix
commonif true a 1-element cov will be treated as the common variance for all diagonal elements

Definition at line 88 of file BeamParameters.cc.

◆ setCovVertex()

void setCovVertex ( const TMatrixDSym &  cov)
inline

Set the covariance matrix of the vertex position.

The upper triangle will be saved.

Definition at line 100 of file BeamParameters.h.

◆ setGenerationFlags()

void setGenerationFlags ( int  flags)
inlineoverridevirtual

Set the generation flags to be used for event generation (ORed combination of EGenerationFlags).

The only difference w.r.t. MCInitialParticles::setGenerationFlags() is that a WARNING is thrown if a flag different from 14 (= all the smearings turned on) is set.

Reimplemented from MCInitialParticles.

Definition at line 74 of file BeamParameters.h.

◆ setHER()

void setHER ( double  energy,
double  angleX,
double  angleY,
const std::vector< double > &  cov 
)

Set the HER FourVector and error matrix from beam energy, angle and covariance entries.

The vector for the covariance matrix can have either 0, 1, 3, 6 or 9 entries:

  • 0 entries means no error.
  • 1 entry will be treated as the variance of the beam energy.
  • 3 entries will be treated as uncorrelated variances for all variables (diagonal of the matrix)
  • 6 entries will be interpreted as the upper triangle of the covariance matrix with the element order being (0, 0), (0, 1), (0, 2), (1, 1), (1, 2), (2, 2)
  • 9 entries for the full covariance matrix in row-major order. The symmetry of the matrix will not be checked but just the lower triangle will be used.
Parameters
energybeam energy
angleXangle between beam direction and z axis (in the x-z plane). Negative values are treated as M_PI - abs(angle)
angleYangle between beam direction and z axis (in the y-z plane). Negative values are treated as M_PI - abs(angle)
coventries of the covariance matrix.

Definition at line 36 of file BeamParameters.cc.

◆ setLER()

void setLER ( double  energy,
double  angleX,
double  angleY,
const std::vector< double > &  cov 
)

Set the LER FourVector and error matrix from beam energy, angle and covariance entries.

The vector for the covariance matrix can have either 0, 1, 3, 6 or 9 entries:

  • 0 entries means no error.
  • 1 entry will be treated as the variance of the beam energy.
  • 3 entries will be treated as uncorrelated variances for all variables (diagonal of the matrix)
  • 6 entries will be interpreted as the upper triangle of the covariance matrix with the element order being (0, 0), (0, 1), (0, 2), (1, 1), (1, 2), (2, 2)
  • 9 entries for the full covariance matrix in row-major order. The symmetry of the matrix will not be checked but just the lower triangle will be used.
Parameters
energybeam energy
angleXangle between beam direction and z axis (in the x-z plane). Negative values are treated as M_PI - abs(angle)
angleYangle between beam direction and z axis (in the y-z plane). Negative values are treated as M_PI - abs(angle)
coventries of the covariance matrix.

Definition at line 29 of file BeamParameters.cc.

◆ setVertex()

void setVertex ( const TVector3 &  vertex,
const std::vector< double > &  cov 
)

Set the vertex position and error matrix.

The vector for the covariance matrix can have either 0, 1, 3, 6 or 9 entries:

  • 0 entries means no error.
  • 1 entry will be treated as the common variance for x, y and z
  • 3 entries will be treated as uncorrelated variances for all variables (diagonal of the matrix)
  • 6 entries will be interpreted as the upper triangle of the covariance matrix with the element order being (0, 0), (0, 1), (0, 2), (1, 1), (1, 2), (2, 2)
  • 9 entries for the full covariance matrix in row-major order. The symmetry of the matrix will not be checked but just the upper triangle will be used.
Parameters
vertexvertex position
coventries of the covariance matrix.

Definition at line 43 of file BeamParameters.cc.

Member Data Documentation

◆ m_CMSToLab

ROOT::Math::LorentzRotation* m_CMSToLab {nullptr}
mutableprivateinherited

transient

Boost from CMS into lab. (calculated on first use, not saved to file)

Definition at line 192 of file MCInitialParticles.h.

◆ m_invariantMass

double m_invariantMass {0.0}
mutableprivateinherited

transient

invariant mass of HER+LER (calculated on first use, not saved to file)

Definition at line 194 of file MCInitialParticles.h.

◆ m_labToCMS

ROOT::Math::LorentzRotation* m_labToCMS {nullptr}
mutableprivateinherited

Boost from Lab into CMS.

(calculated on first use, not saved to file)

Definition at line 190 of file MCInitialParticles.h.

◆ m_validFlag

bool m_validFlag = false
privateinherited

transient

Flag to check if a valid MCInitialParticles object was already generated and filled in an event.

Definition at line 196 of file MCInitialParticles.h.


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