Jpp  master_rocky-40-g5f0272dcd
the software that should make you happy
Functions
JTreeScanner.cc File Reference

Program to test ordered reading using JSUPPORT::JTreeScanner. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <limits>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "JDAQ/JDAQEventIO.hh"
#include "JDAQ/JDAQSummarysliceIO.hh"
#include "JDAQ/JDAQEvaluator.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JTreeScanner.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

Program to test ordered reading using JSUPPORT::JTreeScanner.

Author
mdejong

Definition in file JTreeScanner.cc.

Function Documentation

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 30 of file JTreeScanner.cc.

31 {
32  using namespace std;
33  using namespace JPP;
34  using namespace KM3NETDAQ;
35 
36  JMultipleFileScanner<> inputFile;
37  JLimit_t& numberOfEvents = inputFile.getLimit();
38  string outputFile;
39  int debug;
40 
41  try {
42 
43  JParser<> zap("Program to test ordered reading from a ROOT TTree.");
44 
45  zap['f'] = make_field(inputFile);
46  zap['n'] = make_field(numberOfEvents) = JLimit::max();
47  zap['o'] = make_field(outputFile) = "histogram.root";
48  zap['d'] = make_field(debug) = 2;
49 
50  zap(argc, argv);
51  }
52  catch(const exception& error) {
53  FATAL(error.what() << endl);
54  }
55 
56 
57  TFile out(outputFile.c_str(), "recreate");
58 
59  TH1D h0("h0", NULL, 401, -200.5, +200.5);
60  TH1D h1("h1", NULL, 401, -200.5, +200.5);
61 
62 
64 
66 
67  int frame_index = 0;
68 
69  for (JTreeScanner<JDAQEvent, JDAQEvaluator> in(inputFile); in.hasNext(); ) {
70 
71  STATUS("event: " << setw(10) << in.getCounter() << '\r'); DEBUG(endl);
72 
73  JDAQEvent* event = in.next();
74  Long64_t index = scan.find(*event);
75 
76  JDAQSummaryslice* summary = scan.getEntry(index);
77 
78  h0.Fill(event ->getFrameIndex() - frame_index);
79  h1.Fill(summary->getFrameIndex() - event->getFrameIndex());
80 
81  frame_index = event->getFrameIndex();
82  }
83 
84  out.Write();
85  out.Close();
86 }
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:69
#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
Template definition for direct access of elements in ROOT TChain.
int getFrameIndex() const
Get frame index.
JTriggerCounter_t next()
Increment trigger counter.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
KM3NeT DAQ data structures and auxiliaries.
Definition: DataQueue.cc:39
int getFrameIndex(const double t_ns)
Get frame index for a given time in ns.
Definition: JDAQClock.hh:251
Definition: JSTDTypes.hh:14
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:45