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 <string>
22#include <stdlib.h>
23
24// ROOT
25#include <TRandom.h>
26
27using namespace std;
28using namespace Belle2;
29using namespace he3tube;
30
31
32//-----------------------------------------------------------------
33// Register the Module
34//-----------------------------------------------------------------
35REG_MODULE(He3Digitizer);
36
37//-----------------------------------------------------------------
38// Implementation
39//-----------------------------------------------------------------
40
42{
43 // Set module properties
44 setDescription("He3tube digitizer module");
45
46 addParam("conversionFactor", m_ConversionFactor,
47 "Conversion factor for pulse heights", 0.303132019);
48 addParam("useMCParticles", m_mcpExist, "Use MCParticles to determine if neutrons are present", true);
49
50}
51
53{
54}
55
57{
58 B2INFO("He3Digitizer: Initializing");
59 m_he3tubeHit.registerInDataStore();
60
61 StoreArray<MCParticle> mcParticles;
63 RelationArray relMCSimHit(mcParticles, simHits);
64
65 getXMLData();
66}
67
69{
70}
71
73{
74 B2DEBUG(250, "Beginning of event method");
75
78
79 //skip events with no He3SimHits, but continue the event counter
80 if (He3SimHits.getEntries() == 0) {
81 B2DEBUG(250, "Skipping event #" << Event << " since there are no sim hits");
82 Event++;
83 return;
84 }
85
86 auto peak = new double[numOfTubes]();
87 auto lowTime = new double[numOfTubes](); //earliest deposit in each tube
88 std::fill_n(lowTime, numOfTubes, 9999999999999);
89 auto edepDet = new double[numOfTubes](); //total energy deposited per tube
90 auto NbEle_tot = new double[numOfTubes]();
91
92 auto definiteNeutron = new bool[numOfTubes]();
93 //vector<int> neutronIDs;
94
95 if (m_mcpExist) { //if the mcparticles exist, we know which simhits are from neutron events
96
97 StoreArray<MCParticle> mcParticles;
98 RelationIndex<MCParticle, He3tubeSimHit> relMCSimHit(mcParticles, He3SimHits);
99
100 int nMCParticles = mcParticles.getEntries();
101 for (int i = 0; i < nMCParticles; ++i) {
102 MCParticle& mcp = *mcParticles[i];
103
104 // Find all He3tubeSimHits which point from that MCParticle.
105 typedef RelationIndex<MCParticle, He3tubeSimHit>::Element relMCSimHit_Element;
106 for (const relMCSimHit_Element& relation : relMCSimHit.getElementsFrom(mcp)) {
107
108 //int processNum = mcp.getSecondaryPhysicsProcess();
109
110 MCParticle* mom = mcp.getMother();
111
112 if (mom == NULL) continue;
113
114 //if(processNum==121){
115 if (mom->getPDG() == Const::neutron.getPDGCode()) {
116 He3tubeSimHit* aHit = He3SimHits[relation.indexTo];
117 int detNB = aHit->getdetNb();
118 definiteNeutron[detNB] = true;
119 ProcessHit(aHit, lowTime, edepDet, NbEle_tot);
120 }
121 }
122 }
123
124
125
126 } else { //if mcparticles don't exist, we don't know which simhits are from neutron events
127 int nentries = He3SimHits.getEntries();
128 for (int i = 0; i < nentries; i++) {
129 He3tubeSimHit* aHit = He3SimHits[i];
130 int pdg = aHit->gettkPDG();
131 if (pdg == Const::proton.getPDGCode() || pdg == 1000010030) ProcessHit(aHit, lowTime, edepDet, NbEle_tot);
132 }
133
134 }
135
136
137
138 for (int i = 0; i < numOfTubes; i++) {
139 peak[i] = (double)NbEle_tot[i] * m_ConversionFactor; //convert the number of ions into a value in ADC counts
140 if (peak[i] != 0) {
141
142 //create He3tubeHit
143 He3Hits.appendNew(He3tubeHit(edepDet[i], i, peak[i], lowTime[i], definiteNeutron[i]));
144
145 B2DEBUG(80, "He3Digitizer: " << edepDet[i] << "MeV deposited in tube #" << i << " with waveform peak of " << peak[i]);
146 }
147 }
148
149 Event++;
150
151 //delete arrays
152 delete[] peak;
153 delete[] lowTime;
154 delete[] edepDet;
155 delete[] NbEle_tot;
156 delete[] definiteNeutron;
157
158}
159
160
161void He3DigitizerModule::ProcessHit(He3tubeSimHit* aHit, double* lowTime, double* edepDet, double* NbEle_tot)
162{
163
164 //get various info from this
165 double edep = aHit->getEnergyDep();
166 double niel = aHit->getEnergyNiel();
167 int detNB = aHit->getdetNb();
168 double time = aHit->getGlTime();
169 if (time < lowTime[detNB]) lowTime[detNB] =
170 time; //get the earliest time for this simhit
171
172 if (detNB >= numOfTubes) {
173 B2WARNING("He3Digitizer: Detector number of He3tubeSimHit is greater than number implemented! Ignoring He3tubeSimHit.");
174 return;
175 }
176
177 //ionization energy
178 double ionEn = (edep - niel) * 1e3; //MeV to keV
179
180 if ((ionEn * 1e3) < m_Workfct) return; //keV to eV
181 //if there is enough energy to ionize the gas, get the number if ions produced
182 if ((ionEn * 1e3) > m_Workfct) {
183
184 double meanE1 = ionEn * 1e3 / m_Workfct;
185 double sigma = sqrt(m_Fanofac * meanE1);
186 const int NbEle = (int)gRandom->Gaus(meanE1,
187 sigma); //this is the number of electrons created, and the number which reach the wire, since gain is essentially 1.
188
189 NbEle_tot[detNB] = NbEle_tot[detNB] + NbEle;
190
191 }
192
193 edepDet[detNB] = edepDet[detNB] + edep;
194
195}
196
197
199{
200 GearDir content = GearDir("/Detector/DetectorComponent[@name=\"HE3TUBE\"]/Content/");
201
202 //get the location of the tubes
203 for (const GearDir& activeParams : content.getNodes("Active")) {
204 B2DEBUG(250, "Tubes located at x=" << activeParams.getLength("x_he3tube"));
205 numOfTubes++;
206 }
207
208
209 B2INFO("He3Digitizer: There are " << numOfTubes << " tubes implemented");
210
211}
212
214{
215}
216
218{
219}
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:101
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
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:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:649
double sqrt(double a)
sqrt for double
Definition: beamHelpers.h:28
MCParticle * getMother() const
Returns a pointer to the mother particle.
Definition: MCParticle.h:590
Abstract base class for different kinds of events.
STL namespace.