Belle II Software development
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:
MCInitialParticles

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).
 
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.
 
void setCovVertex (const TMatrixDSym &cov)
 Set the covariance matrix of the vertex position.
 
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.
 
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.
 
void setVertex (const ROOT::Math::XYZVector &vertex, const std::vector< double > &cov)
 Set the vertex position and error matrix.
 
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 ROOT::Math::XYZVector &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 ROOT::Math::XYZVector &vertex)
 Set the initial event values, i.e.
 
void setByLorentzTransformation (double Ecms, double bX, double bY, double bZ, double angleXZ, double angleYZ, const ROOT::Math::XYZVector &vertex)
 Initialize the event values from CMS energy and parameters of the Lorentz transformation between LAB and CMS.
 
void setTime (double time)
 Set collision 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 ROOT::Math::XYZVector & getVertex () const
 Get the position of the collision.
 
double getTime () const
 Get collision time.
 
double getEnergy () const
 Get the actual collision energy (in lab system)
 
double getMass () const
 Get the invariant mass of the collision (= energy in CMS)
 
ROOT::Math::XYZVector getBoostVector () const
 Get the boost vector (velocity of system produced in the collision)
 
ROOT::Math::PxPyPzEVector getBoostedHER () const
 Get the 4-vector of electron beam in CM system obtained by pure boost.
 
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.
 

Static Public Member Functions

static ROOT::Math::PxPyPzEVector getFourVector (double energy, double angleX, double angleY, bool isHER)
 Calculate FourVector of a beam from energy and angles in xz and yz planes.
 
static ROOT::Math::LorentzRotation cmsToLab (double bX, double bY, double bZ, double angleXZ, double angleYZ)
 Return the LorentzRotation from CMS to LAB based on the following parameters.
 

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 Lorentz transformations LAB->CMS & CMS->LAB if necessary.
 
void resetBoost ()
 Reset cached transformations after changing parameters.
 
 ClassDef (MCInitialParticles, 5)
 ROOT Dictionary.
 

Static Private Member Functions

static void setCovMatrix (Double32_t *member, const std::vector< double > &cov, bool common)
 Set covariance matrix from vector of entries.
 
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.
 
ROOT::Math::XYZVector m_vertex
 collision position
 
double m_time = 0
 collision time
 
ROOT::Math::LorentzRotation * m_labToCMS {nullptr}
 Boost from Lab into CMS.
 
ROOT::Math::LorentzRotation * m_CMSToLab {nullptr}
 transient
 
double m_invariantMass {0.0}
 transient
 
ROOT::Math::PxPyPzEVector m_boostedHER
 transient
 
bool m_validFlag = false
 transient
 

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.

36 {
38 c_generateCMS = 1 << 0,
40 c_smearBeamEnergy = 1 << 1,
42 c_smearBeamDirection = 1 << 2,
46 c_smearVertex = 1 << 3,
49 };
@ c_smearBeam
smear the full beam momentum (energy and direction)
@ 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)

Constructor & Destructor Documentation

◆ BeamParameters() [1/2]

BeamParameters ( )
inline

Default constructor.

Definition at line 34 of file BeamParameters.h.

35 {
36 // 14 means all the smearings turned on
38 }
Double32_t m_covLER[6]
Covariance matrix of the low energy beam at the IP.
Double32_t m_covHER[6]
Covariance matrix of the high energy beam at the IP.
Double32_t m_covVertex[6]
Covariance matrix of the vertex position.
int m_generationFlags
Flags to be used when generating events.
MCInitialParticles()
Default constructor.

◆ BeamParameters() [2/2]

BeamParameters ( const BeamParameters b)
inline

Copy constructor.

Definition at line 41 of file BeamParameters.h.

42 {
43 std::copy_n(b.m_covHER, 6, m_covHER);
44 std::copy_n(b.m_covLER, 6, m_covLER);
45 std::copy_n(b.m_covVertex, 6, m_covVertex);
46 }

Member Function Documentation

◆ cmsToLab()

ROOT::Math::LorentzRotation cmsToLab ( double  bX,
double  bY,
double  bZ,
double  angleXZ,
double  angleYZ 
)
staticinherited

Return the LorentzRotation from CMS to LAB based on the following parameters.

