GEMC  1.8
Geant4 Monte-Carlo Framework
bst_strip.cc
Go to the documentation of this file.
1 // %%%%%%%%%%%%%
2 // gemc headers
3 // %%%%%%%%%%%%%
4 #include "bst_strip.h"
5 
6 #include <iostream>
7 #include <cmath>
8 #include <cstdlib>
9 
11 {
12 
13  intersensors = (0.11 + 0.002)*mm; // gap between sensors + microgap
14  alpha = 3.0*deg; // max angle of the strips
15  pitch = 0.156*mm; // pitch of the strips
16 
17  // number of cards by sector for each layer
18  NSensors.push_back(3); NSensors.push_back(3);
19  NSensors.push_back(3); NSensors.push_back(3);
20  NSensors.push_back(3); NSensors.push_back(3);
21  NSensors.push_back(3); NSensors.push_back(3);
22 
23  DZ_inLength = 0.984*mm; // size of the band of dead zones all around in the length of the card
24  DZ_inWidth = 0.835*mm; // size of the band of dead zones all around in the width of the card
25  SensorLength = 111.625*mm; // length of 1 Sensor
26  SensorWidth = 42.02*mm; // width 1 Sensor
27 
28  // Number of strips
29  Nstrips = (int) floor((SensorWidth-2.0*DZ_inLength)/pitch) - 1;
30 
31 }
32 
33 
34 vector<double> bst_strip::FindStrip(int layer, int sector, int isens, G4ThreeVector Lxyz)
35 {
36 
37  // the return vector is always in pairs.
38  // The first number is the ID,
39  // the second number is the sharing percentage
40  // layer and sector are the indexes (i.e. layer is the actual layern-1)
41  vector<double> strip_id;
42 
43  int StripHit = -1;
44  double minDist = 999; // distance between actual x and x of the k-strip at the actual z
45  double dalpha = alpha/((float)Nstrips);
46 
47  // local x position in the sensor: need to add 1/2 active area
48  double lx = Lxyz.x() + SensorWidth/2.0 - DZ_inLength;
49 
50  // local z position in the module: need to add 1/2 active area + 1(2) full sensor lengths if on sensor 2(3)
51  double lz = Lxyz.z() + SensorLength/2.0 - DZ_inWidth + (isens-1)*(SensorLength + intersensors) ;
52 
53  // particle must be in the z active area of the sensor
54  if( fabs(Lxyz.z()) < (SensorLength/2.0 - DZ_inWidth))
55  {
56 
57 
58  // looping over all strips to find the closest
59  // to the point
60  for(int k=0; k<Nstrips; k++)
61  {
62  // angle of the k-strip
63  double alpha_k = k*dalpha;
64 
65  // strip equation is z = mz + intcp
66  double intcp = 0;
67  double m = 0;
68 
69  if(layer%2 == 0)
70  {
71  // intercept: for layer A it starts from the top of the active area of the sensor in sector 1
72  // which is at positive local x
73  intcp = SensorWidth - 2*DZ_inLength - (k+1)*pitch ;
74  m = -tan(alpha_k/rad);
75  }
76 
77  if(layer%2 == 1)
78  {
79  // intercept: for layer B it starts from the bottom of the active area of the sensor in sector 1
80  // which is at negative local x
81  intcp = (k+1)*pitch ;
82  m = tan(alpha_k/rad);
83  }
84 
85  // x position of the strip according to the active local z position
86  double stripx = intcp + m*lz;
87 
88  if(fabs(lx - stripx) < minDist)
89  {
90  minDist = fabs(lx - stripx);
91  StripHit = k+1;
92  }
93 
94  }
95  }
96 
97  if(StripHit != -1)
98  {
99  // correcting for z positioning inside the module
100  double dpitch = pitch + lz*tan(dalpha/rad);
101 
102  // one strip only, all energy to it
103  if(fabs(minDist)<=dpitch/4.0)
104  {
105  strip_id.push_back(StripHit);
106  strip_id.push_back(1);
107  }
108  // two hits, on the right of the strip
109  // 10% loss due to capacitance between strip and backplane
110  if(minDist>dpitch/4.0)
111  {
112  strip_id.push_back(StripHit);
113  strip_id.push_back(0.45);
114  strip_id.push_back(StripHit+1);
115  strip_id.push_back(0.45);
116  }
117  // two hits, on the left of the strip
118  // 10% loss due to capacitance between strip and backplane
119  if(minDist<-dpitch/4.0)
120  {
121  strip_id.push_back(StripHit);
122  strip_id.push_back(0.45);
123  strip_id.push_back(StripHit-1);
124  strip_id.push_back(0.45);
125  }
126  }
127  else
128  {
129  strip_id.push_back(-1);
130  strip_id.push_back(1);
131  }
132  return strip_id;
133 
134 
135 }
136 
137 
void fill_infos()
Definition: bst_strip.cc:10
vector< int > NSensors
Definition: bst_strip.h:15
double alpha
Definition: bst_strip.h:11
double DZ_inLength
Definition: bst_strip.h:17
double pitch
Definition: bst_strip.h:12
double SensorWidth
Definition: bst_strip.h:20
double SensorLength
Definition: bst_strip.h:19
vector< double > FindStrip(int layer, int sector, int isens, G4ThreeVector Lxyz)
Definition: bst_strip.cc:34
double intersensors
Definition: bst_strip.h:14
int Nstrips
Definition: bst_strip.h:21
double DZ_inWidth
Definition: bst_strip.h:18