Class Reference for E1039 Core & Analysis Software
AnaChamPlane.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 <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 "AnaChamPlane.h"
17 using namespace std;
18 
19 AnaChamPlane::AnaChamPlane(const std::string& plane_name)
20  : SubsysReco("AnaChamPlane"+plane_name)
21  , m_plane_name(plane_name)
22  , m_plane_id (0)
23  , m_n_ele (0)
24  , m_n_evt_all (0)
25  , m_n_hit_all (0)
26  , m_evt (0)
27  , m_hit_vec (0)
28  , m_file_out (0)
29  , m_h1_nhit (0)
30  , m_h1_ele (0)
31  , m_h1_pos (0)
32  , m_h1_time (0)
33  , m_h1_dist (0)
34 {
35  ;
36 }
37 
39 {
41 }
42 
44 {
45  m_evt = findNode::getClass<SQEvent >(topNode, "SQEvent");
46  m_hit_vec = findNode::getClass<SQHitVector>(topNode, "SQHitVector");
47  if (!m_evt || !m_hit_vec) return Fun4AllReturnCodes::ABORTEVENT;
48 
49  GeomSvc* geom = GeomSvc::instance();
50  m_plane_id = geom->getDetectorID(m_plane_name);
51  if (m_plane_id <= 0) {
52  cout << "!!ERROR!! " << Name() << " is not given a proper plane name.\n";
54  }
55  m_n_ele = geom->getPlaneNElements(m_plane_id);
56 
57  m_file_out = new TFile(TString::Format("output_%s.root", Name().c_str()).Data(), "RECREATE");
58 
59  m_h1_nhit = new TH1D(
60  TString::Format("h1_nhit_%d_%s", m_plane_id, m_plane_name.c_str()).Data(),
61  TString::Format("#%d: %s;N of hits/plane/event;Hit count", m_plane_id, m_plane_name.c_str()).Data(),
62  m_n_ele+1, -0.5, m_n_ele+0.5);
63 
64  m_h1_ele = new TH1D(
65  TString::Format("h1_ele_%d_%s", m_plane_id, m_plane_name.c_str()).Data(),
66  TString::Format("#%d: %s;Element ID;Hit count", m_plane_id, m_plane_name.c_str()).Data(),
67  m_n_ele, 0.5, m_n_ele+0.5);
68 
69  double plane_width = geom->getPlaneScaleX(m_plane_id);
70  m_h1_pos = new TH1D(
71  TString::Format("h1_pos_%d_%s", m_plane_id, m_plane_name.c_str()).Data(),
72  TString::Format("#%d: %s;Position (cm);Hit count", m_plane_id, m_plane_name.c_str()).Data(),
73  60, -0.6*plane_width, 0.6*plane_width);
74 
75  const double DT = 20/9.0; // 4/9 ns per single count of Taiwan TDC
76  const int NT = 450;
77  const double T0 = 225.5*DT;
78  const double T1 = 675.5*DT;
79  m_h1_time = new TH1D(
80  TString::Format("h1_time_%d_%s", m_plane_id, m_plane_name.c_str()).Data(),
81  TString::Format("#%d: %s;tdcTime;Hit count", m_plane_id, m_plane_name.c_str()).Data(),
82  NT, T0, T1);
83 
84  m_h1_dist = new TH1D(
85  TString::Format("h1_dist_%d_%s", m_plane_id, m_plane_name.c_str()).Data(),
86  TString::Format("#%d: %s;Drift dist;Hit count", m_plane_id, m_plane_name.c_str()).Data(),
87  120, 0.0, 0.6 * geom->getCellWidth(m_plane_id));
88 
90 }
91 
93 {
94  //int run_id = m_evt->get_run_id();
95  //int spill_id = m_evt->get_spill_id();
96  //int event_id = m_evt->get_event_id();
97  m_n_evt_all++;
98 
99  shared_ptr<SQHitVector> hv(UtilSQHit::FindHits(m_hit_vec, m_plane_id));
100  m_n_hit_all += hv->size();
101 
102  m_h1_nhit->Fill(hv->size());
103  for (SQHitVector::ConstIter it = hv->begin(); it != hv->end(); it++) {
104  int ele = (*it)->get_element_id();
105  double pos = (*it)->get_pos();
106  double time = (*it)->get_tdc_time();
107  double dist = (*it)->get_drift_distance();
108  m_h1_ele ->Fill(ele );
109  m_h1_pos ->Fill(pos );
110  m_h1_time->Fill(time);
111  m_h1_dist->Fill(dist);
112  }
113 
115 }
116 
118 {
119  std::ofstream ofs(TString::Format("output_%s.txt", Name().c_str()).Data());
120  ofs << "N of events = " << m_n_evt_all << "\n"
121  << "N of hits = " << m_n_hit_all << "\n";
122  ofs.close();
123 
124  m_file_out->cd();
125  m_file_out->Write();
126  m_file_out->Close();
127 
129 }
int InitRun(PHCompositeNode *topNode)
Definition: AnaChamPlane.cc:43
int process_event(PHCompositeNode *topNode)
Definition: AnaChamPlane.cc:92
AnaChamPlane(const std::string &plane_name)
Definition: AnaChamPlane.cc:19
int End(PHCompositeNode *topNode)
Called at the end of all processing.
int Init(PHCompositeNode *topNode)
Definition: AnaChamPlane.cc:38
virtual const std::string Name() const
Returns the name of this module.
Definition: Fun4AllBase.h:23
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
double getPlaneScaleX(int detectorID) const
Definition: GeomSvc.h:242
double getCellWidth(int detectorID) const
Definition: GeomSvc.h:238
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