Belle II Software  release-06-02-00
ScintillatorSimulator.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 #pragma once
10 
11 /* KLM headers. */
12 #include <klm/dataobjects/bklm/BKLMSimHit.h>
13 #include <klm/dataobjects/eklm/EKLMSimHit.h>
14 #include <klm/dbobjects/KLMScintillatorDigitizationParameters.h>
15 #include <klm/dbobjects/KLMScintillatorFEEData.h>
16 #include <klm/simulation/ScintillatorFirmware.h>
17 #include <klm/time/KLMTime.h>
18 
19 namespace Belle2 {
25  namespace KLM {
26 
31 
32  public:
33 
35  struct Photoelectron {
36  int bin;
37  double expTime;
38  bool isReflected;
39  };
40 
49  ScintillatorFirmware* fitter,
50  double digitizationInitialTime,
51  bool debug);
52 
57 
62 
67 
73  void simulate(
74  const std::multimap<KLMChannelNumber, const BKLMSimHit*>::iterator& firstHit,
75  const std::multimap<KLMChannelNumber, const BKLMSimHit*>::iterator& end);
76 
82  void simulate(
83  const std::multimap<KLMChannelNumber, const EKLMSimHit*>::iterator& firstHit,
84  const std::multimap<KLMChannelNumber, const EKLMSimHit*>::iterator& end);
85 
90 
95  enum ScintillatorFirmwareFitStatus getFitStatus() const;
96 
100  double getNPhotoelectrons();
101 
106 
110  double getEnergy();
111 
115  void setFEEData(const KLMScintillatorFEEData* FEEData);
116 
125  void generatePhotoelectrons(double stripLen, double distSiPM,
126  int nPhotons, double timeShift,
127  bool isReflected);
128 
135  void fillSiPMOutput(float* hist, bool useDirect, bool useReflected);
136 
141  float getMCTime() const
142  {
143  return m_MCTime;
144  }
145 
150  float getSiPMMCTime() const
151  {
152  return m_SiPMMCTime;
153  }
154 
155  private:
156 
161  void reallocPhotoElectronBuffers(int size);
162 
166  void prepareSimulation();
167 
171  void performSimulation();
172 
178  int* sortPhotoelectrons(int nPhotoelectrons);
179 
183  void addRandomSiPMNoise();
184 
188  void simulateADC();
189 
193  void debugOutput();
194 
196  const KLMTime* m_Time;
197 
200 
203 
206 
208  bool m_Debug;
209 
211  double m_histRange;
212 
215 
218 
220  float* m_amplitude;
221 
224 
227 
230 
233 
236 
239 
242 
244  enum ScintillatorFirmwareFitStatus m_FPGAStat;
245 
248 
250  int m_npe;
251 
253  double m_Energy;
254 
256  std::string m_stripName;
257 
259  double m_Pedestal;
260 
263 
266 
268  float m_MCTime;
269 
272 
273  };
274 
275  }
276 
278 }
Class to store KLM scintillator simulation parameters in the database.
KLM time conversion.
Definition: KLMTime.h:27
Digitize EKLMSim2Hits to get EKLM StripHits.
float getSiPMMCTime() const
Get SiPM MC time.
void setFEEData(const KLMScintillatorFEEData *FEEData)
Set FEE data.
int * m_PhotoelectronIndex2
Buffer for photoelectron indices.
float * m_amplitudeReflected
Analog amplitude (reflected).
void debugOutput()
Debug output (signal and fit result histograms).
enum ScintillatorFirmwareFitStatus getFitStatus() const
Get fit status.
double * m_SignalTimeDependence
Buffer for signal time dependence calculation.
void generatePhotoelectrons(double stripLen, double distSiPM, int nPhotons, double timeShift, bool isReflected)
Generate photoelectrons.
double * m_SignalTimeDependenceDiff
Buffer for signal time dependence calculation.
double m_Energy
Total energy deposited in the strip.
float * m_amplitudeDirect
Analog amplitude (direct).
ScintillatorSimulator(const ScintillatorSimulator &)=delete
Copy constructor (disabled).
void addRandomSiPMNoise()
Add random noise to the signal (amplitude-dependend).
void simulateADC()
Simulate ADC (create digital signal from analog),.
void prepareSimulation()
Prepare simulation.
void simulate(const std::multimap< KLMChannelNumber, const BKLMSimHit * >::iterator &firstHit, const std::multimap< KLMChannelNumber, const BKLMSimHit * >::iterator &end)
Simulate BKLM strip.
double getEnergy()
Get total energy deposited in the strip (sum over ssimulation hits).
double m_histRange
Time range, (number of digitizations) * (ADC sampling time).
void reallocPhotoElectronBuffers(int size)
Reallocate photoelectron buffers.
bool m_Debug
Debug mode (generates additional output files with histograms).
double m_DigitizationInitialTime
Initial digitization time.
double getNPhotoelectrons()
Get number of photoelectrons (fit result).
ScintillatorSimulator & operator=(const ScintillatorSimulator &)=delete
Operator = (disabled).
struct Photoelectron * m_Photoelectrons
Buffer for photoelectron data.
ScintillatorFirmware * m_fitter
Fitter.
std::string m_stripName
Name of the strip.
enum ScintillatorFirmwareFitStatus m_FPGAStat
FPGA fit status.
KLMScintillatorFirmwareFitResult * getFPGAFit()
Get fit data.
KLMScintillatorFirmwareFitResult m_FPGAFit
FPGA fit data.
void performSimulation()
Perform common simulation stage.
double m_PhotoelectronAmplitude
Photoelectron amplitude.
int getNGeneratedPhotoelectrons()
Get generated number of photoelectrons.
int * m_PhotoelectronIndex
Buffer for photoelectron indices.
void fillSiPMOutput(float *hist, bool useDirect, bool useReflected)
Fill SiPM output.
int * sortPhotoelectrons(int nPhotoelectrons)
Sort photoelectrons.
int m_npe
Number of photoelectrons (generated).
int m_PhotoelectronBufferSize
Size of photoelectron data buffer.
ScintillatorSimulator(const KLMScintillatorDigitizationParameters *digPar, ScintillatorFirmware *fitter, double digitizationInitialTime, bool debug)
Constructor.
const KLMScintillatorDigitizationParameters * m_DigPar
Parameters.
Abstract base class for different kinds of events.
int bin
Hit time bin in ADC output histogram.
bool isReflected
Direct (false) or reflected (true).
double expTime
exp(-m_DigPar->PEAttenuationFreq * (-time))