Jpp test-rotations-new
the software that should make you happy
Loading...
Searching...
No Matches
JDAQProfile.cc File Reference

Example program to histogram various data profiles. More...

#include <string>
#include <iostream>
#include <iomanip>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TProfile.h"
#include "JDAQ/JDAQSummarysliceIO.hh"
#include "JDAQ/JDAQEventIO.hh"
#include "km3net-dataformat/online/JDAQClock.hh"
#include "JDAQ/JDAQEvaluator.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JTreeScanner.hh"
#include "JSupport/JSupportToolkit.hh"
#include "JSupport/JSupport.hh"
#include "JROOT/JManager.hh"
#include "Jeep/JParser.hh"
#include "Jeep/JMessage.hh"

Go to the source code of this file.

Functions

int main (int argc, char **argv)
 

Detailed Description

Example program to histogram various data profiles.

Author
mdejong

Definition in file JDAQProfile.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 33 of file JDAQProfile.cc.

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 NOTICE("Determine frame index range... " << flush);
63
64 const JFrameIndexRange frame_index = getFrameIndexRange<JDAQSummaryslice>(inputFile);
65
66 NOTICE("= (" << frame_index.first << "," << frame_index.second << ")" << endl);
67
68 const Long64_t NX = frame_index.second - frame_index.first + 1;
69 const Long64_t nx = (NX + prescale - 1) / prescale;
70 const Double_t xmin = frame_index.first;
71 const Double_t xmax = frame_index.first + nx * prescale;
72
73
74 TFile out(outputFile.c_str(), "recreate");
75
76 typedef JManager<int, TH1D> JManager_t;
77
78 JManager_t m_summary(new TH1D("Summary[%]", NULL, nx, xmin, xmax));
79 JManager_t m_trigger(new TH1D("Trigger[%]", NULL, nx, xmin, xmax));
80 JManager_t m_status (new TH1D("Status[%]", NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS - 0.5));
81
82 TH1D hu("#N [all]", NULL, nx, xmin, xmax);
83 TH1D hv("#M [UDP]", NULL, nx, xmin, xmax);
84 TH1D hw("#L [HRV]", NULL, nx, xmin, xmax);
85 TH1D ha("L0 [all]", NULL, nx, xmin, xmax);
86 TH1D hb("L0 [HRV]", NULL, nx, xmin, xmax);
87 TH1D h1("WR", NULL, nx, xmin, xmax);
88 TH1D h2("trigger", NULL, nx, xmin, xmax);
89 TH1D h3("PMT", NULL, nx, xmin, xmax);
90
91 h2.Sumw2();
92
93 for (JMultipleFileScanner<JDAQSummaryslice> in(inputFile, numberOfEvents); in.hasNext(); ) {
94
95 STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
96
97 JDAQSummaryslice* summary = in.next();
98
99 const Double_t x = summary->getFrameIndex();
100
101 Double_t U = summary->size();
102 Double_t V = 0.0;
103 Double_t W = 0.0;
104 Double_t Y[] = { 0.0, 0.0 };
105 Double_t Z = 0.0;
106
107 for (JDAQSummaryslice::const_iterator frame = summary->begin(); frame != summary->end(); ++frame) {
108
109 Double_t y = 0.0;
110
111 if (frame->testDAQStatus()) {
112
113 V += 1.0;
114
115 for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
116
117 if (frame->testHighRateVeto(pmt) || frame->testFIFOStatus(pmt))
118 m_status[frame->getModuleID()]->Fill((Double_t) pmt);
119 else
120 W += 1.0 / NUMBER_OF_PMTS;
121
122 y += frame->getRate(pmt) * 1e-3; // kHz
123 }
124 }
125
126 m_summary[frame->getModuleID()]->Fill(x, y);
127
128 for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
129
130 Y[0] += frame->getRate(pmt) * 1e-3; // kHz
131
132 if (!frame->testHighRateVeto(pmt) && !frame->testFIFOStatus(pmt)) {
133 Y[1] += frame->getRate(pmt) * 1e-3; // kHz
134 }
135 }
136
137 if (frame->testWhiteRabbitStatus()) {
138 Z += 1.0;
139 }
140 }
141
142 Y[0] /= (summary->size() * NUMBER_OF_PMTS);
143 Y[1] /= (summary->size() * NUMBER_OF_PMTS);
144 Z /= summary->size();
145
146 hu.Fill(x, U / prescale);
147 hv.Fill(x, V / prescale);
148 hw.Fill(x, W / prescale);
149 ha.Fill(x, Y[0] / prescale);
150 hb.Fill(x, Y[1] / prescale);
151 h1.Fill(x, Z / prescale);
152 }
153 STATUS(endl);
154
155 for (JMultipleFileScanner<JDAQEvent> in(inputFile, numberOfEvents); in.hasNext(); ) {
156
157 STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
158
159 JDAQEvent* event = in.next();
160
161 const Double_t x = event->getFrameIndex();
162 const Double_t y = 1.0e9 / getFrameTime() / prescale;
163
164 h2.Fill(x, y);
165
166 typedef JDAQTriggeredHit JHit_t;
167 //typedef JDAQSnapshotHit JHit_t;
168
169 for (JDAQEvent::const_iterator<JHit_t> hit = event->begin<JHit_t>(); hit != event->end<JHit_t>(); ++hit) {
170 m_trigger[hit->getModuleID()]->Fill(x);
171 }
172 }
173 STATUS(endl);
174
175 m_summary.Write(out);
176 m_trigger.Write(out);
177 m_status .Write(out);
178
179 out.Write();
180 out.Close();
181}
string outputFile
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#define STATUS(A)
Definition JMessage.hh:63
#define NOTICE(A)
Definition JMessage.hh:64
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:72
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
Utility class to parse command line options.
Definition JParser.hh:1698
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys.
Definition JManager.hh:47
General purpose class for object reading from a list of file names.
virtual bool hasNext() override
Check availability of next element.
JKey_t first
Definition JPair.hh:128
JValue_t second
Definition JPair.hh:129
int getFrameIndex() const
Get frame index.
Template const_iterator.
Definition JDAQEvent.hh:68
JTriggerCounter_t next()
Increment trigger counter.
const double xmax
const double xmin
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
JFrameIndexRange getFrameIndexRange(JTreeScannerInterface< T, KM3NETDAQ::JDAQEvaluator > &in)
Get range of frame indices.
KM3NeT DAQ data structures and auxiliaries.
Definition DataQueue.cc:39
double getFrameTime()
Get frame time duration.
Definition JDAQClock.hh:162
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition JDAQ.hh:26
Auxiliary class to set-up Hit.
Definition JSirene.hh:58
Auxiliary class for defining the range of iterations of objects.
Definition JLimit.hh:45
static counter_type max()
Get maximum counter value.
Definition JLimit.hh:128