Class Reference for E1039 Core & Analysis Software
AnaChamHit.cc
Go to the documentation of this file.
1 #include <iomanip>
2 #include <TSystem.h>
3 #include <TFile.h>
4 #include <TTree.h>
5 #include <TChain.h>
6 #include <TH1D.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 <UtilAna/UtilHist.h>
18 #include "AnaChamHit.h"
19 using namespace std;
20 
21 AnaChamHit::AnaChamHit(const std::string& name)
22  : SubsysReco (name)
23  , m_sq_evt (0)
24  , m_sq_hit_vec(0)
25  , m_file_name ("output.root")
26  , m_file (0)
27  , m_tree (0)
28 {
29  ;
30 }
31 
33 {
34  ;
35 }
36 
38 {
40 }
41 
43 {
44  m_sq_evt = findNode::getClass<SQEvent >(topNode, "SQEvent");
45  m_sq_hit_vec = findNode::getClass<SQHitVector>(topNode, "SQHitVector");
46  if (!m_sq_evt || !m_sq_hit_vec) return Fun4AllReturnCodes::ABORTEVENT;
47 
48  m_file = new TFile(m_file_name.c_str(), "RECREATE");
49  m_tree = new TTree("tree", "Created by AnaChamHit");
50  m_tree->Branch("event" , &m_evt);
51  m_tree->Branch("hit_list", &m_hit_list);
52 
54 }
55 
57 {
58  if (! m_sq_evt->get_trigger(SQEvent::MATRIX1)) {
60  }
61 
62  m_evt.run_id = m_sq_evt->get_run_id();
63  m_evt.spill_id = m_sq_evt->get_spill_id();
64  m_evt.event_id = m_sq_evt->get_event_id();
65  m_evt.fpga_bits = (m_sq_evt->get_trigger() >> SQEvent::MATRIX1) & 0x1f;
66  m_evt.nim_bits = (m_sq_evt->get_trigger() >> SQEvent::NIM1 ) & 0x1f;
67  m_evt.D1 = m_evt.D2 = m_evt.D3p = m_evt.D3m = 0;
68 
69  GeomSvc* geom = GeomSvc::instance();
70  m_hit_list.clear();
71  for (SQHitVector::Iter it = m_sq_hit_vec->begin(); it != m_sq_hit_vec->end(); it++) {
72  SQHit* hit = *it;
73  int det_id = hit->get_detector_id();
74  if (! geom->isChamber(det_id)) continue;
75 
76  if ( 0 < det_id && det_id <= 6) m_evt.D1++;
77  else if (12 < det_id && det_id <= 18) m_evt.D2++;
78  else if (18 < det_id && det_id <= 24) m_evt.D3p++;
79  else if (24 < det_id && det_id <= 30) m_evt.D3m++;
80 
81  HitData hd;
82  hd.det_id = det_id;
83  hd.ele_id = hit->get_element_id();
84  hd.tdc_time = hit->get_tdc_time();
85  hd.drift_dist = hit->get_drift_distance();
86  m_hit_list.push_back(hd);
87  }
88 
89  m_tree->Fill();
91 }
92 
94 {
95  m_file->cd();
96  m_file->Write();
97  m_file->Close();
99 }
100 
101 void AnaChamHit::AnalyzeTree(TChain* tree)
102 {
103  cout << "N of trees = " << tree->GetNtrees() << endl;
104 
105  TH1* h1_ele [31];
106  TH1* h1_time[31];
107  TH1* h1_dist[31];
108  GeomSvc* geom = GeomSvc::instance();
109  ostringstream oss;
110  for (int pl = 1; pl <= 30; pl++) {
111  if (6 < pl && pl <= 12) continue;
112  string name = geom->getDetectorName (pl);
113  int n_ele = geom->getPlaneNElements(pl);
114  double space = geom->getPlaneSpacing (pl);
115  oss.str("");
116  oss << "h1_ele_" << pl;
117  h1_ele[pl] = new TH1D(oss.str().c_str(), "", n_ele, 0.5, n_ele+0.5);
118  oss.str("");
119  oss << setw(2) << pl << ":" << name << ";Element ID;Hit count";
120  h1_ele[pl]->SetTitle(oss.str().c_str());
121 
122  const double DT = 12/9.0; // 4/9 ns per single count of Taiwan TDC
123  oss.str("");
124  oss << "h1_time_" << pl;
125  h1_time[pl] = new TH1D(oss.str().c_str(), "", 1000, 400.5*DT, 1400.5*DT);
126  oss.str("");
127  oss << setw(2) << pl << ":" << name << ";tdcTime (ns);Hit count";
128  h1_time[pl]->SetTitle(oss.str().c_str());
129 
130  oss.str("");
131  oss << "h1_dist_" << pl;
132  h1_dist[pl] = new TH1D(oss.str().c_str(), "", 110, 0.0, 0.55*space);
133  oss.str("");
134  oss << setw(2) << pl << ":" << name << ";Drift distance (cm);Hit count";
135  h1_dist[pl]->SetTitle(oss.str().c_str());
136  }
137 
138  EventData* evt = 0;
139  HitList* hit_list = 0;
140  tree->SetBranchAddress("event" , &evt);
141  tree->SetBranchAddress("hit_list", &hit_list);
142 
143  int n_ent = tree->GetEntries();
144  cout << "N of entries = " << n_ent << endl;
145  for (int i_ent = 0; i_ent < n_ent; i_ent++) {
146  if ((i_ent+1) % (n_ent/10) == 0) cout << " " << 10*(i_ent+1)/(n_ent/10) << "%" << flush;
147  tree->GetEntry(i_ent);
148  for (auto it = hit_list->begin(); it != hit_list->end(); it++) {
149  HitData* hd = &(*it);
150  short det_id = hd->det_id;
151  if (6 < det_id && det_id <= 12) continue;
152  short ele_id = hd->ele_id;
153  double time = hd->tdc_time;
154  double dist = hd->drift_dist;
155  h1_ele [det_id]->Fill(ele_id);
156  h1_time[det_id]->Fill(time);
157  h1_dist[det_id]->Fill(dist);
158  }
159  }
160 
161  gSystem->mkdir("result", true);
162  TCanvas* c1 = new TCanvas("c1", "");
163  c1->SetGrid();
164  //c1->SetLogy(true);
165 
166  oss << setfill('0');
167  for (int pl = 1; pl <= 30; pl++) {
168  if (6 < pl && pl <= 12) continue;
169  h1_ele[pl]->Draw();
170  oss.str("");
171  oss << "result/h1_ele_" << setw(2) << pl << ".png";
172  c1->SaveAs(oss.str().c_str());
173 
174  h1_time[pl]->Draw();
175  UtilHist::AutoSetRange(h1_time[pl]);
176  oss.str("");
177  oss << "result/h1_time_" << setw(2) << pl << ".png";
178  c1->SaveAs(oss.str().c_str());
179 
180  h1_dist[pl]->Draw();
181  oss.str("");
182  oss << "result/h1_dist_" << setw(2) << pl << ".png";
183  c1->SaveAs(oss.str().c_str());
184  }
185 
186  delete c1;
187 }
188 
std::vector< HitData > HitList
Definition: TreeData.h:34
int InitRun(PHCompositeNode *topNode)
Definition: AnaChamHit.cc:42
int End(PHCompositeNode *topNode)
Called at the end of all processing.
Definition: AnaChamHit.cc:93
int Init(PHCompositeNode *topNode)
Definition: AnaChamHit.cc:37
int process_event(PHCompositeNode *topNode)
Definition: AnaChamHit.cc:56
AnaChamHit(const std::string &name="AnaChamHit")
Definition: AnaChamHit.cc:21
virtual ~AnaChamHit()
Definition: AnaChamHit.cc:32
static void AnalyzeTree(TChain *tree)
Definition: AnaChamHit.cc:101
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
bool isChamber(const int detectorID) const
Return "true" for chamber planes.
Definition: GeomSvc.cxx:772
double getPlaneSpacing(int detectorID) const
Definition: GeomSvc.h:236
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 float get_drift_distance() const
Return the drift distance of this hit. Probably the value is not properly set at present....
Definition: SQHit.h:57
virtual short get_element_id() const
Return the element ID of this hit.
Definition: SQHit.h:45
virtual float get_tdc_time() const
Return the TDC time (nsec) of this hit.
Definition: SQHit.h:54
virtual short get_detector_id() const
Return the detector ID of this hit.
Definition: SQHit.h:42
void AutoSetRange(TH1 *h1, const int margin_lo=5, const int margin_hi=5)
Adjust the axis range via "h1->GetXaxis()->SetRange(bin_lo, bin_hi)" to zoom up non-empty bins.
Definition: UtilHist.cc:17
int spill_id
Definition: TreeData.h:9
int event_id
Definition: TreeData.h:10
short nim_bits
Definition: TreeData.h:12
short D1
Definition: TreeData.h:13
short D2
Definition: TreeData.h:14
short fpga_bits
Definition: TreeData.h:11
short run_id
Definition: TreeData.h:8
short D3p
Definition: TreeData.h:15
short D3m
Definition: TreeData.h:16
short det_id
Definition: TreeData.h:24
short ele_id
Definition: TreeData.h:25
double tdc_time
Definition: TreeData.h:26
double drift_dist
Definition: TreeData.h:27