Class Reference for E1039 Core & Analysis Software
AnaTrack.cc
Go to the documentation of this file.
1 #include <fstream>
2 #include <iomanip>
3 #include <TSystem.h>
4 #include <TFile.h>
5 #include <TTree.h>
6 #include <TChain.h>
7 #include <TH1D.h>
8 #include <TH2D.h>
9 #include <THStack.h>
10 #include <TCanvas.h>
11 #include <TLatex.h>
12 #include <interface_main/SQRun.h>
14 #include <interface_main/SQEvent.h>
16 #include <ktracker/SRecEvent.h>
18 #include <phool/PHNodeIterator.h>
19 #include <phool/PHIODataNode.h>
20 #include <phool/getClass.h>
21 #include <geom_svc/GeomSvc.h>
22 //#include <UtilAna/UtilHist.h>
23 #include "AnaTrack.h"
24 using namespace std;
25 
26 AnaTrack::AnaTrack(const std::string& name)
27  : SubsysReco (name)
28  , m_sq_evt (0)
29  , m_sq_hit_vec(0)
30  , m_sq_trk_vec(0)
31  , m_file_name ("output_track.root")
32  , m_file (0)
33  , m_tree (0)
34 {
35  ;
36 }
37 
39 {
40  ;
41 }
42 
43 int AnaTrack::Init(PHCompositeNode* topNode)
44 {
46 }
47 
49 {
50  //GeomSvc* geom = GeomSvc::instance();
51 
52  m_sq_evt = findNode::getClass<SQEvent >(topNode, "SQEvent");
53  m_sq_hit_vec = findNode::getClass<SQHitVector >(topNode, "SQHitVector");
54  m_sq_trk_vec = findNode::getClass<SQTrackVector >(topNode, "SQRecTrackVector");
55  if (!m_sq_evt || !m_sq_hit_vec || !m_sq_trk_vec) return Fun4AllReturnCodes::ABORTEVENT;
56 
57  m_file = new TFile(m_file_name.c_str(), "RECREATE");
58  m_tree = new TTree("tree", "Created by AnaTrack");
59  m_tree->Branch("event" , &m_evt);
60  m_tree->Branch("track_list", &m_trk_list);
61 
62  SQRun* sq_run = findNode::getClass<SQRun>(topNode, "SQRun");
63  if (!sq_run) return Fun4AllReturnCodes::ABORTEVENT;
64 
66 }
67 
69 {
70  m_evt.run_id = m_sq_evt->get_run_id();
71  m_evt.spill_id = m_sq_evt->get_spill_id();
72  m_evt.event_id = m_sq_evt->get_event_id();
73  m_evt.fpga_bits = (m_sq_evt->get_trigger() >> SQEvent::MATRIX1) & 0x1f;
74  m_evt.nim_bits = (m_sq_evt->get_trigger() >> SQEvent::NIM1 ) & 0x1f;
75 
76  m_evt.D1 = m_evt.D2 = m_evt.D3p = m_evt.D3m = 0;
77  for (SQHitVector::Iter it = m_sq_hit_vec->begin(); it != m_sq_hit_vec->end(); it++) {
78  SQHit* hit = *it;
79  int det_id = hit->get_detector_id();
80  if ( 0 < det_id && det_id <= 6) m_evt.D1++;
81  else if (12 < det_id && det_id <= 18) m_evt.D2++;
82  else if (18 < det_id && det_id <= 24) m_evt.D3p++;
83  else if (24 < det_id && det_id <= 30) m_evt.D3m++;
84  }
85 
86  m_trk_list.clear();
87  for (auto it = m_sq_trk_vec->begin(); it != m_sq_trk_vec->end(); it++) {
88  SRecTrack* trk = dynamic_cast<SRecTrack*>(*it);
89  TrackData td;
90  td.charge = trk->get_charge();
91  td.road = trk->getTriggerRoad();
92  td.n_hits = trk->get_num_hits();
93  td.chisq = trk->get_chisq();
94  td.chisq_target = trk->get_chisq_target();
96  td.chisq_dump = trk->get_chisq_dump();
97  td.pos_vtx = trk->get_pos_vtx();
98  td.mom_vtx = trk->get_mom_vtx();
99  td.pos_st1 = trk->get_pos_st1();
100  td.mom_st1 = trk->get_mom_st1();
101  td.pos_st3 = trk->get_pos_st3();
102  td.mom_st3 = trk->get_mom_st3();
103  m_trk_list.push_back(td);
104  }
105 
106  m_tree->Fill();
108 }
109 
110 int AnaTrack::End(PHCompositeNode* topNode)
111 {
112  m_file->cd();
113  m_file->Write();
114  m_file->Close();
116 }
117 
118 void AnaTrack::AnalyzeTree(TChain* tree)
119 {
120  string dir_out = "result/track";
121  cout << "N of trees = " << tree->GetNtrees() << endl;
122  gSystem->mkdir(dir_out.c_str(), true);
123  ofstream ofs(dir_out + "/result.txt");
124 
125  TFile* file_out = new TFile((dir_out+"/result.root").c_str(), "RECREATE");
126 
127  TH2* h2_nhit = new TH2D("h2_nhit", ";N of hits/track;", 6, 12.5, 18.5, 2, -2, 2);
128  TH2* h2_chi2 = new TH2D("h2_chi2", ";Track #chi^{2};", 100, 0, 10 , 2, -2, 2);
129  TH2* h2_chi2_tgt = new TH2D("h2_chi2_tgt", ";Track #chi^{2} at target;" , 100, 0, 10, 2, -2, 2);
130  TH2* h2_chi2_dum = new TH2D("h2_chi2_dum", ";Track #chi^{2} at dump;" , 100, 0, 10, 2, -2, 2);
131  TH2* h2_chi2_ups = new TH2D("h2_chi2_ups", ";Track #chi^{2} at upstream;", 100, 0, 10, 2, -2, 2);
132 
133  TH2* h2_x_vtx = new TH2D("h2_x_vtx" , ";Track x (cm) @ Vertex;Charge;", 100, -1, 1, 2, -2, 2);
134  TH2* h2_y_vtx = new TH2D("h2_y_vtx" , ";Track y (cm) @ Vertex;Charge;", 100, -1, 1, 2, -2, 2);
135  TH2* h2_z_vtx = new TH2D("h2_z_vtx" , ";Track z (cm) @ Vertex;Charge;", 120, -700, 500, 2, -2, 2);
136  TH2* h2_px_vtx = new TH2D("h2_px_vtx", ";Track p_{x} (GeV) @ Vertex;Charge;", 100, -5, 5, 2, -2, 2);
137  TH2* h2_py_vtx = new TH2D("h2_py_vtx", ";Track p_{y} (GeV) @ Vertex;Charge;", 100, -5, 5, 2, -2, 2);
138  TH2* h2_pz_vtx = new TH2D("h2_pz_vtx", ";Track p_{z} (GeV) @ Vertex;Charge;", 100, 0, 100, 2, -2, 2);
139 
140  TH2* h2_x_st1 = new TH2D("h2_x_st1" , ";Track x (cm) @ St 1;Charge;", 100, -50, 50, 2, -2, 2);
141  TH2* h2_y_st1 = new TH2D("h2_y_st1" , ";Track y (cm) @ St 1;Charge;", 100, -50, 50, 2, -2, 2);
142  TH2* h2_z_st1 = new TH2D("h2_z_st1" , ";Track z (cm) @ St 1;Charge;", 100, 550, 650, 2, -2, 2);
143  TH2* h2_px_st1 = new TH2D("h2_px_st1", ";Track p_{x} (GeV) @ St 1;Charge;", 100, -5, 5, 2, -2, 2);
144  TH2* h2_py_st1 = new TH2D("h2_py_st1", ";Track p_{y} (GeV) @ St 1;Charge;", 100, -5, 5, 2, -2, 2);
145  TH2* h2_pz_st1 = new TH2D("h2_pz_st1", ";Track p_{z} (GeV) @ St 1;Charge;", 100, 0, 100, 2, -2, 2);
146 
147  TH2* h2_x_st3 = new TH2D("h2_x_st3" , ";Track x (cm) @ St 3;Charge;", 100, -150, 150, 2, -2, 2);
148  TH2* h2_y_st3 = new TH2D("h2_y_st3" , ";Track y (cm) @ St 3;Charge;", 100, -150, 150, 2, -2, 2);
149  TH2* h2_z_st3 = new TH2D("h2_z_st3" , ";Track z (cm) @ St 3;Charge;", 100, 1900, 2000, 2, -2, 2);
150  TH2* h2_px_st3 = new TH2D("h2_px_st3", ";Track p_{x} (GeV) @ St 3;Charge;", 100, -5, 5, 2, -2, 2);
151  TH2* h2_py_st3 = new TH2D("h2_py_st3", ";Track p_{y} (GeV) @ St 3;Charge;", 100, -5, 5, 2, -2, 2);
152  TH2* h2_pz_st3 = new TH2D("h2_pz_st3", ";Track p_{z} (GeV) @ St 3;Charge;", 100, 0, 100, 2, -2, 2);
153 
154  //GeomSvc* geom = GeomSvc::instance();
155  ostringstream oss;
156 
157  EventData* evt = 0;
158  TrackList* trk_list = 0;
159  tree->SetBranchAddress("event" , &evt);
160  tree->SetBranchAddress("track_list", &trk_list);
161 
162  int n_ent = tree->GetEntries();
163  cout << "N of entries = " << n_ent << endl;
164  for (int i_ent = 0; i_ent < n_ent; i_ent++) {
165  if ((i_ent+1) % (n_ent/10) == 0) cout << " " << 10*(i_ent+1)/(n_ent/10) << "%" << flush;
166  tree->GetEntry(i_ent);
167 
168  //if (! (evt->fpga_bits & 0x1)) continue;
169  //if (! (evt->fpga_bits & 0x4)) continue;
170  if (! (evt->fpga_bits & 0x8)) continue;
171  //if (! (evt->nim_bits & 0x4)) continue;
172 
173  for (auto it = trk_list->begin(); it != trk_list->end(); it++) {
174  TrackData* td = &(*it);
175  int charge = td->charge;
176  if (abs(charge) != 1) continue;
177  //int road = td->road;
178  int nhit = td->n_hits;
179  double chi2 = td->chisq;
180  double chi2_tgt = td->chisq_target;
181  double chi2_dum = td->chisq_dump;
182  double chi2_ups = td->chisq_upstream;
183  TVector3* pos_vtx = &td->pos_vtx;
184  TLorentzVector* mom_vtx = &td->mom_vtx;
185  TVector3* pos_st1 = &td->pos_st1;
186  TLorentzVector* mom_st1 = &td->mom_st1;
187  TVector3* pos_st3 = &td->pos_st3;
188  TLorentzVector* mom_st3 = &td->mom_st3;
189 
190  h2_nhit ->Fill(nhit, charge);
191  h2_chi2 ->Fill(chi2, charge);
192  h2_chi2_tgt->Fill(chi2_tgt, charge);
193  h2_chi2_dum->Fill(chi2_dum, charge);
194  h2_chi2_ups->Fill(chi2_ups, charge);
195 
196  h2_x_vtx ->Fill(pos_vtx->X(), charge);
197  h2_y_vtx ->Fill(pos_vtx->Y(), charge);
198  h2_z_vtx ->Fill(pos_vtx->Z(), charge);
199  h2_px_vtx ->Fill(mom_vtx->X(), charge);
200  h2_py_vtx ->Fill(mom_vtx->Y(), charge);
201  h2_pz_vtx ->Fill(mom_vtx->Z(), charge);
202 
203  h2_x_st1 ->Fill(pos_st1->X(), charge);
204  h2_y_st1 ->Fill(pos_st1->Y(), charge);
205  h2_z_st1 ->Fill(pos_st1->Z(), charge);
206  h2_px_st1 ->Fill(mom_st1->X(), charge);
207  h2_py_st1 ->Fill(mom_st1->Y(), charge);
208  h2_pz_st1 ->Fill(mom_st1->Z(), charge);
209 
210  h2_x_st3 ->Fill(pos_st3->X(), charge);
211  h2_y_st3 ->Fill(pos_st3->Y(), charge);
212  h2_z_st3 ->Fill(pos_st3->Z(), charge);
213  h2_px_st3 ->Fill(mom_st3->X(), charge);
214  h2_py_st3 ->Fill(mom_st3->Y(), charge);
215  h2_pz_st3 ->Fill(mom_st3->Z(), charge);
216  }
217  }
218 
219  DrawHistIn1D(h2_nhit, "nhit", dir_out);
220  DrawHistIn1D(h2_chi2, "chi2", dir_out);
221  DrawHistIn1D(h2_chi2_tgt, "chi2_tgt", dir_out);
222  DrawHistIn1D(h2_chi2_dum, "chi2_dum", dir_out);
223  DrawHistIn1D(h2_chi2_ups, "chi2_ups", dir_out);
224 
225  DrawHistIn1D(h2_x_vtx , "x_vtx", dir_out);
226  DrawHistIn1D(h2_y_vtx , "y_vtx", dir_out);
227  DrawHistIn1D(h2_z_vtx , "z_vtx", dir_out);
228  DrawHistIn1D(h2_px_vtx, "px_vtx", dir_out);
229  DrawHistIn1D(h2_py_vtx, "py_vtx", dir_out);
230  DrawHistIn1D(h2_pz_vtx, "pz_vtx", dir_out);
231 
232  DrawHistIn1D(h2_x_st1 , "x_st1", dir_out);
233  DrawHistIn1D(h2_y_st1 , "y_st1", dir_out);
234  DrawHistIn1D(h2_z_st1 , "z_st1", dir_out);
235  DrawHistIn1D(h2_px_st1, "px_st1", dir_out);
236  DrawHistIn1D(h2_py_st1, "py_st1", dir_out);
237  DrawHistIn1D(h2_pz_st1, "pz_st1", dir_out);
238 
239  DrawHistIn1D(h2_x_st3 , "x_st3", dir_out);
240  DrawHistIn1D(h2_y_st3 , "y_st3", dir_out);
241  DrawHistIn1D(h2_z_st3 , "z_st3", dir_out);
242  DrawHistIn1D(h2_px_st3, "px_st3", dir_out);
243  DrawHistIn1D(h2_py_st3, "py_st3", dir_out);
244  DrawHistIn1D(h2_pz_st3, "pz_st3", dir_out);
245 
246  ofs.close();
247  file_out->Write();
248  file_out->Close();
249 }
250 
251 void AnaTrack::DrawHistIn1D(TH2* h2, const std::string label, const std::string dir_out)
252 {
253  TCanvas* c1 = new TCanvas("c1", "");
254  c1->SetGrid();
255  //c1->SetLogy(true);
256 
257  TH1* h1_neg = h2->ProjectionX("h1_neg", 1, 1);
258  TH1* h1_pos = h2->ProjectionX("h1_pos", 2, 2);
259  h1_pos->SetLineColor(kRed);
260  h1_neg->SetLineColor(kBlue);
261 
262  ostringstream oss;
263  oss << ";" << h2->GetXaxis()->GetTitle() << ";" << h2->GetZaxis()->GetTitle();
264  THStack hs("hs", oss.str().c_str());
265  hs.Add(h1_pos);
266  hs.Add(h1_neg);
267  hs.Draw("nostack");
268 
269  TLatex tex;
270  tex.SetNDC(true);
271  tex.SetTextColor(kRed);
272  tex.DrawLatex(0.25, 0.94, "Red: #mu^{+}");
273  tex.SetTextColor(kBlue);
274  tex.DrawLatex(0.75, 0.94, "Blue: #mu^{-}");
275 
276  oss.str("");
277  oss << dir_out << "/h1_" << label << ".png";
278  c1->SaveAs(oss.str().c_str());
279 
280  delete h1_neg;
281  delete h1_pos;
282  delete c1;
283 }
std::vector< TrackData > TrackList
Definition: TreeData.h:53
AnaTrack(const std::string &name="AnaTrack")
Definition: AnaTrack.cc:24
int InitRun(PHCompositeNode *topNode)
Definition: AnaTrack.cc:46
int Init(PHCompositeNode *topNode)
Definition: AnaTrack.cc:41
virtual ~AnaTrack()
Definition: AnaTrack.h:42
static void AnalyzeTree(TChain *tree)
Definition: AnaTrack.cc:118
static void DrawHistIn1D(TH2 *h2, const std::string label, const std::string dir_out)
Definition: AnaTrack.cc:251
int End(PHCompositeNode *topNode)
Called at the end of all processing.
Definition: AnaTrack.cc:195
int process_event(PHCompositeNode *topNode)
Definition: AnaTrack.cc:104
virtual int get_run_id() const =0
Return the run ID.
@ MATRIX1
Definition: SQEvent.h:27
@ NIM1
Definition: SQEvent.h:22
virtual bool get_trigger(const SQEvent::TriggerMask i) const =0
Return the trigger bit (fired or not) of the selected trigger channel.
virtual int get_spill_id() const =0
Return the spill ID.
virtual int get_event_id() const =0
Return the event ID, which is unique per run.
virtual ConstIter end() const =0
virtual ConstIter begin() const =0
std::vector< SQHit * >::iterator Iter
Definition: SQHitVector.h:38
An SQ interface class to hold one detector hit.
Definition: SQHit.h:20
virtual short get_detector_id() const
Return the detector ID of this hit.
Definition: SQHit.h:42
An SQ interface class to hold the run-level info.
Definition: SQRun.h:18
virtual ConstIter begin() const =0
virtual ConstIter end() const =0
virtual TVector3 get_pos_st3() const
Return the track position at Station 3.
Definition: SRecEvent.h:72
virtual TLorentzVector get_mom_st3() const
Return the track momentum at Station 3.
Definition: SRecEvent.h:81
virtual int get_charge() const
Return the charge, i.e. +1 or -1.
Definition: SRecEvent.h:60
Int_t getTriggerRoad()
Definition: SRecEvent.h:222
virtual double get_chisq_target() const
Definition: SRecEvent.h:85
virtual TVector3 get_pos_st1() const
Return the track position at Station 1.
Definition: SRecEvent.h:69
virtual double get_chisq_upstream() const
Definition: SRecEvent.h:88
virtual TLorentzVector get_mom_st1() const
Return the track momentum at Station 1.
Definition: SRecEvent.h:78
virtual int get_num_hits() const
Return the number of hits associated to this track.
Definition: SRecEvent.h:63
virtual double get_chisq_dump() const
Definition: SRecEvent.h:86
virtual TLorentzVector get_mom_vtx() const
Return the track momentum at vertex.
Definition: SRecEvent.h:75
virtual double get_chisq() const
Definition: SRecEvent.h:84
virtual TVector3 get_pos_vtx() const
Return the track position at vertex.
Definition: SRecEvent.h:66
short fpga_bits
Definition: TreeData.h:11
int road
Definition: TreeData.h:62
double chisq_target
Definition: TreeData.h:65
TVector3 pos_st3
Definition: TreeData.h:72
TLorentzVector mom_st3
Definition: TreeData.h:73
TVector3 pos_st1
Definition: TreeData.h:70
TVector3 pos_vtx
Definition: TreeData.h:24
double chisq
Definition: TreeData.h:64
int charge
Definition: TreeData.h:23
short n_hits
Definition: TreeData.h:63
double chisq_dump
Definition: TreeData.h:66
TLorentzVector mom_vtx
Definition: TreeData.h:25
TLorentzVector mom_st1
Definition: TreeData.h:71
double chisq_upstream
Definition: TreeData.h:67