er  dev
ERNeuRadDigitizer.h
1 /********************************************************************************
2  * Copyright (C) Joint Institute for Nuclear Research *
3  * *
4  * This software is distributed under the terms of the *
5  * GNU Lesser General Public Licence version 3 (LGPL) version 3, *
6  * copied verbatim in the file "LICENSE" *
7  ********************************************************************************/
8 
9 #ifndef ERNEURADDIGITIZER_H
10 #define ERNEURADDIGITIZER_H
11 
12 #include <vector>
13 
14 #include "TStopwatch.h"
15 
16 #include "FairTask.h"
17 
18 #include "ERNeuRadPoint.h"
19 #include "ERNeuRadPhotoElectron.h"
20 #include "ERNeuRadPixelSignal.h"
21 #include "ERNeuRadSetup.h"
22 
23 class TObjectArray;
24 class TF1;
25 
37 class ERNeuRadDigitizer : public FairTask
38 {
39 public:
42 
46  ERNeuRadDigitizer(Int_t verbose);
47 
50 
52  virtual InitStatus Init();
53 
55  virtual void Exec(Option_t* opt);
56 
58  virtual void Finish();
59 
61  virtual void Reset();
62 
64  inline void SetPixelJitter(const Double_t PixelJitter) {fPixelJitter = PixelJitter;}
65  inline void SetPixelDelay(const Double_t PixelDelay) {fPixelDelay = PixelDelay;}
66  inline void SetScincilationTau(const Double_t tau) {fScincilationTau = tau;}
67  inline void SetUseCrosstalks(const Bool_t use) {fUseCrosstalks = use;}
68 
70  Int_t PhotoElectronCount() const;
71  Int_t PixelSignalCount() const;
72 
73 protected:
74 
75  //Digitization parameters
76  ERNeuRadSetup* fNeuRadSetup;
77 
78  //Input arrays
79  TClonesArray *fNeuRadPoints;
80 
81  //Output arrays
82  TClonesArray *fNeuRadPhotoElectron;
83  TClonesArray *fNeuRadPixelSignal;
84 
85  //Event header variable
86  Int_t fPECountF;
87  Int_t fPECountB;
88  Float_t fSumAmplitudeF;
89  Float_t fSumAmplitudeB;
90 
91  // constants
92  static const Double_t cSciFiLightYield; // [photons/MeV]
93  static const Double_t cSpeedOfLight; // [cm/ns]
94  static const Double_t cMaterialSpeedOfLight; // [cm/ns]
95  static const Double_t cLightFractionInTotalIntReflection;
96  // доля света захватываемая файбером в полное внутренне отражение в каждую сторону.
97  static const Double_t cPixelDelay; // [ns]
98  static const Double_t cPixelJitter; // [ns]
99  static const Double_t cScincilationTau; // [ns]
100 
101  // Allow for user params
102  Double_t fPixelJitter; // [ns]
103  Double_t fPixelDelay; // [ns]
104  Double_t fExcessNoiseFactor;
105  Double_t fScincilationTau; // [ns]
106 
107  TStopwatch fPhotoElectronsCreatingTimer;
108  Double_t fPhotoElectronsCreatingTime;
109  TStopwatch fPixelSignalCreatingTimer;
110  Double_t fPixelSignalCreatingTime;
111  Bool_t fUseCrosstalks;
112 
113  TF1* fPEA;
114 
115 protected:
116 
117  ERNeuRadPhotoElectron* AddPhotoElectron(Int_t i_point, Int_t side,
118  Double_t lytime, Double_t cathode_time, Double_t anode_time,
119  Int_t photon_count,Double_t amplitudes);
120 
121  virtual ERNeuRadPixelSignal* AddPixelSignal(Int_t iModule, Int_t iPixel, Int_t fpoints_count, Int_t side);
122 
123  virtual void PhotoElectronsCreating(Int_t i_point, ERNeuRadPoint *point,
124  std::vector<ERNeuRadPhotoElectron* >** pePerPixels,
125  Int_t side, Int_t& sumPECount, Float_t& sumAmplitude);
126 
127  virtual void PixelSignalsCreating(Int_t iModule, Int_t iPixel,
128  std::vector<ERNeuRadPhotoElectron* >** pePerPixels, Int_t side);
129 
130  Int_t Crosstalks(Int_t pointModule, Int_t pointPixel, Int_t& peModule, Int_t& pePixel);
131 
132  static Double_t PeFunc(Double_t *x, Double_t *par);
133 
134 private:
135 
136  virtual void SetParContainers();
137 
138  ClassDef(ERNeuRadDigitizer,1)
139 };
140 
141 #endif // ERNEURADDIGITIZER_H
void SetPixelJitter(const Double_t PixelJitter)
virtual InitStatus Init()
Int_t PhotoElectronCount() const
Class for the NeuRad digital response calculation.
virtual void Exec(Option_t *opt)
The data class for storing pieces of charged tracks in sensitive volumes in NeuRad.
Definition: ERNeuRadPoint.h:22