Belle II Software development
He3DigitizerModule.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/he3tube/modules/He3DigitizerModule.h>
10#include <beast/he3tube/dataobjects/He3tubeSimHit.h>
11
12#include <mdst/dataobjects/MCParticle.h>
13#include <framework/datastore/RelationArray.h>
14#include <framework/datastore/RelationIndex.h>
15#include <framework/logging/Logger.h>
16#include <framework/gearbox/GearDir.h>
17#include <framework/gearbox/Const.h>
18
19//c++
20#include <cmath>
21#include <boost/foreach.hpp>
22#include <string>
23#include <fstream>
24#include <stdlib.h>
25
26// ROOT
27#include <TRandom.h>
28
29using namespace std;
30using namespace Belle2;
31using namespace he3tube;
32
33
34//-----------------------------------------------------------------
35// Register the Module
36//-----------------------------------------------------------------
37REG_MODULE(He3Digitizer);
38
39//-----------------------------------------------------------------
40// Implementation
41//-----------------------------------------------------------------
42
44{
45 // Set module properties
46 setDescription("He3tube digitizer module");
47
48 addParam("conversionFactor", m_ConversionFactor,
49 "Conversion factor for pulse heights", 0.303132019);
50 addParam("useMCParticles", m_mcpExist, "Use MCParticles to determine if neutrons are present", true);
51
52}
53
55{
56}
57
59{
60 B2INFO("He3Digitizer: Initializing");
61 m_he3tubeHit.registerInDataStore();
62
63 StoreArray<MCParticle> mcParticles;
65 RelationArray relMCSimHit(mcParticles, simHits);
66
67 getXMLData();
68}
69
71{
72}
73
75{
76 B2DEBUG(250, "Beginning of event method");
77
80
81 auto peak = new double[numOfTubes]();
82 auto lowTime = new double[numOfTubes](); //earliest deposit in each tube
83 std::fill_n(lowTime, numOfTubes, 9999999999999);
84 auto edepDet = new double[numOfTubes](); //total energy deposited per tube
85 auto NbEle_tot = new double[numOfTubes]();
86
87 //skip events with no He3SimHits, but continue the event counter
88 if (He3SimHits.getEntries() == 0) {
89 B2DEBUG(250, "Skipping event #" << Event << " since there are no sim hits");
90 Event++;
91 return;
92 }
93
94 auto definiteNeutron = new bool[numOfTubes]();
95 //vector<int> neutronIDs;
96
97 if (m_mcpExist) { //if the mcparticles exist, we know which simhits are from neutron events
98
99 StoreArray<MCParticle> mcParticles;
100 RelationIndex<MCParticle, He3tubeSimHit> relMCSimHit(mcParticles, He3SimHits);
101
102 int nMCParticles = mcParticles.getEntries();
103 for (int i = 0; i < nMCParticles; ++i) {
104 MCParticle& mcp = *mcParticles[i];
105
106 //Find all He3tubeSimHits which point from that MCParticle using a typedef and BOOST_FOREACH
107 //The typedef is needed as BOOST_FOREACH is a macro and cannot handle anything including a comma
108 typedef RelationIndex<MCParticle, He3tubeSimHit>::Element relMCSimHit_Element;
109 BOOST_FOREACH(const relMCSimHit_Element & relation, relMCSimHit.getElementsFrom(mcp)) {
110
111 //int processNum = mcp.getSecondaryPhysicsProcess();
112
113 MCParticle* mom = mcp.getMother();
114
115 if (mom == NULL) continue;
116
117 //if(processNum==121){
118 if (mom->getPDG() == Const::neutron.getPDGCode()) {
119 He3tubeSimHit* aHit = He3SimHits[relation.indexTo];
120 int detNB = aHit->getdetNb();
121 definiteNeutron[detNB] = true;
122 ProcessHit(aHit, lowTime, edepDet, NbEle_tot);
123 }
124 }
125 }
126
127
128
129 } else { //if mcparticles don't exist, we don't know which simhits are from neutron events
130 int nentries = He3SimHits.getEntries();
131 for (int i = 0; i < nentries; i++) {
132 He3tubeSimHit* aHit = He3SimHits[i];
133 int pdg = aHit->gettkPDG();
134 if (pdg == Const::proton.getPDGCode() || pdg == 1000010030) ProcessHit(aHit, lowTime, edepDet, NbEle_tot);
135 }
136
137 }
138
139
140
141 for (int i = 0; i < numOfTubes; i++) {
142 peak[i] = (double)NbEle_tot[i] * m_ConversionFactor; //convert the number of ions into a value in ADC counts
143 if (peak[i] != 0) {
144
145 //create He3tubeHit
146 He3Hits.appendNew(He3tubeHit(edepDet[i], i, peak[i], lowTime[i], definiteNeutron[i]));
147
148 B2DEBUG(80, "He3Digitizer: " << edepDet[i] << "MeV deposited in tube #" << i << " with waveform peak of " << peak[i]);
149 }
150 }
151
152 Event++;
153
154 //delete arrays
155 delete[] peak;
156 delete[] lowTime;
157 delete[] edepDet;
158 delete[] NbEle_tot;
159 delete[] definiteNeutron;
160
161}
162
163
164void He3DigitizerModule::ProcessHit(He3tubeSimHit* aHit, double* lowTime, double* edepDet, double* NbEle_tot)
165{
166
167 //get various info from this
168 double edep = aHit->getEnergyDep();
169 double niel = aHit->getEnergyNiel();
170 int detNB = aHit->getdetNb();
171 double time = aHit->getGlTime();
172 if (time < lowTime[detNB]) lowTime[detNB] =
173 time; //get the earliest time for this simhit
174
175 if (detNB >= numOfTubes) {
176 B2WARNING("He3Digitizer: Detector number of He3tubeSimHit is greater than number implemented! Ignoring He3tubeSimHit.");
177 return;
178 }
179
180 //ionization energy
181 double ionEn = (edep - niel) * 1e3; //MeV to keV
182
183 if ((ionEn * 1e3) < m_Workfct) return; //keV to eV
184 //if there is enough energy to ionize the gas, get the number if ions produced
185 if ((ionEn * 1e3) > m_Workfct) {
186
187 double meanE1 = ionEn * 1e3 / m_Workfct;
188 double sigma = sqrt(m_Fanofac * meanE1);
189 const int NbEle = (int)gRandom->Gaus(meanE1,
190 sigma); //this is the number of electrons created, and the number which reach the wire, since gain is essentially 1.
191
192 NbEle_tot[detNB] = NbEle_tot[detNB] + NbEle;
193
194 }
195
196 edepDet[detNB] = edepDet[detNB] + edep;
197
198}
199
200
202{
203 GearDir content = GearDir("/Detector/DetectorComponent[@name=\"HE3TUBE\"]/Content/");
204
205 //get the location of the tubes
206 BOOST_FOREACH(const GearDir & activeParams, content.getNodes("Active")) {
207 B2DEBUG(250, "Tubes located at x=" << activeParams.getLength("x_he3tube"));
208 numOfTubes++;
209 }
210
211
212 B2INFO("He3Digitizer: There are " << numOfTubes << " tubes implemented");
213
214}
215
217{
218}
219
221{
222}
int getPDGCode() const
PDG code.
Definition: Const.h:473
static const ParticleType neutron
neutron particle
Definition: Const.h:675
static const ChargedStable proton
proton particle
Definition: Const.h:663
GearDir is the basic class used for accessing the parameter store.
Definition: GearDir.h:31
ClassHe3Hit - digitization simulated hit for the He3tube detector.
Definition: He3tubeHit.h:26
ClassHe3tubeSimHit - Geant4 simulated hit for the He3tube detector.
Definition: He3tubeSimHit.h:30
float getEnergyDep() const
Return the energy deposition in electrons.
Definition: He3tubeSimHit.h:63
float getdetNb() const
Return the He3tube number.
Definition: He3tubeSimHit.h:71
float getEnergyNiel() const
Return the non-ionization energy in electrons.
Definition: He3tubeSimHit.h:65
int gettkPDG() const
Return the PDG number of the track.
Definition: He3tubeSimHit.h:67
float getGlTime() const
Return the global time.
Definition: He3tubeSimHit.h:73
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:112
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
Low-level class to create/modify relations between StoreArrays.
Definition: RelationArray.h:62
Provides access to fast ( O(log n) ) bi-directional lookups on a specified relation.
Definition: RelationIndex.h:76
range_from getElementsFrom(const FROM *from) const
Return a range of all elements pointing from the given object.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
T * appendNew()
Construct a new T object at the end of the array.
Definition: StoreArray.h:246
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
double getLength(const std::string &path="") const noexcept(false)
Get the parameter path as a double converted to the standard length unit.
Definition: Interface.h:259
StoreArray< He3tubeHit > m_he3tubeHit
Array for He3tubeHit.
void ProcessHit(He3tubeSimHit *aHit, double *lowTime, double *edepDet, double *NbEle_tot)
Process the he3tube simhits.
virtual void initialize() override
Initialize the Module.
virtual void event() override
Event processor.
virtual void endRun() override
End-of-run action.
virtual void getXMLData()
reads data from HE3TUBE.xml: tube location, drift data filename, sigma of impulse response function
bool m_mcpExist
Whether or not mcparticle array exists.
virtual void terminate() override
Termination action.
double m_ConversionFactor
Conversion to ADC counts, set in steering file.
virtual void beginRun() override
Called when entering a new run.
double m_Workfct
ionization energy of He3
He3DigitizerModule()
Constructor: Sets the description, the properties and the parameters of the module.
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
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
MCParticle * getMother() const
Returns a pointer to the mother particle.
Definition: MCParticle.h:600
Abstract base class for different kinds of events.
STL namespace.