Belle II Software  release-08-01-10
ReaderSAD.h
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #pragma once
10 
11 #include <framework/core/FrameworkExceptions.h>
12 
13 #include <mdst/dataobjects/MCParticleGraph.h>
14 
15 #include <TGeoMatrix.h>
16 #include <TFile.h>
17 #include <TTree.h>
18 
19 #include <string>
20 
21 
22 namespace Belle2 {
35  class ReaderSAD {
36 
37  public:
38 
39  //Define exceptions
41  BELLE2_DEFINE_EXCEPTION(SADCouldNotOpenFileError, "Could not open file %1% !");
43  BELLE2_DEFINE_EXCEPTION(SADEndOfFile, "End of the SAD file.");
44 
47  c_HER = 0,
48  c_LER = 1
49  };
50 
54  ReaderSAD();
55 
59  ~ReaderSAD();
60 
68  void initialize(TGeoHMatrix* transMatrix, double sRange, AcceleratorRings accRing, double readoutTime);
69 
74  void open(const std::string& filename);
75 
82  void setMomentumRes(double pxRes, double pyRes) { m_pxRes = pxRes; m_pyRes = pyRes; }
83 
91  double getSADParticle(MCParticleGraph& graph);
92 
100  bool getRealParticle(MCParticleGraph& graph);
101 
107  void addAllSADParticles(MCParticleGraph& graph);
108 
115  TGeoHMatrix SADtoGeant(ReaderSAD::AcceleratorRings accRing, double s);
116 
117 
118  protected:
119 
120  TFile* m_file;
121  TTree* m_tree;
123  TGeoHMatrix* m_transMatrix;
124  double m_sRange;
126  double m_pxRes;
127  double m_pyRes;
130  double m_readoutTime;
132  unsigned int m_realPartNum;
133  unsigned int m_realPartEntry;
136  double m_lostX;
137  double m_lostY;
138  double m_lostS;
139  double m_lostPx;
140  double m_lostPy;
141  double m_lostRate;
142  double m_lostE;
146  double m_inputSAD_ss;
147  double m_inputSAD_s;
148  double m_inputSAD_Lss;
150  double m_inputSAD_x;
151  double m_inputSAD_y;
152  double m_inputSAD_px;
153  double m_inputSAD_py;
156  double m_inputSAD_r;
157  double m_inputSAD_rr;
159  double m_inputSAD_E;
163  private:
164 
167 
173  void addParticleToMCParticles(MCParticleGraph& graph, bool gaussSmearing = false);
174 
180  int calculateRealParticleNumber(double rate);
181 
183 // TGeoHMatrix SADtoGeant(ReaderSAD::AcceleratorRings accRing, double s);
184 
186  double x0;
187  double z0;
188  double l;
189  double phi;
190  };
191 
193  struct bendingElement {
194  double rt;
195  double x0;
196  double z0;
197  double sphi;
198  double dphi;
199  };
200 
201  };
202 
204 }
205 
Class to build, validate and sort a particle decay chain.
Class to read files that have been created by SAD and store their content in a MCParticle graph.
Definition: ReaderSAD.h:35
double m_pyRes
The resolution for the y momentum component of the SAD real particle.
Definition: ReaderSAD.h:127
double m_lostPx
x momentum at lost position [m].
Definition: ReaderSAD.h:139
double getSADParticle(MCParticleGraph &graph)
Reads one SAD particle from the file and creates one event per SAD particle.
Definition: ReaderSAD.cc:122
double m_inputSAD_x
x at lost position [m].
Definition: ReaderSAD.h:150
double m_lostS
lost position [m] along ring.
Definition: ReaderSAD.h:138
TTree * m_tree
The input root tree.
Definition: ReaderSAD.h:121
double m_inputSAD_ssraw
scattered position [m]
Definition: ReaderSAD.h:144
~ReaderSAD()
Destructor.
Definition: ReaderSAD.cc:64
void addAllSADParticles(MCParticleGraph &graph)
Reads all SAD particles from the file into the MCParticles collection which are inside the specified ...
Definition: ReaderSAD.cc:212
AcceleratorRings m_accRing
The accelerator ring from which the particles originate.
Definition: ReaderSAD.h:125
double m_inputSAD_r
sqrt(x*x+y*y) [m]
Definition: ReaderSAD.h:156
BELLE2_DEFINE_EXCEPTION(SADCouldNotOpenFileError, "Could not open file %1% !")
Exception is thrown if the SAD file could not be opened.
double m_lostX
x at lost position [m].
Definition: ReaderSAD.h:136
int calculateRealParticleNumber(double rate)
Calculates the number of real particles for a SAD particle.
Definition: ReaderSAD.cc:403
TGeoHMatrix * m_transMatrix
Transformation matrix from local SAD to global geant4 space.
Definition: ReaderSAD.h:123
unsigned int m_realPartEntry
The current number of the created real particles.
Definition: ReaderSAD.h:133
void convertParamsToSADUnits()
Convert the parameters from the SAD units to the basf2 units.
Definition: ReaderSAD.cc:236
void open(const std::string &filename)
Opens a root file and prepares it for reading.
Definition: ReaderSAD.cc:80
double m_inputSAD_rr
sqrt(x*x+y*y) [m] before matching onto beam pipe inner surface
Definition: ReaderSAD.h:157
void initialize(TGeoHMatrix *transMatrix, double sRange, AcceleratorRings accRing, double readoutTime)
Initializes the reader, sets the particle parameters and calculates important values.
Definition: ReaderSAD.cc:71
double m_SADToRealFactor
The factor to calculate the number of real particles from a SAD particle.
Definition: ReaderSAD.h:129
ReaderSAD()
Constructor of the ReaderSAD class.
Definition: ReaderSAD.cc:28
double m_inputSAD_watt
loss wattage [W]
Definition: ReaderSAD.h:161
int m_readEntry
The current number of the SAD entry that is read.
Definition: ReaderSAD.h:134
TFile * m_file
The input root file.
Definition: ReaderSAD.h:120
void addParticleToMCParticles(MCParticleGraph &graph, bool gaussSmearing=false)
Adds the current particle described by the member variables to the MCParticles collection.
Definition: ReaderSAD.cc:266
TGeoHMatrix SADtoGeant(ReaderSAD::AcceleratorRings accRing, double s)
Transformation matrix.
Definition: ReaderSAD.cc:416
double m_lostPy
y momentum at lost position [m].
Definition: ReaderSAD.h:140
double m_inputSAD_E
energy at lost position [m].
Definition: ReaderSAD.h:159
double m_inputSAD_yraw
y at lost position [m] before matching onto beam pipe inner surface
Definition: ReaderSAD.h:155
double m_inputSAD_dp_over_p0
dp_over_p0
Definition: ReaderSAD.h:158
double m_inputSAD_xraw
x at lost position [m] before matching onto beam pipe inner surface
Definition: ReaderSAD.h:154
double m_inputSAD_y
y at lost position [m].
Definition: ReaderSAD.h:151
double m_inputSAD_sraw
lost position [m
Definition: ReaderSAD.h:145
int m_inputSAD_nturn
number of turns from scattered to lost
Definition: ReaderSAD.h:149
void setMomentumRes(double pxRes, double pyRes)
Sets the resolution of the momentum for the real particles.
Definition: ReaderSAD.h:82
double m_inputSAD_Lss
length of element in which scattered [m]
Definition: ReaderSAD.h:148
double m_lostY
y at lost position [m].
Definition: ReaderSAD.h:137
unsigned int m_realPartNum
The current number of the created real particles.
Definition: ReaderSAD.h:132
AcceleratorRings
The both accelerator rings.
Definition: ReaderSAD.h:46
@ c_HER
High Energy Ring (electrons)
Definition: ReaderSAD.h:47
@ c_LER
Low Energy Ring (positrons)
Definition: ReaderSAD.h:48
bool getRealParticle(MCParticleGraph &graph)
Reads one SAD particle from the file, calculates the number of real particles which are represented b...
Definition: ReaderSAD.cc:165
double m_pxRes
The resolution for the x momentum component of the SAD real particle.
Definition: ReaderSAD.h:126
double m_readoutTime
The readout time.
Definition: ReaderSAD.h:130
double m_inputSAD_rate
loss rate [Hz]
Definition: ReaderSAD.h:160
double m_inputSAD_px
x momentum at lost position [m].
Definition: ReaderSAD.h:152
double m_inputSAD_py
y momentum at lost position [m].
Definition: ReaderSAD.h:153
double m_lostRate
loss rate [Hz]>
Definition: ReaderSAD.h:141
double m_inputSAD_s
lost position (|s|<Ltot/2) [m]
Definition: ReaderSAD.h:147
double m_inputSAD_ss
scattered position (|s|<Ltot/2) [m]
Definition: ReaderSAD.h:146
BELLE2_DEFINE_EXCEPTION(SADEndOfFile, "End of the SAD file.")
Exception is thrown if the end of the SAD file has been reached.
double m_sRange
The +- range for the s value for which particles are loaded.
Definition: ReaderSAD.h:124
double m_lostE
energy at lost position [m].
Definition: ReaderSAD.h:142
Abstract base class for different kinds of events.
double z0
Initial position in Z
Definition: ReaderSAD.h:196
double sphi
Bending Angle.
Definition: ReaderSAD.h:197
double dphi
Bending Angle for torus in phi.
Definition: ReaderSAD.h:198
double x0
Initial position in X
Definition: ReaderSAD.h:195
Calculates the transformation matrix from local SAD to global geant4 space.
Definition: ReaderSAD.h:185
double z0
Initial position in Z.
Definition: ReaderSAD.h:187
double x0
Initial position in X.
Definition: ReaderSAD.h:186