Belle II Software development
ECLDigiStudyModule.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/**************************************************************************
10 * *
11 * Description: This module write ECL digi information in a root tuple *
12 * to study amplitude and time info *
13 * *
14 **************************************************************************/
15
16/* Own header. */
17#include <ecl/modules/eclDataAnalysis/ECLDigiStudyModule.h>
18
19/* ECL headers. */
20#include <ecl/dataobjects/ECLDigit.h>
21#include <ecl/dataobjects/ECLDsp.h>
22#include <ecl/dataobjects/ECLHit.h>
23#include <ecl/dataobjects/ECLTrig.h>
24#include <ecl/geometry/ECLGeometryPar.h>
25
26/* ROOT headers. */
27#include <TFile.h>
28#include <TTree.h>
29
30/* C++ headers. */
31#include <algorithm>
32
33using namespace std;
34using namespace Belle2;
35using namespace ECL;
36//-----------------------------------------------------------------
37// Register the Module
38//-----------------------------------------------------------------
39REG_MODULE(ECLDigiStudy);
40
41//-----------------------------------------------------------------
42// Implementation
43//-----------------------------------------------------------------
45{
46 //Set module properties
47 setDescription("EclCovMatrixNtuple: write ECL waveform and fitted time and amplitude in a root file");
49 //Parameters definition
50 addParam("outputFileName", m_dataOutFileName,
51 "Output root file name of this module", string("digistudy"));
52 addParam("dspArrayName1", m_dspArrayName1, "name of input ECLDsp Array 1", string("ECLDsps"));
53 addParam("digiArrayName1", m_digiArrayName1, "name of input ECLDigit Array 1", string("ECLDigits"));
54 addParam("dspArrayName2", m_dspArrayName2, "name of input ECLDsp Array 2", string("ECLDspsPureCsI"));
55 addParam("digiArrayName2", m_digiArrayName2, "name of input ECLDigit Array 2", string("ECLDigitsPureCsI"));
56}
57
59{
60
61 m_eclDspArray1.registerInDataStore(m_dspArrayName1);
62 m_eclDspArray2.registerInDataStore(m_dspArrayName2);
63 m_eclDigiArray1.registerInDataStore(m_digiArrayName1);
64 m_eclDigiArray2.registerInDataStore(m_digiArrayName2);
65
66 B2INFO("[EclCovMatrixNtuple Module]: Starting initialization of EclCovMatrixNtuple Module.");
67
68 // Initializing the output root file
69 string dataFileName = m_dataOutFileName + ".root";
70
71 m_rootFile = new TFile(dataFileName.c_str(), "RECREATE");
72 m_tree = new TTree("m_tree", "EclCovMatrixNtuple tree");
74
75 m_tree->Branch("nhits", &m_nhits, "nhits/I");
76 m_tree->Branch("cellid", m_cellId, "cellid[nhits]/I");
77 m_tree->Branch("energy", m_energy, "energy[nhits]/D");
78 m_tree->Branch("alle", m_allenergy, "alle[nhits]/D");
79 m_tree->Branch("time", m_time, "time[nhits]/D");
80 m_tree->Branch("theta", m_theta, "theta[nhits]/I");
81 m_tree->Branch("phi", m_phi, "phi[nhits]/I");
82 m_tree->Branch("a1", m_DspHit1, "a1[nhits][31]/I");
83 m_tree->Branch("a2", m_DspHit2, "a2[nhits][31]/I");
84 m_tree->Branch("base1", m_baseline1, "base1[nhits][16]/I");
85 m_tree->Branch("base2", m_baseline2, "base2[nhits][16]/I");
86 m_tree->Branch("baseavg1", m_baselineAvg1, "baseavg1[nhits]/D");
87 m_tree->Branch("baseavg2", m_baselineAvg2, "baseavg2[nhits]/D");
88 m_tree->Branch("mv1", m_maxVal1, "mv1[nhits]/I");
89 m_tree->Branch("mv2", m_maxVal2, "mv2[nhits]/I");
90 m_tree->Branch("t1", m_digiTime1, "t1[nhits]/D");
91 m_tree->Branch("t2", m_digiTime2, "t2[nhits]/D");
92 m_tree->Branch("e1", m_digiE1, "e1[nhits]/D");
93 m_tree->Branch("e2", m_digiE2, "e2[nhits]/D");
94 m_tree->Branch("q1", m_digiQual1, "q1[nhits]/I");
95 m_tree->Branch("q2", m_digiQual2, "q2[nhits]/I");
96
97 m_tree->Branch("trig1", &m_trig1, "trig1/D");
98 m_tree->Branch("trig2", &m_trig2, "trig2/D");
99 m_tree->Branch("neclhits", &m_neclhits, "neclhits/I");
100
101 B2INFO("[ECLDigiStudyModule Module]: Finished initialising");
102}
103
105{
106 m_rootFile->cd();
107 m_tree->Write();
108 m_rootFile->Close();
109}
110
112{
113
115
116 int i;
117 m_neclhits = m_eclHitsArray.getEntries();
124
125 for (int j = 0; j < ECLElementNumbers::c_NCrystals; j++) {
126 fill_n(m_DspHit1[j], 31, 0);
127 fill_n(m_DspHit2[j], 31, 0);
128 fill_n(m_baseline1[j], 16, 0);
129 fill_n(m_baseline2[j], 16, 0);
130 }
131
136
143
144 m_trig1 = m_eclTrigArray[0]->getTimeTrig();
145 m_trig2 = 0;
146
147 for (const auto& hit : m_eclHitsArray) {
148
149 i = hit.getCellId() - 1;
150 m_cellId[i] = i;
151 double e = hit.getEnergyDep();
152 if (hit.getBackgroundTag() == BackgroundMetaData::bg_none) {
153 if (m_energy[i] < e)
154 m_time[i] = hit.getTimeAve();
155 m_energy[i] += e;
156 }
157 m_allenergy[i] += e;
158 eclgeo->Mapping(i);
159 m_theta[i] = eclgeo->GetThetaID();
160 m_phi[i] = eclgeo->GetPhiID();
161
162 m_digiTime1[i] = m_digiTime2[i] = -10000;
163 m_digiE1[i] = m_digiE2[i] = -10.;
164 }
165
166 for (const auto& digit1 : m_eclDigiArray1) {
167 i = digit1.getCellId() - 1;
168 m_digiTime1[i] = digit1.getTimeFit();
169 m_digiE1[i] = digit1.getAmp() / 20000.;
170 m_digiQual1[i] = digit1.getQuality();
171 }
172
173 for (const auto& digit2 : m_eclDigiArray2) {
174 i = digit2.getCellId() - 1;
175 m_digiTime2[i] = digit2.getTimeFit() / 1000.;
176 m_digiE2[i] = digit2.getAmp() / 20000.;
177 m_digiQual2[i] = digit2.getQuality();
178 }
179
180 for (const auto& dsp : m_eclDspArray1) {
181 i = dsp.getCellId() - 1;
182 dsp.getDspA(m_DspHit1[i]);
183 for (int j = 0; j < 16; j++) {
184 m_baseline1[i][j] = m_DspHit1[i][j] - 3000;
185 m_baselineAvg1[i] += m_baseline1[i][j];
186 }
187 m_baselineAvg1[i] /= 16;
188 assert(m_DspHit1[i] == (&m_DspHit1[i][0]));
189 m_maxVal1[i] = * (max_element(& m_DspHit1[i][16], (& m_DspHit1[i][16]) + 15));
190 }
191
192 for (const auto& dsp : m_eclDspArray2) {
193 i = dsp.getCellId() - 1;
194 dsp.getDspA(m_DspHit2[i]);
195 for (int j = 0; j < 16; j++) {
196 m_baseline2[i][j] = m_DspHit2[i][j] - 3000;
197 m_baselineAvg2[i] += m_baseline2[i][j];
198 }
199 m_baselineAvg2[i] /= 16;
200 m_maxVal2[i] = * (max_element(& m_DspHit2[i][16], (& m_DspHit2[i][16]) + 15));
201 }
202
203 m_tree->Fill();
204}
int m_theta[ECLElementNumbers::c_NCrystals]
Array oh ThetaID.
TTree * m_tree
Root tree for saving the output.
int m_digiQual2[ECLElementNumbers::c_NCrystals]
Digit quality for second digit array.
double m_trig2
Trigger time for array 2.
std::string m_dataOutFileName
output root file name (given as Module parameter)
StoreArray< ECLDigit > m_eclDigiArray2
Store array: ECLDigit.
virtual void initialize() override
Initializes the module.
int m_DspHit1[ECLElementNumbers::c_NCrystals][31]
WF sampling points for first digit array.
virtual void event() override
Method is called for each event.
double m_time[ECLElementNumbers::c_NCrystals]
Array of digit time.
TFile * m_rootFile
Root file for saving the output.
int m_neclhits
Actual number of hits.
virtual void terminate() override
terminate
int m_digiQual1[ECLElementNumbers::c_NCrystals]
Digit quality for first digit array.
std::string m_dspArrayName1
Name of first DSP array.
int m_maxVal1[ECLElementNumbers::c_NCrystals]
WF maximum for first digit array.
std::string m_digiArrayName2
Name of second digit array.
int m_cellId[ECLElementNumbers::c_NCrystals]
Array of cellIDs.
std::string m_digiArrayName1
Name of first digit array.
double m_digiTime2[ECLElementNumbers::c_NCrystals]
Digit time for second digit array.
int m_nhits
Maximum number of hits.
StoreArray< ECLTrig > m_eclTrigArray
Store array: ECLTrig.
std::string m_dspArrayName2
Name of second DSP array.
double m_energy[ECLElementNumbers::c_NCrystals]
Array of deposited MC energy.
double m_digiTime1[ECLElementNumbers::c_NCrystals]
Digit time for first digit array.
StoreArray< ECLDigit > m_eclDigiArray1
Store array: ECLDigit.
double m_digiE2[ECLElementNumbers::c_NCrystals]
Deposited energy for second digit array.
double m_baselineAvg1[ECLElementNumbers::c_NCrystals]
Baseline energy for first digit array.
int m_maxVal2[ECLElementNumbers::c_NCrystals]
WF maximum for second digit array.
int m_baseline2[ECLElementNumbers::c_NCrystals][16]
Baseline sampling points for second digit array.
StoreArray< ECLHit > m_eclHitsArray
Store array: ECLHit.
int m_baseline1[ECLElementNumbers::c_NCrystals][16]
Baseline sampling points for first digit array.
double m_digiE1[ECLElementNumbers::c_NCrystals]
Deposited energy for first digit array.
StoreArray< ECLDsp > m_eclDspArray1
Store array: ECLDsp.
double m_trig1
Trigger time for array 1.
double m_allenergy[ECLElementNumbers::c_NCrystals]
Array of deposited energy.
StoreArray< ECLDsp > m_eclDspArray2
Store array: ECLDsp.
int m_phi[ECLElementNumbers::c_NCrystals]
Array of PhiID.
int m_DspHit2[ECLElementNumbers::c_NCrystals][31]
WF sampling points for second digit array.
double m_baselineAvg2[ECLElementNumbers::c_NCrystals]
Baseline energy for second digit array.
The Class for ECL Geometry Parameters.
static ECLGeometryPar * Instance()
Static method to get a reference to the ECLGeometryPar instance.
void Mapping(int cid)
Mapping theta, phi Id.
int GetThetaID()
Get Theta Id.
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
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
const int c_NCrystals
Number of crystals.
Abstract base class for different kinds of events.
STL namespace.