GEMC  1.8
Geant4 Monte-Carlo Framework
FT_hitprocess.cc
Go to the documentation of this file.
1 // %%%%%%%%%%
2 // G4 headers
3 // %%%%%%%%%%
4 #include "G4UnitsTable.hh"
5 #include "G4Poisson.hh"
6 #include "Randomize.hh"
7 
8 // %%%%%%%%%%%%%
9 // gemc headers
10 // %%%%%%%%%%%%%
11 #include "FT_hitprocess.h"
12 
14 {
15  PH_output out;
16  out.identity = aHit->GetId();
17  HCname = "FT Hit Process";
18 
19  // %%%%%%%%%%%%%%%%%%%
20  // Raw hit information
21  // %%%%%%%%%%%%%%%%%%%
22  int nsteps = aHit->GetPos().size();
23 
24  // Get Total Energy deposited
25  double Etot = 0;
26  vector<G4double> Edep = aHit->GetEdep();
27  for(int s=0; s<nsteps; s++) Etot = Etot + Edep[s];
28 
29  // average global, local positions of the hit
30  double x, y, z;
31  double lx, ly, lz;
32  x = y = z = lx = ly = lz = 0;
33  vector<G4ThreeVector> pos = aHit->GetPos();
34  vector<G4ThreeVector> Lpos = aHit->GetLPos();
35 
36  if(Etot>0)
37  for(int s=0; s<nsteps; s++)
38  {
39  x = x + pos[s].x()*Edep[s]/Etot;
40  y = y + pos[s].y()*Edep[s]/Etot;
41  z = z + pos[s].z()*Edep[s]/Etot;
42  lx = lx + Lpos[s].x()*Edep[s]/Etot;
43  ly = ly + Lpos[s].y()*Edep[s]/Etot;
44  lz = lz + Lpos[s].z()*Edep[s]/Etot;
45  }
46  else
47  {
48  x = pos[0].x();
49  y = pos[0].y();
50  z = pos[0].z();
51  lx = Lpos[0].x();
52  ly = Lpos[0].y();
53  lz = Lpos[0].z();
54  }
55 
56  // average time
57  double time = 0;
58  vector<G4double> times = aHit->GetTime();
59  for(int s=0; s<nsteps; s++) time = time + times[s]/nsteps;
60 
61  // Energy of the track
62  double Ene = aHit->GetE();
63 
64  out.raws.push_back(Etot);
65  out.raws.push_back(x);
66  out.raws.push_back(y);
67  out.raws.push_back(z);
68  out.raws.push_back(lx);
69  out.raws.push_back(ly);
70  out.raws.push_back(lz);
71  out.raws.push_back(time);
72  out.raws.push_back((double) aHit->GetPID());
73  out.raws.push_back(aHit->GetVert().getX());
74  out.raws.push_back(aHit->GetVert().getY());
75  out.raws.push_back(aHit->GetVert().getZ());
76  out.raws.push_back(Ene);
77  out.raws.push_back((double) aHit->GetmPID());
78  out.raws.push_back(aHit->GetmVert().getX());
79  out.raws.push_back(aHit->GetmVert().getY());
80  out.raws.push_back(aHit->GetmVert().getZ());
81 
82  // %%%%%%%%%%%%
83  // Digitization
84  // %%%%%%%%%%%%
85 
86  // R.De Vita (April 2009)
87 
88  // relevant parameter for digitization (in the future should be read from database)
89  double tdc_time_to_channel=20; // conversion factor from time(ns) to TDC channels)
90  double tdc_max=8191; // TDC range
91  double time_res=0.2; // time resolution
92  double adc_charge_tochannel=20; // conversion factor from charge(pC) to ADC channels
93  double PbWO4_light_yield =240/MeV; // Lead Tungsten Light Yield (APD have similar QE for
94  // fast component, lambda=420nm-ly=120ph/MeV, and slow component,
95  // lambda=560nm-ly=20ph/MeV, taking fast component only)
96 // double PbWO4_light_yield =672/MeV; // LY at -25 deg=2.8 x LY at +18 deg
97  double APD_qe = 0.70; // APD Quantum Efficiency (Hamamatsu S8664-55)
98  double APD_size = 100*mm*mm; // APD size ( 10 mm x 10 mm)
99  double APD_gain = 150; // based on FT note
100  double APD_noise = 0.0033; // relative noise based on a Voltage and Temperature stability of 10 mV (3.9%/V) and 0.1 C (3.3%/C)
101  double AMP_input_noise = 9000; // preamplifier input noise in number of electrons
102  double AMP_gain = 1800; // preamplifier gain = 5V/pC x 25 ns (tipical signal duration)
103  double light_speed =15;
104 
105  // Get the crystal length: in the FT crystal are BOXes and the half-length is the 3rd element
106  double length = 2 * aHit->GetDetector().dimensions[2];
107  // Get the crystal width (rear face): in the FT crystal are BOXes and the half-length is the 2th element
108  double width = 2 * aHit->GetDetector().dimensions[1];
109 
110  // use Crystal ID to define IDX and IDY
111  int IDX = out.identity[0].id;
112  int IDY = out.identity[1].id;
113 
114  // initialize ADC and TDC
115  int ADC = 0;
116  int TDC = 8191;
117 
118  if(Etot>0)
119  {
120 /*
121 // commented out to use average time instead of minimum time in TDC calculation
122  for(int s=0; s<nsteps; s++)
123  {
124  double dRight = length/2 - Lpos[s].z(); // distance along z between the hit position and the end of the crystal
125  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)
126  if(Edep[s]>1*MeV) Tmin=min(Tmin,timeR); // signal time is set to first hit time with energy above 1 MeV
127  }
128  TDC=int(Tmin*tdc_time_to_channel);
129  if(TDC>tdc_max) TDC=(int)tdc_max;
130 */
131  double dRight = length/2 - lz; // distance along z between the hit position and the end of the crystal
132  double timeR = time + dRight/cm/light_speed; // arrival time of the signal at the end of the crystal (speed of light in the crystal=15 cm/ns)
133 // adding spread on time
134  timeR=timeR+G4RandGauss::shoot(0.,time_res);
135 
136  TDC=int(timeR*tdc_time_to_channel);
137  if(TDC>tdc_max) TDC=(int)tdc_max;
138 
139  // calculate number of photoelectrons detected by the APD considering the light yield, the q.e., and the size of the sensor
140  double npe=G4Poisson(Etot*PbWO4_light_yield*0.5*APD_qe*APD_size/width/width);
141 // double npe=(Etot*PbWO4_light_yield*0.5*APD_qe*APD_size/width/width);
142 
143  // for PMT, an addition factor of 0.5625 is needed to reproduce the 13.5 photoelectrons with a 20% QE
144  // double npe=G4Poisson(Etot*PbWO4_light_yield*0.5*0.5625*APD_qe*APD_size/width/width);
145 
146  // calculating APD output charge (in number of electrons) and adding noise
147  double nel=npe*APD_gain;
148  nel=nel*G4RandGauss::shoot(1.,APD_noise);
149  if(nel<0) nel=0;
150  // adding preamplifier input noise
151  nel=nel+AMP_input_noise*G4RandGauss::shoot(0.,1.);
152  if(nel<0) nel=0;
153  // converting to charge (in picoCoulomb)
154  double crg=nel*AMP_gain*1.6e-7;
155  // converting to ADC channels
156  ADC= (int) (crg*adc_charge_tochannel);
157 
158  }
159  out.dgtz.push_back(IDX);
160  out.dgtz.push_back(IDY);
161  out.dgtz.push_back(ADC);
162  out.dgtz.push_back(TDC);
163 
164 
165  return out;
166 }
167 
168 vector<identifier> FT_HitProcess :: ProcessID(vector<identifier> id, G4Step* aStep, detector Detector, gemc_opts Opt)
169 {
170  id[id.size()-1].id_sharing = 1;
171  return id;
172 }
173 
174 
175 
176 
177 
178 
179 
180 
181 
182 
183 
PH_output ProcessHit(MHit *, gemc_opts)
Method to process the hit.
G4ThreeVector GetVert()
Definition: MHit.h:72
vector< double > raws
Raw information.
Definition: MPHBaseClass.h:26
int GetPID()
Definition: MHit.h:107
int GetmPID()
Definition: MHit.h:122
vector< identifier > GetId()
Definition: MHit.h:96
string HCname
Hit Collection name.
Definition: MPHBaseClass.h:41
G4ThreeVector GetmVert()
Definition: MHit.h:127
vector< identifier > identity
Identifier.
Definition: MPHBaseClass.h:28
vector< G4ThreeVector > GetLPos()
Definition: MHit.h:69
double GetE()
Definition: MHit.h:89
vector< G4ThreeVector > GetPos()
Definition: MHit.h:65
Definition: MHit.h:29
vector< double > GetTime()
Definition: MHit.h:82
vector< identifier > ProcessID(vector< identifier >, G4Step *, detector, gemc_opts)
Method to calculate new identifier.
vector< double > dimensions
vector of dimensions. Size, units depends on solid type
Definition: detector.h:63
vector< double > GetEdep()
Definition: MHit.h:76
detector GetDetector()
Definition: MHit.h:101
vector< int > dgtz
Digitized information.
Definition: MPHBaseClass.h:27