Class Reference for E1039 Core & Analysis Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
OnlMonH4.cc
Go to the documentation of this file.
1 #include <iomanip>
3 #include <TH1D.h>
4 #include <TH2D.h>
5 #include <interface_main/SQRun.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/UtilHist.h>
14 #include <UtilAna/UtilSQHit.h>
15 #include "OnlMonServer.h"
16 #include "OnlMonH4.h"
17 using namespace std;
18 
19 OnlMonH4::OnlMonH4(const HodoType_t type) : m_type(type)
20 {
21  NumCanvases(2);
22  switch (m_type) {
23  case H4T:
24  Name("OnlMonH4T");
25  Title("H4T PMT");
26  m_det_name = "H4T";
27  m_list_det_name.push_back("H4Tu");
28  m_list_det_name.push_back("H4Td");
29  break;
30  case H4B:
31  Name("OnlMonH4B");
32  Title("H4B PMT");
33  m_det_name = "H4B";
34  m_list_det_name.push_back("H4Bu");
35  m_list_det_name.push_back("H4Bd");
36  break;
37  case H4Y1L:
38  Name("OnlMonH4Y1L");
39  Title("H4Y1L PMT");
40  m_det_name = "H4Y1L";
41  m_list_det_name.push_back("H4Y1Ll");
42  m_list_det_name.push_back("H4Y1Lr");
43  break;
44  case H4Y1R:
45  Name("OnlMonH4Y1R");
46  Title("H4Y1R PMT");
47  m_det_name = "H4Y1R";
48  m_list_det_name.push_back("H4Y1Rl");
49  m_list_det_name.push_back("H4Y1Rr");
50  break;
51  case H4Y2L:
52  Name("OnlMonH4Y2L");
53  Title("H4Y2L PMT");
54  m_det_name = "H4Y2L";
55  m_list_det_name.push_back("H4Y2Ll");
56  m_list_det_name.push_back("H4Y2Lr");
57  break;
58  case H4Y2R:
59  Name("OnlMonH4Y2R");
60  Title("H4Y2R PMT");
61  m_det_name = "H4Y2R";
62  m_list_det_name.push_back("H4Y2Rl");
63  m_list_det_name.push_back("H4Y2Rr");
64  break;
65  }
66  // m_list_det_name.size() must be N_PL (=2)
67 }
68 
70 {
72 }
73 
75 {
76  const double DT = 40/9.0; // 4/9 ns per single count of Taiwan TDC
77  const int NT = 350;
78  const double T0 = 100.5*DT;
79  const double T1 = 450.5*DT;
80 
81  GeomSvc* geom = GeomSvc::instance();
82  int n_ele = geom->getPlaneNElements( geom->getDetectorID(m_det_name) );
83 
84  string det_name[N_PL];
85  ostringstream oss;
86  m_list_det.clear();
87  for (int i_pl = 0; i_pl < N_PL; i_pl++) {
88  det_name[i_pl] = m_list_det_name[i_pl];
89  short det_id = geom->getDetectorID(det_name[i_pl]);
90  m_list_det.push_back(det_id);
91 
92  oss.str("");
93  oss << "h1_ele_" << i_pl;
94  h1_ele[i_pl] = new TH1D(oss.str().c_str(), "", n_ele, 0.5, n_ele+0.5);
95  oss.str("");
96  oss << det_name[i_pl] << ";Element ID;Hit count";
97  h1_ele[i_pl]->SetTitle(oss.str().c_str());
98 
99 
100  oss.str("");
101  oss << "h1_time_" << i_pl;
102  h1_time[i_pl] = new TH1D(oss.str().c_str(), "", NT, T0, T1);
103  oss.str("");
104  oss << det_name[i_pl] << ";tdcTime;Hit count";
105  h1_time[i_pl]->SetTitle(oss.str().c_str());
106 
107  RegisterHist(h1_ele [i_pl]);
108  RegisterHist(h1_time[i_pl]);
109  }
110 
111  h2_ele = new TH2D("h2_ele", "", n_ele, 0.5, n_ele+0.5, n_ele, 0.5, n_ele+0.5);
112  oss.str("");
113  oss << ";Element ID of " << det_name[0] << ";Element ID of " << det_name[1] << ";Hit count";
114  h2_ele->SetTitle(oss.str().c_str());
115 
116  h2_time = new TH2D("h2_time", "", NT, T0, T1, NT, T0, T1);
117  oss.str("");
118  oss << ";tdcTime of " << det_name[0] << ";tdcTime ID of " << det_name[1] << ";Hit count";
119  h2_time->SetTitle(oss.str().c_str());
120 
121  RegisterHist(h2_ele );
122  RegisterHist(h2_time);
123 
125 }
126 
128 {
129  SQEvent* event_header = findNode::getClass<SQEvent >(topNode, "SQEvent");
130  SQHitVector* hit_vec = findNode::getClass<SQHitVector>(topNode, "SQHitVector");
131  if (!event_header || !hit_vec) return Fun4AllReturnCodes::ABORTEVENT;
132 
133  SQHitVector* hv[N_PL];
134  for (int i_pl = 0; i_pl < N_PL; i_pl++) {
135  hv[i_pl] = UtilSQHit::FindHits(hit_vec, m_list_det[i_pl]);
136  for (SQHitVector::ConstIter it = hv[i_pl]->begin(); it != hv[i_pl]->end(); it++) {
137  h1_ele [i_pl]->Fill((*it)->get_element_id());
138  h1_time[i_pl]->Fill((*it)->get_tdc_time ());
139  }
140  }
141  if (hv[0]->size() == 1 && hv[1]->size() == 1) {
142  h2_ele ->Fill( hv[0]->at(0)->get_element_id(), hv[1]->at(0)->get_element_id() );
143  h2_time->Fill( hv[0]->at(0)->get_tdc_time (), hv[1]->at(0)->get_tdc_time () );
144  }
145  for (int i_pl = 0; i_pl < N_PL; i_pl++) delete hv[i_pl];
146 
148 }
149 
151 {
153 }
154 
156 {
157  ostringstream oss;
158  for (int pl = 0; pl < N_PL; pl++) {
159  oss.str("");
160  oss << "h1_ele_" << pl;
161  h1_ele[pl] = FindMonHist(oss.str().c_str());
162  if (! h1_ele[pl]) return 1;
163  oss.str("");
164  oss << "h1_time_" << pl;
165  h1_time[pl] = FindMonHist(oss.str().c_str());
166  if (! h1_time[pl]) return 1;
167  }
168  h2_ele = FindMonHist("h2_ele" );
169  h2_time = FindMonHist("h2_time");
170  if (! h2_ele || ! h2_time) return 1;
171  return 0;
172 }
173 
175 {
176  OnlMonCanvas* can0 = GetCanvas(0);
178  TPad* pad0 = can0->GetMainPad();
179  pad0->SetGrid();
180  pad0->Divide(1, 2);
181  TVirtualPad* pad00 = pad0->cd(1); // for 1D histograms
182  pad00->Divide(2, 1);
183  pad00->cd(1); h1_ele[0]->Draw();
184  pad00->cd(2); h1_ele[1]->Draw();
185  pad0->cd(2); // for 2D histogram
186  h2_ele->Draw("colz");
187  if (h1_ele[0]->Integral() + h1_ele[1]->Integral() == 0) {
189  can0->AddMessage("Empty.");
190  }
191 
192  UtilHist::AutoSetRange (h1_time[0]);
193  UtilHist::AutoSetRange (h1_time[1]);
194  UtilHist::AutoSetRangeX((TH2*)h2_time);
195  UtilHist::AutoSetRangeY((TH2*)h2_time);
196  OnlMonCanvas* can1 = GetCanvas(1);
198  TPad* pad1 = can1->GetMainPad();
199  pad1->SetGrid();
200  pad1->Divide(1, 2);
201  TVirtualPad* pad10 = pad1->cd(1); // for 1D histograms
202  pad10->Divide(2, 1);
203  pad10->cd(1); h1_time[0]->Draw();
204  pad10->cd(2); h1_time[1]->Draw();
205  pad1->cd(2); // for 2D histogram
206  h2_time->Draw("colz");
207  if (h1_time[0]->Integral() + h1_time[1]->Integral() == 0) {
209  can1->AddMessage("Empty.");
210  }
211 
212  return 0;
213 }
void AutoSetRangeX(TH2 *h2, const int margin_lo=5, const int margin_hi=5)
Definition: UtilHist.cc:29
int NumCanvases()
Definition: OnlMonClient.h:101
int DrawMonitor()
Definition: OnlMonH4.cc:174
void AutoSetRange(TH1 *h1, const int margin_lo=5, const int margin_hi=5)
Adjust the axis range via &quot;h1-&gt;GetXaxis()-&gt;SetRange(bin_lo, bin_hi)&quot; to zoom up non-empty bins...
Definition: UtilHist.cc:17
HodoType_t
Definition: OnlMonH4.h:7
OnlMonCanvas * GetCanvas(const int num=0)
std::vector< SQHit * >::const_iterator ConstIter
Definition: SQHitVector.h:37
int InitOnlMon(PHCompositeNode *topNode)
Definition: OnlMonH4.cc:69
virtual ConstIter end() const
Definition: SQHitVector.h:59
void SetStatus(const MonStatus_t stat)
Definition: OnlMonCanvas.h:42
int EndOnlMon(PHCompositeNode *topNode)
Definition: OnlMonH4.cc:150
void RegisterHist(TH1 *h1, const HistMode_t mode=MODE_ADD)
std::string Title()
Definition: OnlMonClient.h:73
int InitRunOnlMon(PHCompositeNode *topNode)
Definition: OnlMonH4.cc:74
virtual const std::string Name() const
Returns the name of this module.
Definition: Fun4AllBase.h:23
An SQ interface class to hold a list of SQHit objects.
Definition: SQHitVector.h:32
An SQ interface class to hold one event header.
Definition: SQEvent.h:17
void AutoSetRangeY(TH2 *h2, const int margin_lo=5, const int margin_hi=5)
Definition: UtilHist.cc:43
int ProcessEventOnlMon(PHCompositeNode *topNode)
Definition: OnlMonH4.cc:127
TPad * GetMainPad()
Definition: OnlMonCanvas.cc:74
static GeomSvc * instance()
singlton instance
Definition: GeomSvc.cxx:211
virtual size_t size() const
Definition: SQHitVector.h:50
int getDetectorID(const std::string &detectorName) const
Get the plane position.
Definition: GeomSvc.h:184
int FindAllMonHist()
Definition: OnlMonH4.cc:155
TH1 * FindMonHist(const std::string name, const bool non_null=true)
int getPlaneNElements(int detectorID)
Definition: GeomSvc.h:208
void AddMessage(const char *msg)
Definition: OnlMonCanvas.cc:50
OnlMonH4(const HodoType_t type)
Definition: OnlMonH4.cc:19
SQHitVector * FindHits(const SQHitVector *vec_in, const std::string det_name)
Extract a set of hits that are of the given detector (det_name).
Definition: UtilSQHit.cc:15