GEMC  1.8
Geant4 Monte-Carlo Framework
ECwithG4strips_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 // %%%%%%%%%%%%%
12 
13 // Process the ID and hit for the EC using individual EC scintillator strips.
14 
16 {
17  PH_output out;
18  out.identity = aHit->GetId();
19  // get sector, stack (inner or outer), view (U, V, W), and strip.
20  int sector = out.identity[0].id;
21  int stack = out.identity[1].id;
22  int view = out.identity[2].id;
23  int strip = out.identity[3].id;
24  HCname = "ECwithG4strips Hit Process";
25 
26  //cout << "from ECwithG4strips ProcessHit strip = " << strip << endl;
27 
28  // %%%%%%%%%%%%%%%%%%%
29  // Raw hit information
30  // %%%%%%%%%%%%%%%%%%%
31  int nsteps = aHit->GetPos().size();
32 
33  // Get Total Energy deposited
34  double Etot = 0;
35  vector<G4double> Edep = aHit->GetEdep();
36  for(int s=0; s<nsteps; s++) Etot = Etot + Edep[s];
37 
38  // average global, local positions of the hit
39  double x, y, z;
40  double lx, ly, lz;
41  x = y = z = lx = ly = lz = 0;
42  vector<G4ThreeVector> pos = aHit->GetPos();
43  vector<G4ThreeVector> Lpos = aHit->GetLPos();
44 
45  if(Etot>0)
46  for(int s=0; s<nsteps; s++)
47  {
48  x = x + pos[s].x()*Edep[s]/Etot;
49  y = y + pos[s].y()*Edep[s]/Etot;
50  z = z + pos[s].z()*Edep[s]/Etot;
51  lx = lx + Lpos[s].x()*Edep[s]/Etot;
52  ly = ly + Lpos[s].y()*Edep[s]/Etot;
53  lz = lz + Lpos[s].z()*Edep[s]/Etot;
54  }
55  else
56  {
57  x = pos[0].x();
58  y = pos[0].y();
59  z = pos[0].z();
60  lx = Lpos[0].x();
61  ly = Lpos[0].y();
62  lz = Lpos[0].z();
63  }
64 
65  // average time
66  double time = 0;
67  vector<G4double> times = aHit->GetTime();
68  for(int s=0; s<nsteps; s++) time = time + times[s]/nsteps;
69 
70  // Energy of the track
71  double Ene = aHit->GetE();
72 
73  out.raws.push_back(Etot);
74  out.raws.push_back(x);
75  out.raws.push_back(y);
76  out.raws.push_back(z);
77  out.raws.push_back(lx);
78  out.raws.push_back(ly);
79  out.raws.push_back(lz);
80  out.raws.push_back(time);
81  out.raws.push_back((double) aHit->GetPID());
82  out.raws.push_back(aHit->GetVert().getX());
83  out.raws.push_back(aHit->GetVert().getY());
84  out.raws.push_back(aHit->GetVert().getZ());
85  out.raws.push_back(Ene);
86  out.raws.push_back((double) aHit->GetmPID());
87  out.raws.push_back(aHit->GetmVert().getX());
88  out.raws.push_back(aHit->GetmVert().getY());
89  out.raws.push_back(aHit->GetmVert().getZ());
90 
91 
92  // %%%%%%%%%%%%
93  // Digitization
94  // %%%%%%%%%%%%
95 
96  // Jerry Gilfoyle, Feb, 2010
97  // parameters: changed from hardwiring (see below) to gpars array - gilfoyle 8/31/12
98  int ec_tdc_time_to_channel = (int) gpars["EC/ec_tdc_time_to_channel"]; // conversion from time (ns) to TDC channels.
99  double ECfactor = (double) gpars["EC/ECfactor"]; // number of photons divided by the energy deposited in MeV; value taken from gsim. see EC NIM paper table 1.
100  int EC_TDC_MAX = (int) gpars["EC/EC_TDC_MAX"]; // max value for EC tdc.
101  int ec_charge_to_channel = (int) gpars["EC/ec_charge_to_channel"]; // conversion from charge (pC) to ADC channels; taken from IC code.
102  double ec_npe_to_charge = (double) gpars["EC/ec_npe_to_charge"]; // unjustified (!!) conversion from number of photoelectrons to charge in pC.
103 
104  // parameters were initially hardwired into the code with the values shown below.
105  // Moved into trunk/database_io/clas/geo/ec/parameters.txt - gilfoyle 8/31/12
106  //
107  // double ECfactor = 35; // number of photons divided by the energy deposited in MeV; value taken from gsim. see EC NIM paper table 1.
108  // int EC_TDC_MAX = 4095; // max value for EC tdc.
109  // int ec_charge_to_channel = 20; // conversion from charge (pC) to ADC channels; taken from IC code.
110  // double ec_npe_to_charge = 0.0005; // unjustified (!!) conversion from number of photoelectrons to charge in pC.
111 
112  // initialize ADC and TDC
113  int EC_ADC = 0;
114  int EC_TDC = EC_TDC_MAX;
115 
116  // simulate the adc value.
117  if (Etot > 0) {
118  // number of photoelectrons.
119  double EC_npe = G4Poisson(Etot*ECfactor);
120  // Fluctuations in PMT gain distributed using Gaussian with sigma SNR = sqrt(ngamma)/sqrt(del/del-1) del = dynode gain = 3 (From RCA PMT Handbook) p. 169)
121  // algorithm, values, and comment above taken from gsim.
122  double sigma = sqrt(EC_npe)*1.15;
123  double EC_charge = G4RandGauss::shoot(EC_npe,sigma)*ec_npe_to_charge*ec_charge_to_channel;
124  if (EC_charge <= 0) EC_charge=0.0; // guard against weird, rare events.
125  EC_ADC = (int) EC_charge;
126  }
127 
128  // simulate the tdc.
129  EC_TDC = (int) (time*ec_tdc_time_to_channel);
130  if (EC_TDC > EC_TDC_MAX) EC_TDC = EC_TDC_MAX;
131 
132  out.dgtz.push_back(sector);
133  out.dgtz.push_back(stack);
134  out.dgtz.push_back(view);
135  out.dgtz.push_back(strip);
136  out.dgtz.push_back(EC_ADC);
137  out.dgtz.push_back(EC_TDC);
138 
139  //cout << "sector = " << sector << " stack = " << stack << " view = " << view << " strip = " << strip << " EC_ADC = " << EC_ADC << " EC_TDC = " << EC_TDC << " Edep = " << Etot << endl;
140 
141  return out;
142 }
143 
144 vector<identifier> ECwithG4strips_HitProcess :: ProcessID(vector<identifier> id, G4Step* aStep, detector Detector, gemc_opts Opt)
145 {
146  vector<identifier> yid = id;
147 
148  //int view = yid[2].id; // get the view: U->1, V->2, W->3
149 
150  G4StepPoint *prestep = aStep->GetPreStepPoint();
151  G4StepPoint *poststep = aStep->GetPostStepPoint();
152  G4VTouchable* TH = (G4VTouchable*) aStep->GetPreStepPoint()->GetTouchable();
153 
154  string name = TH->GetVolume(0)->GetName();
155  G4ThreeVector xyz = poststep->GetPosition();
156  G4ThreeVector Lxyz = prestep->GetTouchableHandle()->GetHistory()
157  ->GetTopTransform().TransformPoint(xyz);
158 
159  //double xlocal = Lxyz.x();
160  //double ylocal = Lxyz.y();
161 
162  //double pDy1 = Detector.dimensions[3]; ///< G4Trap Semilength. these points are needed to get the strip number.
163  //double pDx2 = Detector.dimensions[5]; ///< G4Trap Semilength.
164 
165  int strip_number = yid[3].id; // scintillator strip from geometry.
166  string detName = Detector.name;
167 
168  //cout << "strip from the geometry = " << strip_number << " view=" << view << " name=" << detName << endl;
169 
170  //
171  // Some EC geometry:
172  //
173  // B ---------------------- C
174  // \ /
175  // \ /
176  // \ /
177  // \ / face of EC sector looking from the target; z-axis is out
178  // \ x /
179  // \ / | y
180  // \ / | coordinates: origin is at the
181  // \ / _______| geometric center of the triangle.
182  // \ / x
183  // \/
184  // A
185  //
186  // U - strips parallel to BC, start numbering at A
187  // V - strips parallel to AB, start numbering at C
188  // W - strips parallel tp CA, start numbering at B
189  //
190  // other points: D - point where line from C crosses AB at right angle.
191  // F - CF is the component of the hit along CD.
192  // G - point where line from B crosses AC at right angle.
193  // H - BH is the component of the hit position along BG.
194  //
195 
196 
197  // ******************************************************************************************************************************
198  //
199  // replace the strip number calculated from the layer position with the strip number from the revised geometry. - gilfoyle 7/11/12
200  //
201  // ******************************************************************************************************************************
202 
203  // if (view == 1) {
204  // cout << "test = " << setw(10) << testStrip << " strip calculated from slab =" << setw(2) << strip << "; strip from G4 strip geometry = " << setw(2) << strip_number << "; view=" <<
205  // view << " name=" << setw(20) << detName << " ylocal = " << setw(8) << ylocal <<
206  // " pDy1 = " << setw(5) << pDy1 << endl;
207  // }
208 
209  yid[3].id = strip_number;
210  yid[3].id_sharing = 1;
211 
212  //if (stack==1 && view == 3) cout << "sector = " << sector << " stack = " << stack << " view = " << view << " xlocal = " << xlocal << " ylocal = " << ylocal << " strip = " << strip << " pDx2 = " << pDx2 << " pDy1 = " << pDy1 << endl;
213 
214  return yid;
215 }
216 
217 
218 
219 
220 
221 
222 
223 
224 
225 
226 
vector< identifier > ProcessID(vector< identifier >, G4Step *, detector, gemc_opts)
Method to calculate new identifier.
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
string name
Name of the volume. Since this is the key of the STL map, it has to be unique.
Definition: detector.h:53
map< string, double > gpars
Definition: MPHBaseClass.h:42
PH_output ProcessHit(MHit *, gemc_opts)
Method to process the hit.
Definition: MHit.h:29
vector< double > GetTime()
Definition: MHit.h:82
vector< double > GetEdep()
Definition: MHit.h:76
vector< int > dgtz
Digitized information.
Definition: MPHBaseClass.h:27