Parameters
bXx-component of the boost vector, i.e. of (pHER + pLER) / (eHER + eLER), where pHER & pLER are momentum 3-vectors
bYy-component of the boost vector, i.e. of (pHER + pLER) / (eHER + eLER), where pHER & pLER are momentum 3-vectors
bZz-component of the boost vector, i.e. of (pHER + pLER) / (eHER + eLER), where pHER & pLER are momentum 3-vectors
angleXZangle in the XZ plane of the collision axis in the CM system obtained by pure boost
angleYZangle in the YZ plane of the collision axis in the CM system obtained by pure boost

Definition at line 31 of file MCInitialParticles.cc.

32{
33 //boost to CM frame
34 ROOT::Math::LorentzRotation boost(ROOT::Math::Boost(-1.*ROOT::Math::XYZVector(bX, bY, bZ)));
35
36
37 //rotation such that the collision axis is aligned with z-axis
38 ROOT::Math::XYZVector zaxis(0., 0., 1.); //target collision axis
39
40 double tanAngleXZ = tan(angleXZ);
41 double tanAngleYZ = tan(angleYZ);
42 double Norm = 1 / sqrt(1 + pow(tanAngleXZ, 2) + pow(tanAngleYZ, 2));
43 ROOT::Math::XYZVector electronCMS(Norm * tanAngleXZ, Norm * tanAngleYZ, Norm); //current collision axis
44
45 ROOT::Math::XYZVector rotAxis = zaxis.Cross(electronCMS);
46 double rotangle = asin(rotAxis.R());
47
48 ROOT::Math::LorentzRotation rotation(ROOT::Math::AxisAngle(rotAxis, -rotangle));
49
50
51 ROOT::Math::LorentzRotation trans = rotation * boost;
52 ROOT::Math::LorentzRotation transI = trans.Inverse();
53 return transI;
54}
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28

◆ getBoostedHER()

ROOT::Math::PxPyPzEVector getBoostedHER ( ) const
inlineinherited

Get the 4-vector of electron beam in CM system obtained by pure boost.

Definition at line 185 of file MCInitialParticles.h.

185{ calculateBoost(); return m_boostedHER; }
ROOT::Math::PxPyPzEVector m_boostedHER
transient
void calculateBoost() const
Calculate the Lorentz transformations LAB->CMS & CMS->LAB if necessary.

◆ getBoostVector()

ROOT::Math::XYZVector getBoostVector ( ) const
inlineinherited

Get the boost vector (velocity of system produced in the collision)

Definition at line 182 of file MCInitialParticles.h.

182{ calculateBoost(); return -(m_her + m_ler).BoostToCM(); }
ROOT::Math::PxPyPzEVector m_ler
LER 4vector.
ROOT::Math::PxPyPzEVector m_her
HER 4vector.

◆ getCMSToLab()

const ROOT::Math::LorentzRotation & getCMSToLab ( ) const
inlineinherited

Return the LorentzRotation to convert from CMS to lab frame.

Definition at line 194 of file MCInitialParticles.h.

195 {
196 calculateBoost(); return *m_CMSToLab;
197 }
ROOT::Math::LorentzRotation * m_CMSToLab
transient

◆ getCovHER()

TMatrixDSym getCovHER ( ) const
inline

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.

Definition at line 157 of file BeamParameters.h.

157{ return getCovMatrix(m_covHER); }
static TMatrixDSym getCovMatrix(const Double32_t *member)
Obtain covariance matrix from a given float array.

◆ getCovLER()

TMatrixDSym getCovLER ( ) const
inline

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.

Definition at line 162 of file BeamParameters.h.

162{ return getCovMatrix(m_covLER); }

◆ getCovVertex()

TMatrixDSym getCovVertex ( ) const
inline

Get the covariance matrix of the vertex position.

Definition at line 165 of file BeamParameters.h.

165{ return getCovMatrix(m_covVertex); }

◆ getEnergy()

double getEnergy ( ) const
inlineinherited

Get the actual collision energy (in lab system)

Definition at line 176 of file MCInitialParticles.h.

176{ return (m_her + m_ler).E(); }

◆ getFourVector()

