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

Example program to histogram KM3NETDAQ::JDAQSummaryslice. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <limits>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "TProfile.h"
#include "JDAQ/JDAQSummarysliceIO.hh"
#include "JDetector/JDetectorSupportkit.hh"
#include "JTrigger/JTriggerToolkit.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JSupport.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 KM3NETDAQ::JDAQSummaryslice.

Author
mdejong

Definition in file JDAQSummaryslice.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 30 of file JDAQSummaryslice.cc.

31{
32 using namespace std;
33 using namespace JPP;
34 using namespace KM3NETDAQ;
35
37 JLimit_t& numberOfEvents = inputFile.getLimit();
38 string outputFile;
39 bool correct;
40 int debug;
41
42 try {
43
44 JParser<> zap("Example program to histogram summary data.");
45
46 zap['f'] = make_field(inputFile);
47 zap['o'] = make_field(outputFile) = "summaryslice.root";
48 zap['n'] = make_field(numberOfEvents) = JLimit::max();
49 zap['c'] = make_field(correct);
50 zap['d'] = make_field(debug) = 2;
51
52 zap(argc, argv);
53 }
54 catch(const exception& error) {
55 FATAL(error.what() << endl);
56 }
57
58
59 const double factor = 1.0e-3; // [kHz]
60
61
62 TFile out(outputFile.c_str(), "recreate");
63
64 TH1D h0("h0", NULL, JDAQRate::getN(), JDAQRate::getData(factor));
65 TProfile h1("h1", NULL, 32,-0.5, 31.5);
66 TProfile h2("h2", NULL, 32,-0.5, 31.5);
67
68 TH2D hu("hu", NULL, 201, -0.5, +200.5, 201, -0.5, +200.5);
69
70 TH2D hv("hv", NULL, 51, -0.5, +50.5, 200, 0.0, 50.0);
71 TH2D hw("hw", NULL, 51, -0.5, +50.5, 200, 0.0, 50.0);
72
73 TH2D hx("hx", NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS - 0.5, 200, 0.0, 50.0);
74 TH2D hy("hy", NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS - 0.5, 200, 0.0, 50.0);
75
76
77 while (inputFile.hasNext()) {
78
79 STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
80
81 JDAQSummaryslice* summaryslice = inputFile.next();
82
83 DEBUG("Summary: "
84 << setw(8) << inputFile.getCounter() << ' '
85 << setw(8) << summaryslice->getRunNumber() << ' '
86 << setw(8) << summaryslice->getFrameIndex() << ' '
87 << setw(6) << summaryslice->size() << endl);
88
89 const JDetectorAddressMap& demo = getDetectorAddressMap(summaryslice->getDetectorID());
90
91 for (JDAQSummaryslice::const_iterator frame = summaryslice->begin(); frame != summaryslice->end(); ++frame) {
92
93 const JModuleAddressMap& memo = demo.get(frame->getModuleID());
94
95 int N[2] = { 0 };
96 double R[2] = { 0.0 };
97
98 int lower[2] = { 0 };
99 int upper[2] = { 0 };
100
101 for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
102
103 const int index = (!frame->testHighRateVeto(pmt) && !frame->testFIFOStatus(pmt) ? 0 : 1);
104
105 N[index] += 1;
106 R[index] += correct ? getRate(*frame, pmt, factor) : frame->getRate(pmt, factor);
107
108 if (memo.getAddressTranslator(pmt).ring <= 'D')
109 lower[index] += 1;
110 else
111 upper[index] += 1;
112
113 h0.Fill(correct ? getRate(*frame, pmt, factor) : frame->getRate(pmt, factor), frame->getWeight(pmt, factor));
114
115 h1.Fill((Double_t) pmt, (frame->testHighRateVeto(pmt) ? 1.0 : 0.0));
116 h2.Fill((Double_t) pmt, (frame->testFIFOStatus (pmt) ? 1.0 : 0.0));
117 }
118
119 hu.Fill((double) frame->getUDPMaximalSequenceNumber(),
120 (double) frame->getUDPNumberOfReceivedPackets());
121
122 if (N[0] != 0) { hv.Fill((double) frame->getUDPMaximalSequenceNumber(), R[0] / N[0]); }
123 if (N[1] != 0) { hw.Fill((double) frame->getUDPMaximalSequenceNumber(), R[1] / N[1]); }
124
125 if (N[0] != 0) {
126
127 hx.Fill((double) N[1], R[0] / N[0]);
128
129 if (lower[0] != 0 && upper[0] != 0) {
130 hy.Fill((double) N[1], R[0] / N[0]);
131 }
132 }
133
134 const bool status = (frame->getUDPNumberOfReceivedPackets() == frame->getUDPMaximalSequenceNumber() + 1 &&
135 frame->hasUDPTrailer());
136
137 h1.Fill((Double_t) 31, (frame->testWhiteRabbitStatus() ? 1.0 : 0.0));
138 h2.Fill((Double_t) 31, (status ? 1.0 : 0.0));
139 }
140 }
141 STATUS(endl);
142
143 out.Write();
144 out.Close();
145}
string outputFile
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#define STATUS(A)
Definition JMessage.hh:63
#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
Lookup table for PMT addresses in detector.
const JModuleAddressMap & get(const int id) const
Get module address map.
Lookup table for PMT addresses in optical module.
const JPMTAddressTranslator & getAddressTranslator(const int tdc) const
Get PMT address translator.
Utility class to parse command line options.
Definition JParser.hh:1698
General purpose class for object reading from a list of file names.
virtual bool hasNext() override
Check availability of next element.
counter_type getCounter() const
Get counter.
virtual const pointer_type & next() override
Get next element.
int getDetectorID() const
Get detector identifier.
int getRunNumber() const
Get run number.
int getFrameIndex() const
Get frame index.
static int getN()
Get number of bins.
static const double * getData(const double factor=1.0)
Get abscissa values.
JDetectorAddressMap & getDetectorAddressMap()
Get detector address map.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
double getRate(const JDAQSummaryFrame &frame, const int pmt, const double factor=1.0)
Get corrected rate of PMT.
KM3NeT DAQ data structures and auxiliaries.
Definition DataQueue.cc:39
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition JDAQ.hh:26
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