Belle II Software development
SADInputModule.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#include <generators/modules/SADInputModule.h>
10#include <mdst/dataobjects/MCParticleGraph.h>
11
12#include <framework/logging/Logger.h>
13#include <framework/utilities/FileSystem.h>
14
15#include <framework/gearbox/Unit.h>
16#include <framework/gearbox/GearDir.h>
17
18#include <framework/datastore/StoreArray.h>
19
20#include <TGeoMatrix.h>
21
22using namespace std;
23using namespace Belle2;
24
25//-----------------------------------------------------------------
26// Register the Module
27//-----------------------------------------------------------------
28REG_MODULE(SADInput);
29
30//-----------------------------------------------------------------
31// Implementation
32//-----------------------------------------------------------------
33
35{
36 //Set module properties
37 setDescription("Reads the SAD data from a root file and stores it into the MCParticle collection.");
39
40 //Parameter definition
41 addParam("AccRing", m_accRing, "The accelerator ring: 0 = LER, 1 = HER");
42 addParam("ReadoutTime", m_readoutTime, "The readout time of the detector [ns]", 20 * Unit::us);
43 addParam("ReadMode", m_readMode,
44 "The read mode: 0 = one SAD particle per event, 1 = one real particle per event, 2 = all SAD particles per event", 0);
45 addParam("Filename", m_filename, "The filename of the SAD input file.");
46 addParam("Range", m_range, "All particles within the range around the IP are loaded [cm].", 300.0 * Unit::cm);
47 addParam("PxResolution", m_pxRes, "The resolution for the x momentum component of the SAD real particle.", 0.01);
48 addParam("PyResolution", m_pyRes, "The resolution for the y momentum component of the SAD real particle.", 0.01);
49 addParam("RotateParticles", m_rotateParticles,
50 "Rotate the SAD particles around the nominal beam axis [deg] (just for unphysical tests !!!).", 0.0);
51
52 m_PipePartMatrix = NULL;
53}
54
55
57{
58 //Register inputs/outputs
59 m_eventMetaDataPtr.isRequired();
60 StoreArray<MCParticle> mcParticles;
61 mcParticles.registerInDataStore();
62
63 //Check parameters
65 B2ERROR("Parameter <Filename>: Could not open the file. The filename " << m_filename << " does not exist !");
66 } else m_reader.open(m_filename);
67
68 //Initialize the SAD reader.
69 //Set the transformation from local SAD plane to global geant4 space.
70 m_PipePartMatrix = new TGeoHMatrix("SADTrafo");
72
73 GearDir ler("/Detector/SuperKEKB/LER/");
74 GearDir her("/Detector/SuperKEKB/HER/");
75
76 switch (m_accRing) {
77 case 0 : m_PipePartMatrix->RotateY(ler.getDouble("angle") / Unit::deg);
79 break;
80 case 1 : m_PipePartMatrix->RotateY(her.getDouble("angle") / Unit::deg);
82 break;
83 default: B2FATAL("Please specify a valid number for the accRing parameter (0 = LER, 1 = HER) !");
84 break;
85 }
86
88}
89
90
92{
93 try {
95
96 try {
97 //----------------------------------
98 // Read the SAD data
99 //----------------------------------
100 switch (m_readMode) {
101 case 0: readSADParticle(m_reader, mpg);
102 break;
103 case 1: readRealParticle(m_reader, mpg);
104 break;
105 case 2: m_reader.addAllSADParticles(mpg);
106 break;
107 default: B2FATAL("Please specify a valid number for the readMode parameter !");
108 break;
109 }
110
111 //----------------------------------
112 // Generate MCParticles collection
113 //----------------------------------
115
116 } catch (ReaderSAD::SADEndOfFile& exc) {
117 B2DEBUG(10, exc.what());
118 m_eventMetaDataPtr->setEndOfData();
119 return;
120 }
121 } catch (runtime_error& exc) {
122 B2ERROR(exc.what());
123 }
124}
125
126
127//====================================================================
128// Private methods
129//====================================================================
130
132{
133 double rate = reader.getSADParticle(mpg);
134 if (rate < 0) return;
135 m_eventMetaDataPtr->setGeneratedWeight(rate);
136}
137
139{
140 if (!reader.getRealParticle(mpg)) return;
141 m_eventMetaDataPtr->setGeneratedWeight(1.0);
142}
static bool fileExists(const std::string &filename)
Check if the file with given filename exists.
Definition: FileSystem.cc:32
GearDir is the basic class used for accessing the parameter store.
Definition: GearDir.h:31
Class to build, validate and sort a particle decay chain.
@ c_checkCyclic
Check for cyclic dependencies.
@ c_setDecayInfo
Set decay time and vertex.
void generateList(const std::string &name="", int options=c_setNothing)
Generates the MCParticle list and stores it in the StoreArray with the given name.
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_Input
This module is an input module (reads data).
Definition: Module.h:78
Class to read files that have been created by SAD and store their content in a MCParticle graph.
Definition: ReaderSAD.h:35
void addAllSADParticles(MCParticleGraph &graph)
Reads all SAD particles from the file into the MCParticles collection which are inside the specified ...
Definition: ReaderSAD.cc:212
void open(const std::string &filename)
Opens a root file and prepares it for reading.
Definition: ReaderSAD.cc:80
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
void setMomentumRes(double pxRes, double pyRes)
Sets the resolution of the momentum for the real particles.
Definition: ReaderSAD.h:82
@ c_HER
High Energy Ring (electrons)
Definition: ReaderSAD.h:47
@ c_LER
Low Energy Ring (positrons)
Definition: ReaderSAD.h:48
double m_pyRes
The resolution for the y momentum component of the SAD real particle.
TGeoHMatrix * m_PipePartMatrix
Transformation matrix from SAD space into geant4 space.
virtual void initialize() override
Checks the validity of the module parameters.
double m_range
All particles within the range around the IP are loaded.
virtual void event() override
Reads the data and stores it into the MCParticle collection.
int m_readMode
The read mode: 0 = one SAD particle per event, 1 = one real particle per event, 2 = all SAD particles...
void readRealParticle(ReaderSAD &reader, MCParticleGraph &mpg)
Reads one real particle per event into the MCParticle graph.
double m_rotateParticles
Rotate the SAD particles around the nominal beam axis (just for unphysical tests !...
int m_accRing
The accelerator ring: 0 = LER, 1 = HER.
SADInputModule()
Constructor.
std::string m_filename
The filename of the SAD file.
StoreObjPtr< EventMetaData > m_eventMetaDataPtr
EventMetaData (for setting event weights).
double m_pxRes
The resolution for the x momentum component of the SAD real particle.
double m_readoutTime
The readout time of the detector [ns].
ReaderSAD m_reader
The SAD reader object for the SAD data.
void readSADParticle(ReaderSAD &reader, MCParticleGraph &mpg)
Reads one SAD particle per event into the MCParticle graph.
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
static const double deg
degree to radians
Definition: Unit.h:109
static const double us
[microsecond]
Definition: Unit.h:97
static const double cm
Standard units with the value = 1.
Definition: Unit.h:47
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:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.
STL namespace.