Class Reference for E1039 Core & Analysis Software
AnaTriggerHit.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 <TH1D.h>
6 #include <TCanvas.h>
7 #include <interface_main/SQRun.h>
11 #include <phool/PHNodeIterator.h>
12 #include <phool/PHIODataNode.h>
13 #include <phool/getClass.h>
14 #include <geom_svc/GeomSvc.h>
15 #include <UtilAna/UtilSQHit.h>
16 #include "AnaTriggerHit.h"
17 using namespace std;
18 
20 const vector<string> AnaTriggerHit::m_list_det_name = {
21  "H1T", "H1B", "H2T", "H2B", "H3T", "H3B", "H4T", "H4B"
22 };
23 
24 AnaTriggerHit::AnaTriggerHit(const std::string& name)
25  : SubsysReco(name)
26  , m_intime(false)
27  , m_level(1)
28  , m_file(0)
29  , m_tree(0)
30 {
31  ;
32 }
33 
35 {
37 }
38 
40 {
44  m_evt = findNode::getClass<SQEvent >(topNode, "SQEvent");
45  m_hit_vec = findNode::getClass<SQHitVector>(topNode, "SQTriggerHitVector");
46  if (!m_evt || !m_hit_vec) return Fun4AllReturnCodes::ABORTEVENT;
47 
51  gSystem->mkdir("result", true);
52 
53  m_file = new TFile("result/output.root", "RECREATE");
54  m_tree = new TTree("tree", "Created by AnaTriggerHit");
55  m_tree->Branch("det_name", &b_det_name, "det_name/C");
56  m_tree->Branch("det_id" , &b_det_id , "det_id/I");
57  m_tree->Branch("ele_id" , &b_ele_id , "ele_id/I");
58  m_tree->Branch("time" , &b_time , "time/D");
59 
60  ostringstream oss;
61  GeomSvc* geom = GeomSvc::instance();
62  m_list_det_id.clear();
63  for (unsigned int i_det = 0; i_det < m_list_det_name.size(); i_det++) {
64  string name = m_list_det_name[i_det];
65  int id = geom->getDetectorID(name);
66  if (id <= 0) {
67  cerr << "!ERROR! AnaTriggerHit::InitRun(): Invalid ID (" << id << "). Probably the detector name that you specified in 'm_list_det_name' (" << name << ") is not valid. Abort." << endl;
68  exit(1);
69  }
70  m_list_det_id.push_back(id);
71  int n_ele = geom->getPlaneNElements(id);
72  cout << " " << setw(6) << name << " = " << id << endl;
73 
74  oss.str("");
75  oss << "h1_ele_id_" << name;
76  m_h1_ele[i_det] = new TH1D(oss.str().c_str(), "", n_ele, 0.5, n_ele+0.5);
77  oss.str("");
78  oss << "V1495 Level-" << m_level << " " << name << ";Element ID;Hit count";
79  m_h1_ele[i_det]->SetTitle(oss.str().c_str());
80 
81  oss.str("");
82  oss << "h1_nhit_" << name;
83  m_h1_nhit[i_det] = new TH1D(oss.str().c_str(), "", 20, -0.5, 19.5);
84  oss.str("");
85  oss << "V1495 Level-" << m_level << " " << name << ";N of hits/plane/event;Hit count";
86  m_h1_nhit[i_det]->SetTitle(oss.str().c_str());
87 
88  oss.str("");
89  oss << "h1_time_" << name;
90  m_h1_time[i_det] = new TH1D(oss.str().c_str(), "", 2000, -0.5, 1999.5);
91  oss.str("");
92  oss << "V1495 Level-" << m_level << " " << name << ";TDC time (ns);Hit count";
93  m_h1_time[i_det]->SetTitle(oss.str().c_str());
94  }
95 
97 }
98 
100 {
101  //int spill_id = m_evt->get_spill_id();
102  //int event_id = m_evt->get_event_id();
103 
107 // if (! m_evt->get_trigger(SQEvent::NIM2)) {
108 // return Fun4AllReturnCodes::EVENT_OK;
109 // }
110 
114  for (unsigned int i_det = 0; i_det < m_list_det_name.size(); i_det++) {
115  strncpy(b_det_name, m_list_det_name[i_det].c_str(), sizeof(b_det_name));
116  b_det_id = m_list_det_id[i_det];
117  shared_ptr<SQHitVector> hv(UtilSQHit::FindHits(m_hit_vec, b_det_id));
118  for (SQHitVector::ConstIter it = hv->begin(); it != hv->end(); it++) {
119  if ((*it)->get_level() != m_level) continue;
120  if (m_intime && ! (*it)->is_in_time()) continue;
121 
122  b_ele_id = (*it)->get_element_id();
123  b_time = (*it)->get_tdc_time ();
124  m_tree->Fill();
125 
126  m_h1_ele [i_det]->Fill(b_ele_id);
127  m_h1_time[i_det]->Fill(b_time);
128  }
129  m_h1_nhit[i_det]->Fill(hv->size());
130  }
131 
133 }
134 
136 {
137  ostringstream oss;
138  TCanvas* c1 = new TCanvas("c1", "");
139  c1->SetGrid();
140  for (unsigned int i_det = 0; i_det < m_list_det_id.size(); i_det++) {
141  m_h1_ele[i_det]->Draw();
142  oss.str("");
143  oss << "result/" << m_h1_ele[i_det]->GetName() << ".png";
144  c1->SaveAs(oss.str().c_str());
145 
146  m_h1_nhit[i_det]->Draw();
147  oss.str("");
148  oss << "result/" << m_h1_nhit[i_det]->GetName() << ".png";
149  c1->SaveAs(oss.str().c_str());
150 
151  m_h1_time[i_det]->Draw();
152  oss.str("");
153  oss << "result/" << m_h1_time[i_det]->GetName() << ".png";
154  c1->SaveAs(oss.str().c_str());
155  }
156  delete c1;
157 
158  m_file->cd();
159  m_file->Write();
160  m_file->Close();
161 
163 }
int End(PHCompositeNode *topNode)
Called at the end of all processing.
int InitRun(PHCompositeNode *topNode)
int Init(PHCompositeNode *topNode)
AnaTriggerHit(const std::string &name="AnaTriggerHit")
int process_event(PHCompositeNode *topNode)
User interface class about the geometry of detector planes.
Definition: GeomSvc.h:164
int getDetectorID(const std::string &detectorName) const
Get the plane position.
Definition: GeomSvc.h:219
static GeomSvc * instance()
singlton instance
Definition: GeomSvc.cxx:212
int getPlaneNElements(int detectorID) const
Definition: GeomSvc.h:260
std::vector< SQHit * >::const_iterator ConstIter
Definition: SQHitVector.h:37
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