GEMC  2.3
Geant4 Monte-Carlo Framework
FMT_hitprocess.cc
Go to the documentation of this file.
1 // gemc headers
2 #include "FMT_hitprocess.h"
3 #include "fmt_strip.h"
4 
5 // CLHEP units
6 #include "CLHEP/Units/PhysicalConstants.h"
7 using namespace CLHEP;
8 
9 map<string, double>FMT_HitProcess :: integrateDgt(MHit* aHit, int hitn)
10 {
11  map<string, double> dgtz;
12  vector<identifier> identity = aHit->GetId();
13 
14  // FMT ID:
15  // layer, type, sector, strip
16 
17  int layer = 1*identity[0].id + identity[1].id - 1 ; // modified on 7/27/2015 to match new geometry (Frederic Georges)
18  //int layer = 2*identity[0].id + identity[1].id - 2 ;
19  int sector = identity[2].id;
20  int strip = identity[3].id;
21 
22  if(verbosity>4)
23  {
24  trueInfos tInfos(aHit);
25  cout << log_msg << " layer: " << layer << " sector: " << sector << " Strip: " << strip
26  << " x=" << tInfos.x << " y=" << tInfos.y << " z=" << tInfos.z << endl;
27  }
28 
29  dgtz["hitn"] = hitn;
30  dgtz["layer"] = layer;
31  dgtz["sector"] = sector;
32  dgtz["strip"] = strip;
33 
34  return dgtz;
35 }
36 
37 
38 
39 vector<identifier> FMT_HitProcess :: processID(vector<identifier> id, G4Step* aStep, detector Detector)
40 {
41  double x, y, z;
42  G4ThreeVector xyz = aStep->GetPostStepPoint()->GetPosition();
43  x = xyz.x()/mm;
44  y = xyz.y()/mm;
45  z = xyz.z()/mm;
46 
47  vector<identifier> yid = id;
48  class fmt_strip fmts;
49  fmts.fill_infos();
50 
51  int layer = 1*yid[0].id + yid[1].id - 1 ; // modified on 7/27/2015 to match new geometry (Frederic Georges)
52  //int layer = 2*yid[0].id + yid[1].id - 2 ;
53  int sector = yid[2].id;
54 
55  //yid[3].id = fmts.FindStrip(layer-1, sector-1, x, y, z);
56  double depe = aStep->GetTotalEnergyDeposit();
57  //cout << "resolMM " << layer << " " << x << " " << y << " " << z << " " << depe << " " << aStep->GetTrack()->GetTrackID() << endl;
58  vector<double> multi_hit = fmts.FindStrip(layer-1, sector-1, x, y, z, depe);
59 
60  int n_multi_hits = multi_hit.size()/2;
61 
62  // closest strip
63  //yid[4].id = (int) multi_hit[0];
64  yid[3].id = (int) multi_hit[0];
65 
66  yid[0].id_sharing = multi_hit[1];
67  yid[1].id_sharing = multi_hit[1];
68  yid[2].id_sharing = multi_hit[1];
69  yid[3].id_sharing = multi_hit[1];
70  // yid[4].id_sharing = multi_hit[1];
71 
72  // additional strip
73  for(int h=1; h<n_multi_hits; h++)
74  {
75  for(int j=0; j<3; j++)
76  {
77  identifier this_id;
78  this_id.name = yid[j].name;
79  this_id.rule = yid[j].rule;
80  this_id.id = yid[j].id;
81  this_id.time = yid[j].time;
82  this_id.TimeWindow = yid[j].TimeWindow;
83  this_id.TrackId = yid[j].TrackId;
84  this_id.id_sharing = multi_hit[3];
85  yid.push_back(this_id);
86  }
87  // last id is strip
88  identifier this_id;
89  this_id.name = yid[3].name;
90  this_id.rule = yid[3].rule;
91  this_id.id = (int) multi_hit[2];
92  this_id.time = yid[3].time;
93  this_id.TimeWindow = yid[3].TimeWindow;
94  this_id.TrackId = yid[3].TrackId;
95  this_id.id_sharing = multi_hit[3];
96  yid.push_back(this_id);
97  }
98 
99  return yid;
100 }
101 
102 
103 
104 
105 map< string, vector <int> > FMT_HitProcess :: multiDgt(MHit* aHit, int hitn)
106 {
107  map< string, vector <int> > MH;
108 
109  return MH;
110 }
111 
112 
113 
114 
115 
116 
117 
118 
119 
120 
121 
122 
void fill_infos()
Definition: fmt_strip.cc:7
double verbosity
vector< identifier > processID(vector< identifier >, G4Step *, detector)
vector< identifier > GetId()
Definition: Hit.h:103
int id
manually assing ID. 0 if "ncopy" (will be set at hit processing time)
Definition: identifier.h:39
string name
Name of the detector.
Definition: identifier.h:37
map< string, double > integrateDgt(MHit *, int)
double z
Definition: HitProcess.h:41
double x
Definition: HitProcess.h:41
double TimeWindow
Time Window. If abs(steptime - time) is smaller than TimeWindow, it&#39;s the same hit.
Definition: identifier.h:41
int TrackId
If Time Window is 0, it&#39;s a "flux" detector: if it&#39;s the same track, it&#39;s the same hit...
Definition: identifier.h:42
vector< double > FindStrip(int layer, int sector, double x, double y, double z, double Edep)
Definition: fmt_strip.cc:69
Definition: Hit.h:22
double y
Definition: HitProcess.h:41
double id_sharing
A single step can generate multiple identifiers. This variable represent the percentage sharing of th...
Definition: identifier.h:43
double time
Time of the first step.
Definition: identifier.h:40
string rule
"manual" or "ncopy"
Definition: identifier.h:38
map< string, vector< int > > multiDgt(MHit *, int)