Jpp  master_rocky
the software that should make you happy
JCalibrateToT.cc
Go to the documentation of this file.
1 
2 #include <string>
3 #include <iostream>
4 #include <iomanip>
5 #include <limits>
6 
7 #include "TROOT.h"
8 #include "TFile.h"
9 #include "TH2D.h"
10 
13 
14 #include "JDAQ/JDAQTimesliceIO.hh"
15 #include "JDetector/JDetector.hh"
18 #include "JTrigger/JHitR0.hh"
19 #include "JTrigger/JMatchL0.hh"
24 #include "JSupport/JSupport.hh"
25 #include "JSupport/JMeta.hh"
26 #include "JROOT/JManager.hh"
27 #include "JMath/JMathToolkit.hh"
28 
30 #include "JROOT/JRootToolkit.hh"
31 #include "JROOT/JRootFileWriter.hh"
36 
37 #include "Jeep/JPrint.hh"
38 #include "Jeep/JParser.hh"
39 #include "Jeep/JMessage.hh"
40 
41 /**
42  * \file
43  *
44  * Monitoring of PMT time-over-threshold distributions.
45  * \author mkarel
46  */
47 int main(int argc, char **argv)
48 {
49  using namespace std;
50  using namespace JPP;
51  using namespace KM3NETDAQ;
52 
54  JLimit_t& numberOfEvents = inputFile.getLimit();
55  string outputFile;
56  string detectorFile;
57  JRange<double> T_ns;
58  double ctMax;
59  JRange<int> multiplicity;
60  double deadTime_us;
61  JROOTClassSelector selector;
62  JPreprocessor option;
63  int debug;
64 
65  try {
66 
67  JParser<> zap("Monitoring of PMT time-over-threshold distributions.");
68 
69  zap['f'] = make_field(inputFile, "input file.");
70  zap['o'] = make_field(outputFile, "output file.") = "calibrate_tot.root";
71  zap['n'] = make_field(numberOfEvents) = JLimit::max();
72  zap['a'] = make_field(detectorFile, "detector file.");
73  zap['T'] = make_field(T_ns, "time window [ns].") = JRange<double>(10.0, 25.0);
74  zap['c'] = make_field(ctMax, "maximal cosine space angle between PMT axes.") = 0.0;
75  zap['M'] = make_field(multiplicity, "multiplicity range of hits on DOM.") = JRange<int>(1, 2);
76  zap['D'] = make_field(deadTime_us, "L0 dead time (us)") = 10;
77  zap['C'] = make_field(selector, "timeslice selector, e.g. JDAQTimesliceL1.") = getROOTClassSelection<JDAQTimesliceTypes_t>();
78  zap['O'] = make_field(option, "hit pre-processing option.") = JPreprocessor::getOptions();
79  zap['d'] = make_field(debug, "debug flag.") = 1;
80 
81  zap(argc, argv);
82  }
83  catch(const exception &error) {
84  FATAL(error.what() << endl);
85  }
86 
87 
88  //-----------------------------------------------------------
89  // check the input parameters
90  //-----------------------------------------------------------
91 
92  if (!T_ns.is_valid()) {
93  FATAL("Invalid time window [ns] " << T_ns << endl);
94  }
95 
96  if (selector == JDAQTimeslice ::Class_Name() ||
97  selector == JDAQTimesliceL1::Class_Name()) {
98 
99  JTriggerParameters parameters;
100 
101  try {
102  parameters = getTriggerParameters(inputFile);
103  }
104  catch(const JException& error) {
105  FATAL("No trigger parameters from input." << endl);
106  }
107 
108  if ((selector == JDAQTimeslice ::Class_Name() && parameters.writeL1.prescale > 0) ||
109  (selector == JDAQTimesliceL1::Class_Name())) {
110 
111  if (parameters.TMaxLocal_ns < T_ns.getUpperLimit()) {
112  FATAL("Option -T <T_ns> = " << T_ns.getUpperLimit() << " is larger than in the trigger " << parameters.TMaxLocal_ns << endl);
113  }
114  }
115  }
116 
117  if (!multiplicity.is_valid()) { FATAL("Invalid multiplicity " << multiplicity << endl); }
118  if ( multiplicity.getLowerLimit() < 1) { FATAL("Invalid multiplicity " << multiplicity << endl); }
119 
120  //-----------------------------------------------------------
121  // load detector file
122  //-----------------------------------------------------------
123 
125 
126  try {
127  load(detectorFile, detector);
128  }
129  catch(const JException& error) {
130  FATAL(error);
131  }
132 
133  if (detector.empty()) {
134  FATAL("Empty detector." << endl);
135  }
136 
137  const JModuleRouter router(detector);
138 
139  //-----------------------------------------------------------
140  // initialise histograms
141  //-----------------------------------------------------------
142 
143  const double zmin = -0.5; // [ns]
144  const double zmax = 256.5; // [ns]
145  const int nz = (int) ((zmax-zmin) / 1.0);
146 
147  JManager<int, TH2D> manager(new TH2D(MAKE_CSTRING("%" << _2SToT), NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS-0.5, nz, zmin, zmax));
148 
149  typedef JHitR0 hit_type;
150  typedef JSuperFrame2D<hit_type> JSuperFrame2D_t;
151  typedef JSuperFrame1D<hit_type> JSuperFrame1D_t;
152 
153  const JMatchL0<hit_type> match(T_ns.getUpperLimit()); // time window self-coincidences [ns]
154 
155  const double deadTime_ns = deadTime_us * 1e3;
156 
158 
159 
160  TFile out(outputFile.c_str(), "recreate");
161 
162  putObject(out, JMeta(argc, argv));
163 
164  JDAQHeader header;
165 
166  for (counter_type counter = 0; in.hasNext() && counter != inputFile.getLimit(); ++counter) {
167 
168  STATUS("event: " << setw(10) << counter << '\r'); DEBUG(endl);
169 
170  const JDAQTimeslice* timeslice = in.next();
171 
172  if (header.getDetectorID() != timeslice->getDetectorID() ||
173  header.getRunNumber () != timeslice->getRunNumber ()) {
174 
175  header = timeslice->getDAQHeader();
176 
177  putObject(out, header);
178  }
179 
180  for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
181 
182  if (router.hasModule(frame->getModuleID())) {
183 
184  TH2D* h2 = manager[frame->getModuleID()];
185  const JModule& module = router.getModule(frame->getModuleID());
186 
187  JSuperFrame2D_t& buffer = JSuperFrame2D_t::demultiplex(*frame, module);
188 
189  buffer.preprocess(option, match);
190 
191  JSuperFrame1D_t& data = JSuperFrame1D_t::multiplex(buffer);
192 
193  vector<double> t0(NUMBER_OF_PMTS, numeric_limits<double>::lowest());
194 
195  for (JSuperFrame1D_t::const_iterator p = data.begin(); p != data.end(); ) {
196 
197  JSuperFrame1D_t::const_iterator q = p;
198 
199  // reject genuine coincidences
200 
201  double ct_max = -1.0;
202  double dt_min = numeric_limits<double>::max();
203 
204  while (++q != data.end() && q->getT() - p->getT() < T_ns.getUpperLimit()) {
205 
206  const double ct = getDot(module.getPMT(p->getPMT()), module.getPMT(q->getPMT()));
207  const double dt = q->getT() - p->getT();
208 
209  if (ct > ct_max) {
210  ct_max = ct;
211  }
212 
213  if (dt < dt_min) {
214  dt_min = dt;
215  }
216  }
217 
218  if (multiplicity(distance(p,q)) && ct_max < ctMax && dt_min > T_ns.getLowerLimit()) {
219 
220  for (JSuperFrame1D_t::const_iterator i = p; i != q; ++i) {
221 
222  if (i->getT() > t0[i->getPMT()] + deadTime_ns) {
223  h2->Fill(i->getPMT(), i->getToT());
224  }
225 
226  t0[i->getPMT()] = i->getT();
227  }
228  }
229 
230  p = q;
231  }
232  }
233  }
234  }
235 
236  STATUS(endl);
237 
238  out << manager;
239 
240  out.Close();
241  out.Write();
242 }
int main(int argc, char **argv)
string outputFile
KM3NeT DAQ constants, bit handling, etc.
Data structure for detector geometry and calibration.
Basic data structure for L0 hit.
Dynamic ROOT object management.
Match operator for consecutive hits.
Auxiliary methods for geometrical methods.
General purpose messaging.
#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
ROOT I/O of application specific meta data.
Direct access to module in detector data structure.
Scanning of objects from multiple files according a format that follows from the extension of each fi...
Utility class to parse command line options.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:2142
Auxiliaries for pre-processing of hits.
I/O formatting auxiliaries.
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:72
ROOT TTree parameter settings of various packages.
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
Detector data structure.
Definition: JDetector.hh:96
Router for direct addressing of module data in detector data structure.
const JModule & getModule(const JObjectID &id) const
Get module parameters.
bool hasModule(const JObjectID &id) const
Has module.
Data structure for a composite optical module.
Definition: JModule.hh:75
const JPMT & getPMT(const int index) const
Get PMT.
Definition: JModule.hh:172
General exception.
Definition: JException.hh:24
Auxiliary class for multiplexing object iterators.
virtual bool hasNext() override
Check availability of next element.
virtual const pointer_type & next() override
Get next element.
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
void Write(TDirectory &out, const bool wm=false)
Write objects to file.
Definition: JManager.hh:304
General purpose class for object reading from a list of file names.
bool is_valid() const
Check validity of range.
Definition: JRange.hh:311
T getLowerLimit() const
Get lower limit.
Definition: JRange.hh:202
T getUpperLimit() const
Get upper limit.
Definition: JRange.hh:213
Reduced data structure for L0 hit.
Definition: JHitR0.hh:27
L0 match criterion.
Definition: JMatchL0.hh:29
1-dimensional frame with time calibrated data from one optical module.
2-dimensional frame with time calibrated data from one optical module.
int getDetectorID() const
Get detector identifier.
int getRunNumber() const
Get run number.
const JDAQHeader & getDAQHeader() const
Get DAQ header.
Definition: JDAQHeader.hh:49
double getDot(const JNeutrinoDirection &first, const JNeutrinoDirection &second)
Dot product.
Definition: JAstronomy.hh:676
static const char *const _2SToT
Histogram naming.
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
bool putObject(TDirectory &dir, const TObject &object)
Write object to ROOT directory.
Long64_t counter_type
Type definition for counter.
JTriggerParameters getTriggerParameters(const JMultipleFileScanner_t &file_list)
Get trigger parameters.
KM3NeT DAQ data structures and auxiliaries.
Definition: DataQueue.cc:39
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
Definition: JSTDTypes.hh:14
Detector file.
Definition: JHead.hh:227
Auxiliary class to select ROOT class based on class name.
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:45
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:72
Auxiliary class for specifying the way of pre-processing of hits.