Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JDAQProfile.cc
Go to the documentation of this file.
1 
2 #include <string>
3 #include <iostream>
4 #include <iomanip>
5 
6 #include "TROOT.h"
7 #include "TFile.h"
8 #include "TH1D.h"
9 #include "TProfile.h"
10 
11 #include "JDAQ/JDAQSummaryslice.hh"
12 #include "JDAQ/JDAQEvent.hh"
13 #include "JDAQ/JDAQClock.hh"
14 #include "JDAQ/JDAQEvaluator.hh"
15 
17 #include "JSupport/JTreeScanner.hh"
19 #include "JSupport/JSupport.hh"
20 
21 #include "JGizmo/JManager.hh"
22 
23 #include "Jeep/JParser.hh"
24 #include "Jeep/JMessage.hh"
25 
26 
27 /**
28  * \file
29  *
30  * Example program to histogram various data profiles.
31  * \author mdejong
32  */
33 int main(int argc, char **argv)
34 {
35  using namespace std;
36  using namespace JPP;
37  using namespace KM3NETDAQ;
38 
39  string inputFile;
40  JLimit_t numberOfEvents;
41  string outputFile;
42  Long64_t prescale;
43  int debug;
44 
45  try {
46 
47  JParser<> zap("Example program to histogram various data profiles.");
48 
49  zap['f'] = make_field(inputFile);
50  zap['o'] = make_field(outputFile) = "profile.root";
51  zap['n'] = make_field(numberOfEvents) = JLimit::max();
52  zap['P'] = make_field(prescale) = 1;
53  zap['d'] = make_field(debug) = 1;
54 
55  zap(argc, argv);
56  }
57  catch(const exception& error) {
58  FATAL(error.what() << endl);
59  }
60 
61 
62  cout.tie(&cerr);
63 
64 
65  NOTICE("Determine frame index range... " << flush);
66 
67  JFrameIndexRange frame_index = getFrameIndexRange<JDAQSummaryslice>(inputFile);
68 
69  NOTICE("= (" << frame_index.first << "," << frame_index.second << ")" << endl);
70 
71  const Long64_t NX = frame_index.second - frame_index.first + 1;
72  const Long64_t nx = (NX + prescale - 1) / prescale;
73  const Double_t xmin = (Double_t) frame_index.first - 0.5;
74  const Double_t xmax = (Double_t) frame_index.second + 0.5;
75 
76 
77  TFile out(outputFile.c_str(), "recreate");
78 
79  typedef JManager<int, TH1D> JManager_t;
80 
81  JManager_t map_summary(new TH1D("Summary[%]", NULL, NX, xmin, xmax));
82  JManager_t map_trigger(new TH1D("Trigger[%]", NULL, nx, xmin, xmax));
83  JManager_t map_status (new TH1D("Status[%]", NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS - 0.5));
84 
85  TH1D h0("L0", NULL, NX, xmin, xmax);
86  TH1D h1("WR", NULL, NX, xmin, xmax);
87  TH1D h2("trigger", NULL, nx, xmin, xmax);
88 
89  h2.Sumw2();
90 
91  for (JMultipleFileScanner<JDAQSummaryslice> in(inputFile, numberOfEvents); in.hasNext(); ) {
92 
93  STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
94 
95  JDAQSummaryslice* summaryslice = in.next();
96 
97  const Double_t x = summaryslice->getFrameIndex();
98 
99 
100  for (JDAQSummaryslice::const_iterator frame = summaryslice->begin(); frame != summaryslice->end(); ++frame) {
101 
102  if (frame->testHighRateVeto()) {
103 
104  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
105 
106  if (frame->testHighRateVeto(pmt)) {
107  map_status[frame->getModuleID()]->Fill((Double_t) pmt);
108  }
109  }
110  }
111 
112  Double_t y = 0.0;
113 
114  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
115  if (!frame->testHighRateVeto(pmt)) {
116  y += frame->getRate(pmt) * 1e-3; // kHz
117  }
118  }
119 
120  map_summary[frame->getModuleID()]->Fill(x, y);
121  }
122 
123 
124  Double_t y = 0.0;
125 
126  for (JDAQSummaryslice::const_iterator frame = summaryslice->begin(); frame != summaryslice->end(); ++frame) {
127 
128  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
129  y += frame->getRate(pmt) * 1e-3; // kHz
130  }
131 
132  h1.Fill(x, (frame->testWhiteRabbitStatus() ? 1.0 : 0.0));
133  }
134 
135  h0.Fill(x, y);
136  }
137  STATUS(endl);
138 
139 
140  for (JMultipleFileScanner<JDAQEvent> in(inputFile, numberOfEvents); in.hasNext(); ) {
141 
142  STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
143 
144  JDAQEvent* event = in.next();
145 
146  const Double_t x = event->getFrameIndex();
147  const Double_t y = 1.0e9 / getFrameTime() / prescale;
148 
149  h2.Fill(x, y);
150 
151  typedef JDAQTriggeredHit JHit_t;
152  //typedef JDAQSnapshotHit JHit_t;
153 
154  for (JDAQEvent::const_iterator<JHit_t> hit = event->begin<JHit_t>(); hit != event->end<JHit_t>(); ++hit) {
155  map_trigger[hit->getModuleID()]->Fill(x);
156  }
157  }
158  STATUS(endl);
159 
160  map_summary.Write(out);
161  map_trigger.Write(out);
162  map_status .Write(out);
163 
164  out.Write();
165  out.Close();
166 }
Utility class to parse command line options.
Definition: JParser.hh:1410
JTOOLS::JRange< int > JFrameIndexRange
Type definition for frame index range.
Auxiliary class to manage set of histograms.
#define STATUS(A)
Definition: JMessage.hh:61
Template const_iterator.
Definition: JDAQEvent.hh:69
Structure to store the ToT mean and standard deviation of the hits produced by a nanobeacon in a sour...
Dynamic ROOT object management.
string outputFile
int getFrameIndex() const
Get frame index.
JLimit JLimit_t
Type definition of limit.
Definition: JLimit.hh:214
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1836
double getFrameTime()
Get frame time duration.
Definition: JDAQClock.hh:162
#define NOTICE(A)
Definition: JMessage.hh:62
Support methods.
int debug
debug level
Definition: JSirene.cc:59
General purpose messaging.
#define FATAL(A)
Definition: JMessage.hh:65
Scanning of objects from multiple files according a format that follows from the extension of each fi...
Utility class to parse command line options.
ROOT TTree parameter settings.
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
JTriggerCounter_t next()
Increment trigger counter.
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:60
int main(int argc, char *argv[])
Definition: Main.cpp:15