Class Reference for E1039 Core & Analysis Software
AnaEffHodo.cc
Go to the documentation of this file.
1 #include <iomanip>
2 #include <TFile.h>
3 #include <TTree.h>
4 #include <TH1D.h>
5 #include <THStack.h>
6 #include <TEfficiency.h>
7 #include <TCanvas.h>
8 #include <interface_main/SQRun.h>
12 #include <phool/PHNodeIterator.h>
13 #include <phool/PHIODataNode.h>
14 #include <phool/getClass.h>
15 #include <geom_svc/GeomSvc.h>
16 #include <UtilAna/UtilSQHit.h>
17 #include "AnaEffHodo.h"
18 using namespace std;
19 
21 {
22  ;
23 }
24 
26 {
27  n_evt_all = n_evt_trig = n_evt_nhit = 0;
28 
29  ostringstream oss;
30  ofs.open("output.txt");
31  f_out = new TFile("output.root", "RECREATE");
32 
33  h1_eff_all = new TH1D("h1_eff_all", "", 1, 0, 1);
34  h1_eff_ok = new TH1D("h1_eff_ok" , "", 1, 0, 1);
35 
36  h1_nhit = new TH1D("h1_nhit", ";N of hits/plane/event;Hit count", 10, -0.5, 9.5);
37  h1_ele = new TH1D("h1_ele", ";Element ID;Hit count", 23, 0.5, 23.5);
38 
39  const double DT = 20/9.0; // 4/9 ns per single count of Taiwan TDC
40  const int NT = 1000;
41  const double T0 = 0.5*DT;
42  const double T1 = (NT+0.5)*DT;
43  h1_time = new TH1D("h1_time", ";tdcTime;Hit count", NT, T0, T1);
44 
46 }
47 
49 {
51 }
52 
54 {
55  SQEvent* event = findNode::getClass<SQEvent >(topNode, "SQEvent");
56  SQHitVector* hit_vec = findNode::getClass<SQHitVector>(topNode, "SQHitVector");
57  if (!event || !hit_vec) return Fun4AllReturnCodes::ABORTEVENT;
58  //int run_id = event->get_run_id();
59  //int spill_id = event->get_spill_id();
60  //int event_id = event->get_event_id();
61  n_evt_all++;
62 
66  if (! event->get_trigger(SQEvent::NIM1)) { // depends on the setting of run(s) you analyze
68  }
69  n_evt_trig++;
70 
71  shared_ptr<SQHitVector> hv_h1t(UtilSQHit::FindFirstHits(hit_vec, "H1T"));
72  shared_ptr<SQHitVector> hv_h1b(UtilSQHit::FindFirstHits(hit_vec, "H1B"));
73  shared_ptr<SQHitVector> hv_h2t(UtilSQHit::FindFirstHits(hit_vec, "H2T"));
74  shared_ptr<SQHitVector> hv_h2b(UtilSQHit::FindFirstHits(hit_vec, "H2B"));
75  shared_ptr<SQHitVector> hv_h3t(UtilSQHit::FindFirstHits(hit_vec, "H3T"));
76  shared_ptr<SQHitVector> hv_h3b(UtilSQHit::FindFirstHits(hit_vec, "H3B"));
77  shared_ptr<SQHitVector> hv_h4t(UtilSQHit::FindFirstHits(hit_vec, "H4T"));
78  shared_ptr<SQHitVector> hv_h4b(UtilSQHit::FindFirstHits(hit_vec, "H4B"));
79  if (hv_h1t->size() + hv_h1b->size() != 1 ||
80  hv_h2t->size() + hv_h2b->size() != 1 ||
81  hv_h3t->size() + hv_h3b->size() != 1 ||
82  hv_h4t->size() + hv_h4b->size() != 1
84 
85  n_evt_nhit++;
86 
90  shared_ptr<SQHitVector> hv_h2l(UtilSQHit::FindFirstHits(hit_vec, "H2L"));
91  shared_ptr<SQHitVector> hv_h2r(UtilSQHit::FindFirstHits(hit_vec, "H2R"));
92 
93  int nhit = hv_h2l->size() + hv_h2r->size();
94  h1_nhit->Fill(nhit);
95  for (SQHitVector::ConstIter it = hv_h2l->begin(); it != hv_h2l->end(); it++) {
96  h1_ele ->Fill((*it)->get_element_id());
97  h1_time->Fill((*it)->get_tdc_time ());
98  }
99  for (SQHitVector::ConstIter it = hv_h2r->begin(); it != hv_h2r->end(); it++) {
100  h1_ele ->Fill((*it)->get_element_id());
101  h1_time->Fill((*it)->get_tdc_time ());
102  }
103 
104  bool is_eff = nhit > 0; // very simple judgment for now
105  if (is_eff) h1_eff_ok->Fill(0);
106  h1_eff_all->Fill(0);
107 
109 }
110 
112 {
113  ostringstream oss;
114 
115  ofs << "N of events:\n"
116  << " All = " << n_evt_all << "\n"
117  << " Trigger OK = " << n_evt_trig << "\n"
118  << " n_hit OK = " << n_evt_nhit << "\n"
119  << endl;
120 
121  TCanvas* c1 = new TCanvas("c1", "");
122  c1->SetGrid();
123 
124  h1_eff_all->SetLineColor(kRed);
125  h1_eff_all->SetLineWidth(2);
126  THStack hs("hs", ";;N of hits");
127  hs.Add(h1_eff_all);
128  hs.Add(h1_eff_ok );
129  hs.Draw("nostack");
130  c1->SaveAs("h1_eff.png");
131 
132  TEfficiency* teff_pl = new TEfficiency(*h1_eff_ok, *h1_eff_all);
133  teff_pl->SetMarkerStyle(21);
134  teff_pl->SetMarkerColor(kRed);
135  teff_pl->SetLineColor (kRed);
136  teff_pl->Draw("APE1");
137  teff_pl->SetTitle(";Plane;Plane efficiency");
138  c1->SaveAs("teff_pl.png");
139 
140  h1_nhit->Draw();
141  c1->SaveAs("h1_nhit.png");
142 
143  h1_ele->Draw();
144  c1->SaveAs("h1_ele.png");
145 
146  h1_time->Draw();
147  c1->SaveAs("h1_time.png");
148 
149  delete c1;
150 
151  f_out->cd();
152  f_out->Write();
153  f_out->Close();
154  ofs.close();
155 
157 }
int InitRun(PHCompositeNode *topNode)
Definition: AnaEffHodo.cc:48
int End(PHCompositeNode *topNode)
Called at the end of all processing.
Definition: AnaEffHodo.cc:111
int Init(PHCompositeNode *topNode)
Definition: AnaEffHodo.cc:25
int process_event(PHCompositeNode *topNode)
Definition: AnaEffHodo.cc:53
An SQ interface class to hold one event header.
Definition: SQEvent.h:17
@ NIM1
Definition: SQEvent.h:22
An SQ interface class to hold a list of SQHit objects.
Definition: SQHitVector.h:32
std::vector< SQHit * >::const_iterator ConstIter
Definition: SQHitVector.h:37
SQHitVector * FindFirstHits(const SQHitVector *vec_in, const std::string det_name, const bool in_time=false)
Extract a set of first hits that are of the given detector (det_name), where "first" means the earlie...
Definition: UtilSQHit.cc:40