Belle II Software development
DosiStudyModule.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 <beast/dosi/modules/DosiStudyModule.h>
10#include <beast/dosi/dataobjects/DosiSimHit.h>
11#include <beast/dosi/dataobjects/DosiHit.h>
12
13#include <mdst/dataobjects/MCParticle.h>
14#include <framework/datastore/StoreArray.h>
15#include <framework/datastore/RelationIndex.h>
16#include <framework/logging/Logger.h>
17#include <framework/gearbox/Const.h>
18
19//c++
20#include <cmath>
21#include <string>
22
23// ROOT
24#include <TString.h>
25
26#include <generators/SAD/dataobjects/SADMetaHit.h>
27
28using namespace std;
29using namespace Belle2;
30using namespace dosi;
31
32//-----------------------------------------------------------------
33// Register the Module
34//-----------------------------------------------------------------
35REG_MODULE(DosiStudy);
36
37//-----------------------------------------------------------------
38// Implementation
39//-----------------------------------------------------------------
40
42{
43 // Set module properties
44 setDescription("Study module for Dosis (BEAST)");
45}
46
48{
49}
50
51//This module is a histomodule. Any histogram created here will be saved by the HistoManager module
53{
54 for (int i = 0; i < 5; i++) {
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72 h_dosi_edep0[i] = new TH1F(TString::Format("dosi_edep0_%d", i), "Energy deposited [MeV]", 50000, 0., 400.);
73 h_dosi_edep1[i] = new TH1F(TString::Format("dosi_edep1_%d", i), "Energy deposited [MeV]", 50000, 0., 400.);
74 h_dosi_edep2[i] = new TH1F(TString::Format("dosi_edep2_%d", i), "Energy deposited [MeV]", 50000, 0., 400.);
75 h_dosi_edep3[i] = new TH1F(TString::Format("dosi_edep3_%d", i), "Energy deposited [MeV]", 50000, 0., 400.);
76 h_dosi_edep4[i] = new TH1F(TString::Format("dosi_edep4_%d", i), "Energy deposited [MeV]", 50000, 0., 400.);
77 h_dosi_edep5[i] = new TH1F(TString::Format("dosi_edep5_%d", i), "Energy deposited [MeV]", 50000, 0., 400.);
78 h_dosi_edep6[i] = new TH1F(TString::Format("dosi_edep6_%d", i), "Energy deposited [MeV]", 50000, 0., 400.);
79 h_dosi_edep7[i] = new TH1F(TString::Format("dosi_edep7_%d", i), "Energy deposited [MeV]", 50000, 0., 400.);
80 h_dosi_edep8[i] = new TH1F(TString::Format("dosi_edep8_%d", i), "Energy deposited [MeV]", 50000, 0., 400.);
81
82 h_dosi_edep0[i]->Sumw2();
83 h_dosi_edep1[i]->Sumw2();
84 h_dosi_edep2[i]->Sumw2();
85 h_dosi_edep3[i]->Sumw2();
86 h_dosi_edep4[i]->Sumw2();
87 h_dosi_edep5[i]->Sumw2();
88 h_dosi_edep6[i]->Sumw2();
89 h_dosi_edep7[i]->Sumw2();
90 h_dosi_edep8[i]->Sumw2();
91
92 h_dosi_rs_edep0[i] = new TH2F(TString::Format("dosi_rs_edep0_%d", i), "Energy deposited [MeV]", 50000, 0., 400., 12, 0., 12.);
93 h_dosi_rs_edep1[i] = new TH2F(TString::Format("dosi_rs_edep1_%d", i), "Energy deposited [MeV]", 50000, 0., 400., 12, 0., 12.);
94 h_dosi_rs_edep2[i] = new TH2F(TString::Format("dosi_rs_edep2_%d", i), "Energy deposited [MeV]", 50000, 0., 400., 12, 0., 12.);
95 h_dosi_rs_edep3[i] = new TH2F(TString::Format("dosi_rs_edep3_%d", i), "Energy deposited [MeV]", 50000, 0., 400., 12, 0., 12.);
96 h_dosi_rs_edep4[i] = new TH2F(TString::Format("dosi_rs_edep4_%d", i), "Energy deposited [MeV]", 50000, 0., 400., 12, 0., 12.);
97 h_dosi_rs_edep5[i] = new TH2F(TString::Format("dosi_rs_edep5_%d", i), "Energy deposited [MeV]", 50000, 0., 400., 12, 0., 12.);
98 h_dosi_rs_edep6[i] = new TH2F(TString::Format("dosi_rs_edep6_%d", i), "Energy deposited [MeV]", 50000, 0., 400., 12, 0., 12.);
99 h_dosi_rs_edep7[i] = new TH2F(TString::Format("dosi_rs_edep7_%d", i), "Energy deposited [MeV]", 50000, 0., 400., 12, 0., 12.);
100 h_dosi_rs_edep8[i] = new TH2F(TString::Format("dosi_rs_edep8_%d", i), "Energy deposited [MeV]", 50000, 0., 400., 12, 0., 12.);
101
102 h_dosi_rs_edep0[i]->Sumw2();
103 h_dosi_rs_edep1[i]->Sumw2();
104 h_dosi_rs_edep2[i]->Sumw2();
105 h_dosi_rs_edep3[i]->Sumw2();
106 h_dosi_rs_edep4[i]->Sumw2();
107 h_dosi_rs_edep5[i]->Sumw2();
108 h_dosi_rs_edep6[i]->Sumw2();
109 h_dosi_rs_edep7[i]->Sumw2();
110 h_dosi_rs_edep8[i]->Sumw2();
111 }
112
113}
114
115
117{
118 B2INFO("DosiStudyModule: Initialize");
119
120 //read dosi xml file
121 //getXMLData();
122 //StoreArray<MCParticle> mcParticles;
123 //StoreArray<DosiSimHit> SimHits;
124 //RelationArray relMCSimHit(mcParticles, SimHits);
125
126 REG_HISTOGRAM
127
128}
129
131{
132}
133
135{
136 //Here comes the actual event processing
137 StoreArray<MCParticle> mcParticles;
140 StoreArray<SADMetaHit> MetaHits;
141
142 //Look at the meta data to extract IR rate and scattering ring section
143
144 //double rate = 0;
145 int ring_section = -1;
146 const int section_ordering[12] = {1, 12, 11, 10, 9, 8, 7, 6, 5, 4, 3, 2};
147 for (const auto& MetaHit : MetaHits) {
148 //rate = MetaHit.getrate();
149 double sad_ssraw = MetaHit.getssraw();
150 double ssraw = 0;
151 if (sad_ssraw >= 0) ssraw = sad_ssraw / 100.;
152 else ssraw = 3000. + sad_ssraw / 100.;
153 //else if (sad_ssraw < 0) ssraw = 3000. + sad_ssraw / 100.;
154 ring_section = section_ordering[(int)((ssraw) / 250.)] - 1;
155 //ring_section = MetaHit.getring_section() - 1;
156 }
157
160
161 for (const auto& SimHit : SimHits) {
162 const int detNB = SimHit.getCellId();
163 const double Edep = SimHit.getEnergyDep() * 1e3; //GeV -> MeV
164 h_dosi_edep0[detNB]->Fill(Edep);
165 h_dosi_rs_edep0[detNB]->Fill(Edep, ring_section);
166 for (const RelationElement& rel : relMCSimHit.getElementsTo(SimHit)) {
167 const MCParticle* particle = rel.from;
168 const float Mass = particle->getMass();
169 const float Energy = particle->getEnergy();
170 const float Kinetic = (Energy - Mass) * 1e3; //GeV to MeV
171 const int pdg = particle->getPDG();
172 h_dosi_edep1[detNB]->Fill(Edep);
173 if (0.005 <= Kinetic && Kinetic <= 10.0)
174 h_dosi_edep2[detNB]->Fill(Edep);
175 if (0.005 <= Kinetic && Kinetic <= 10.0 && (fabs(pdg) == Const::electron.getPDGCode()))
176 h_dosi_edep3[detNB]->Fill(Edep);
177 if (0.1 <= Kinetic && Kinetic <= 10.0 && (fabs(pdg) == Const::electron.getPDGCode()))
178 h_dosi_edep4[detNB]->Fill(Edep);
179 if (0.001 <= Kinetic && Kinetic <= 0.050)
180 h_dosi_edep5[detNB]->Fill(Edep);
181 if (pdg == Const::neutron.getPDGCode())
182 h_dosi_edep6[detNB]->Fill(Edep);
183 if (pdg == Const::neutron.getPDGCode() && (0.001 <= Kinetic && Kinetic <= 10.0))
184 h_dosi_edep7[detNB]->Fill(Edep);
185
186
187 h_dosi_rs_edep1[detNB]->Fill(Edep, ring_section);
188 if (0.005 <= Kinetic && Kinetic <= 10.0)
189 h_dosi_rs_edep2[detNB]->Fill(Edep, ring_section);
190 if (0.005 <= Kinetic && Kinetic <= 10.0 && (fabs(pdg) == Const::electron.getPDGCode()))
191 h_dosi_rs_edep3[detNB]->Fill(Edep, ring_section);
192 if (0.1 <= Kinetic && Kinetic <= 10.0 && (fabs(pdg) == Const::electron.getPDGCode()))
193 h_dosi_rs_edep4[detNB]->Fill(Edep, ring_section);
194 if (0.001 <= Kinetic && Kinetic <= 0.050)
195 h_dosi_rs_edep5[detNB]->Fill(Edep, ring_section);
196 if (pdg == Const::neutron.getPDGCode())
197 h_dosi_rs_edep6[detNB]->Fill(Edep, ring_section);
198 if (pdg == Const::neutron.getPDGCode() && (0.001 <= Kinetic && Kinetic <= 10.0))
199 h_dosi_rs_edep7[detNB]->Fill(Edep, ring_section);
200 }
201 }
202}
203/*
204//read tube centers, impulse response, and garfield drift data filename from DOSI.xml
205void DosiStudyModule::getXMLData()
206{
207 GearDir content = GearDir("/Detector/DetectorComponent[@name=\"DOSI\"]/Content/");
208 m_Ethres = content.getDouble("Ethres");
209
210 B2INFO("DosiStudy");
211
212}
213*/
215{
216
217
218}
219
221{
222}
223
224
int getPDGCode() const
PDG code.
Definition: Const.h:473
static const ParticleType neutron
neutron particle
Definition: Const.h:675
static const ChargedStable electron
electron particle
Definition: Const.h:659
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
Class to store a single element of a relation.
Provides access to fast ( O(log n) ) bi-directional lookups on a specified relation.
Definition: RelationIndex.h:76
range_to getElementsTo(const TO *to) const
Return a range of all elements pointing to the given object.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
TH1F * h_dosi_edep4[18]
Energy deposited.
TH1F * h_dosi_edep6[18]
Energy deposited.
TH1F * h_dosi_edep1[18]
Energy deposited.
virtual ~DosiStudyModule()
Destructor.
TH2F * h_dosi_rs_edep6[18]
Energy deposited.
TH2F * h_dosi_rs_edep3[18]
Energy deposited.
TH1F * h_dosi_edep7[18]
Energy deposited.
virtual void initialize() override
Initialize the Module.
TH1F * h_dosi_edep8[18]
Energy deposited.
virtual void event() override
Event processor.
TH2F * h_dosi_rs_edep5[18]
Energy deposited.
TH2F * h_dosi_rs_edep2[18]
Energy deposited.
TH1F * h_dosi_edep3[18]
Energy deposited.
TH2F * h_dosi_rs_edep7[18]
Energy deposited.
virtual void endRun() override
End-of-run action.
TH2F * h_dosi_rs_edep0[18]
Energy deposited.
virtual void terminate() override
Termination action.
virtual void beginRun() override
Called when entering a new run.
TH1F * h_dosi_edep2[18]
Energy deposited.
TH2F * h_dosi_rs_edep4[18]
Energy deposited.
DosiStudyModule()
Constructor: Sets the description, the properties and the parameters of the module.
TH2F * h_dosi_rs_edep1[18]
Energy deposited.
TH1F * h_dosi_edep5[18]
Energy deposited.
TH2F * h_dosi_rs_edep8[18]
Energy deposited.
TH1F * h_dosi_edep0[18]
reads data from DOSI.xml: tube location, drift data filename, sigma of impulse response function
virtual void defineHisto() override
Defines the histograms.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:649
Abstract base class for different kinds of events.
STL namespace.