Class Reference for E1039 Core & Analysis Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
SQDigitizer.cc
Go to the documentation of this file.
1 #include "SQDigitizer.h"
2 
6 #include <g4main/PHG4Hitv1.h>
9 #include <phool/PHNodeIterator.h>
10 #include <phool/PHCompositeNode.h>
11 #include <phool/PHIODataNode.h>
12 #include <phool/getClass.h>
13 
14 #include <geom_svc/GeomSvc.h>
15 
16 #include <iomanip>
17 #include <cmath>
18 #include <fstream>
19 #include <sstream>
20 #include <string>
21 #include <vector>
22 
23 #include <Geant4/G4SystemOfUnits.hh>
24 #include <Geant4/Randomize.hh>
25 
26 #include <TMath.h>
27 #include <TSQLServer.h>
28 #include <TSQLResult.h>
29 #include <TSQLRow.h>
30 
31 #ifndef __CINT__
32 #include <boost/lexical_cast.hpp>
33 #endif
34 
35 #ifndef __CINT__
37 {
38  if(Verbosity() > 2) std::cout << "SQDigitizer: Init " << detIDByName.size() << std::endl;
39 
40  p_geomSvc = GeomSvc::instance();
42  {
43  if(!enableDC1 && (i >= 7 && i <= 12)) continue;
44  if(!enableDPHodo && (i > nChamberPlanes+nHodoPlanes+nPropPlanes && i <= nChamberPlanes+nHodoPlanes+nPropPlanes+nDarkPhotonPlanes)) continue;
45 
46  std::string detName = p_geomSvc->getDetectorName(i);
47  detIDByName[detName] = i;
48  }
49 
50  //TODO: Load the detector realization setup
51 
53 }
54 #endif
55 
56 SQDigitizer::SQDigitizer(const std::string& name, const int verbose):
57  SubsysReco(name),
58  p_geomSvc(nullptr),
59  enableDC1(false),
60  enableDPHodo(true)
61 {
62  detIDByName.clear();
63  Verbosity(0);
64 }
65 
67 {}
68 
70 {
71  if(Verbosity() > 2) std::cout << "SQDigitizer: InitRun " << detIDByName.size() << std::endl;
72 
73  PHNodeIterator iter(topNode);
74 
75  // Looking for the DST node
76  PHCompositeNode *dstNode;
77  dstNode = dynamic_cast<PHCompositeNode*>(iter.findFirst("PHCompositeNode", "DST"));
78  if(!dstNode)
79  {
80  std::cerr << Name() << " DST Node missing, abort." << std::endl;
82  }
83 
84  // Book input g4hits - everything should be present
85  hitContainerByName.clear();
86  for(auto it = detIDByName.begin(); it != detIDByName.end(); ++it)
87  {
88 
89  std::string detName = it->first;
90  std::string g4hitNodeName = "G4HIT_" + detName;
91  if(Verbosity() > 2)
92  {
93  std::cout << Name() << ": booking input G4HIT node " << g4hitNodeName << std::endl;
94  }
95 
96  PHG4HitContainer* hits = findNode::getClass<PHG4HitContainer>(topNode, g4hitNodeName.c_str());
97  if(!hits)
98  {
99  if(Verbosity() > 2) std::cout << Name() << ": Could not locate g4 hit node " << g4hitNodeName << std::endl;
100  //return Fun4AllReturnCodes::ABORTRUN;
101  }
102  else
103  {
104  hitContainerByName[detName] = hits;
105  }
106  }
107 
108  // Exit if there is nothing to work with
109  if(hitContainerByName.empty())
110  {
111  std::cerr << Name() << ": no g4hit node foud, abort " << std::endl;
113  }
114 
115  //Book output SQHits
116  digits = findNode::getClass<SQHitVector>(topNode, "SQHitVector");
117  if(!digits)
118  {
119  if(Verbosity() > 2) std::cout << Name() << ": booking output node SQHitVector" << std::endl;
120 
121  digits = new SQHitVector_v1();
122  PHIODataNode<PHObject>* newNode = new PHIODataNode<PHObject>(digits, "SQHitVector", "PHObject");
123  dstNode->addNode(newNode);
124  }
125 
127 }
128 
130 {
131  for(auto it = hitContainerByName.begin(); it != hitContainerByName.end(); ++it)
132  {
134  {
135  std::cout << Name() << ": Digitizing hits for " << it->first << std::endl;
136  }
137 
138  int detID = detIDByName[it->first];
140  {
141  digitizePlane(it->first);
142  }
143  else
144  {
145  digitizeEMCal(it->first);
146  }
147  }
148 
150 }
151 
152 void SQDigitizer::digitizePlane(const std::string& detName)
153 {
154  PHG4HitContainer* g4hits = hitContainerByName[detName];
155  if(g4hits->size() < 1) return;
156 
157  int detID = detIDByName[detName];
158  for(PHG4HitContainer::ConstIterator it = g4hits->getHits().first; it != g4hits->getHits().second; ++it)
159  {
160  const PHG4Hit& g4hit = *(it->second);
161 
162  int track_id = g4hit.get_trkid();
163  if(track_id < 0) continue; //only save primary track hits
164 
165  //get average momentum and position
166  double x = 0.5*(g4hit.get_x(0) + g4hit.get_x(1));
167  double y = 0.5*(g4hit.get_y(0) + g4hit.get_y(1));
168  if(!p_geomSvc->isInPlane(detID, x, y)) continue; //only save in-acceptance hits
169 
170  double z = 0.5*(g4hit.get_z(0) + g4hit.get_z(1));
171  double px = g4hit.get_px(0);
172  double py = g4hit.get_py(0);
173  double pz = g4hit.get_pz(0);
174 
175  double tx = px/pz;
176  double ty = py/pz;
177  double x0 = x - tx*z;
178  double y0 = y - ty*z;
179 
180  double z_ref = p_geomSvc->getPlanePosition(detID);
181  double x_ref = tx*z_ref + x0;
182  double y_ref = ty*z_ref + y0;
183 
184  double w = p_geomSvc->getInterception(detID, tx, ty, x0, y0);
185  int eleID = p_geomSvc->getExpElementID(detID, w);
186  if(eleID < 1 || eleID > p_geomSvc->getPlaneNElements(detID)) continue; //only save hits within active region
187 
188  SQMCHit_v1 digiHit;
189  digiHit.set_track_id(track_id);
190  digiHit.set_g4hit_id(g4hit.get_hit_id());
191  digiHit.set_truth_x(x_ref);
192  digiHit.set_truth_y(y_ref);
193  digiHit.set_truth_z(z_ref);
194  digiHit.set_truth_px(px);
195  digiHit.set_truth_py(py);
196  digiHit.set_truth_pz(pz);
197  digiHit.set_hit_id(digits->size());
198  digiHit.set_in_time(1);
199  digiHit.set_hodo_mask(0);
200  digiHit.set_detector_id(detID);
201  digiHit.set_element_id(eleID);
202  digiHit.set_tdc_time(0.);
203  digiHit.set_pos(p_geomSvc->getMeasurement(detID, eleID));
204 
205  //Drift distance is calculated differently for chamber and hodos
206  if(detID <= nChamberPlanes || (detID > nChamberPlanes+nHodoPlanes && detID <= nChamberPlanes+nHodoPlanes+nPropPlanes))
207  {
208  digiHit.set_drift_distance(p_geomSvc->getDCA(detID, eleID, tx, ty, x0, y0));
209  }
210  else
211  {
212  digiHit.set_drift_distance(0.);
213  }
214 
215  //push the hit to vector
216  digits->push_back(&digiHit);
217 
218  //special treatment for hodoscopes
219  double dw = w - digiHit.get_pos();
220  double hodoWidth = p_geomSvc->getCellWidth(detID)/2.;
221  double hodoOverlap = p_geomSvc->getPlaneOverlap(detID);
222  if(fabs(dw) > hodoWidth - hodoOverlap && fabs(dw) < hodoWidth) //hit happens in the overlap region
223  {
224  if(dw < 0. && eleID != 1)
225  {
226  digiHit.set_element_id(eleID-1);
227  digiHit.set_pos(p_geomSvc->getMeasurement(detID, eleID-1));
228  digiHit.set_hit_id(digits->size());
229 
230  digits->push_back(&digiHit);
231  }
232  else if(dw > 0. && eleID != p_geomSvc->getPlaneNElements(detID))
233  {
234  digiHit.set_element_id(eleID+1);
235  digiHit.set_pos(p_geomSvc->getMeasurement(detID, eleID+1));
236  digiHit.set_hit_id(digits->size());
237 
238  digits->push_back(&digiHit);
239  }
240  }
241  }
242 }
243 
244 void SQDigitizer::digitizeEMCal(const std::string& detName)
245 {
246  PHG4HitContainer* g4hits = hitContainerByName[detName];
247  if(g4hits->size() < 1) return;
248 
249  int detID = detIDByName[detName];
250  std::map<int, SQCalHit_v1> digiHits; //key -> towerID, val -> calHit
251  for(PHG4HitContainer::ConstIterator it = g4hits->getHits().first; it != g4hits->getHits().second; ++it)
252  {
253  const PHG4Hit& g4hit = *(it->second);
254 
255  int towerID = g4hit.get_scint_id();
256  if(digiHits.find(towerID) == digiHits.end())
257  {
258  SQCalMCHit_v1 digiHit;
259  digiHit.set_detector_id(detID);
260  digiHit.set_element_id(towerID);
261  digiHit.set_shower_id(g4hit.get_shower_id());
262  digiHit.set_track_id(g4hit.get_trkid());
263  digiHit.set_g4hit_id(g4hit.get_hit_id());
264 
265  digiHit.set_in_time(1);
266  digiHit.set_hodo_mask(0);
267  digiHit.set_tdc_time(0.);
268  digiHit.set_drift_distance(0.);
269  digiHit.set_pos(0.);
270 
271  digiHit.set_truth_x(g4hit.get_x(0));
272  digiHit.set_truth_y(g4hit.get_y(0));
273  digiHit.set_truth_z(g4hit.get_z(0));
274  digiHit.set_truth_px(g4hit.get_px(0));
275  digiHit.set_truth_py(g4hit.get_py(0));
276  digiHit.set_truth_pz(g4hit.get_pz(0));
277 
278  digiHits[towerID] = digiHit;
279  }
280 
281  digiHits[towerID].add_cell(g4hit.get_index_l(), g4hit.get_edep());
282  }
283 
284  for(auto iter = digiHits.begin(); iter != digiHits.end(); ++iter)
285  {
286  //iter->second.identify();
287  iter->second.set_hit_id(digits->size());
288  digits->push_back(&(iter->second));
289  }
290 }
291 
293 {
294  //TODO: implement efficiency and calibration smearing
295  return true;
296 }
virtual void set_element_id(const short id)
Definition: SQHit_v1.h:40
virtual void set_truth_px(const float a)
Definition: SQMCHit_v1.h:47
virtual void set_pos(const float a)
Definition: SQHit.h:61
virtual float get_x(const int i) const
Definition: PHG4Hit.h:21
virtual void set_hodo_mask(const bool a)
Definition: SQHit.h:94
std::string getDetectorName(const int &detectorID) const
Definition: GeomSvc.h:188
virtual float get_pz(const int i) const
Definition: PHG4Hit.h:26
virtual int get_shower_id() const
Definition: PHG4Hit.h:38
virtual void set_truth_px(const float a)
Definition: SQCalMCHit_v1.h:43
virtual float get_pos() const
Return the absolute position of this hit. Probably the value is not properly set at present...
Definition: SQHit_v1.h:51
PHNode * findFirst(const std::string &, const std::string &)
double getCellWidth(int detectorID)
Definition: GeomSvc.h:200
#define nPropPlanes
Definition: GlobalConsts.h:8
virtual void set_tdc_time(const float a)
Definition: SQHit.h:55
virtual void set_truth_py(const float a)
Definition: SQMCHit_v1.h:50
PHBoolean addNode(PHNode *)
virtual void set_detector_id(const short a)
Definition: SQHit_v1.h:37
void digitizePlane(const std::string &detName)
main external call, fill the digi hit vector
Definition: SQDigitizer.cc:152
virtual int get_trkid() const
Definition: PHG4Hit.h:41
virtual void set_hit_id(const int a)
Definition: SQHit_v1.h:34
An SQ interface class to hold one detector hit.
Definition: SQHit.h:20
void digitizeEMCal(const std::string &detName)
digitize the emcal hits
Definition: SQDigitizer.cc:244
SQDigitizer(const std::string &name="SQDigitizer", const int verbose=0)
Definition: SQDigitizer.cc:56
virtual void set_g4hit_id(const PHG4HitDefs::keytype a)
Definition: SQCalMCHit_v1.h:31
virtual int get_index_l() const
Definition: PHG4Hit.h:49
virtual void set_g4hit_id(const PHG4HitDefs::keytype a)
Definition: SQMCHit_v1.h:35
virtual void set_truth_pz(const float a)
Definition: SQMCHit_v1.h:53
#define nHodoPlanes
Definition: GlobalConsts.h:7
virtual void set_detector_id(const short a)
Definition: SQCalHit_v1.h:29
virtual PHG4HitDefs::keytype get_hit_id() const
Definition: PHG4Hit.h:36
virtual void set_track_id(const int a)
Definition: SQCalMCHit_v1.h:25
double getPlaneOverlap(int detectorID) const
Definition: GeomSvc.h:199
virtual void set_drift_distance(const float a)
Definition: SQHit.h:58
virtual const std::string Name() const
Returns the name of this module.
Definition: Fun4AllBase.h:23
virtual float get_edep() const
Definition: PHG4Hit.h:31
ConstRange getHits(const unsigned int detid) const
return all hits matching a given detid
virtual float get_y(const int i) const
Definition: PHG4Hit.h:22
virtual void set_truth_x(const float a)
Definition: SQMCHit_v1.h:38
double getInterception(int detectorID, double tx, double ty, double x0, double y0) const
Get the interception of a line an a plane.
Definition: GeomSvc.h:232
virtual void set_truth_y(const float a)
Definition: SQMCHit_v1.h:41
virtual int Verbosity() const
Gets the verbosity of this module.
Definition: Fun4AllBase.h:64
virtual void push_back(const SQHit *hit)
Definition: SQHitVector.h:55
Output a lot of messages.
Definition: Fun4AllBase.h:48
double getDCA(int detectorID, int elementID, double tx, double ty, double x0, double y0)
Definition: GeomSvc.cxx:940
Definition: PHG4Hit.h:9
Map::const_iterator ConstIterator
virtual void set_pos(const float a)
Definition: SQHit_v1.h:52
virtual void set_drift_distance(const float a)
Definition: SQHit_v1.h:49
int InitRun(PHCompositeNode *topNode)
module initialization
Definition: SQDigitizer.cc:69
virtual float get_z(const int i) const
Definition: PHG4Hit.h:23
int process_event(PHCompositeNode *topNode)
event processing
Definition: SQDigitizer.cc:129
int Init(PHCompositeNode *topNode)
Definition: SQDigitizer.cc:36
virtual void set_in_time(const bool a)
Definition: SQHit.h:91
#define nChamberPlanes
Definition: GlobalConsts.h:6
static GeomSvc * instance()
singlton instance
Definition: GeomSvc.cxx:211
virtual size_t size() const
Definition: SQHitVector.h:50
virtual float get_px(const int i) const
Definition: PHG4Hit.h:24
virtual void set_hodo_mask(const bool a)
Definition: SQHit_v1.h:58
virtual void set_truth_pz(const float a)
Definition: SQCalMCHit_v1.h:49
virtual void set_shower_id(const int a)
Definition: SQCalMCHit_v1.h:28
int getExpElementID(int detectorID, double pos_exp)
Definition: GeomSvc.cxx:807
virtual void set_truth_y(const float a)
Definition: SQCalMCHit_v1.h:37
int getPlaneNElements(int detectorID)
Definition: GeomSvc.h:208
virtual void set_truth_z(const float a)
Definition: SQMCHit_v1.h:44
unsigned int size(void) const
virtual void set_truth_x(const float a)
Definition: SQCalMCHit_v1.h:34
void getMeasurement(int detectorID, int elementID, double &measurement, double &dmeasurement)
Convert the detectorID and elementID to the actual hit position.
Definition: GeomSvc.cxx:781
#define nDarkPhotonPlanes
Definition: GlobalConsts.h:9
bool isInPlane(int detectorID, double x, double y)
See if a point is in a plane.
Definition: GeomSvc.cxx:752
virtual void set_element_id(const short id)
Definition: SQCalHit_v1.h:32
virtual void set_in_time(const bool a)
Definition: SQHit_v1.h:55
virtual void set_track_id(const int a)
Definition: SQMCHit_v1.h:32
virtual float get_py(const int i) const
Definition: PHG4Hit.h:25
virtual ~SQDigitizer()
Definition: SQDigitizer.cc:66
virtual void set_truth_py(const float a)
Definition: SQCalMCHit_v1.h:46
virtual void set_tdc_time(const float a)
Definition: SQHit_v1.h:46
virtual int get_scint_id() const
Definition: PHG4Hit.h:39
virtual void set_truth_z(const float a)
Definition: SQCalMCHit_v1.h:40
double getPlanePosition(int detectorID) const
Definition: GeomSvc.h:197
bool realize(SQHit &dHit)
realization process
Definition: SQDigitizer.cc:292