Jpp  18.0.1-rc.2
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
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

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  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 = (Double_t) frame_index.first - 0.5;
71  const Double_t xmax = (Double_t) frame_index.second + 0.5;
72 
73 
74  TFile out(outputFile.c_str(), "recreate");
75 
76  typedef JManager<int, TH1D> JManager_t;
77 
78  JManager_t map_summary(new TH1D("Summary[%]", NULL, NX, xmin, xmax));
79  JManager_t map_trigger(new TH1D("Trigger[%]", NULL, nx, xmin, xmax));
80  JManager_t map_status (new TH1D("Status[%]", NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS - 0.5));
81 
82  TH1D hn("N", NULL, NX, xmin, xmax);
83  TH1D h0("L0", NULL, NX, xmin, xmax);
84  TH1D h1("WR", NULL, NX, xmin, xmax);
85  TH1D h2("trigger", NULL, nx, xmin, xmax);
86 
87  h2.Sumw2();
88 
89  for (JMultipleFileScanner<JDAQSummaryslice> in(inputFile, numberOfEvents); in.hasNext(); ) {
90 
91  STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
92 
93  JDAQSummaryslice* summaryslice = in.next();
94 
95  const Double_t x = summaryslice->getFrameIndex();
96 
97 
98  for (JDAQSummaryslice::const_iterator frame = summaryslice->begin(); frame != summaryslice->end(); ++frame) {
99 
100  if (frame->testHighRateVeto()) {
101 
102  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
103 
104  if (frame->testHighRateVeto(pmt)) {
105  map_status[frame->getModuleID()]->Fill((Double_t) pmt);
106  }
107  }
108  }
109 
110  Double_t y = 0.0;
111 
112  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
113  if (!frame->testHighRateVeto(pmt)) {
114  y += frame->getRate(pmt) * 1e-3; // kHz
115  }
116  }
117 
118  map_summary[frame->getModuleID()]->Fill(x, y);
119  }
120 
121 
122  Double_t y = 0.0;
123 
124  for (JDAQSummaryslice::const_iterator frame = summaryslice->begin(); frame != summaryslice->end(); ++frame) {
125 
126  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
127  y += frame->getRate(pmt) * 1e-3; // kHz
128  }
129 
130  h1.Fill(x, (frame->testWhiteRabbitStatus() ? 1.0 : 0.0));
131  }
132 
133  hn.Fill(x, summaryslice->size());
134  h0.Fill(x, y);
135  }
136  STATUS(endl);
137 
138 
139  for (JMultipleFileScanner<JDAQEvent> in(inputFile, numberOfEvents); in.hasNext(); ) {
140 
141  STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
142 
143  JDAQEvent* event = in.next();
144 
145  const Double_t x = event->getFrameIndex();
146  const Double_t y = 1.0e9 / getFrameTime() / prescale;
147 
148  h2.Fill(x, y);
149 
150  typedef JDAQTriggeredHit JHit_t;
151  //typedef JDAQSnapshotHit JHit_t;
152 
153  for (JDAQEvent::const_iterator<JHit_t> hit = event->begin<JHit_t>(); hit != event->end<JHit_t>(); ++hit) {
154  map_trigger[hit->getModuleID()]->Fill(x);
155  }
156  }
157  STATUS(endl);
158 
159  map_summary.Write(out);
160  map_trigger.Write(out);
161  map_status .Write(out);
162 
163  out.Write();
164  out.Close();
165 }
const double xmax
Definition: JQuadrature.cc:24
Utility class to parse command line options.
Definition: JParser.hh:1514
#define STATUS(A)
Definition: JMessage.hh:63
Template const_iterator.
Definition: JDAQEvent.hh:68
string outputFile
int getFrameIndex() const
Get frame index.
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys...
Definition: JManager.hh:43
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1989
double getFrameTime()
Get frame time duration.
Definition: JDAQClock.hh:162
#define NOTICE(A)
Definition: JMessage.hh:64
#define FATAL(A)
Definition: JMessage.hh:67
const double xmin
Definition: JQuadrature.cc:23
Auxiliary class to set-up Hit.
Definition: JSirene.hh:57
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable INPUT_FILE $argv[2] eval JPrintDetector a $DETECTOR O IDENTIFIER eval JPrintDetector a $DETECTOR O SUMMARY JAcoustics sh $DETECTOR_ID source JAcousticsToolkit sh CHECK_EXIT_CODE typeset A EMITTERS get_tripods $WORKDIR tripod txt EMITTERS get_transmitters $WORKDIR transmitter txt EMITTERS for EMITTER in
Definition: JCanberra.sh:46
int debug
debug level
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62