Class Reference for E1039 Core & Analysis Software
AnaTrkQA.h
Go to the documentation of this file.
1 
9 #ifndef _H_AnaTrkQA_H_
10 #define _H_AnaTrkQA_H_
11 
12 // ROOT
13 #include <TSQLServer.h>
14 #include <TSQLResult.h>
15 #include <TSQLRow.h>
16 #include <TMatrixD.h>
17 #include <TVector3.h>
18 
19 // Fun4All includes
20 #include <fun4all/SubsysReco.h>
21 
22 // STL includes
23 #include <vector>
24 #include <string>
25 #include <iostream>
26 #include <list>
27 #include <map>
28 
29 class TVector3;
30 class TLorentzVector;
31 
32 class SQRun;
33 class SQSpillMap;
34 class SQEvent;
35 class SQHitMap;
36 class SQHitVector;
37 class SQHit;
38 
40 class PHG4HitContainer;
41 class SRecEvent;
42 class SRecTrack;
43 class GeomSvc;
44 
45 class TFile;
46 class TTree;
47 
48 class AnaTrkQA: public SubsysReco {
49 
50  public:
51 
52  AnaTrkQA(const std::string &name = "AnaTrkQA.root");
53  virtual ~AnaTrkQA() {
54  }
55 
56  int Init(PHCompositeNode *topNode);
57  int InitRun(PHCompositeNode *topNode);
58  int process_event(PHCompositeNode *topNode);
59  int End(PHCompositeNode *topNode);
60 
61  int InitEvalTree();
62  int ResetEvalVars();
63 
64  const std::string& get_hit_container_choice() const {
65  return _hit_container_type;
66  }
67 
68  void set_hit_container_choice(const std::string& hitContainerChoice) {
69  _hit_container_type = hitContainerChoice;
70  }
71 
72  const std::string& get_out_name() const {
73  return _out_name;
74  }
75 
76  void set_out_name(const std::string& outName) {
77  _out_name = outName;
78  }
79 
80  private:
81 
82  int GetNodes(PHCompositeNode *topNode);
83 
84  int TruthRecoEval(PHCompositeNode *topNode);
85 
86  bool FindG4HitAtStation(const int trk_id, const PHG4HitContainer* g4hc, TVector3* pos, TLorentzVector* mom);
87  int FindCommonHitIDs(std::vector<int>& hitidvec1, std::vector<int>& hitidvec2);
88  SRecTrack* FindBestMomRecTrack(SRecEvent *recEvent, const float true_P);
89  bool FindG4HitAtHodo(const int trk_id, const PHG4HitContainer* g4hc);
90  bool FindG4HitAtProp(const int trk_id, const PHG4HitContainer* g4hc);
91 
92  std::string _hit_container_type;
93 
94  size_t _event;
95  SQRun* _run_header;
96  SQSpillMap * _spill_map;
97 
98  SQEvent * _event_header;
99  SQHitMap *_hit_map;
100  SQHitVector *_hit_vector;
101 
102  PHG4TruthInfoContainer* _truth;
103  SRecEvent* _recEvent;
104 
105  PHG4HitContainer *g4hc_d1x;
106  PHG4HitContainer *g4hc_d2xp;
107  PHG4HitContainer *g4hc_d3px;
108  PHG4HitContainer *g4hc_d3mx;
109 
110  PHG4HitContainer *g4hc_h1t;
111  PHG4HitContainer *g4hc_h1b;
112  PHG4HitContainer *g4hc_h2t;
113  PHG4HitContainer *g4hc_h2b;
114  PHG4HitContainer *g4hc_h3t;
115  PHG4HitContainer *g4hc_h3b;
116  PHG4HitContainer *g4hc_h4t;
117  PHG4HitContainer *g4hc_h4b;
118 
119  PHG4HitContainer *g4hc_p1y1;
120  PHG4HitContainer *g4hc_p1y2;
121  PHG4HitContainer *g4hc_p1x1;
122  PHG4HitContainer *g4hc_p1x2;
123  PHG4HitContainer *g4hc_p2x1;
124  PHG4HitContainer *g4hc_p2x2;
125  PHG4HitContainer *g4hc_p2y1;
126  PHG4HitContainer *g4hc_p2y2;
127 
128  std::string _out_name;
129 
130  TTree* _tout_reco;
131  TTree* _qa_tree;
132  TFile *file;
133 
134  int run_id;
135  int spill_id;
136  float target_pos;
137  int event_id;
138  int krecstat;
139  int kalman_stat;
140  unsigned short emu_trigger;
141 
142  int n_hits;
143  int nhits_st1[100];
144  int nhits_st2[100];
145  int nhits_st3[100];
146  int hit_id[100];
147  int detector_id[100];
148  int element_id[100];
149  int hodo_mask[100];
150  float drift_distance[100];
151  float pos[100];
152  float detector_z[100];
153 
154  float truth_x[100];
155  float truth_y[100];
156  float truth_z[100];
157  float truth_pos[100];
158 
159  int n_tracks;
160  int n_recTracks;
161  int rec_id[100];
162  int par_id[100];
163  int pid[100];
164 
165  float gvx[100];
166  float gvy[100];
167  float gvz[100];
168  float gpx[100];
169  float gpy[100];
170  float gpz[100];
171  float gx_st1[100];
172  float gy_st1[100];
173  float gz_st1[100];
174 
175  float ac_gpx[100];
176  float ac_gpy[100];
177  float ac_gpz[100];
178 
179  float sq_pos_st1[100];
180  float sq_drift_st1[100];
181  float sq_decID[100];
182  float sq_x_st1[100];
183  float sq_y_st1[100];
184  float sq_z_st1[100];
185  float sq_px_st1[100];
186  float sq_py_st1[100];
187  float sq_pz_st1[100];
188 
189 
190  float sq_pos_st2[100];
191  float sq_drift_st2[100];
192  float sq_x_st2[100];
193  float sq_y_st2[100];
194  float sq_z_st2[100];
195  float sq_px_st2[100];
196  float sq_py_st2[100];
197  float sq_pz_st2[100];
198 
199 
200  float sq_pos_st3[100];
201  float sq_drift_st3[100];
202  float sq_x_st3[100];
203  float sq_y_st3[100];
204  float sq_z_st3[100];
205  float sq_px_st3[100];
206  float sq_py_st3[100];
207  float sq_pz_st3[100];
208 
209 
210  float chisq_st1[100];
211  float prob_st1[100];
212  float quality[100];
213  float fit_time[100];
214 
215  float pull_q2p_st1[100];
216  float pull_q2p_st2[100];
217  float pull_q2p_st3[100];
218 
219 
220  float rec_drift_st1[100];
221  float rec_px_st1[100];
222  float rec_py_st1[100];
223  float rec_pz_st1[100];
224  float rec_p_st1[100];
225  float rec_x_st1[100];
226  float rec_y_st1[100];
227  float rec_z_st1[100];
228 
229  float rec_drift_st2[100];
230  float rec_px_st2[100];
231  float rec_py_st2[100];
232  float rec_pz_st2[100];
233  float rec_p_st2[100];
234  float rec_x_st2[100];
235  float rec_y_st2[100];
236  float rec_z_st2[100];
237 
238 
239  float rec_drift_st3[100];
240  float rec_px_st3[100];
241  float rec_py_st3[100];
242  float rec_pz_st3[100];
243  float rec_p_st3[100];
244  float rec_x_st3[100];
245  float rec_y_st3[100];
246  float rec_z_st3[100];
247 
248  float gpx_st1[100];
249  float gpy_st1[100];
250  float gpz_st1[100];
251  float gpt[100];
252  float geta[100];
253  float gphi[100];
254  int gnhits[100];
255  int gndc[100];
256  int gnhodo[100];
257  int gnprop[100];
258  int gndp[100];
259  int ntruhits[100];
260  int nhits[100];
261  int charge[100];
262  float rec_vx[100];
263  float rec_vy[100];
264  float rec_vz[100];
265  float rec_px[100];
266  float rec_py[100];
267  float rec_pz[100];
268  float rec_pt[100];
269  float rec_eta[100];
270  float rec_phi[100];
271  float pull_state00[100];
272 
273  int gelmid[1000][128];
274 
275 
276 
277  GeomSvc *p_geomSvc;
278 };
279 
280 
281 #endif /* _H_AnaTrkQA_H_ */
int ResetEvalVars()
Definition: AnaTrkQA.cxx:598
int process_event(PHCompositeNode *topNode)
Definition: AnaTrkQA.cxx:88
AnaTrkQA(const std::string &name="AnaTrkQA.root")
Definition: AnaTrkQA.cxx:57
virtual ~AnaTrkQA()
Definition: AnaTrkQA.h:53
const std::string & get_out_name() const
Definition: AnaTrkQA.h:72
int InitEvalTree()
Definition: AnaTrkQA.cxx:472
const std::string & get_hit_container_choice() const
Definition: AnaTrkQA.h:64
int End(PHCompositeNode *topNode)
===========================
Definition: AnaTrkQA.cxx:462
void set_out_name(const std::string &outName)
Definition: AnaTrkQA.h:76
int Init(PHCompositeNode *topNode)
Definition: AnaTrkQA.cxx:71
void set_hit_container_choice(const std::string &hitContainerChoice)
Definition: AnaTrkQA.h:68
int InitRun(PHCompositeNode *topNode)
Definition: AnaTrkQA.cxx:75
User interface class about the geometry of detector planes.
Definition: GeomSvc.h:164
An SQ interface class to hold one event header.
Definition: SQEvent.h:17
An SQ interface class to hold a list of SQHit objects as std::map.
Definition: SQHitMap.h:23
An SQ interface class to hold a list of SQHit objects.
Definition: SQHitVector.h:32
An SQ interface class to hold one detector hit.
Definition: SQHit.h:20
An SQ interface class to hold the run-level info.
Definition: SQRun.h:18
An SQ interface class to hold a list of SQSpill objects.
Definition: SQSpillMap.h:19