Class Reference for E1039 Core & Analysis Software
ExtractTdcDist.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 <TSystem.h>
6 #include <TCanvas.h>
9 #include <phool/PHNodeIterator.h>
10 #include <phool/PHIODataNode.h>
11 #include <phool/getClass.h>
12 #include <geom_svc/GeomSvc.h>
13 #include <UtilAna/UtilSQHit.h>
14 #include "ExtractTdcDist.h"
15 using namespace std;
16 
17 ExtractTdcDist::ExtractTdcDist(const std::string& name)
18  : SubsysReco("ExtractTdcDist"+name)
19  , m_name (name)
20  , m_n_evt_all (0)
21  , m_n_evt_trig(0)
22  , m_n_evt_ana (0)
23 {
24  if (name == "d0" ) m_type = D0 ;
25  else if (name == "d1" ) m_type = D1 ;
26  else if (name == "d2" ) m_type = D2 ;
27  else if (name == "d3p") m_type = D3p;
28  else if (name == "d3m") m_type = D3m;
29  else {
30  cout << "ExtractTdcDist: Bad name (" << name << "). Abort." << endl;
31  exit(1);
32  }
33 }
34 
42 {
43  double t0_lo;
44  double t0_hi;
45  switch (m_type) {
46  case D0 : m_pl0 = 1; m_trig_mask = SQEvent::NIM2; t0_lo = 1000; t0_hi = 1400; break;
47  case D1 : m_pl0 = 7; m_trig_mask = SQEvent::NIM2; t0_lo = 850; t0_hi = 1250; break;
48  case D2 : m_pl0 = 13; m_trig_mask = SQEvent::NIM4; t0_lo = 850; t0_hi = 1250; break;
49  case D3p: m_pl0 = 19; m_trig_mask = SQEvent::NIM4; t0_lo = 740; t0_hi = 1140; break;
50  case D3m: m_pl0 = 25; m_trig_mask = SQEvent::NIM4; t0_lo = 850; t0_hi = 1250; break;
51  }
52 
53  m_dir_out = "tdc_dist_" + m_name;
54  gSystem->mkdir(m_dir_out.c_str(), true);
55 
56  GeomSvc* geom = GeomSvc::instance();
57  m_ofs.open( (m_dir_out+"/output.txt").c_str() );
58  m_file = new TFile( (m_dir_out+"/output.root").c_str(), "RECREATE");
59 
60  ostringstream oss;
61  oss << setfill('0');
62 
63  for (int i_pl = 0; i_pl < N_PL; i_pl++) {
64  int det_id = m_pl0 + i_pl;
65  string name = geom->getDetectorName(det_id);
66  int n_ele = geom->getPlaneNElements(det_id);
67  m_list_pl_id .push_back(det_id);
68  m_list_pl_name .push_back(name);
69  m_list_pl_n_ele.push_back(n_ele);
70  cout << " " << setw(6) << name << " " << det_id << " " << n_ele << endl;
71 
72  oss.str("");
73  oss << "h1_nhit_" << setw(2) << det_id;
74  m_h1_nhit[i_pl] = new TH1D(oss.str().c_str(), "", 50, -0.5, 49.5);
75  oss.str("");
76  oss << name << ";N of hits/plane/event;Hit count";
77  m_h1_nhit[i_pl]->SetTitle(oss.str().c_str());
78 
79  oss.str("");
80  oss << "h1_ele_" << setw(2) << det_id;
81  m_h1_ele[i_pl] = new TH1D(oss.str().c_str(), "", n_ele, 0.5, n_ele+0.5);
82  oss.str("");
83  oss << name << ";Element ID;Hit count";
84  m_h1_ele[i_pl]->SetTitle(oss.str().c_str());
85 
86  const double DT = 20/9.0; // 4/9 ns per single count of Taiwan TDC
87  int i_t_lo = (int)(t0_lo / DT);
88  int i_t_hi = (int)(t0_hi / DT) + 1;
89  double t_lo = i_t_lo * DT;
90  double t_hi = i_t_hi * DT;
91  int t_n = i_t_hi - i_t_lo;
92  oss.str("");
93  oss << "h1_time_" << setw(2) << det_id;
94  m_h1_time[i_pl] = new TH1D(oss.str().c_str(), "", t_n, t_lo, t_hi);
95  oss.str("");
96  oss << name << ";tdcTime;Hit count";
97  m_h1_time[i_pl]->SetTitle(oss.str().c_str());
98 
99  t_lo = 0 * DT;
100  t_hi = 1000 * DT;
101  t_n = 1000;
102  oss.str("");
103  oss << "h1_time_wide_" << setw(2) << det_id;
104  m_h1_time_wide[i_pl] = new TH1D(oss.str().c_str(), "", t_n, t_lo, t_hi);
105  oss.str("");
106  oss << name << ";tdcTime;Hit count";
107  m_h1_time_wide[i_pl]->SetTitle(oss.str().c_str());
108  }
109 
111 }
112 
114 {
116 }
117 
119 {
120  SQEvent* event = findNode::getClass<SQEvent >(topNode, "SQEvent");
121  SQHitVector* hit_vec = findNode::getClass<SQHitVector>(topNode, "SQHitVector");
122  if (!event || !hit_vec) return Fun4AllReturnCodes::ABORTEVENT;
123  m_n_evt_all++;
124 
125  if (! event->get_trigger(m_trig_mask)) return Fun4AllReturnCodes::EVENT_OK;
126  m_n_evt_trig++;
127 
128  bool good_event = true;
129  int n_hit_tot = 0;
130  shared_ptr<SQHitVector> list_hv[N_PL];
131  for (int i_pl = 0; i_pl < N_PL; i_pl++) {
132  list_hv[i_pl] = shared_ptr<SQHitVector>(UtilSQHit::FindHits(hit_vec, m_list_pl_id[i_pl]));
133 
134  int n_hit = list_hv[i_pl]->size();
135  m_h1_nhit[i_pl]->Fill(n_hit);
136  if (n_hit > 1) {
137  good_event = false;
138  break;
139  }
140  n_hit_tot += n_hit;
141  }
142  if (!good_event || n_hit_tot < 5) return Fun4AllReturnCodes::EVENT_OK;
143 
144  m_n_evt_ana++;
145 
149  for (int i_pl = 0; i_pl < N_PL; i_pl++) {
150  for (auto it = list_hv[i_pl]->begin(); it != list_hv[i_pl]->end(); it++) {
151  int ele = (*it)->get_element_id();
152  double time = (*it)->get_tdc_time ();
153  m_h1_ele [i_pl]->Fill(ele );
154  m_h1_time [i_pl]->Fill(time);
155  m_h1_time_wide[i_pl]->Fill(time);
156  }
157  }
158 
160 }
161 
163 {
164  ostringstream oss;
165 
166  m_ofs << "N of events:\n"
167  << " All = " << m_n_evt_all << "\n"
168  << " Trigger OK = " << m_n_evt_trig << "\n"
169  << " Analyzed = " << m_n_evt_ana << "\n"
170  << endl;
171 
172  TCanvas* c1 = new TCanvas("c1", "");
173  c1->SetGrid();
174  for (int i_pl = 0; i_pl < N_PL; i_pl++) {
175  m_h1_nhit[i_pl]->Draw();
176  oss.str("");
177  oss << m_dir_out << "/" << m_h1_nhit[i_pl]->GetName() << ".png";
178  c1->SaveAs(oss.str().c_str());
179 
180  m_h1_ele[i_pl]->Draw();
181  oss.str("");
182  oss << m_dir_out << "/" << m_h1_ele[i_pl]->GetName() << ".png";
183  c1->SaveAs(oss.str().c_str());
184 
185  m_h1_time[i_pl]->Draw();
186  oss.str("");
187  oss << m_dir_out << "/" << m_h1_time[i_pl]->GetName() << ".png";
188  c1->SaveAs(oss.str().c_str());
189 
190  m_h1_time_wide[i_pl]->Draw();
191  oss.str("");
192  oss << m_dir_out << "/" << m_h1_time_wide[i_pl]->GetName() << ".png";
193  c1->SaveAs(oss.str().c_str());
194  }
195  delete c1;
196 
197  m_file->cd();
198  m_file->Write();
199  m_file->Close();
200  m_ofs.close();
201 
203 }
int process_event(PHCompositeNode *topNode)
static const int N_PL
ExtractTdcDist(const std::string &name="d0")
int End(PHCompositeNode *topNode)
Called at the end of all processing.
int Init(PHCompositeNode *topNode)
int InitRun(PHCompositeNode *topNode)
User interface class about the geometry of detector planes.
Definition: GeomSvc.h:164
static GeomSvc * instance()
singlton instance
Definition: GeomSvc.cxx:212
int getPlaneNElements(int detectorID) const
Definition: GeomSvc.h:260
std::string getDetectorName(const int &detectorID) const
Definition: GeomSvc.h:223
An SQ interface class to hold one event header.
Definition: SQEvent.h:17
@ NIM4
Definition: SQEvent.h:25
@ NIM2
Definition: SQEvent.h:23
An SQ interface class to hold a list of SQHit objects.
Definition: SQHitVector.h:32
SQHitVector * FindHits(const SQHitVector *vec_in, const std::string det_name, const bool in_time=false)
Extract a set of hits that are of the given detector (det_name).
Definition: UtilSQHit.cc:16