static ROOT::Math::PxPyPzEVector getFourVector ( double  energy,
double  angleX,
double  angleY,
bool  isHER 
)
static

Calculate FourVector of a beam from energy and angles in xz and yz planes.

if isHER=true, the angles are measured wrt +p if isHER=false, the angles are measured wrt -p, as is the standard convention for LER

Parameters
energybeam energy
angleXhorizontal angle wrt z-axis, i.e. angle measured in xz plane
angleYvertical angle wrt z-axis, i.e. angle measured in yz plane
isHERisHER=true for HER, isHER=false for LER

◆ getGenerationFlags()

int getGenerationFlags ( ) const
inlineinherited

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

Definition at line 203 of file MCInitialParticles.h.

203{ return m_generationFlags; }

◆ 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.

14{
15 std::string flags = "";
16 const std::vector<std::pair<int, std::string>> flagvalues{
17 {c_generateCMS, "generateCMS"},
18 {c_smearBeamEnergy, "smearBeamEnergy"},
19 {c_smearBeamDirection, "smearBeamDirection"},
20 {c_smearVertex, "smearVertex"}
21 };
22 for (auto& i : flagvalues) {
23 if (hasGenerationFlags(i.first)) {
24 if (flags.size() > 0) flags += separator;
25 flags += i.second;
26 }
27 }
28 return flags;
29}
bool hasGenerationFlags(int flags) const
Check if a certain set of EGenerationFlags is set.

◆ getHER()

const ROOT::Math::PxPyPzEVector & getHER ( ) const
inlineinherited

Get 4vector of the high energy beam.

Definition at line 164 of file MCInitialParticles.h.

164{ return m_her; }

◆ getLabToCMS()

const ROOT::Math::LorentzRotation & getLabToCMS ( ) const
inlineinherited

Return the LorentzRotation to convert from lab to CMS frame.

Definition at line 188 of file MCInitialParticles.h.

189 {
190 calculateBoost(); return *m_labToCMS;
191 }
ROOT::Math::LorentzRotation * m_labToCMS
Boost from Lab into CMS.

◆ getLER()

const ROOT::Math::PxPyPzEVector & getLER ( ) const
inlineinherited

Get 4vector of the low energy beam.

Definition at line 167 of file MCInitialParticles.h.

167{ return m_ler; }

◆ getMass()

double getMass ( ) const
inlineinherited

Get the invariant mass of the collision (= energy in CMS)

Definition at line 179 of file MCInitialParticles.h.

179{ calculateBoost(); return m_invariantMass; }

◆ getTime()

double getTime ( ) const
inlineinherited

Get collision time.

Definition at line 173 of file MCInitialParticles.h.

173{return m_time;}
double m_time
collision time

◆ getValidFlag()

bool getValidFlag ( ) const
inlineinherited

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

Definition at line 200 of file MCInitialParticles.h.

200{ return m_validFlag; }

◆ getVertex()

const ROOT::Math::XYZVector & getVertex ( ) const
inlineinherited

Get the position of the collision.

Definition at line 170 of file MCInitialParticles.h.

170{ return m_vertex; }
ROOT::Math::XYZVector m_vertex
collision position

◆ hasGenerationFlags()

bool hasGenerationFlags ( int  flags) const
inlineinherited

Check if a certain set of EGenerationFlags is set.

Definition at line 206 of file MCInitialParticles.h.

206{ return (m_generationFlags & flags) == flags; }

◆ operator=()

BeamParameters & operator= ( const BeamParameters b)
inline

Assignment operator.

Definition at line 49 of file BeamParameters.h.

50 {
52 std::copy_n(b.m_covHER, 6, m_covHER);
53 std::copy_n(b.m_covLER, 6, m_covLER);
54 std::copy_n(b.m_covVertex, 6, m_covVertex);
55 return *this;
56 }
MCInitialParticles & operator=(const MCInitialParticles &b)
Assignment operator.

◆ operator==() [1/2]

bool operator== ( const BeamParameters b) const
inline

Equality operator.

Definition at line 59 of file BeamParameters.h.

