Jpp  19.1.0
the software that should make you happy
Functions
JDAQHitRouter.cc File Reference

Example program to histogram string and floor hits. More...

#include <string>
#include <iostream>
#include <iomanip>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TH2D.h"
#include "JTools/JRange.hh"
#include "JROOT/JManager.hh"
#include "JROOT/JRootToolkit.hh"
#include "JDAQ/JDAQEventIO.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JDAQHitRouter.hh"
#include "JDetector/JStringRouter.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 string and floor hits.

Author
mdejong

Definition in file JDAQHitRouter.cc.

Function Documentation

◆ main()

int main ( int  argc,
char **  argv 
)

Definition at line 35 of file JDAQHitRouter.cc.

36 {
37  using namespace std;
38  using namespace JPP;
39  using namespace KM3NETDAQ;
40 
42  JLimit_t& numberOfEvents = inputFile.getLimit();
43  string outputFile;
44  string detectorFile;
45  int debug;
46 
47  try {
48 
49  JParser<> zap("Example program to histogram string and floor hits.");
50 
51  zap['f'] = make_field(inputFile);
52  zap['o'] = make_field(outputFile) = "router.root";
53  zap['a'] = make_field(detectorFile);
54  zap['n'] = make_field(numberOfEvents) = JLimit::max();
55  zap['d'] = make_field(debug) = 2;
56 
57  zap(argc, argv);
58  }
59  catch(const exception& error) {
60  FATAL(error.what() << endl);
61  }
62 
63 
65 
66  try {
67  load(detectorFile, detector);
68  }
69  catch(const JException& error) {
70  FATAL(error);
71  }
72 
73  const floor_range range = getRangeOfFloors(detector);
74 
75  const JDAQHitRouter router(detector);
76  const JStringRouter string(detector);
77 
78 
79  TH1D h1("h1", NULL, 100, 0.0, 1.0e1);
80 
81  JManager<int, TH2D> H2(new TH2D("h2[%]", NULL,
82  string.size(), -0.5, string.size() - 0.5,
83  range.getUpperLimit(), 1 - 0.5, range.getUpperLimit() + 0.5));
84 
85  for (Int_t i = 1; i <= H2->GetXaxis()->GetNbins(); ++i) {
86  H2->GetXaxis()->SetBinLabel(i, MAKE_CSTRING(string.at(i-1)));
87  }
88  for (Int_t i = 1; i <= H2->GetYaxis()->GetNbins(); ++i) {
89  H2->GetYaxis()->SetBinLabel(i, MAKE_CSTRING(i));
90  }
91 
92  TH2D* h3 = (TH2D*) H2->Clone("h3");
93 
94  while (inputFile.hasNext()) {
95 
96  STATUS("event: " << setw(10) << inputFile.getCounter() << '\r'); DEBUG(endl);
97 
98  JDAQEvent* event = inputFile.next();
99 
100  JTimeRange tx = JTimeRange::DEFAULT_RANGE();
101 
102  {
103  typedef JDAQTriggeredHit JHit_t;
104 
105  for (JDAQEvent::const_iterator<JHit_t> hit = event->begin<JHit_t>(); hit != event->end<JHit_t>(); ++hit) {
106 
107  const JPMTChannel& channel = router.getPMTChannel(*hit);
108 
109  H2->Fill((double) string.getIndex(channel.getString()), (double) channel.getFloor());
110 
111  for (int i = 0; i != NUMBER_OF_TRIGGER_BITS; ++i) {
112  if (hit->hasTriggerBit(i)) {
113  H2[i]->Fill((double) string.getIndex(channel.getString()), (double) channel.getFloor());
114  }
115  }
116 
117  tx.include(getTime(*hit, router.getPMT(*hit)));
118  }
119  }
120 
121  h1.Fill(tx.getLength() * 1.0e-3);
122 
123  {
124  typedef JDAQSnapshotHit JHit_t;
125 
126  for (JDAQEvent::const_iterator<JHit_t> hit = event->begin<JHit_t>(); hit != event->end<JHit_t>(); ++hit) {
127 
128  const JPMTChannel& channel = router.getPMTChannel(*hit);
129 
130  h3->Fill((double) string.getIndex(channel.getString()), (double) channel.getFloor());
131  }
132  }
133  }
134  STATUS(endl);
135 
136  TFile out(outputFile.c_str(), "recreate");
137 
138  out << h1;
139  out << *H2;
140  out << H2;
141  out << *h3;
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:69
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:2142
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:72
Simple wrapper around JModuleRouter class for direct addressing of PMT data in detector data structur...
Detector data structure.
Definition: JDetector.hh:96
int getFloor() const
Get floor number.
Definition: JLocation.hh:146
int getString() const
Get string number.
Definition: JLocation.hh:135
Auxiliary class to uniquely identify PMT readout channel.
Definition: JPMTChannel.hh:35
General exception.
Definition: JException.hh:24
Utility class to parse command line options.
Definition: JParser.hh:1698
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys.
Definition: JManager.hh:47
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.
range_type & include(argument_type x)
Include given value to range.
Definition: JRange.hh:397
T getLength() const
Get length (difference between upper and lower limit).
Definition: JRange.hh:289
T getUpperLimit() const
Get upper limit.
Definition: JRange.hh:213
Template const_iterator.
Definition: JDAQEvent.hh:68
double getTime(const Hit &hit)
Get true time of hit.
floor_range getRangeOfFloors(const JDetector &detector)
Get range of floors.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
int getIndex()
Get index for user I/O manipulation.
Definition: JManip.hh:26
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
KM3NeT DAQ data structures and auxiliaries.
Definition: DataQueue.cc:39
static const unsigned int NUMBER_OF_TRIGGER_BITS
Number of trigger bits.
Definition: JSTDTypes.hh:14
Detector file.
Definition: JHead.hh:227
Router for mapping of string identifier to index.
Auxiliary class to set-up Hit.
Definition: JSirene.hh:58
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:45