Belle II Software development
KLMScintillatorSimulatorModule.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 <klm/modules/KLMScintillatorSimulator/KLMScintillatorSimulatorModule.h>
11
12/* KLM headers. */
13#include <klm/eklm/geometry/GeometryData.h>
14#include <klm/simulation/ScintillatorSimulator.h>
15
16/* Basf2 headers. */
17#include <framework/gearbox/Unit.h>
18
19/* ROOT headers. */
20#include <TH1F.h>
21
22using namespace Belle2;
23
24REG_MODULE(KLMScintillatorSimulator)
25
26static const char MemErr[] = "Memory allocation error.";
27
29 Module(),
30 m_fout(nullptr),
31 m_SciSimPar(nullptr),
32 m_hDir(nullptr),
33 m_hRef(nullptr)
34{
35 setDescription("Standalone generation and studies of ADC output.");
37 addParam("OutputFile", m_out, "Output file.",
38 std::string("KLMScintillatorSimulator.root"));
39 addParam("Mode", m_mode, "Mode (\"Shape\" or \"Strips\").",
40 std::string("Strips"));
41}
42
44{
45}
46
48 const char* name, double l, double d, int npe)
49{
50 int j;
51 double t, s;
52 KLM::ScintillatorSimulator fe(m_SciSimPar, nullptr, 0, false);
53 KLMTime* klmTime = &(KLMTime::Instance());
54 klmTime->updateConstants();
55 TH1F* h = nullptr;
58 try {
59 h = new TH1F(name, "", m_SciSimPar->getNDigitizations(), 0, t);
60 } catch (std::bad_alloc& ba) {
61 B2FATAL(MemErr);
62 }
63 for (j = 0; j < m_SciSimPar->getNDigitizations(); j++) {
64 m_hDir[j] = 0;
65 m_hRef[j] = 0;
66 }
67 fe.generatePhotoelectrons(l, d, npe, 0, false);
68 fe.generatePhotoelectrons(l, d, npe, 0, true);
69 fe.fillSiPMOutput(m_hDir, true, false);
70 fe.fillSiPMOutput(m_hRef, false, true);
71 s = 0;
72 for (j = 0; j < m_SciSimPar->getNDigitizations(); j++)
73 s = s + m_hDir[j] + m_hRef[j];
74 for (j = 1; j <= m_SciSimPar->getNDigitizations(); j++)
75 h->SetBinContent(j, (m_hDir[j - 1] + m_hRef[j - 1]) / s);
76 h->Write();
77 delete h;
78}
79
81{
82 /* cppcheck-suppress variableScope */
83 char str[32];
84 /* cppcheck-suppress variableScope */
85 int i;
86 /* cppcheck-suppress variableScope */
87 double l;
88 if (!m_SciSimParDatabase.isValid())
89 B2FATAL("EKLM digitization parameters are not available.");
92 try {
93 m_fout = new TFile(m_out.c_str(), "recreate");
94 } catch (std::bad_alloc& ba) {
95 B2FATAL(MemErr);
96 }
97 m_hDir = (float*)malloc(m_SciSimPar->getNDigitizations() * sizeof(float));
98 if (m_hDir == nullptr)
99 B2FATAL(MemErr);
100 m_hRef = (float*)malloc(m_SciSimPar->getNDigitizations() * sizeof(float));
101 if (m_hRef == nullptr)
102 B2FATAL(MemErr);
103 if (m_mode.compare("Strips") == 0) {
104 for (i = 1; i <= geoDat->getNStrips(); i++) {
105 l = geoDat->getStripLength(i) / CLHEP::mm * Unit::mm;
106 snprintf(str, 32, "h%d_near", i);
107 generateHistogram(str, l, 0, 10000);
108 snprintf(str, 32, "h%d_far", i);
109 generateHistogram(str, l, l, 10000);
110 }
111 } else if (m_mode.compare("Shape") == 0) {
113 generateHistogram("FitShape", 0, 0, 1000000);
114 } else
115 B2FATAL("Unknown operation mode.");
116 free(m_hDir);
117 free(m_hRef);
118 m_fout->Close();
119 delete m_fout;
120}
121
123{
124}
125
127{
128}
129
131{
132}
133
135{
136 delete m_SciSimPar;
137}
138
int getNStrips() const
Get number of strips.
EKLM geometry data.
Definition: GeometryData.h:38
static const GeometryData & Instance(enum DataSource dataSource=c_Database, const GearDir *gearDir=nullptr)
Instantiation.
Definition: GeometryData.cc:33
double getStripLength(int strip) const
Get strip length.
Definition: GeometryData.h:71
Class to store KLM scintillator simulation parameters in the database.
void setMirrorReflectiveIndex(float reflectiveIndex)
Set mirror reflective index.
int getADCSamplingTDCPeriods() const
Get ADC sampling time in TDC periods.
int getNDigitizations() const
Get number of digitizations (points) in one sample.
void generateHistogram(const char *name, double l, double d, int npe)
Generate output histogram.
void event() override
This method is called for each event.
void endRun() override
This method is called if the current run ends.
void terminate() override
This method is called at the end of the event processing.
KLMScintillatorDigitizationParameters * m_SciSimPar
Scintillator simulation parameters.
void beginRun() override
Called when entering a new run.
DBObjPtr< KLMScintillatorDigitizationParameters > m_SciSimParDatabase
Scintillator simulation parameters.
KLM time conversion.
Definition: KLMTime.h:27
double getTDCPeriod() const
Get TDC period.
Definition: KLMTime.h:45
void updateConstants()
Update constants from database objects.
Definition: KLMTime.cc:20
static KLMTime & Instance()
Instantiation.
Definition: KLMTime.cc:14
Digitize EKLMSim2Hits to get EKLM StripHits.
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_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
static const double mm
[millimeters]
Definition: Unit.h:70
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.