60 {
61 // since we only save the covariance matrices with float precision we
62 // need to also do the comparison with float precision.
63 auto floatcmp = [](double dbl_a, double dbl_b) { return (float)dbl_a == (float)dbl_b; };
65 std::equal(m_covHER, m_covHER + 6, b.m_covHER, floatcmp) &&
66 std::equal(m_covLER, m_covLER + 6, b.m_covLER, floatcmp) &&
67 std::equal(m_covVertex, m_covVertex + 6, b.m_covVertex, floatcmp);
68 }
bool operator==(const MCInitialParticles &b) const
Equality operator.

◆ operator==() [2/2]

bool operator== ( const MCInitialParticles b) const
inlineinherited

Equality operator.

Definition at line 68 of file MCInitialParticles.h.

69 {
70 // FIXME: sin(x) returns slightly different values on
71 // different platforms in some cases so we cannot just do an equality
72 // comparison. We need to do this more elegantly, this is just for
73 // testing if it solves all problems
74#if defined(MCP_DBL_CMP) || defined(MCP_VEC3_CMP) || defined(MCP_VEC4_CMP)
75#error Macro already defined, cannot continue
76#endif
77#define MCP_DBL_CMP(a,b,x) ((a.X()==b.X())||(std::abs(a.X()-b.X())<1e-10))
78#define MCP_VEC3_CMP(a,b) (MCP_DBL_CMP(a,b,X) && MCP_DBL_CMP(a,b,Y) && MCP_DBL_CMP(a,b,Z))
79#define MCP_VEC4_CMP(a,b) (MCP_VEC3_CMP(a,b) && MCP_DBL_CMP(a,b,E))
80 return MCP_VEC4_CMP(m_her, b.m_her) && MCP_VEC4_CMP(m_ler, b.m_ler) && MCP_VEC3_CMP(m_vertex, b.m_vertex)
81 && m_validFlag == b.m_validFlag &&
82 m_generationFlags == b.m_generationFlags;
83#undef MCP_DBL_CMP
84#undef MCP_VEC3_CMP
85#undef MCP_VEC4_CMP
86 }

◆ set()

void set ( const ROOT::Math::PxPyPzEVector &  her,
const ROOT::Math::PxPyPzEVector &  ler,
const ROOT::Math::XYZVector &  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.

98 {
99 m_her = her;
100 m_ler = ler;
101 m_vertex = vertex;
102 m_validFlag = true;
103 resetBoost();
104 }
void resetBoost()
Reset cached transformations after changing parameters.

◆ setByLorentzTransformation()

void setByLorentzTransformation ( double  Ecms,
double  bX,
double  bY,
double  bZ,
double  angleXZ,
double  angleYZ,
const ROOT::Math::XYZVector &  vertex 
)
inlineinherited

Initialize the event values from CMS energy and parameters of the Lorentz transformation between LAB and CMS.

In addition the vertex is also initialized.

Parameters
Ecmscentre-of-mass energy of the collision
bXx-component of the boost vector, i.e. of (pHER + pLER) / (eHER + eLER), where pHER & pLER are momentum 3-vectors
bYy-component of the boost vector, i.e. of (pHER + pLER) / (eHER + eLER), where pHER & pLER are momentum 3-vectors
bZz-component of the boost vector, i.e. of (pHER + pLER) / (eHER + eLER), where pHER & pLER are momentum 3-vectors
angleXZangle in the XZ plane of the collision axis in the CM system obtained by pure boost
angleYZangle in the YZ plane of the collision axis in the CM system obtained by pure boost
vertexposition of the actual collision vertex

Definition at line 116 of file MCInitialParticles.h.

118 {
119 if (m_labToCMS) delete m_labToCMS;
120 if (m_CMSToLab) delete m_CMSToLab;
121
122 m_invariantMass = Ecms;
123 m_CMSToLab = new ROOT::Math::LorentzRotation();
124 m_labToCMS = new ROOT::Math::LorentzRotation();
125 *m_CMSToLab = cmsToLab(bX, bY, bZ, angleXZ, angleYZ);
126 *m_labToCMS = m_CMSToLab->Inverse();
127
128 const double me = Const::electron.getMass();
129 double p = sqrt(Ecms * Ecms / 4 - me * me);
130 m_her = (*m_CMSToLab) * ROOT::Math::PxPyPzEVector(0.0, 0.0, p, Ecms / 2);
131 m_ler = (*m_CMSToLab) * ROOT::Math::PxPyPzEVector(0.0, 0.0, -p, Ecms / 2);
132
133 m_vertex = vertex;
134 m_validFlag = true;
135 }
double getMass() const
Particle mass.
Definition: UnitConst.cc:356
static const ChargedStable electron
electron particle
Definition: Const.h:659
static ROOT::Math::LorentzRotation cmsToLab(double bX, double bY, double bZ, double angleXZ, double angleYZ)
Return the LorentzRotation from CMS to LAB based on the following parameters.
static const double me
electron mass
Definition: beamHelpers.h:177

◆ setCovHER()

void setCovHER ( const TMatrixDSym &  cov)
inline

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.

Definition at line 91 of file BeamParameters.h.

91{ setCovMatrix(m_covHER, cov); }
static void setCovMatrix(Double32_t *member, const std::vector< double > &cov, bool common)
Set covariance matrix from vector of entries.

◆ 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); }

