Class Reference for E1039 Core & Analysis Software
AnaHodoHit.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 "AnaHodoHit.h"
17 using namespace std;
18 
20 const vector<string> AnaHodoHit::m_list_det_name = { "H1T", "H1B", "H2T", "H2B" };
21 
22 AnaHodoHit::AnaHodoHit(const std::string& name)
23  : SubsysReco(name)
24  , m_file(0)
25  , m_tree(0)
26 {
27  ;
28 }
29 
31 {
33 }
34 
36 {
40  m_evt = findNode::getClass<SQEvent >(topNode, "SQEvent");
41  m_hit_vec = findNode::getClass<SQHitVector>(topNode, "SQTriggerHitVector");
42  if (!m_evt || !m_hit_vec) return Fun4AllReturnCodes::ABORTEVENT;
43 
47  gSystem->mkdir("result", true);
48 
49  m_file = new TFile("result/output.root", "RECREATE");
50  m_tree = new TTree("tree", "Created by AnaHodoHit");
51  m_tree->Branch("det_name", &b_det_name, "det_name/C");
52  m_tree->Branch("det" , &b_det , "det/I");
53  m_tree->Branch("ele" , &b_ele , "ele/I");
54  m_tree->Branch("time" , &b_time , "time/D");
55 
56  ostringstream oss;
57  GeomSvc* geom = GeomSvc::instance();
58  for (unsigned int i_det = 0; i_det < m_list_det_name.size(); i_det++) {
59  string name = m_list_det_name[i_det];
60  int id = geom->getDetectorID(name);
61  if (id <= 0) {
62  cerr << "!ERROR! AnaHodoHit::InitRun(): Invalid ID (" << id << "). Probably the detector name that you specified in 'list_det_name' (" << name << ") is not valid. Abort." << endl;
63  exit(1);
64  }
65  m_list_det_id.push_back(id);
66  int n_ele = geom->getPlaneNElements(id);
67  cout << " " << setw(6) << name << " = " << id << endl;
68 
69  oss.str("");
70  oss << "h1_ele_" << name;
71  m_h1_ele[i_det] = new TH1D(oss.str().c_str(), "", n_ele, 0.5, n_ele+0.5);
72  oss.str("");
73  oss << name << ";Element ID;Hit count";
74  m_h1_ele[i_det]->SetTitle(oss.str().c_str());
75 
76  oss.str("");
77  oss << "h1_nhit_" << name;
78  m_h1_nhit[i_det] = new TH1D(oss.str().c_str(), "", 10, -0.5, 9.5);
79  oss.str("");
80  oss << name << ";N of hits/plane/event;Hit count";
81  m_h1_nhit[i_det]->SetTitle(oss.str().c_str());
82  }
83 
85 }
86 
88 {
89  //int spill_id = m_evt->get_spill_id();
90  //int event_id = m_evt->get_event_id();
91 
95  if (! m_evt->get_trigger(SQEvent::NIM2)) {
97  }
98 
100  static int did_h3t = 0;
101  static int did_h3b = 0;
102  if (did_h3t == 0) { // Resolve IDs only once.
103  GeomSvc* geom = GeomSvc::instance();
104  did_h3t = geom->getDetectorID("H3T");
105  did_h3b = geom->getDetectorID("H3B");
106  cout << "H3T = " << did_h3t << ", H3B = " << did_h3b << endl;
107  }
108  shared_ptr<SQHitVector> hv_h3t(UtilSQHit::FindHits(m_hit_vec, did_h3t));
109  shared_ptr<SQHitVector> hv_h3b(UtilSQHit::FindHits(m_hit_vec, did_h3b));
110  if (hv_h3t->size() + hv_h3b->size() != 1) return Fun4AllReturnCodes::EVENT_OK;
111 
115  for (unsigned int i_det = 0; i_det < m_list_det_name.size(); i_det++) {
116  strncpy(b_det_name, m_list_det_name[i_det].c_str(), sizeof(b_det_name));
117  b_det = m_list_det_id[i_det];
118  shared_ptr<SQHitVector> hv(UtilSQHit::FindHits(m_hit_vec, b_det));
119  for (SQHitVector::ConstIter it = hv->begin(); it != hv->end(); it++) {
120  b_ele = (*it)->get_element_id();
121  b_time = (*it)->get_tdc_time ();
122  m_tree->Fill();
123 
124  m_h1_ele[i_det]->Fill(b_ele);
125  }
126  m_h1_nhit[i_det]->Fill(hv->size());
127  }
128 
130 }
131 
133 {
134  ostringstream oss;
135  TCanvas* c1 = new TCanvas("c1", "");
136  c1->SetGrid();
137  for (unsigned int i_det = 0; i_det < m_list_det_id.size(); i_det++) {
138  m_h1_ele[i_det]->Draw();
139  oss.str("");
140  oss << "result/" << m_h1_ele[i_det]->GetName() << ".png";
141  c1->SaveAs(oss.str().c_str());
142 
143  m_h1_nhit[i_det]->Draw();
144  oss.str("");
145  oss << "result/" << m_h1_nhit[i_det]->GetName() << ".png";
146  c1->SaveAs(oss.str().c_str());
147  }
148  delete c1;
149 
150  m_file->cd();
151  m_file->Write();
152  m_file->Close();
153 
155 }
int process_event(PHCompositeNode *topNode)
Definition: AnaHodoHit.cc:87
int Init(PHCompositeNode *topNode)
Definition: AnaHodoHit.cc:30
AnaHodoHit(const std::string &name="AnaHodoHit")
Definition: AnaHodoHit.cc:22
int InitRun(PHCompositeNode *topNode)
Definition: AnaHodoHit.cc:35
int End(PHCompositeNode *topNode)
Called at the end of all processing.
Definition: AnaHodoHit.cc:132
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
@ NIM2
Definition: SQEvent.h:23
virtual bool get_trigger(const SQEvent::TriggerMask i) const =0
Return the trigger bit (fired or not) of the selected trigger channel.
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