Belle II Software development
SVDDigitizerModule.h
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#ifndef SVDDigitizerModule_H
10#define SVDDigitizerModule_H
11
12#include <framework/core/Module.h>
13#include <svd/dataobjects/SVDSimHit.h>
14#include <svd/simulation/SVDWaveform.h>
15#include <svd/geometry/SensorInfo.h>
16#include <svd/calibration/SVDNoiseCalibrations.h>
17#include <svd/calibration/SVDPulseShapeCalibrations.h>
18#include <svd/calibration/SVDChargeSimulationCalibrations.h>
19#include <svd/calibration/SVDFADCMaskedStrips.h>
20#include <svd/online/SVDOnlineToOfflineMap.h>
21#include <framework/database/PayloadFile.h>
22#include <svd/dataobjects/SVDEventInfo.h>
23#include <framework/dbobjects/HardwareClockSettings.h>
24
25#include <string>
26
27#include <Math/Vector3D.h>
28#include <root/TFile.h>
29#include <root/TTree.h>
30#include <root/TH1D.h>
31#include <root/TH2F.h>
32
33namespace Belle2 {
38 namespace SVD {
39
41 typedef std::map<short int, SVDWaveform> StripWaveforms;
42
44 typedef std::pair<StripWaveforms, StripWaveforms> SensorWaveforms;
45
47 typedef std::map<VxdID, SensorWaveforms> Waveforms;
48
54 class SVDDigitizerModule : public Module {
55 public:
58
60 void processHit();
61
68 void driftCharge(const ROOT::Math::XYZVector& position, double carriers, SVD::SensorInfo::CarrierType carrierType);
69
75 double addNoise(double charge, double noise);
76
80 void saveDigits();
81
85 void saveWaveforms();
86
90 void saveSignals();
91
93 virtual void initialize() override;
95 virtual void beginRun() override;
97 virtual void event() override;
99 virtual void terminate() override;
100
101 protected:
102
103 // Members holding module parameters:
104
105 // 1. Collections
107 std::string m_storeMCParticlesName = "";
109 std::string m_storeSimHitsName = "";
111 std::string m_storeTrueHitsName = "";
115 std::string m_relTrueHitSimHitName = "";
117 std::string m_storeShaperDigitsName = "";
123 std::string m_svdEventInfoName = "SVDEventInfoSim";
124
125 // 2. Physics
127 double m_segmentLength = 0.020;
128
129 // 3. Noise
131 bool m_applyPoisson = true;
133 double m_SNAdjacent = 3.0;
135 bool m_roundZS = true;
139 double m_noiseFraction = 0.01;
140
141 // 4. Timing
145 double m_betaPrimeDecayTimeU = 250.0;
147 double m_betaPrimeDecayTimeV = 250.0;
149 double m_samplingTime = -1;
153 double m_startSampling = -2;
157 double m_initTime = 0;
158
165 float m_minTimeFrame = -300;
167 float m_maxTimeFrame = 150;
172
173 // 5. 3-mixed-6 and 3-sample daqMode
175 bool m_is3sampleEvent = false;
178
179 // 6. Reporting
181 std::string m_rootFilename = "";
183 bool m_storeWaveforms = false;
185 std::string m_signalsList = "";
186
187
188 // Other data members:
189
192
194 const SVDSimHit* m_currentHit = nullptr;
204 double m_currentTime = 0;
206 double m_sensorThickness = 0.03;
207
212
214 int getFirstSample(int triggerBin, int relativShift);
215
216 //payloads:
218 static std::string m_xmlFileName ;
220 std::unique_ptr<SVDOnlineToOfflineMap> m_map;
226 // ROOT stuff:
228 TFile* m_rootFile = nullptr;
230 TH1D* m_histChargeSharing_u = nullptr;
232 TH1D* m_histChargeSharing_v = nullptr;
233
235 TH1D* m_histMobility_e = nullptr;
237 TH1D* m_histMobility_h = nullptr;
239 TH1D* m_histVelocity_e = nullptr;
241 TH1D* m_histVelocity_h = nullptr;
243 TH1D* m_histDistanceToPlane_e = nullptr;
245 TH1D* m_histDistanceToPlane_h = nullptr;
247 TH1D* m_histDriftTime_e = nullptr;
249 TH1D* m_histDriftTime_h = nullptr;
251 TH1D* m_histHitTime = nullptr;
253 TH2F* m_histHitTimeTB = nullptr;
254
256 TH1D* m_histLorentz_u = nullptr;
258 TH1D* m_histLorentz_v = nullptr;
260 TH1D* m_signalDist_u = nullptr;
262 TH1D* m_signalDist_v = nullptr;
264 TTree* m_waveTree = nullptr;
265
266 };//end class declaration
267
268
269 } // end namespace SVD
271} // end namespace Belle2
272
273#endif // SVDDigitizerModule_H
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
Base class for Modules.
Definition: Module.h:72
This class defines the dbobject and the methods to access SVD simulation calibrations; coupling const...
This class defines the dbobject and the method to access strips which are masked at FADC level.
This class defines the dbobject and the method to access SVD calibrations from the noise local runs.
This class defines the dbobject and the methods to access the SVD calibrations from the local runs pr...
Class SVDSimHit - Geant4 simulated hit for the SVD.
Definition: SVDSimHit.h:26
The SVD Digitizer module.
double m_startSampling
Time window start, excluding trigger bin effect.
TH1D * m_histHitTime
Histogram showing the hit time.
TH1D * m_histDistanceToPlane_e
Histogram showing the distance to plane for e.
TH1D * m_histChargeSharing_u
Histogram showing the charge sharing + diffusion in u (r-phi).
double m_sensorThickness
Thickness of current sensor (read from m_currentSensorInfo)
void processHit()
Process one SVDSimHit by dividing the step in smaller steps and drifting the charge.
TH1D * m_histVelocity_h
Histogram showing the velocity of h.
SVDFADCMaskedStrips m_MaskedStr
FADC masked strip payload.
TH1D * m_histVelocity_e
Histogram showing the velocity of e-.
Waveforms m_waveforms
Structure containing waveforms in all existing sensors.
double m_samplingTime
Interval between two waveform samples, by default taken from HardwareClockSettings.
std::string m_relTrueHitSimHitName
Name of the relation between SVDTrueHits and SVDSimHits.
TH1D * m_histLorentz_v
Histogram showing the Lorentz angles in v (z).
void saveDigits()
Save digits to the DataStore Saves samples of generated waveforms.
bool m_randomizeEventTimes
Randomize event times? If set to true, event times will be randomized uniformly from m_minTimeFrame t...
int getFirstSample(int triggerBin, int relativShift)
return the starting sample
int m_startingSample
Starting sample for the selection of 3 samples in 3-mixed-6.
int m_currentParticle
Index of the particle which caused the current hit.
std::string m_relShaperDigitMCParticleName
Name of the relation between SVDShaperDigits and MCParticles.
int m_nSamplesOverZS
Keep digit if at least m_nSamplesOverZS are over threshold.
TH1D * m_histChargeSharing_v
Histogram showing the charge sharing + diffusion in v (z).
TH1D * m_histMobility_h
Histogram showing the mobility of h.
virtual void initialize() override
Initialize the module and check module parameters.
bool m_roundZS
Round ZS cut to nearest ADU.
std::string m_storeShaperDigitsName
Name of the collection for the SVDShaperDigits.
virtual void event() override
Digitize one event.
bool m_is3sampleEvent
True if the event should be simulated with 3 sample.
double m_SNAdjacent
Zero-suppression cut.
double m_currentTime
Time of the current SimHit.
TFile * m_rootFile
Pointer to the ROOT filename for statistics.
SVDNoiseCalibrations m_NoiseCal
SVDNoise calibrations db object.
SensorWaveforms * m_currentSensorWaveforms
Pointer to the sensor in which the current hit occurred.
std::string m_storeTrueHitsName
Name of the collection for the SVDTrueHits.
TH1D * m_signalDist_u
Histogram showing the distribution of digit signals in u (r-phi).
DBObjPtr< PayloadFile > m_mapping
channel mapping payload
virtual void terminate() override
Terminate the module.
void saveSignals()
Save signals to a root-delimited file (to be analyzed in Python).
std::string m_storeMCParticlesName
Name of the collection for the MCParticles.
TH1D * m_signalDist_v
Histogram showing the distribution of digit signals in v (z).
double m_noiseFraction
(derived from SNAdjacent) Fraction of noisy strips per sensor.
static std::string m_xmlFileName
< channel mapping xml filename
double m_betaPrimeDecayTimeU
Decay time of betaprime waveform U-side.
SVDPulseShapeCalibrations m_PulseShapeCal
SVDPulseShapeCalibrations calibrations db object.
TH1D * m_histMobility_e
Histogram showing the mobility of e-.
std::string m_relShaperDigitTrueHitName
Name of the relation between SVDShaperDigits and SVDTrueHits.
SVDChargeSimulationCalibrations m_ChargeSimCal
SVDChargeSimulationCalibrations calibrations db object.
int m_nAPV25Samples
number of digitized samples read from SVDEventInfo
TTree * m_waveTree
Tree for waveform storage.
bool m_applyPoisson
Whether or not to apply poisson fluctuation of charge (Fano factor)
virtual void beginRun() override
Initialize the list of existing SVD Sensors.
TH1D * m_histDistanceToPlane_h
Histogram showing the distance to plane for h.
TH1D * m_histDriftTime_e
Histogram showing the drift time of e.
DBObjPtr< HardwareClockSettings > m_hwClock
Hardware Clocks.
std::string m_svdEventInfoName
Name of the SVDEventInfo object.
std::string m_signalsList
Name of the tab-delimited listing of waveforms.
double addNoise(double charge, double noise)
Calculate the noise contribution to one strip with given charge.
const SensorInfo * m_currentSensorInfo
Pointer to the SensorInfo of the current sensor.
TH1D * m_histLorentz_u
Histogram showing the Lorentz angles in u (r-phi).
int m_relativeShift
relative shift in SVDEventInfo obj
TH2F * m_histHitTimeTB
Histogram showing the hit time vs TB.
bool m_storeWaveforms
Store waveform data in the reporting file?
void driftCharge(const ROOT::Math::XYZVector &position, double carriers, SVD::SensorInfo::CarrierType carrierType)
Drift the charge inside the silicon.
void saveWaveforms()
Save waveforms to the statistics file.
std::unique_ptr< SVDOnlineToOfflineMap > m_map
channel mapping map
float m_maxTimeFrame
High edge of randomization time frame.
TH1D * m_histDriftTime_h
Histogram showing the drift time of h.
double m_initTime
Time window start, including the triggerBin effect.
int m_currentTrueHit
Index of the TrueHit the current hit belongs to.
float m_minTimeFrame
Low edge of randomization time frame.
const SVDSimHit * m_currentHit
Pointer to the SVDSimhit currently digitized.
std::string m_relMCParticleSimHitName
Name of the relation between MCParticles and SVDSimHits.
double m_betaPrimeDecayTimeV
Decay time of betaprime waveform V-side.
std::string m_rootFilename
Name of the ROOT filename to output statistics.
std::string m_storeSimHitsName
Name of the collection for the SVDSimhits.
float m_currentEventTime
Current event time.
Specific implementation of SensorInfo for SVD Sensors which provides additional sensor specific infor...
Definition: SensorInfo.h:25
CarrierType
Enum to flag charge carriers.
Definition: SensorInfo.h:38
std::map< VxdID, SensorWaveforms > Waveforms
Map of all waveforms in all sensors.
std::pair< StripWaveforms, StripWaveforms > SensorWaveforms
Waveforms of u- and v- channels in one sensor.
std::map< short int, SVDWaveform > StripWaveforms
Map of all channels' waveforms in one sensor side.
Abstract base class for different kinds of events.