Belle II Software  release-08-01-10
ReaderSAD Class Reference

Class to read files that have been created by SAD and store their content in a MCParticle graph. More...

#include <ReaderSAD.h>

Collaboration diagram for ReaderSAD:

Classes

struct  bendingElement
 Sensitive Element
More...
 
struct  straightElement
 Calculates the transformation matrix from local SAD to global geant4 space. More...
 

Public Types

enum  AcceleratorRings {
  c_HER = 0 ,
  c_LER = 1
}
 The both accelerator rings. More...
 

Public Member Functions

 BELLE2_DEFINE_EXCEPTION (SADCouldNotOpenFileError, "Could not open file %1% !")
 Exception is thrown if the SAD file could not be opened.
 
 BELLE2_DEFINE_EXCEPTION (SADEndOfFile, "End of the SAD file.")
 Exception is thrown if the end of the SAD file has been reached.
 
 ReaderSAD ()
 Constructor of the ReaderSAD class.
 
 ~ReaderSAD ()
 Destructor.
 
void initialize (TGeoHMatrix *transMatrix, double sRange, AcceleratorRings accRing, double readoutTime)
 Initializes the reader, sets the particle parameters and calculates important values. More...
 
void open (const std::string &filename)
 Opens a root file and prepares it for reading. More...
 
void setMomentumRes (double pxRes, double pyRes)
 Sets the resolution of the momentum for the real particles. More...
 
double getSADParticle (MCParticleGraph &graph)
 Reads one SAD particle from the file and creates one event per SAD particle. More...
 
bool getRealParticle (MCParticleGraph &graph)
 Reads one SAD particle from the file, calculates the number of real particles which are represented by the SAD particle and creates one event per real particle. More...
 
void addAllSADParticles (MCParticleGraph &graph)
 Reads all SAD particles from the file into the MCParticles collection which are inside the specified s range. More...
 
TGeoHMatrix SADtoGeant (ReaderSAD::AcceleratorRings accRing, double s)
 Transformation matrix. More...
 

Protected Attributes

TFile * m_file
 The input root file.
 
TTree * m_tree
 The input root tree.
 
TGeoHMatrix * m_transMatrix
 Transformation matrix from local SAD to global geant4 space.
 
double m_sRange
 The +- range for the s value for which particles are loaded.
 
AcceleratorRings m_accRing
 The accelerator ring from which the particles originate.
 
double m_pxRes
 The resolution for the x momentum component of the SAD real particle.
 
double m_pyRes
 The resolution for the y momentum component of the SAD real particle.
 
double m_SADToRealFactor
 The factor to calculate the number of real particles from a SAD particle.
 
double m_readoutTime
 The readout time.
 
unsigned int m_realPartNum
 The current number of the created real particles.
 
unsigned int m_realPartEntry
 The current number of the created real particles.
 
int m_readEntry
 The current number of the SAD entry that is read.
 
double m_lostX
 x at lost position [m].
 
double m_lostY
 y at lost position [m].
 
double m_lostS
 lost position [m] along ring. More...
 
double m_lostPx
 x momentum at lost position [m].
 
double m_lostPy
 y momentum at lost position [m].
 
double m_lostRate
 loss rate [Hz]>
 
double m_lostE
 energy at lost position [m].
 
double m_inputSAD_ssraw
 scattered position [m]
 
