GEMC  1.8
Geant4 Monte-Carlo Framework
CTOF_hitprocess.cc
Go to the documentation of this file.
1 // %%%%%%%%%%
2 // G4 headers
3 // %%%%%%%%%%
4 #include "G4UnitsTable.hh"
5 
6 
7 // %%%%%%%%%%%%%
8 // gemc headers
9 // %%%%%%%%%%%%%
10 #include "CTOF_hitprocess.h"
11 
12 
14 {
15  string hd_msg = Opt.args["LOG_MSG"].args + " CTOF Hit Process " ;
16  double HIT_VERBOSITY = Opt.args["HIT_VERBOSITY"].arg;
17  PH_output out;
18  out.identity = aHit->GetId();
19  HCname = "CTOF Hit Process";
20 
21  // %%%%%%%%%%%%%%%%%%%
22  // Raw hit information
23  // %%%%%%%%%%%%%%%%%%%
24  int nsteps = aHit->GetPos().size();
25 
26  // Get Total Energy deposited
27  double Etot = 0;
28  vector<G4double> Edep = aHit->GetEdep();
29  for(int s=0; s<nsteps; s++) Etot = Etot + Edep[s];
30 
31  // average global, local positions of the hit
32  double x, y, z;
33  double lx, ly, lz;
34  x = y = z = lx = ly = lz = 0;
35  vector<G4ThreeVector> pos = aHit->GetPos();
36  vector<G4ThreeVector> Lpos = aHit->GetLPos();
37 
38  if(Etot>0)
39  for(int s=0; s<nsteps; s++)
40  {
41  x = x + pos[s].x()*Edep[s]/Etot;
42  y = y + pos[s].y()*Edep[s]/Etot;
43  z = z + pos[s].z()*Edep[s]/Etot;
44  lx = lx + Lpos[s].x()*Edep[s]/Etot;
45  ly = ly + Lpos[s].y()*Edep[s]/Etot;
46  lz = lz + Lpos[s].z()*Edep[s]/Etot;
47  }
48  else
49  {
50  x = pos[0].x();
51  y = pos[0].y();
52  z = pos[0].z();
53  lx = Lpos[0].x();
54  ly = Lpos[0].y();
55  lz = Lpos[0].z();
56  }
57 
58  // average time
59  double time = 0;
60  vector<G4double> times = aHit->GetTime();
61  for(int s=0; s<nsteps; s++) time = time + times[s]/nsteps;
62 
63  // Energy of the track
64  double Ene = aHit->GetE();
65 
66  out.raws.push_back(Etot);
67  out.raws.push_back(x);
68  out.raws.push_back(y);
69  out.raws.push_back(z);
70  out.raws.push_back(lx);
71  out.raws.push_back(ly);
72  out.raws.push_back(lz);
73  out.raws.push_back(time);
74  out.raws.push_back((double) aHit->GetPID());
75  out.raws.push_back(aHit->GetVert().getX());
76  out.raws.push_back(aHit->GetVert().getY());
77  out.raws.push_back(aHit->GetVert().getZ());
78  out.raws.push_back(Ene);
79  out.raws.push_back((double) aHit->GetmPID());
80  out.raws.push_back(aHit->GetmVert().getX());
81  out.raws.push_back(aHit->GetmVert().getY());
82  out.raws.push_back(aHit->GetmVert().getZ());
83 
84 
85  // %%%%%%%%%%%%
86  // Digitization
87  // %%%%%%%%%%%%
88 
89  int paddle = out.identity[0].id;
90 
91  if(HIT_VERBOSITY>4)
92  cout << hd_msg << " paddle: " << paddle << " x=" << x/cm << " y=" << y/cm << " z=" << z/cm << endl;
93 
94  out.dgtz.push_back(paddle);
95 
96  return out;
97 }
98 
99 
100 vector<identifier> CTOF_HitProcess :: ProcessID(vector<identifier> id, G4Step *step, detector Detector, gemc_opts Opt)
101 {
102  id[id.size()-1].id_sharing = 1;
103  return id;
104 }
105 
106 
107 
108 
109 
110 
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
vector< identifier > ProcessID(vector< identifier >, G4Step *, detector, gemc_opts)
Method to calculate new identifier.
Definition: MHit.h:29
vector< double > GetTime()
Definition: MHit.h:82
map< string, opts > args
Options map.
Definition: usage.h:68
vector< double > GetEdep()
Definition: MHit.h:76
PH_output ProcessHit(MHit *, gemc_opts)
Method to process the hit.
vector< int > dgtz
Digitized information.
Definition: MPHBaseClass.h:27