GEMC  2.3
Geant4 Monte-Carlo Framework
IC_hitprocess.cc
Go to the documentation of this file.
1 // G4 headers
2 #include "G4Poisson.hh"
3 #include "Randomize.hh"
4 
5 // gemc headers
6 #include "IC_hitprocess.h"
7 
8 map<string, double> IC_HitProcess :: integrateDgt(MHit* aHit, int hitn)
9 {
10  map<string, double> dgtz;
11  vector<identifier> identity = aHit->GetId();
12 
13 
14  // R.De Vita (April 2009)
15 
16  // relevant parameter for digitization (in the future should be read from database)
17  double tdc_time_to_channel=20; // conversion factor from time(ns) to TDC channels)
18  double tdc_max=4095; // TDC range
19  double adc_charge_tochannel=20; // conversion factor from charge(pC) to ADC channels
20  double PbWO4_light_yield =120/MeV; // Lead Tungsten Light Yield (PAD have similar QE for
21  // fast component, lambda=420nm-ly=120ph/MeV, and slow component,
22  // lambda=560nm-ly=20ph/MeV, taking fast component only)
23  double APD_qe = 0.75; // APD Quantum Efficiency (Hamamatsu S8664-55)
24  double APD_size = 25*mm*mm; // APD size ( 5 mm x 5 mm)
25  double APD_gain = 250; // based on IC note
26  double APD_noise = 0.006; // relative noise based on a Voltage and Temperature stability of 30 mV (10%/V) and 0.1 C (-5%/C)
27  double AMP_input_noise = 4500; // preamplifier input noise in number of electrons
28  double AMP_gain = 2500; // preamplifier gain = 5V/pC x 25 ns (tipical signal duration)
29  double light_speed =15;
30 
31  // Get the crystal length: in the IC crystal are trapezoid (TRD) and the half-length is the 5th element
32  double length = 2 * aHit->GetDetector().dimensions[4];
33  // Get the crystal width (rear face): in the IC crystal are trapezoid (TRD) and the half-length is the 2th element
34  double width = 2 * aHit->GetDetector().dimensions[1];
35 
36  // use Crystal ID to define IDX and IDY
37  int IDX = identity[0].id;
38  int IDY = identity[1].id;
39 
40  // initialize ADC and TDC
41  int ADC = 0;
42  int TDC = 8192;
43 
44 
45  double Tmin = 99999.;
46  vector<G4double> times = aHit->GetTime();
47  vector<G4ThreeVector> Lpos = aHit->GetLPos();
48  vector<G4double> Edep = aHit->GetEdep();
49  if(Etot>0)
50  {
51  for(unsigned int s=0; s<nsteps; s++)
52  {
53  double dRight = length/2 - Lpos[s].z(); // distance along z between the hit position and the end of the crystal
54  double timeR = times[s] + dRight/cm/light_speed; // arrival time of the signal at the end of the crystal (speed of light in the crystal=15 cm/ns)
55  if(Edep[s]>1*MeV) Tmin=min(Tmin,timeR); // signal time is set to first hit time with energy above 1 MeV
56  }
57  TDC=int(Tmin*tdc_time_to_channel);
58  if(TDC>tdc_max) TDC=(int)tdc_max;
59  // calculate number of photoelectrons detected by the APD considering the light yield, the q.e., and the size of the sensor
60  double npe=G4Poisson(Etot*PbWO4_light_yield/2*APD_qe*APD_size/width/width);
61  // calculating APD output charge (in number of electrons) and adding noise
62  double nel=npe*APD_gain;
63  nel=nel*G4RandGauss::shoot(1.,APD_noise);
64  if(nel<0) nel=0;
65  // adding preamplifier input noise
66  nel=nel+AMP_input_noise*G4RandGauss::shoot(0.,1.);
67  if(nel<0) nel=0;
68  // converting to charge (in picoCoulomb)
69  double crg=nel*AMP_gain*1.6e-7;
70  // converting to ADC channels
71  ADC= (int) (crg*adc_charge_tochannel);
72 
73  }
74 
75  dgtz["hitn"] = hitn;
76  dgtz["IDX"] = IDX;
77  dgtz["IDY"] = IDY;
78  dgtz["ADC"] = ADC;
79  dgtz["TDC"] = TDC;
80 
81  return dgtz;
82 }
83 
84 vector<identifier> IC_HitProcess :: processID(vector<identifier> id, G4Step* aStep, detector Detector)
85 {
86  id[id.size()-1].id_sharing = 1;
87  return id;
88 }
89 
90 
91 map< string, vector <int> > IC_HitProcess :: multiDgt(MHit* aHit, int hitn)
92 {
93  map< string, vector <int> > MH;
94 
95  return MH;
96 }
97 
98 
99 
100 
101 
102 
103 
104 
105 
106 
107 
108 
109 
110 
vector< identifier > processID(vector< identifier >, G4Step *, detector)
vector< identifier > GetId()
Definition: Hit.h:103
map< string, double > integrateDgt(MHit *, int)
Definition: IC_hitprocess.cc:8
map< string, vector< int > > multiDgt(MHit *, int)
vector< G4ThreeVector > GetLPos()
Definition: Hit.h:76
Definition: Hit.h:22
vector< double > GetTime()
Definition: Hit.h:89
vector< double > dimensions
vector of dimensions. Size, units depends on solid type
Definition: detector.h:77
vector< double > GetEdep()
Definition: Hit.h:83
detector GetDetector()
Definition: Hit.h:108