Belle II Software development
BeamBkgGeneratorModule.cc
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// Own header.
10#include <background/modules/BeamBkgGenerator/BeamBkgGeneratorModule.h>
11
12// framework - DataStore
13#include <framework/datastore/StoreArray.h>
14#include <framework/datastore/StoreObjPtr.h>
15
16// framework aux
17#include <framework/gearbox/Unit.h>
18#include <framework/gearbox/Const.h>
19#include <framework/logging/Logger.h>
20#include <framework/gearbox/GearDir.h>
21
22// data objects
23#include <framework/dataobjects/EventMetaData.h>
24#include <mdst/dataobjects/MCParticle.h>
25#include <generators/SAD/dataobjects/SADMetaHit.h>
26
27// random generator
28#include <TRandom.h>
29
30// coordinates translation
31#include <TGeoMatrix.h>
32#include <generators/SAD/ReaderSAD.h>
33
34using namespace std;
35using namespace Belle2;
36
37//-----------------------------------------------------------------
38// Register module
39//-----------------------------------------------------------------
40
41REG_MODULE(BeamBkgGenerator);
42
43//-----------------------------------------------------------------
44// Implementation
45//-----------------------------------------------------------------
46
48{
49 // set module description
50 setDescription("Beam background generator based on SAD files. "
51 "The generator picks up particles from the SAD file randomly "
52 "according to their rates. "
53 "Number of events is determined from 'realTime' and overall rate, and "
54 "the generator terminates the execution when this number is reached.");
55
56 // Add parameters
57 addParam("fileName", m_fileName, "name of the SAD file converted to root");
58 addParam("treeName", m_treeName, "name of the TTree in the SAD file", string("sad"));
59 addParam("ringName", m_ringName, "name of the superKEKB ring (LER or HER)");
60 addParam("realTime", m_realTime,
61 "equivalent superKEKB running time to generate sample [ns].");
62}
63
65{
66}
67
69{
70 // register MCParticles
71
72 StoreArray<MCParticle> mcParticles;
73 mcParticles.registerInDataStore();
74
76 sadHits.registerInDataStore();
77
78
79 // check steering parameters
80
81 if (m_ringName != "LER" and m_ringName != "HER") {
82 B2ERROR("ringName can only be LER or HER");
83 }
84 if (m_realTime <= 0) {
85 B2ERROR("realTime must be positive");
86 }
87
88 // open SAD file and get TTree
89
90 m_file = TFile::Open(m_fileName.c_str());
91 if (!m_file) {
92 B2ERROR(m_fileName << ": can't open file");
93 return;
94 }
95 m_tree = (TTree*) m_file->Get(m_treeName.c_str());
96 if (!m_tree) {
97 B2ERROR(m_treeName << ": no such TTree in the SAD file");
98 return;
99 }
100
101 m_tree->SetBranchAddress("s", &m_sad.s);
102 m_tree->SetBranchAddress("x", &m_sad.x);
103 m_tree->SetBranchAddress("px", &m_sad.px);
104 m_tree->SetBranchAddress("y", &m_sad.y);
105 m_tree->SetBranchAddress("py", &m_sad.py);
106 m_tree->SetBranchAddress("E", &m_sad.E);
107 m_tree->SetBranchAddress("rate", &m_sad.rate);
108
109 // for SADMetaHit
110 m_tree->SetBranchAddress("ss", &m_sad.ss);
111 m_tree->SetBranchAddress("sraw", &m_sad.sraw);
112 m_tree->SetBranchAddress("ssraw", &m_sad.ssraw);
113 m_tree->SetBranchAddress("nturn", &m_sad.nturn);
114 m_tree->SetBranchAddress("xraw", &m_sad.xraw);
115 m_tree->SetBranchAddress("yraw", &m_sad.yraw);
116 m_tree->SetBranchAddress("r", &m_sad.r);
117 m_tree->SetBranchAddress("rr", &m_sad.rr);
118 m_tree->SetBranchAddress("dp_over_p0", &m_sad.dp_over_p0);
119 m_tree->SetBranchAddress("watt", &m_sad.watt);
120
121 int numEntries = m_tree->GetEntries();
122 if (numEntries <= 0) {
123 B2ERROR("SAD tree is empty");
124 return;
125 }
126
127 // construct vector of cumulative rates and determine number of events to generate
128
129 m_rates.push_back(0);
130 for (int i = 0; i < numEntries; i++) {
131 m_tree->GetEntry(i);
132 m_rates.push_back(m_sad.rate + m_rates.back());
133 m_counters.push_back(0);
134 }
135
136 m_numEvents = gRandom->Poisson(m_rates.back() * m_realTime / Unit::s);
137
138 // set rotation from SAD to Belle II
139
140 if (m_ringName == "LER") {
141 GearDir ring("/Detector/SuperKEKB/LER/");
142 m_rotation.SetAngle(ring.getAngle("angle"));
143 m_ring = 2;
144 } else {
145 GearDir ring("/Detector/SuperKEKB/HER/");
146 m_rotation.SetAngle(ring.getAngle("angle"));
147 m_ring = 1;
148 }
149
150 m_sectionOrdering.insert(m_sectionOrdering.end(), {1, 12, 11, 10, 9, 8, 7, 6, 5, 4, 3, 2});
151
152 B2RESULT("BG rate: " << m_rates.back() / 1e6 << " MHz, events to generate: "
153 << m_numEvents);
154
155}
156
158{
159}
160
162{
163
164 // check event counter
165
167 StoreObjPtr<EventMetaData> eventMetaData;
168 eventMetaData->setEndOfData(); // stop event processing
169 return;
170 }
172
173 // get SAD entry
174
175 int i = generateEntry();
176 m_tree->GetEntry(i);
177 m_counters[i]++;
178
179 // make SADMetaHit
180 int ring_section = -1;
181 double ssraw = m_sad.ss;
182 if (m_sad.ss < 0) ssraw += 3016.;
183 int section = (int)(ssraw * 12 / 3016.);
184 if ((unsigned)section < m_sectionOrdering.size()) ring_section = m_sectionOrdering[section];
185
186 StoreArray<SADMetaHit> SADMetaHits;
188 0., m_sad.nturn,
191 m_sad.watt, m_ring, ring_section));
192
193
194 // transform to Belle II (flip sign of y and s, rotate)
195
196 ROOT::Math::XYZVector position(m_sad.x * Unit::m, -m_sad.y * Unit::m, -m_sad.s * Unit::m);
197 position = m_rotation * position;
198
199 double pz = sqrt(m_sad.E * m_sad.E - m_sad.px * m_sad.px - m_sad.py * m_sad.py);
200 if (m_ringName == "LER") pz = -pz;
201 ROOT::Math::XYZVector momentum(m_sad.px, -m_sad.py, pz);
202 momentum = m_rotation * momentum;
203
204 // PDG code and mass
205
206 int pdgCode = Const::electron.getPDGCode();
207 double mass = Const::electron.getMass();
208 if (m_ringName == "LER") pdgCode = -pdgCode;
209
210 // append SAD particle to MCParticles
211
212 StoreArray<MCParticle> MCParticles;
213 MCParticle* part = MCParticles.appendNew();
214 part->setPDG(pdgCode);
215 part->setMass(mass);
216 part->setProductionVertex(position);
217 part->setProductionTime(0);
218 part->setMomentum(momentum);
219 part->setEnergy(sqrt(momentum.Mag2() + mass * mass));
220 part->setValidVertex(true);
223
224 // FarBeamLine region transformation
225 if (abs(m_sad.s * Unit::m) > 4.0 * Unit::m) {
226 // initial coordinates in SAD space
227 double particlePosSADfar[] = {m_sad.x* Unit::m, -m_sad.y* Unit::m, 0.0 * Unit::m};
228 double particleMomSADfar[] = {m_sad.px* Unit::GeV, -m_sad.py* Unit::GeV, pz* Unit::GeV};
229 // final coordinates in Geant4 space
230 double particlePosGeant4[] = {0.0, 0.0, 0.0};
231 double particleMomGeant4[] = {0.0, 0.0, 0.0};
232
233 // create a transformation matrix for a given ring
234 TGeoHMatrix transMatrix;
235 if (m_ringName == "LER") {
237 } else {
239 }
240
241 // calculate a new set of coordinates in Geant4 space
242 transMatrix.LocalToMaster(particlePosSADfar, particlePosGeant4); // position
243 transMatrix.LocalToMasterVect(particleMomSADfar, particleMomGeant4); // momentum
244 // apply a new set of coordinates
245 part->setMomentum(ROOT::Math::XYZVector(particleMomGeant4[0], particleMomGeant4[1], particleMomGeant4[2]));
246 part->setProductionVertex(ROOT::Math::XYZVector(particlePosGeant4[0], particlePosGeant4[1], particlePosGeant4[2]));
247 }
248}
249
251{
252}
253
255{
256 // close SAD file
257 m_file->Close();
258
259 B2RESULT("BG rate: " << m_rates.back() / 1e6 << " MHz, equivalent time: "
260 << m_realTime / Unit::us << " us, "
261 << m_eventCounter << " events generated");
263 B2ERROR("Number of generated events does not match the equivalent running time: "
264 << m_numEvents << " events needed, but "
265 << m_eventCounter << " generated");
266
267 int imax = 0;
268 double average = 0;
269 for (unsigned i = 0; i < m_counters.size(); i++) {
270 if (m_counters[i] > m_counters[imax]) imax = i;
271 average += m_counters[i];
272 }
273 average /= m_counters.size();
274 B2RESULT("SAD particle usage: on average " << average << " times, max "
275 << m_counters[imax] << " times (entry " << imax << ")");
276
277}
278
279
281{
282 double rate = gRandom->Uniform(m_rates.back());
283 int i1 = 0;
284 int i2 = m_rates.size() - 1;
285 while (i2 - i1 > 1) {
286 int i = (i1 + i2) / 2;
287 if (rate > m_rates[i]) {
288 i1 = i;
289 } else {
290 i2 = i;
291 }
292 }
293 return i1;
294}
std::vector< int > m_counters
counters: how many times SAD particles are used
virtual void initialize() override
Initialize the Module.
virtual ~BeamBkgGeneratorModule()
Destructor.
virtual void event() override
Event processor.
virtual void endRun() override
End-of-run action.
ReaderSAD m_readerSAD
the transformation from SAD to Belle II system for the far beamline
std::vector< double > m_rates
cumulative rates of SAD particles [Hz]
virtual void terminate() override
Termination action.
ROOT::Math::RotationY m_rotation
rotation from SAD to Belle II frame
std::string m_fileName
name of the SAD file converted to root
int generateEntry() const
Pick up particle randomly from the SAD file according to its rate.
int m_ring
ring number, 1-HER, 2-LER
int m_numEvents
number of events to generate
double m_realTime
equivalent SuperKEKB running time in [ns]
virtual void beginRun() override
Called when entering a new run.
std::vector< int > m_sectionOrdering
superKEKB section ordering
std::string m_ringName
name of the superKEKB ring (LER or HER)
std::string m_treeName
name of the TTree in the SAD file
int getPDGCode() const
PDG code.
Definition: Const.h:473
double getMass() const
Particle mass.
Definition: UnitConst.cc:353
static const ChargedStable electron
electron particle
Definition: Const.h:659
GearDir is the basic class used for accessing the parameter store.
Definition: GearDir.h:31
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
@ c_PrimaryParticle
bit 0: Particle is primary particle.
Definition: MCParticle.h:47
@ c_StableInGenerator
bit 1: Particle is stable, i.e., not decaying in the generator.
Definition: MCParticle.h:49
void setMass(float mass)
Set particle mass.
Definition: MCParticle.h:355
void addStatus(unsigned short int bitmask)
Add bitmask to current status.
Definition: MCParticle.h:342
void setEnergy(float energy)
Set energy.
Definition: MCParticle.h:361
void setValidVertex(bool valid)
Set indication whether vertex and time information is valid or just default.
Definition: MCParticle.h:367
void setProductionVertex(const ROOT::Math::XYZVector &vertex)
Set production vertex position.
Definition: MCParticle.h:385
void setPDG(int pdg)
Set PDG code of the particle.
Definition: MCParticle.h:324
void setMomentum(const ROOT::Math::XYZVector &momentum)
Set particle momentum.
Definition: MCParticle.h:406
void setStatus(unsigned short int status)
Set Status code for the particle.
Definition: MCParticle.h:335
void setProductionTime(float time)
Set production time.
Definition: MCParticle.h:373
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
TGeoHMatrix SADtoGeant(ReaderSAD::AcceleratorRings accRing, double s)
Transformation matrix.
Definition: ReaderSAD.cc:416
@ c_HER
High Energy Ring (electrons)
Definition: ReaderSAD.h:47
@ c_LER
Low Energy Ring (positrons)
Definition: ReaderSAD.h:48
ClassSADMetaHit - digitization simulated metahit for the SAD detector.
Definition: SADMetaHit.h:25
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
T * appendNew()
Construct a new T object at the end of the array.
Definition: StoreArray.h:246
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:95
static const double us
[microsecond]
Definition: Unit.h:97
static const double m
[meters]
Definition: Unit.h:69
static const double s
[second]
Definition: Unit.h:95
static const double GeV
Standard of [energy, momentum, mass].
Definition: Unit.h:51
double getAngle(const std::string &path="") const noexcept(false)
Get the parameter path as a double converted to the standard angle unit.
Definition: Interface.h:299
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:649
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
Abstract base class for different kinds of events.
STL namespace.
int nturn
number of turns from scattered to lost
double E
E at lost position [GeV] (in fact momentum magnitude!)
double sraw
s at lost position [m] before matching G4 beam pipe inner surface
double ss
scattered position (|s|<Ltot/2) [m]
double yraw
y at lost position [m] before matching G4 beam pipe inner surface
double rr
sqrt(xraw*xraw+yraw*yraw) [m]
double dp_over_p0
momentum deviation of the lost particle
double xraw
x at lost position [m] before matching G4 beam pipe inner surface
double s
lost position measured from IP along the ring [m]