GEMC  2.3
Geant4 Monte-Carlo Framework
rich_hitprocess.cc
Go to the documentation of this file.
1 // G4 headers
2 #include "G4Poisson.hh"
3 #include "Randomize.hh"
4 
5 
6 // CLHEP units
7 #include "CLHEP/Units/PhysicalConstants.h"
8 using namespace CLHEP;
9 
10 // gemc headers
11 #include "rich_hitprocess.h"
12 
13 map<string, double> rich_HitProcess :: integrateDgt(MHit* aHit, int hitn)
14 {
15  map<string, double> dgtz;
16  vector<identifier> identity = aHit->GetId();
17 
18  int sector = identity[0].id;
19  int pad = identity[1].id;
20 
21 
22  // Ahmed El Alaoui, May, 2010
23  dgtz["hitn"] = hitn;
24  dgtz["sector"] = sector;
25  dgtz["pad"] = pad;
26 
27  return dgtz;
28 }
29 
30 
31 // this routine needs to be modified
32 // no point drawing should be made here, but in MHit
33 // finding the PMT should be in a different function,
34 // with parameters coming from DB
35 
36 #include "G4VVisManager.hh"
37 #include "G4Circle.hh"
38 #include "G4VisAttributes.hh"
39 #include "G4ParticleTable.hh"
40 
41 vector<identifier> rich_HitProcess :: processID(vector<identifier> id, G4Step* aStep, detector Detector)
42 {
43  vector<identifier> id2 = id;
44 
45  G4StepPoint *prestep = aStep->GetPreStepPoint();
46  G4ThreeVector xyz = aStep->GetPostStepPoint()->GetPosition();
47  G4ThreeVector Lxyz = prestep->GetTouchableHandle()->GetHistory()->GetTopTransform().TransformPoint(xyz);
48  const double xpos = Lxyz.x();
49  const double ypos = Lxyz.y();
50  const double Energy = aStep->GetTrack()->GetTotalEnergy()/eV;
51 
52  // H8500 specs:
53  static const int nbins=15;
54  static const double p[nbins]={1.86, 1.96, 2.03, 2.13, 2.28, 2.37, 2.48, 2.59, 2.71, 3.11, 3.37, 3.73, 4.13, 4.39, 4.64};
55  static const double q[nbins]={0.0002, 0.002, 0.007, 0.02, 0.04, 0.08, 0.13, 0.17, 0.20, 0.26, 0.27, 0.28, 0.22, 0.18, 0.10};
56  static const int NIPXL = 8;
57  static const int NJPXL = 8;
58  static const double PhCsize_x(49);
59  static const double PhCsize_y(49);
60  static const double PXLDeadSpace_x(0.15);
61  static const double PXLDeadSpace_y(0.15);
62  static const double PXLsize_x((PhCsize_x-7*PXLDeadSpace_x)/8);
63  static const double PXLsize_y((PhCsize_y-7*PXLDeadSpace_y)/8);
64 
65  //shift the local reference center to down left corner
66  const double new_xpos = xpos+PhCsize_x/2.;
67  const double new_ypos = ypos+PhCsize_y/2.;
68 
69  //find the pixel row
70  const int iPMT = (int) floor(new_xpos/(PXLDeadSpace_x+PXLsize_x))+1;
71  const double lx = new_xpos -(iPMT-1)*(PXLDeadSpace_x+PXLsize_x);
72 
73  //find the pixel column
74  const int jPMT = (int) floor(new_ypos/(PXLDeadSpace_y+PXLsize_y))+1;
75  const double ly = new_ypos -(jPMT-1)*(PXLDeadSpace_y+PXLsize_y);
76 
77  // find the pixel number
78  int ijPMT = NIPXL*(NJPXL-jPMT) + iPMT;
79 
80  // photon hit the dead space:
81  if(lx>PXLsize_x || ly>PXLsize_y)
82  {
83  ijPMT = 0;
84  }
85  else
86  {
87  // pixel# += 100 if fails quauntum efficiency:
88  double QE = -1;
89  if(Energy > p[0] || Energy <= p[nbins-1])
90  {
91  for(int ie=0; ie<nbins-1; ie++)
92  {
93  if(Energy > p[ie] && Energy <= p[ie+1])
94  {
95  QE = q[ie] + (Energy-p[ie])*(q[ie+1]-q[ie])/(p[ie+1]-p[ie]);
96  if(G4UniformRand()>QE)
97  {
98  ijPMT += 100;
99  }
100  break;
101  }
102  }
103  }
104  }
105 
106  G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
107  if(pVVisManager)
108  {
109  if (ijPMT<100)
110  {
111  // All hits are already colored green by MHit.cc. Here,
112  // draw white circle for photons that hit dead area, and
113  // red circle for photons that pass QE.
114  G4Circle circle(xyz);
115  circle.SetFillStyle(G4Circle::filled);
116  circle.SetScreenSize(8);
117  if (ijPMT>0)
118  {
119  G4Colour colour_touch (1.0, 0.0, 0.0);
120  circle.SetVisAttributes(G4VisAttributes(colour_touch));
121  }
122  pVVisManager->Draw(circle);
123  }
124  }
125 
126  id2[2].id = ijPMT;
127  id2[2].id_sharing = 1;
128  return id2;
129 }
130 
131 
132 
133 map< string, vector <int> > rich_HitProcess :: multiDgt(MHit* aHit, int hitn)
134 {
135  map< string, vector <int> > MH;
136 
137  return MH;
138 }
139 
140 
141 
142 
143 
144 
145 
146 
147 
148 
149 
vector< identifier > processID(vector< identifier >, G4Step *, detector)
vector< identifier > GetId()
Definition: Hit.h:103
map< string, double > integrateDgt(MHit *, int)
map< string, vector< int > > multiDgt(MHit *, int)
Definition: Hit.h:22