Belle II Software development
EclCovMatrixNtupleModule.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/EclCovMatrixNtupleModule.h>
18
19/* ECL headers. */
20#include <ecl/dataobjects/ECLDigit.h>
21#include <ecl/dataobjects/ECLDsp.h>
22#include <ecl/dataobjects/ECLTrig.h>
23#include <ecl/geometry/ECLGeometryPar.h>
24
25/* ROOT headers. */
26#include <TFile.h>
27#include <TTree.h>
28
29using namespace std;
30using namespace Belle2;
31using namespace ECL;
32//-----------------------------------------------------------------
33// Register the Module
34//-----------------------------------------------------------------
35REG_MODULE(EclCovMatrixNtuple);
36
37//-----------------------------------------------------------------
38// Implementation
39//-----------------------------------------------------------------
41{
42 //Set module properties
43 setDescription("EclCovMatrixNtuple: write ECL waveform and fitted time and amplitude in a root file");
45 //Parameters definition
46 addParam("outputFileName", m_dataOutFileName,
47 "Output root file name of this module", string("EclCovMatrixNtuple"));
48 addParam("dspArrayName", m_dspArrayName, "name of input ECLDsp Array", string("ECLDsps"));
49 addParam("digiArrayName", m_digiArrayName, "name of input ECLDigit Array", string("ECLDigits"));
50}
51
53{
54 m_eclDspArray.registerInDataStore(m_dspArrayName);
55 m_eclDigiArray.registerInDataStore(m_digiArrayName);
56 B2INFO("[EclCovMatrixNtuple Module]: Starting initialization of EclCovMatrixNtuple Module.");
57
58 // Initializing the output root file
59 string dataFileName = m_dataOutFileName + ".root";
60
61 m_rootFile = new TFile(dataFileName.c_str(), "RECREATE");
62 m_tree = new TTree("m_tree", "EclCovMatrixNtuple tree");
63
64 m_tree->Branch("energy", &m_energy, "energy/D");
65 m_tree->Branch("nhits", &m_nhits, "nhits/I");
66 m_tree->Branch("cellID", m_cellID, "cellID[nhits]/I");
67 m_tree->Branch("theta", m_theta, "theta[nhits]/I");
68 m_tree->Branch("phi", m_phi, "phi[nhits]/I");
69 m_tree->Branch("hitA", m_DspHit, "hitA[nhits][31]/I");
70 m_tree->Branch("hitT", m_hitTime, "hitT[nhits]/D");
71 m_tree->Branch("deltaT", m_DeltaT, "deltaT[nhits]/D");
72
73 m_tree->Branch("digiE", m_hitE, "hitE[nhits]/D");
74 m_tree->Branch("digiT", m_DigiTime, "digiT[nhits]/I");
75
76
77 B2INFO("[EclCovMatrixNtuple Module]: Finished initialising the Time Information Study in Gamma Reconstruction Module.");
78}
79
81{
82 m_rootFile->cd();
83 m_tree->Write();
84 m_rootFile->Close();
85}
86
88{
89
90
91 m_nevt++;
92 m_energy = 0;
93
94
95
96 m_nhits = 0;
97
98 for (int i = 0; i < ECLElementNumbers::c_NCrystals; i++) {
99 m_cellID[i] = 0;
100 m_theta[i] = 0;
101 m_phi[i] = 0;
102 m_hitTime[i] = 0;
103 m_hitE[i] = 0.0;
104 m_hitTime[i] = 0.0;
105 m_DigiTime[i] = 0;
106 m_DeltaT[i] = 0.0;
107 for (int j = 0; j < 31; j++) {
108 m_DspHit[i][j] = 0;
109 }
110 }
111
112
113 m_nhits = m_eclDigiArray.getEntries();
114
115 // There is only 1 ECLTrig per event
116 assert(m_eclTrigArray.getEntries() == 1);
118 for (const auto& adigit : m_eclDigiArray) {
119 size_t cellIndex = static_cast<size_t>(adigit.getCellId() - 1);
120 eclgeo->Mapping(cellIndex);
121 m_theta[cellIndex] = eclgeo->GetThetaID();
122 m_phi[cellIndex] = eclgeo->GetPhiID();
123 /*
124 cout << "cid : " << cellIndex
125 << " theta: " << m_theta[cellIndex]
126 << " phi: " << m_phi[cellIndex]
127 << endl;
128 */
129 m_cellID[cellIndex] = cellIndex;
130 m_hitE[cellIndex] = adigit.getAmp();
131 m_DigiTime[cellIndex] = adigit.getTimeFit();
132 // The following DeltaT IS DIFFERENT DeltaT in igitizer by the last factor 12.!
133 double deltaT = m_eclTrigArray[0]->getTimeTrig() * 508.0 / 12.0;
134 m_DeltaT[cellIndex] = deltaT;
135 m_hitTime[cellIndex] = 1520 - adigit.getTimeFit() - 64 * deltaT * 12.0 * 24.0 / 508.0 / 1536.0;
136 }
137
138 //The following works only because position of an ECLDigit in the ECLDigiArray
139 //is the same of position of the corresponding ECLDsp in the ECLDspArray
140 //Since the two a tightly related this must be enforced in a safer way
141 //in the ecl data objects model (Guglielmo De Nardo)
142 assert(m_eclDspArray.getEntries() == m_eclDigiArray.getEntries());
143 for (const auto& eclDsp : m_eclDspArray) {
144 size_t cellIndex = static_cast<size_t>(eclDsp.getCellId() - 1);
145 eclDsp.getDspA(m_DspHit[cellIndex]);
146 }
147 m_tree->Fill();
148}
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.
StoreArray< ECLDsp > m_eclDspArray
Store array: ECLDsp.
double m_hitTime[ECLElementNumbers::c_NCrystals]
eclHit Time
int m_theta[ECLElementNumbers::c_NCrystals]
Crystal Theta ID.
std::string m_dataOutFileName
output root file name (given as Module parameter)
virtual void initialize() override
Initializes the module.
virtual void event() override
Method is called for each event.
double m_hitE[ECLElementNumbers::c_NCrystals]
eclHit Energy
TFile * m_rootFile
Root file for saving the output.
double m_DeltaT[ECLElementNumbers::c_NCrystals]
eclTrig Time
virtual void terminate() override
Terminates the module.
std::string m_digiArrayName
eclDigit array name
StoreArray< ECLTrig > m_eclTrigArray
Store array: ECLTrig.
int m_cellID[ECLElementNumbers::c_NCrystals]
Crystal ID.
std::string m_dspArrayName
eclDSPs array name
int m_DspHit[ECLElementNumbers::c_NCrystals][31]
eclDsp sample Array
int m_phi[ECLElementNumbers::c_NCrystals]
Crystal Phi ID.
StoreArray< ECLDigit > m_eclDigiArray
Store array: ECLDigit.
int m_DigiTime[ECLElementNumbers::c_NCrystals]
eclDigit Time
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.