Jpp  18.1.0
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
examples/JDAQ/JDAQTimeslice.cc File Reference

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

#include <string>
#include <iostream>
#include <iomanip>
#include <limits>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "km3net-dataformat/online/JDAQPMTIdentifier.hh"
#include "JDAQ/JDAQTimesliceIO.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JSupport.hh"
#include "JROOT/JManager.hh"
#include "Jeep/JeepToolkit.hh"
#include "JLang/JObjectMultiplexer.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::JDAQTimeslice.

Author
mdejong

Definition in file examples/JDAQ/JDAQTimeslice.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 31 of file examples/JDAQ/JDAQTimeslice.cc.

32 {
33  using namespace std;
34  using namespace JPP;
35  using namespace KM3NETDAQ;
36 
38  JLimit_t& numberOfEvents = inputFile.getLimit();
39  string outputFile;
41  int debug;
42 
43  try {
44 
45  JParser<> zap("Example program to histogram timeslice data.");
46 
47  zap['f'] = make_field(inputFile);
48  zap['o'] = make_field(outputFile) = "timeslice.root";
49  zap['n'] = make_field(numberOfEvents) = JLimit::max();
50  zap['P'] = make_field(PMT) = JDAQPMTIdentifier(-1, -1);
51  zap['d'] = make_field(debug) = 2;
52 
53  zap(argc, argv);
54  }
55  catch(const exception& error) {
56  FATAL(error.what() << endl);
57  }
58 
59 
60  const bool zoom = (PMT.getModuleID() != -1 &&
61  PMT.getPMTAddress() != -1);
62 
63  JManager<string, TH1D> H0(new TH1D("h0[%]", NULL, numeric_limits<JDAQHit::JPMT_t>::max(), -0.5, numeric_limits<JDAQHit::JPMT_t>::max() - 0.5));
64  JManager<string, TH1D> H1(new TH1D("h1[%]", NULL, numeric_limits<JDAQHit::JTOT_t>::max(), -0.5, numeric_limits<JDAQHit::JTOT_t>::max() - 0.5));
65  JManager<string, TH1D> H2(new TH1D("h2[%]", NULL, (zoom ? 20000000 : 1000), -0.5, getFrameTime() + 0.5));
66  JManager<string, TH1D> H3(new TH1D("h3[%]", NULL, 100, 0.0, 7.0));
67 
68  H0->Sumw2();
69  H1->Sumw2();
70  H2->Sumw2();
71 
73 
74  for (counter_type counter = 0; in.hasNext(); ++counter) {
75 
76  STATUS("event: " << setw(10) << counter << '\r'); DEBUG(endl);
77 
78  JDAQTimeslice* timeslice = in.next();
79 
80  TH1D* h0 = H0[getClassname(timeslice->GetName())];
81  TH1D* h1 = H1[getClassname(timeslice->GetName())];
82  TH1D* h2 = H2[getClassname(timeslice->GetName())];
83  TH1D* h3 = H3[getClassname(timeslice->GetName())];
84 
85  for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
86 
87  h3->Fill(log10((double) frame->size()));
88 
89  for (JDAQSuperFrame::const_iterator hit = frame->begin(); hit != frame->end(); ++hit) {
90 
91  if (JDAQPMTIdentifier::compare(PMT, JDAQPMTIdentifier(frame->getModuleID(), hit->getPMT()))) {
92  h0->Fill(hit->getPMT());
93  h1->Fill(hit->getToT());
94  h2->Fill(hit->getT());
95  }
96  }
97  }
98  }
99  STATUS(endl);
100 
101  TFile out(outputFile.c_str(), "recreate");
102 
103  out << H0 << H1 << H2 << H3;
104 
105  out.Write();
106  out.Close();
107 }
Utility class to parse command line options.
Definition: JParser.hh:1514
static const JPBS_t PMT(3, 4, 2, 3)
PBS of photo-multiplier tube (PMT)
#define STATUS(A)
Definition: JMessage.hh:63
Long64_t counter_type
Type definition for counter.
Auxiliary class for multiplexing object iterators.
string outputFile
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
std::string getClassname(const std::string &type_name)
Get type name, i.e. part after JEEP::TYPENAME_SEPARATOR.
Definition: JeepToolkit.hh:284
Hit data structure.
Definition: JDAQHit.hh:34
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1989
set_variable E_E log10(E_{fit}/E_{#mu})"
double getFrameTime()
Get frame time duration.
Definition: JDAQClock.hh:162
Data time slice.
#define FATAL(A)
Definition: JMessage.hh:67
General purpose class for object reading from a list of file names.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:84
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