double m_inputSAD_sraw
 lost position [m
 
double m_inputSAD_ss
 scattered position (|s|<Ltot/2) [m]
 
double m_inputSAD_s
 lost position (|s|<Ltot/2) [m]
 
double m_inputSAD_Lss
 length of element in which scattered [m]
 
int m_inputSAD_nturn
 number of turns from scattered to lost
 
double m_inputSAD_x
 x at lost position [m].
 
double m_inputSAD_y
 y at lost position [m].
 
double m_inputSAD_px
 x momentum at lost position [m].
 
double m_inputSAD_py
 y momentum at lost position [m].
 
double m_inputSAD_xraw
 x at lost position [m] before matching onto beam pipe inner surface
 
double m_inputSAD_yraw
 y at lost position [m] before matching onto beam pipe inner surface
 
double m_inputSAD_r
 sqrt(x*x+y*y) [m]
 
double m_inputSAD_rr
 sqrt(x*x+y*y) [m] before matching onto beam pipe inner surface
 
double m_inputSAD_dp_over_p0
 dp_over_p0
 
double m_inputSAD_E
 energy at lost position [m].
 
double m_inputSAD_rate
 loss rate [Hz]
 
double m_inputSAD_watt
 loss wattage [W]
 

Private Member Functions

void convertParamsToSADUnits ()
 Convert the parameters from the SAD units to the basf2 units.
 
void addParticleToMCParticles (MCParticleGraph &graph, bool gaussSmearing=false)
 Adds the current particle described by the member variables to the MCParticles collection. More...
 
int calculateRealParticleNumber (double rate)
 Calculates the number of real particles for a SAD particle. More...
 

Detailed Description

Class to read files that have been created by SAD and store their content in a MCParticle graph.

The input data is stored in a root file and contains the particles together with their lost rate. The reader reads one particle from the file, calculates the number of 'real' particles and creates a new event for each of them.

Definition at line 35 of file ReaderSAD.h.

Member Enumeration Documentation

◆ AcceleratorRings

The both accelerator rings.

Enumerator
c_HER 

High Energy Ring (electrons)

c_LER 

Low Energy Ring (positrons)

Definition at line 46 of file ReaderSAD.h.

46  {
47  c_HER = 0,
48  c_LER = 1
49  };
@ c_HER
High Energy Ring (electrons)
Definition: ReaderSAD.h:47
@ c_LER
Low Energy Ring (positrons)
Definition: ReaderSAD.h:48

Member Function Documentation

◆ addAllSADParticles()

void addAllSADParticles ( MCParticleGraph graph)

Reads all SAD particles from the file into the MCParticles collection which are inside the specified s range.

Parameters
graphReference to the graph which should be filled with the information from the SAD file.

Definition at line 212 of file ReaderSAD.cc.

213 {
214  if (m_tree == NULL) {
215  B2ERROR("The SAD tree doesn't exist !");
216  return ;
217  }
218 
219  int nPart = m_tree->GetEntries();
220 
221  for (int iPart = 0; iPart < nPart; ++iPart) {
222  m_tree->GetEntry(iPart);
224  if (fabs(m_lostS) <= m_sRange) addParticleToMCParticles(graph);
225 
226  //std::cout << " in sad " << m_inputSAD_sraw << std::endl;
227 
228  }
229 }
double m_lostS
lost position [m] along ring.
Definition: ReaderSAD.h:138
TTree * m_tree
The input root tree.
Definition: ReaderSAD.h:121
void convertParamsToSADUnits()
Convert the parameters from the SAD units to the basf2 units.
Definition: ReaderSAD.cc:236
void addParticleToMCParticles(MCParticleGraph &graph, bool gaussSmearing=false)
Adds the current particle described by the member variables to the MCParticles collection.
Definition: ReaderSAD.cc:266
double m_sRange
The +- range for the s value for which particles are loaded.
Definition: ReaderSAD.h:124

◆ addParticleToMCParticles()

void addParticleToMCParticles ( MCParticleGraph graph,
bool  gaussSmearing = false 
)
private

Adds the current particle described by the member variables to the MCParticles collection.

Parameters
graphReference to the graph which should be filled with the information from the SAD file.
gaussSmearingIf set to true the particle momentum is smeared using a Gaussian.

Definition at line 266 of file ReaderSAD.cc.

◆ calculateRealParticleNumber()

int calculateRealParticleNumber ( double  rate)
private

Calculates the number of real particles for a SAD particle.

Parameters
rateThe loss rate of the SAD particle.
Returns
The number of real particles for the given loss rate.

Definition at line 403 of file ReaderSAD.cc.

◆ getRealParticle()

bool getRealParticle ( MCParticleGraph graph)

Reads one SAD particle from the file, calculates the number of real particles which are represented by the SAD particle and creates one event per real particle.

Parameters
graphReference to the graph which should be filled with the information from the SAD file.
Returns
True if the particle could be read.

Definition at line 165 of file ReaderSAD.cc.

◆ getSADParticle()

double getSADParticle ( MCParticleGraph graph)

Reads one SAD particle from the file and creates one event per SAD particle.

The loss rate of the SAD particle is stored in the weight attribute of the event meta info.

Parameters
graphReference to the graph which should be filled with the information from the SAD file.
Returns
The loss rate of the SAD particle which was read. Returns -1 if an error occurred.

Definition at line 122 of file ReaderSAD.cc.

◆ initialize()

void initialize ( TGeoHMatrix *  transMatrix,
double  sRange,
ReaderSAD::AcceleratorRings  accRing,
double  readoutTime 
)

Initializes the reader, sets the particle parameters and calculates important values.

Parameters
transMatrixPointer to the matrix which transforms the particles from the local SAD to the global geant4 coordinate system.
sRangeThe +- range for the s value for which particles are loaded.
accRingThe accelerator ring from which the particles originate.
readoutTimeThe readout time of the detector in [ns].

Definition at line 71 of file ReaderSAD.cc.

◆ open()

void open ( const std::string &  filename)

Opens a root file and prepares it for reading.

Parameters
filenameThe filename of the SAD root file which should be read.

Definition at line 80 of file ReaderSAD.cc.

◆ SADtoGeant()

TGeoHMatrix SADtoGeant ( ReaderSAD::AcceleratorRings  accRing,
double  s 
)

Transformation matrix.

Parameters
accRingThe accelerator ring from which the particles originate.
ss value.

Definition at line 416 of file ReaderSAD.cc.

◆ setMomentumRes()

void setMomentumRes ( double  pxRes,
double  pyRes 
)
inline

Sets the resolution of the momentum for the real particles.

Allows each real particle to be smeared according to the specified momentum resolution.

Parameters
pxResThe resolution for the x momentum component.
pyResThe resolution for the y momentum component.

Definition at line 82 of file ReaderSAD.h.

Member Data Documentation

◆ m_lostS

double m_lostS
protected

lost position [m] along ring.

0 is the IP. Range goes from -L/2 to L/2, where L is the ring circumference.

Definition at line 138 of file ReaderSAD.h.


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