◆ setCovMatrix()

static 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

◆ 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.

100{ setCovMatrix(m_covVertex, cov); }

◆ 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.

75 {
76 if (flags != 14)
77 B2WARNING(R"RAW(The generation flags are not set to the default value.
78
79 The default value for the generation flags is 14: a different value means
80 that some smearings (beam energy, vertex position, etc.) may be turned off.
81
82 This is fine only for local tests or for special MC productions, but not
83 for default MC productions.)RAW");
84 m_generationFlags = flags;
85 }

◆ setHER() [1/2]

void setHER ( const ROOT::Math::PxPyPzEVector &  her)
inline

Set the High Energy Beam 4-momentum.

Definition at line 138 of file MCInitialParticles.h.

139 {
140 m_her = her;
141 resetBoost();
142 }

◆ setHER() [2/2]

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.

◆ setLER() [1/2]

void setLER ( const ROOT::Math::PxPyPzEVector &  ler)
inline

Set the Low Energy Beam 4-momentum.

Definition at line 145 of file MCInitialParticles.h.

146 {
147 m_ler = ler;
148 resetBoost();
149 }

◆ setLER() [2/2]

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.

◆ setTime()

void setTime ( double  time)
inlineinherited

Set collision time.

Definition at line 158 of file MCInitialParticles.h.

158{m_time = time;}

◆ setVertex() [1/2]

void setVertex ( const ROOT::Math::XYZVector &  vertex)
inline

Set the vertex position.

Definition at line 152 of file MCInitialParticles.h.

153 {
154 m_vertex = vertex;
155 }

◆ setVertex() [2/2]

void setVertex ( const ROOT::Math::XYZVector &  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.

Member Data Documentation

◆ m_boostedHER

ROOT::Math::PxPyPzEVector m_boostedHER
mutableprivateinherited

transient

HER 4-momentum in CM frame obtained by pure boost (calculated on first use, not saved to file)

Definition at line 243 of file MCInitialParticles.h.

◆ 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 239 of file MCInitialParticles.h.

◆ m_covHER

Double32_t m_covHER[6]
private

Covariance matrix of the high energy beam at the IP.

Definition at line 207 of file BeamParameters.h.

◆ m_covLER

Double32_t m_covLER[6]
private

Covariance matrix of the low energy beam at the IP.

Definition at line 210 of file BeamParameters.h.

◆ m_covVertex

Double32_t m_covVertex[6]
private

Covariance matrix of the vertex position.

Definition at line 213 of file BeamParameters.h.

◆ m_generationFlags

int m_generationFlags {0}
protectedinherited

Flags to be used when generating events.

Definition at line 250 of file MCInitialParticles.h.

◆ m_her

ROOT::Math::PxPyPzEVector m_her
privateinherited

HER 4vector.

Definition at line 229 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 241 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 237 of file MCInitialParticles.h.

◆ m_ler

ROOT::Math::PxPyPzEVector m_ler
privateinherited

LER 4vector.

Definition at line 231 of file MCInitialParticles.h.

◆ m_time

double m_time = 0
privateinherited

collision time

Definition at line 235 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 245 of file MCInitialParticles.h.

◆ m_vertex

ROOT::Math::XYZVector m_vertex
privateinherited

collision position

Definition at line 233 of file MCInitialParticles.h.


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