Class Reference for E1039 Core & Analysis Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
OnlMonMainDaq.cc
Go to the documentation of this file.
1 #include <iomanip>
3 #include <TH1D.h>
4 #include <TCanvas.h>
5 #include <TPaveText.h>
6 #include <interface_main/SQRun.h>
9 #include <phool/PHNodeIterator.h>
10 #include <phool/PHIODataNode.h>
11 #include <phool/getClass.h>
12 #include <UtilAna/UtilHist.h>
13 #include "OnlMonServer.h"
14 #include "OnlMonMainDaq.h"
15 using namespace std;
16 
18 {
19  Name("OnlMonMainDaq");
20  Title("Main DAQ");
21  NumCanvases(1);
22 }
23 
25 {
27 }
28 
30 {
31  h1_trig = new TH1D("h1_trig", "Trigger Status;Trigger;N of events", 10, 0.5, 10.5);
32  h1_n_taiwan = new TH1D("h1_n_taiwan", "Taiwan TDC Status;N of boards/event;N of events", 200, -0.5, 199.5);
33  h1_evt_qual = new TH1D("h1_evt_qual", "Event Status;Event-quality bit;N of events", 33, -0.5, 32.5);
34  h1_flag_v1495 = new TH1D("h1_flag_v1495", "V1495 Status;v1495 status bit; N of v1495 events", 5, -0.5, 4.5);
35  h1_cnt = new TH1D("h1_cnt", ";Type;Count", 15, 0.5, 15.5);
36 
37  h1_trig->GetXaxis()->SetBinLabel( 1, "FPGA1");
38  h1_trig->GetXaxis()->SetBinLabel( 2, "FPGA2");
39  h1_trig->GetXaxis()->SetBinLabel( 3, "FPGA3");
40  h1_trig->GetXaxis()->SetBinLabel( 4, "FPGA4");
41  h1_trig->GetXaxis()->SetBinLabel( 5, "FPGA5");
42  h1_trig->GetXaxis()->SetBinLabel( 6, "NIM1");
43  h1_trig->GetXaxis()->SetBinLabel( 7, "NIM2");
44  h1_trig->GetXaxis()->SetBinLabel( 8, "NIM3");
45  h1_trig->GetXaxis()->SetBinLabel( 9, "NIM4");
46  h1_trig->GetXaxis()->SetBinLabel(10, "NIM5");
47 
48  ostringstream oss;
49  h1_evt_qual->GetXaxis()->SetBinLabel(1, "OK");
50  for (int ii = 1; ii <= 32; ii++) {
51  oss.str(""); oss << ii;
52  h1_evt_qual->GetXaxis()->SetBinLabel(ii+1, oss.str().c_str());
53  }
54 
55  h1_flag_v1495->GetXaxis()->SetBinLabel(1, "OK");
56  h1_flag_v1495->GetXaxis()->SetBinLabel(2, "d1ad");
57  h1_flag_v1495->GetXaxis()->SetBinLabel(3, "d2ad");
58  h1_flag_v1495->GetXaxis()->SetBinLabel(4, "d3ad");
59  h1_flag_v1495->GetXaxis()->SetBinLabel(5, "Other");
60 
61  RegisterHist(h1_trig);
62  RegisterHist(h1_n_taiwan);
63  RegisterHist(h1_evt_qual);
64  RegisterHist(h1_flag_v1495);
65  RegisterHist(h1_cnt, OnlMonClient::MODE_UPDATE);
66 
68 }
69 
71 {
72  SQRun* run = findNode::getClass<SQRun >(topNode, "SQRun");
73  SQEvent* event = findNode::getClass<SQEvent>(topNode, "SQEvent");
74  if (! run || ! event) return Fun4AllReturnCodes::ABORTEVENT;
75 
76  if (event->get_trigger(SQEvent::MATRIX1)) h1_trig->Fill( 1);
77  if (event->get_trigger(SQEvent::MATRIX2)) h1_trig->Fill( 2);
78  if (event->get_trigger(SQEvent::MATRIX3)) h1_trig->Fill( 3);
79  if (event->get_trigger(SQEvent::MATRIX4)) h1_trig->Fill( 4);
80  if (event->get_trigger(SQEvent::MATRIX5)) h1_trig->Fill( 5);
81  if (event->get_trigger(SQEvent::NIM1 )) h1_trig->Fill( 6);
82  if (event->get_trigger(SQEvent::NIM2 )) h1_trig->Fill( 7);
83  if (event->get_trigger(SQEvent::NIM3 )) h1_trig->Fill( 8);
84  if (event->get_trigger(SQEvent::NIM4 )) h1_trig->Fill( 9);
85  if (event->get_trigger(SQEvent::NIM5 )) h1_trig->Fill(10);
86 
87  h1_cnt->SetBinContent( 1, run->get_n_spill ());
88  h1_cnt->SetBinContent( 2, run->get_n_evt_all ());
89  h1_cnt->SetBinContent( 3, run->get_n_evt_dec ());
90  h1_cnt->SetBinContent( 4, run->get_n_phys_evt ());
91  h1_cnt->SetBinContent( 5, run->get_n_phys_evt_bad ());
92  h1_cnt->SetBinContent( 6, run->get_n_flush_evt ());
93  h1_cnt->SetBinContent( 7, run->get_n_flush_evt_bad());
94  h1_cnt->SetBinContent( 8, run->get_n_hit ());
95  h1_cnt->SetBinContent( 9, run->get_n_t_hit ());
96  h1_cnt->SetBinContent(10, run->get_n_hit_bad ());
97  h1_cnt->SetBinContent(11, run->get_n_t_hit_bad ());
98  h1_cnt->SetBinContent(12, run->get_n_v1495 ());
99  h1_cnt->SetBinContent(13, run->get_n_v1495_d1ad ());
100  h1_cnt->SetBinContent(14, run->get_n_v1495_d2ad ());
101  h1_cnt->SetBinContent(15, run->get_n_v1495_d3ad ());
102 
103  h1_n_taiwan->Fill(event->get_n_board_taiwan());
104 
105  int dq = event->get_data_quality();
106  if (dq == 0) h1_evt_qual->Fill(0);
107  for (int bit = 0; bit < 32; bit++) {
108  if ((dq >> bit) & 0x1) h1_evt_qual->Fill(bit + 1);
109  }
110 
111  int v1495 = event->get_flag_v1495();
112  if (v1495 == 0) h1_flag_v1495->Fill(0);
113  for (int bit = 0; bit < 4; bit++) {
114  if ((v1495 >> bit) & 0x1) h1_flag_v1495->Fill(bit + 1);
115  }
116 
118 }
119 
121 {
123 }
124 
126 {
127  h1_trig = FindMonHist("h1_trig");
128  h1_n_taiwan = FindMonHist("h1_n_taiwan");
129  h1_evt_qual = FindMonHist("h1_evt_qual");
130  h1_flag_v1495 = FindMonHist("h1_flag_v1495");
131  h1_cnt = FindMonHist("h1_cnt");
132  return (h1_trig && h1_evt_qual && h1_flag_v1495 && h1_cnt ? 0 : 1);
133 }
134 
136 {
137  UtilHist::AutoSetRange(h1_n_taiwan);
138 
139  OnlMonCanvas* can = GetCanvas();
141 
142  TPad* pad = can->GetMainPad();
143  pad->SetGrid();
144  pad->Divide(1, 3);
145 
146  pad->cd(1);
147  TPad* pad11 = new TPad("pad11", "", 0.0, 0.0, 0.6, 1.0);
148  TPad* pad12 = new TPad("pad12", "", 0.6, 0.0, 1.0, 1.0);
149  pad11->Draw();
150  pad12->Draw();
151  pad11->cd(); h1_trig->Draw();
152  pad12->cd(); h1_flag_v1495->Draw();
153 
154  pad->cd(2);
155  TPad* pad21 = new TPad("pad21", "", 0.0, 0.0, 0.7, 1.0);
156  TPad* pad22 = new TPad("pad22", "", 0.7, 0.0, 1.0, 1.0);
157  pad21->Draw();
158  pad22->Draw();
159  pad21->cd(); h1_evt_qual->Draw();
160  pad22->cd(); h1_n_taiwan->Draw();
161 
162  pad->cd(3);
163  TPaveText* pate = new TPaveText(.02, .02, .98, .98);
164  ostringstream oss;
165  oss << "N of spill-counter events = " << h1_cnt->GetBinContent(1);
166  pate->AddText(oss.str().c_str());
167  oss.str("");
168  oss << "N of triggered events = " << h1_cnt->GetBinContent(2) << " (all), " << h1_cnt->GetBinContent(3) << " (decoded)";
169  pate->AddText(oss.str().c_str());
170 
171  double n_all = h1_cnt->GetBinContent(4);
172  double n_bad = h1_cnt->GetBinContent(5);
173  oss.str("");
174  oss << "N of Coda physics events = " << n_all << " (all), " << n_bad << " (bad)";
175  pate->AddText(oss.str().c_str());
176 
177  n_all = h1_cnt->GetBinContent(6);
178  n_bad = h1_cnt->GetBinContent(7);
179  oss.str("");
180  oss << "N of Coda flush events = " << n_all << " (all), " << n_bad << " (bad)";
181  pate->AddText(oss.str().c_str());
182  if (n_bad / n_all > 0.05) {
184  can->AddMessage("Errors in >5% of Coda flush events.");
185  } else if (n_bad / n_all > 0.01) {
187  can->AddMessage("Errors in >1% of Coda flush events.");
188  }
189 
190  n_all = h1_cnt->GetBinContent(8);
191  n_bad = h1_cnt->GetBinContent(10);
192  oss.str("");
193  oss << "N of Taiwan TDC hits = " << n_all << " (all), " << n_bad << " (bad)";
194  pate->AddText(oss.str().c_str());
195  if (n_bad / n_all > 0.05) {
197  can->AddMessage("Errors in >5% of Taiwan TDC hits.");
198  } else if (n_bad / n_all > 0.01) {
200  can->AddMessage("Errors in >1% of Taiwan TDC hits.");
201  }
202 
203  n_all = h1_cnt->GetBinContent(9);
204  n_bad = h1_cnt->GetBinContent(11);
205  oss.str("");
206  oss << "N of v1495 TDC hits = " << n_all << " (all), " << n_bad << " (bad)";
207  pate->AddText(oss.str().c_str());
208  if (n_bad / n_all > 0.05) {
210  can->AddMessage("Errors in >5% of v1495 TDC hits.");
211  } else if (n_bad / n_all > 0.01) {
213  can->AddMessage("Errors in >1% of v1495 TDC hits.");
214  }
215 
216  n_all = h1_cnt->GetBinContent(12);
217  double n_d1ad = h1_cnt->GetBinContent(13);
218  double n_d2ad = h1_cnt->GetBinContent(14);
219  double n_d3ad = h1_cnt->GetBinContent(15);
220  oss.str("");
221  oss << "N of v1495 events = " << n_all << " (all), " << n_d1ad << " (d1ad), " << n_d2ad << " (d2ad), " << n_d3ad << " (d3ad)";
222  pate->AddText(oss.str().c_str());
223  pate->Draw();
224  if ((n_d1ad + n_d2ad + n_d3ad) / n_all > 0.05) {
226  can->AddMessage("Errors in >5% of v1495 events.");
227  } else if ((n_d1ad + n_d2ad + n_d3ad) / n_all > 0.01) {
229  can->AddMessage("Errors in >1% of v1495 events.");
230  }
231 
232  return 0;
233 }
virtual int get_n_flush_evt_bad() const
Return the number of bad FLUSH events.
Definition: SQRun.h:83
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
virtual int get_n_evt_dec() const
Return the number of decoded events. The online decoding usually skips a part of events in order to f...
Definition: SQRun.h:71
virtual int get_n_t_hit_bad() const
Return the number of bad V1495-TDC (i.e. trigger) hits.
Definition: SQRun.h:95
virtual int get_n_flush_evt() const
Return the number of all FLUSH events.
Definition: SQRun.h:80
int InitRunOnlMon(PHCompositeNode *topNode)
virtual int get_n_t_hit() const
Return the number of all V1495-TDC (i.e. trigger) hits.
Definition: SQRun.h:89
virtual int get_n_phys_evt() const
Return the number of all PHYSICS events.
Definition: SQRun.h:74
int run(const int nEvents=1)
Definition: run.C:10
void SetWorseStatus(const MonStatus_t stat)
Definition: OnlMonCanvas.cc:56
virtual int get_n_spill() const
Return the number of spill events.
Definition: SQRun.h:65
virtual int get_n_v1495_d3ad() const
Return the number of V1495 events having &#39;d3ad&#39;.
Definition: SQRun.h:107
void SetStatus(const MonStatus_t stat)
Definition: OnlMonCanvas.h:42
virtual int get_n_hit_bad() const
Return the number of bad Taiwan-TDC hits.
Definition: SQRun.h:92
virtual int get_n_v1495() const
Return the number of all V1495 events.
Definition: SQRun.h:98
int EndOnlMon(PHCompositeNode *topNode)
virtual int get_n_phys_evt_bad() const
Return the number of bad PHYSICS events.
Definition: SQRun.h:77
virtual int get_n_v1495_d1ad() const
Return the number of V1495 events having &#39;d1ad&#39;.
Definition: SQRun.h:101
virtual int get_n_evt_all() const
Return the number of all recorded events.
Definition: SQRun.h:68
An SQ interface class to hold one event header.
Definition: SQEvent.h:17
virtual int get_n_v1495_d2ad() const
Return the number of V1495 events having &#39;d2ad&#39;.
Definition: SQRun.h:104
TPad * GetMainPad()
Definition: OnlMonCanvas.cc:74
An SQ interface class to hold the run-level info.
Definition: SQRun.h:18
void AddMessage(const char *msg)
Definition: OnlMonCanvas.cc:50
int ProcessEventOnlMon(PHCompositeNode *topNode)
int InitOnlMon(PHCompositeNode *topNode)
virtual int get_n_hit() const
Return the number of all Taiwan-TDC hits.
Definition: SQRun.h:86