Belle II Software  release-08-01-10
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 
33 namespace 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 
57  class SVDDigitizerModule : public Module {
58  public:
61 
63  void processHit();
64 
71  void driftCharge(const ROOT::Math::XYZVector& position, double carriers, SVD::SensorInfo::CarrierType carrierType);
72 
78  double addNoise(double charge, double noise);
79 
83  void saveDigits();
84 
88  void saveWaveforms();
89 
93  void saveSignals();
94 
96  virtual void initialize() override;
98  virtual void beginRun() override;
100  virtual void event() override;
102  virtual void terminate() override;
103 
104  protected:
105 
106  // Members holding module parameters:
107 
108  // 1. Collections
110  std::string m_storeMCParticlesName = "";
112  std::string m_storeSimHitsName = "";
114  std::string m_storeTrueHitsName = "";
116  std::string m_relMCParticleSimHitName = "";
118  std::string m_relTrueHitSimHitName = "";
120  std::string m_storeShaperDigitsName = "";
124  std::string m_relShaperDigitTrueHitName = "";
126  std::string m_svdEventInfoName = "SVDEventInfoSim";
127 
128  // 2. Physics
130  double m_segmentLength = 0.020;
131 
132  // 3. Noise
134  bool m_applyPoisson = true;
136  double m_SNAdjacent = 3.0;
138  bool m_roundZS = true;
142  double m_noiseFraction = 0.01;
143 
144  // 4. Timing
148  double m_betaPrimeDecayTimeU = 250.0;
150  double m_betaPrimeDecayTimeV = 250.0;
152  double m_samplingTime = -1;
156  double m_startSampling = -2;
160  double m_initTime = 0;
161 
166  bool m_randomizeEventTimes = false;
168  float m_minTimeFrame = -300;
170  float m_maxTimeFrame = 150;
174  float m_currentEventTime = 0.0;
175 
176  // 5. 3-mixed-6 and 3-sample daqMode
178  bool m_is3sampleEvent = false;
181 
182  // 6. Reporting
184  std::string m_rootFilename = "";
186  bool m_storeWaveforms = false;
188  std::string m_signalsList = "";
189 
190 
191  // Other data members:
192 
195 
197  const SVDSimHit* m_currentHit = nullptr;
205  const SensorInfo* m_currentSensorInfo = nullptr;
207  double m_currentTime = 0;
209  double m_sensorThickness = 0.03;
210 
215 
217  int getFirstSample(int triggerBin, int relativShift);
218 
219  //payloads:
221  static std::string m_xmlFileName ;
223  std::unique_ptr<SVDOnlineToOfflineMap> m_map;
229  // ROOT stuff:
231  TFile* m_rootFile = nullptr;
233  TH1D* m_histChargeSharing_u = nullptr;
235  TH1D* m_histChargeSharing_v = nullptr;
236 
238  TH1D* m_histMobility_e = nullptr;
240  TH1D* m_histMobility_h = nullptr;
242  TH1D* m_histVelocity_e = nullptr;
244  TH1D* m_histVelocity_h = nullptr;
246  TH1D* m_histDistanceToPlane_e = nullptr;
248  TH1D* m_histDistanceToPlane_h = nullptr;
250  TH1D* m_histDriftTime_e = nullptr;
252  TH1D* m_histDriftTime_h = nullptr;
254  TH1D* m_histHitTime = nullptr;
256  TH2F* m_histHitTimeTB = nullptr;
257 
259  TH1D* m_histLorentz_u = nullptr;
261  TH1D* m_histLorentz_v = nullptr;
263  TH1D* m_signalDist_u = nullptr;
265  TH1D* m_signalDist_v = nullptr;
267  TTree* m_waveTree = nullptr;
268 
269  };//end class declaration
270 
271 
272  } // end namespace SVD
274 } // end namespace Belle2
275 
276 #endif // SVDDigitizerModule_H
Specialization of DBObjPtr in case of PayloadFiles.
Definition: PayloadFile.h:54
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.