Jpp
 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 "JDAQ/JDAQPMTIdentifier.hh"
#include "JDAQ/JDAQTimeslice.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JSupport.hh"
#include "JGizmo/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 
37  JMultipleFileScanner<JDAQTimesliceTypes_t> inputFile;
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) = 1;
52 
53  zap(argc, argv);
54  }
55  catch(const exception& error) {
56  FATAL(error.what() << endl);
57  }
58 
59 
60  cout.tie(&cerr);
61 
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, 1000, -0.5, getFrameTime() + 0.5));
66 
67  h0->Sumw2();
68  h1->Sumw2();
69  h2->Sumw2();
70 
71  JObjectMultiplexer<JDAQTimesliceTypes_t> in(inputFile);
72 
73  for (counter_type counter = 0; in.hasNext(); ++counter) {
74 
75  STATUS("event: " << setw(10) << counter << '\r'); DEBUG(endl);
76 
77  JDAQTimeslice* timeslice = in.next();
78 
79  for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
80  for (JDAQSuperFrame::const_iterator hit = frame->begin(); hit != frame->end(); ++hit) {
81 
82  if (JDAQPMTIdentifier::compare(PMT, JDAQPMTIdentifier(frame->getModuleID(), hit->getPMT()))) {
83  h0[getClassname(timeslice->GetName())]->Fill(hit->getPMT());
84  h1[getClassname(timeslice->GetName())]->Fill(hit->getToT());
85  h2[getClassname(timeslice->GetName())]->Fill(hit->getT());
86  }
87  }
88  }
89  }
90  STATUS(endl);
91 
92  TFile out(outputFile.c_str(), "recreate");
93 
94  h0.Write(out);
95  h1.Write(out);
96  h2.Write(out);
97 
98  out.Write();
99  out.Close();
100 }
Utility class to parse command line options.
Definition: JParser.hh:1410
#define STATUS(A)
Definition: JMessage.hh:61
Long64_t counter_type
Type definition for counter.
Definition: JCounter.hh:24
string outputFile
JLimit JLimit_t
Type definition of limit.
Definition: JLimit.hh:214
std::string getClassname(const std::string &type_name)
Get type name, i.e.
Definition: JeepToolkit.hh:226
Hit data structure.
Definition: JDAQHit.hh:36
#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
Data time slice.
int debug
debug level
Definition: JSirene.cc:59
#define FATAL(A)
Definition: JMessage.hh:65
JDAQPMTIdentifier PMT
Command line options.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:72
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:60