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 // register MCParticles
67
68 StoreArray<MCParticle> mcParticles;
69 mcParticles.registerInDataStore();
70
72 sadHits.registerInDataStore();
73
74 // check steering parameters
75
76 if (m_ringName != "LER" and m_ringName != "HER") {
77 B2ERROR("ringName can only be LER or HER");
78 }
79 if (m_realTime <= 0) {
80 B2ERROR("realTime must be positive");
81 }
82
83 // open SAD file and get TTree
84
85 m_file = TFile::Open(m_fileName.c_str());
86 if (!m_file) {
87 B2ERROR(m_fileName << ": can't open file");
88 return;
89 }
90 m_tree = static_cast<TTree*>(m_file->Get(m_treeName.c_str()));
91 if (!m_tree) {
92 B2ERROR(m_treeName << ": no such TTree in the SAD file");
93 return;
94 }
95
96 m_tree->SetBranchAddress("s", &m_sad.s);
97 m_tree->SetBranchAddress("x", &m_sad.x);
98 m_tree->SetBranchAddress("px", &m_sad.px);
99 m_tree->SetBranchAddress("y", &m_sad.y);
100 m_tree->SetBranchAddress("py", &m_sad.py);
101 m_tree->SetBranchAddress("E", &m_sad.E);
102 m_tree->SetBranchAddress("rate", &m_sad.rate);
103
104 // for SADMetaHit
105 m_tree->SetBranchAddress("ss", &m_sad.ss);
106 m_tree->SetBranchAddress("sraw", &m_sad.sraw);
107 m_tree->SetBranchAddress("ssraw", &m_sad.ssraw);
108 m_tree->SetBranchAddress("nturn", &m_sad.nturn);
109 m_tree->SetBranchAddress("xraw", &m_sad.xraw);
110 m_tree->SetBranchAddress("yraw", &m_sad.yraw);
111 m_tree->SetBranchAddress("r", &m_sad.r);
112 m_tree->SetBranchAddress("rr", &m_sad.rr);
113 m_tree->SetBranchAddress("dp_over_p0", &m_sad.dp_over_p0);
114 m_tree->SetBranchAddress("watt", &m_sad.watt);
115
116 int numEntries = m_tree->GetEntries();
117 if (numEntries <= 0) {
118 B2ERROR("SAD tree is empty");
119 return;
120 }
121
122 // construct vector of cumulative rates and determine number of events to generate
123
124 m_rates.push_back(0);
125 for (int i = 0; i < numEntries; i++) {
126 m_tree->GetEntry(i);
127 m_rates.push_back(m_sad.rate + m_rates.back());
128 m_counters.push_back(0);
129 }
130
131 m_numEvents = gRandom->Poisson(m_rates.back() * m_realTime / Unit::s);
132
133 // set rotation from SAD to Belle II
134
135 if (m_ringName == "LER") {
136 GearDir ring("/Detector/SuperKEKB/LER/");
137 m_rotation.SetAngle(ring.getAngle("angle"));
138 m_ring = 2;
139 } else {
140 GearDir ring("/Detector/SuperKEKB/HER/");
141 m_rotation.SetAngle(ring.getAngle("angle"));
142 m_ring = 1;
143 }
144
145 m_sectionOrdering.insert(m_sectionOrdering.end(), {1, 12, 11, 10, 9, 8, 7, 6, 5, 4, 3, 2});
146
147 B2RESULT("BG rate: " << m_rates.back() / 1e6 << " MHz, events to generate: "
148 << m_numEvents);
149
150}
151
153{
154
155 // check event counter
156
158 StoreObjPtr<EventMetaData> eventMetaData;
159 eventMetaData->setEndOfData(); // stop event processing
160 return;
161 }
163
164 // get SAD entry
165
166 int i = generateEntry();
167 m_tree->GetEntry(i);
168 m_counters[i]++;
169
170 // make SADMetaHit
171 int ring_section = -1;
172 double ssraw = m_sad.ss;
173 if (m_sad.ss < 0) ssraw += 3016.;
174 int section = (int)(ssraw * 12 / 3016.);
175 if ((unsigned)section < m_sectionOrdering.size()) ring_section = m_sectionOrdering[section];
176
177 StoreArray<SADMetaHit> SADMetaHits;
178 SADMetaHits.appendNew(SADMetaHit(m_sad.ssraw, m_sad.sraw, m_sad.ss, m_sad.s,
179 0., m_sad.nturn,
180 m_sad.x, m_sad.y, m_sad.px, m_sad.py, m_sad.xraw, m_sad.yraw,
181 m_sad.r, m_sad.rr, m_sad.dp_over_p0, m_sad.E, m_sad.rate,
182 m_sad.watt, m_ring, ring_section));
183
184
185 // transform to Belle II (flip sign of y and s, rotate)
186
187 ROOT::Math::XYZVector position(m_sad.x * Unit::m, -m_sad.y * Unit::m, -m_sad.s * Unit::m);
188 position = m_rotation * position;
189
190 double pz = sqrt(m_sad.E * m_sad.E - m_sad.px * m_sad.px - m_sad.py * m_sad.py);
191 if (m_ringName == "LER") pz = -pz;
192 ROOT::Math::XYZVector momentum(m_sad.px, -m_sad.py, pz);
193 momentum = m_rotation * momentum;
194
195 // PDG code and mass
196
197 int pdgCode = Const::electron.getPDGCode();
198 double mass = Const::electron.getMass();
199 if (m_ringName == "LER") pdgCode = -pdgCode;
200
201 // append SAD particle to MCParticles
202
203 StoreArray<MCParticle> MCParticles;
204 MCParticle* part = MCParticles.appendNew();
205 part->setPDG(pdgCode);
206 part->setMass(mass);
207 part->setProductionVertex(position);
208 part->setProductionTime(0);
209 part->setMomentum(momentum);
210 part->setEnergy(sqrt(momentum.Mag2() + mass * mass));
211 part->setValidVertex(true);
214
215 // FarBeamLine region transformation
216 if (abs(m_sad.s * Unit::m) > 4.0 * Unit::m) {
217 // initial coordinates in SAD space
218 double particlePosSADfar[] = {m_sad.x* Unit::m, -m_sad.y* Unit::m, 0.0 * Unit::m};
219 double particleMomSADfar[] = {m_sad.px* Unit::GeV, -m_sad.py* Unit::GeV, pz* Unit::GeV};
220 // final coordinates in Geant4 space
221 double particlePosGeant4[] = {0.0, 0.0, 0.0};
222 double particleMomGeant4[] = {0.0, 0.0, 0.0};
223
224 // create a transformation matrix for a given ring
225 TGeoHMatrix transMatrix;
226 if (m_ringName == "LER") {
227 transMatrix = m_readerSAD.SADtoGeant(ReaderSAD::c_LER, m_sad.s * Unit::m);
228 } else {
229 transMatrix = m_readerSAD.SADtoGeant(ReaderSAD::c_HER, m_sad.s * Unit::m);
230 }
231
232 // calculate a new set of coordinates in Geant4 space
233 transMatrix.LocalToMaster(particlePosSADfar, particlePosGeant4); // position
234 transMatrix.LocalToMasterVect(particleMomSADfar, particleMomGeant4); // momentum
235 // apply a new set of coordinates
236 part->setMomentum(ROOT::Math::XYZVector(particleMomGeant4[0], particleMomGeant4[1], particleMomGeant4[2]));
237 part->setProductionVertex(ROOT::Math::XYZVector(particlePosGeant4[0], particlePosGeant4[1], particlePosGeant4[2]));
238 }
239}
240
242{
243 // close SAD file
244 m_file->Close();
245
246 B2RESULT("BG rate: " << m_rates.back() / 1e6 << " MHz, equivalent time: "
247 << m_realTime / Unit::us << " us, "
248 << m_eventCounter << " events generated");
250 B2ERROR("Number of generated events does not match the equivalent running time: "
251 << m_numEvents << " events needed, but "
252 << m_eventCounter << " generated");
253
254 int imax = 0;
255 double average = 0;
256 for (unsigned i = 0; i < m_counters.size(); i++) {
257 if (m_counters[i] > m_counters[imax]) imax = i;
258 average += m_counters[i];
259 }
260 average /= m_counters.size();
261 B2RESULT("SAD particle usage: on average " << average << " times, max "
262 << m_counters[imax] << " times (entry " << imax << ")");
263
264}
265
266
268{
269 double rate = gRandom->Uniform(m_rates.back());
270 int i1 = 0;
271 int i2 = m_rates.size() - 1;
272 while (i2 - i1 > 1) {
273 int i = (i1 + i2) / 2;
274 if (rate > m_rates[i]) {
275 i1 = i;
276 } else {
277 i2 = i;
278 }
279 }
280 return i1;
281}
std::vector< int > m_counters
counters: how many times SAD particles are used
virtual void initialize() override
Initialize the Module.
virtual void event() override
Event processor.
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]
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
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
void setDescription(const std::string &description)
Sets the description of the module.
Definition Module.cc:214
Module()
Constructor.
Definition Module.cc:30
@ 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:96
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.