GEMC  2.3
Geant4 Monte-Carlo Framework
fmt_strip.cc
Go to the documentation of this file.
1 // gemc headers
2 #include "fmt_strip.h"
3 #include "Randomize.hh"
4 #include <iostream>
5 #include <cmath>
6 
8 {
9  // all dimensions are in mm
10 
11  Pi = 3.14159265358;
12  // interlayer = (0.5+0.1+0.015+0.128+0.030+2.500)*2.; // distance between 2 layers of a superlayer // not usefull anymore, 07/27/2015 (Frederic Georges)
13  // intersuperlayer = 20.0; // distance between 2 superlayers
14  // intersuperlayer = 21.0; // distance between 2 superlayers // modified on 5/20/2015 to match new geometry (R.De Vita)
15  intersuperlayer = 10.5; // distance between 2 superlayers // modified on 7/27/2015 to match new geometry (Frederic Georges)
16  // pitch = 0.500; // pitch of the strips
17  pitch = 0.525; // pitch of the strips // modified on 7/27/2015 to match new geometry (Frederic Georges)
18  hDrift = 5.0;
19  //hStrip2Det = hDrift/2.;
20  sigma_td_max = 0.01; // very small transverse diffusion because of B field
21  w_i = 25.0;
22 
23  //R_min = 12.900;
24  //R_max = 215.0; // outer radius of strip part
25  R_min = 45.500; // modified on 7/27/2015 to match new geometry (Frederic Georges)
26  R_max = 185.400; // outer radius of strip part // modified on 7/27/2015 to match new geometry (Frederic Georges)
27  // Z_1stlayer = 295.-0.5-0.1-0.015-0.128-0.030-2.500; // z position of the 1st layer (middle of the drift gap)
28  //Z_1stlayer = 305.25-0.5-0.1-0.015-0.128-0.030-2.500; // z position of the 1st layer (middle of the drift gap) // modified on 5/20/2015 to match new geometry (R.De Vita)
29  Z_1stlayer = 305.250-0.005-0.200-2.000-0.200-0.005-0.050-0.020-0.128-0.030-2.500; // z position of the 1st layer (middle of the drift gap) // modified on 7/27/2015 to match new geometry (Frederic Georges)
30 
31  // z of the upstream part of the layer
32  /*
33  Z0.push_back(Z_1stlayer);
34  Z0.push_back(Z0[0]+interlayer);
35  Z0.push_back(Z_1stlayer+intersuperlayer);
36  Z0.push_back(Z0[2]+interlayer);
37  Z0.push_back(Z_1stlayer+2.*intersuperlayer);
38  Z0.push_back(Z0[4]+interlayer);
39  */
40  // modified on 7/27/2015 to match new geometry (Frederic Georges)
41  Z0.push_back(Z_1stlayer);
42  Z0.push_back(Z0[0]+intersuperlayer);
43  Z0.push_back(Z0[1]+intersuperlayer);
44  Z0.push_back(Z0[2]+intersuperlayer);
45  Z0.push_back(Z0[3]+intersuperlayer);
46  Z0.push_back(Z0[4]+intersuperlayer);
47 
48  // angles of each layer
49  /*
50  alpha.push_back(0);
51  alpha.push_back(Pi/2.);
52  alpha.push_back(Pi/3.);
53  alpha.push_back(Pi/2+Pi/3.);
54  alpha.push_back(2.*Pi/3.);
55  alpha.push_back(2.*Pi/3+Pi/2.);
56  */
57  // modified on 7/27/2015 to match new geometry (Frederic Georges)
58  alpha.push_back(19.*Pi/180.);
59  alpha.push_back(alpha[0]+Pi/3.);
60  alpha.push_back(alpha[0]+2.*Pi/3.);
61  alpha.push_back(alpha[0]+Pi);
62  alpha.push_back(alpha[0]+4.*Pi/3.);
63  alpha.push_back(alpha[0]+5.*Pi/3.);
64 
65  // Number of strips and pixels
66  N_str = 1024; //16 connectors * 64 strips = 1024 strips for each fmt
67 }
68 
69 vector<double> fmt_strip::FindStrip(int layer, int sector, double x, double y, double z, double Edep)
70 {
71  // the return vector is always in pairs.
72  // The first number is the ID,
73  // the second number is the sharing percentage
74  vector<double> strip_id;
75  // number of electrons (Nt)
76  Nel = (int) (1e6*Edep/w_i);
77  if(z-Z0[layer]>0.5001*hDrift || z-Z0[layer]<-0.5001*hDrift) cout << "Warning! z position of the FMT hit is not in the sensitive volume: " << z-Z0[layer]<< endl;
78  sigma_td = sigma_td_max* sqrt((z-Z0[layer])/hDrift); // expression without Lorentz angle
79 
80  int ClosestStrip=0;
81  if(Nel>0)
82  {
83  for(int iel=0;iel<Nel;iel++)
84  { // loop over (total) electrons
85  x_real = (double) (G4RandGauss::shoot(x*cos(alpha[layer])+y*sin(alpha[layer]),sigma_td));
86  y_real = (double) (G4RandGauss::shoot(y*cos(alpha[layer])-x*sin(alpha[layer]),sigma_td));
87 
88  if(y_real > -84.1 && y_real < 84.1 && x_real < 0){ // R_max - 3*64*0.525 - 0.5 = 84.1;
89  ClosestStrip = (int) (floor((84.1-y_real)/pitch)+1);
90  }
91  else if(y_real < -84.1 && y_real > -184.9){ // R_max - 3*64*0.525 - 0.5 = 84.1; R_max - 0.5 = 184.9
92  ClosestStrip = (int) (floor((-y_real-84.1)/pitch)+1) + 320; // 5*64 = 320
93  }
94  else if(y_real > -84.1 && y_real < 84.1 && x_real > 0){ // R_max - 3*64*0.525 - 0.5 = 84.1;
95  ClosestStrip = (int) (floor((y_real+84.1)/pitch)+1) + 512; // (5+3)*64 = 512
96  }
97  else if(y_real > 84.1 && y_real < 184.9){ // R_max - 3*64*0.525 - 0.5 = 84.1; R_max - 0.5 = 184.9
98  ClosestStrip = (int) (floor((y_real-84.1)/pitch)+1) + 832; // (5+3+5)*64 = 832
99  }
100 
101  if(sqrt(x*x+y*y)<R_max && sqrt(x*x+y*y)>R_min && ClosestStrip>=1 && ClosestStrip<=N_str)
102  { // strip is in the acceptance (~)
103  for(int istrip=0;istrip< (int) (strip_id.size()/2);istrip++)
104  {
105  if(strip_id[2*istrip]==ClosestStrip)
106  {// already hit strip - add Edep
107  strip_id[2*istrip+1]=strip_id[2*istrip+1]+1./((double) Nel); // no gain fluctuation yet
108  ClosestStrip=-1; // not to use it anymore
109  }
110  }
111  if(ClosestStrip>-1)
112  { // this is a new strip
113  strip_id.push_back(ClosestStrip);
114  strip_id.push_back(1./((double) Nel)); // no gain fluctuation yet
115  }
116  }
117  else
118  {// not in the acceptance
119  strip_id.push_back(-1);
120  strip_id.push_back(1);
121  }
122  }
123  }
124  else
125  { // Nel=0, consider the Edep is 0
126  strip_id.push_back(-1);
127  strip_id.push_back(1);
128  }
129  return strip_id;
130 }
131 
132 
133 
134 
135 
136 
137 
138 
139 
140 
141 
double Pi
Definition: fmt_strip.h:8
double y_real
Definition: fmt_strip.h:26
void fill_infos()
Definition: fmt_strip.cc:7
int N_str
Definition: fmt_strip.h:19
double R_max
Definition: fmt_strip.h:12
double x_real
Definition: fmt_strip.h:27
double sigma_td_max
Definition: fmt_strip.h:22
int Nel
Definition: fmt_strip.h:25
double Z_1stlayer
Definition: fmt_strip.h:15
double sigma_td
Definition: fmt_strip.h:23
double w_i
Definition: fmt_strip.h:24
double R_min
Definition: fmt_strip.h:13
vector< double > alpha
Definition: fmt_strip.h:18
double pitch
Definition: fmt_strip.h:7
vector< double > FindStrip(int layer, int sector, double x, double y, double z, double Edep)
Definition: fmt_strip.cc:69
vector< double > Z0
Definition: fmt_strip.h:17
double hDrift
Definition: fmt_strip.h:20
double intersuperlayer
Definition: fmt_strip.h:11