Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JCalibrateToT.cc File Reference

Monitoring of PMT time-over-threshold distributions. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <limits>
#include "TROOT.h"
#include "TFile.h"
#include "TH2D.h"
#include "JDAQ/JDAQ.hh"
#include "JDAQ/JDAQTimeslice.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.hh"
#include "JTrigger/JHitR0.hh"
#include "JTrigger/JSuperFrame1D.hh"
#include "JTrigger/JSuperFrame2D.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JMeta.hh"
#include "JGizmo/JManager.hh"
#include "JMath/JMathToolkit.hh"
#include "JROOT/JROOTClassSelector.hh"
#include "JROOT/JRootToolkit.hh"
#include "JLang/JObjectMultiplexer.hh"
#include "JTrigger/JTriggerParameters.hh"
#include "JSupport/JTriggerParametersSupportkit.hh"
#include "JCalibrate/JCalibrateToT.hh"
#include "Jeep/JPrint.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

Monitoring of PMT time-over-threshold distributions.

Author
mkarel

Definition in file JCalibrateToT.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 49 of file JCalibrateToT.cc.

50 {
51  using namespace std;
52  using namespace JPP;
53  using namespace KM3NETDAQ;
54 
55  JMultipleFileScanner<JDAQTimesliceTypes_t> inputFile;
56  JLimit_t& numberOfEvents = inputFile.getLimit();
57  string outputFile;
58  string detectorFile;
59  JRange<double> T_ns;
60  double ctMax;
61  JRange<int> multiplicity;
62  double deadTime_us;
63  JROOTClassSelector selector;
64  string option;
65  int debug;
66 
67  try {
68 
69  JParser<> zap("Monitoring of PMT time-over-threshold distributions.");
70 
71  zap['f'] = make_field(inputFile, "input file.");
72  zap['o'] = make_field(outputFile, "output file.") = "calibrate_tot.root";
73  zap['n'] = make_field(numberOfEvents) = JLimit::max();
74  zap['a'] = make_field(detectorFile, "detector file.");
75  zap['T'] = make_field(T_ns, "time window [ns].") = JRange<double>(10.0, 25.0);
76  zap['c'] = make_field(ctMax, "maximal cosine space angle between PMT axes.") = 0.0;
77  zap['M'] = make_field(multiplicity, "multiplicity range of hits on DOM.") = JRange<int>(1, 2);
78  zap['D'] = make_field(deadTime_us, "L0 dead time (us)") = 10;
79  zap['C'] = make_field(selector, "timeslice selector, e.g. JDAQTimesliceL1.") = getROOTClassSelection<JDAQTimesliceTypes_t>();
80  zap['O'] = make_field(option, "hit processing option.") = none_t, join_t, filter_t;
81  zap['d'] = make_field(debug, "debug flag.") = 1;
82 
83  zap(argc, argv);
84  }
85  catch(const exception &error) {
86  FATAL(error.what() << endl);
87  }
88 
89  cout.tie(&cerr);
90 
91  //-----------------------------------------------------------
92  // check the input parameters
93  //-----------------------------------------------------------
94 
95  JTriggerParameters parameters;
96 
97  try {
98  parameters = getTriggerParameters(inputFile);
99  }
100  catch(const JException& error) {
101  FATAL("No trigger parameters from input." << endl);
102  }
103 
104  if (!T_ns.is_valid()) {
105  FATAL("Invalid time window [ns] " << T_ns << 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  if (!multiplicity.is_valid()) { FATAL("Invalid multiplicity " << multiplicity << endl); }
117  if ( multiplicity.getLowerLimit() < 1) { FATAL("Invalid multiplicity " << multiplicity << endl); }
118 
119  //-----------------------------------------------------------
120  // load detector file
121  //-----------------------------------------------------------
122 
123  JDetector detector;
124 
125  try {
126  load(detectorFile, detector);
127  }
128  catch(const JException& error) {
129  FATAL(error);
130  }
131 
132  if (detector.empty()) {
133  FATAL("Empty detector." << endl);
134  }
135 
136  const JModuleRouter router(detector);
137 
138  //-----------------------------------------------------------
139  // initialise histograms
140  //-----------------------------------------------------------
141 
142  const double zmin = -0.5; // [ns]
143  const double zmax = 256.5; // [ns]
144  const int nz = (int) ((zmax-zmin) / 1.0);
145 
146  JManager<int, TH2D> manager(new TH2D(MAKE_CSTRING("%" << _2SToT), NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS-0.5, nz, zmin, zmax));
147 
148  typedef JHitR0 hit_type;
149  typedef JSuperFrame2D<hit_type> JSuperFrame2D_t;
150  typedef JSuperFrame1D<hit_type> JSuperFrame1D_t;
151 
152  const JMatch_t match(T_ns.getUpperLimit()); // time window self-coincidences [ns]
153 
154  const double deadTime_ns = deadTime_us * 1e3;
155 
156  JObjectMultiplexer<JDAQTimesliceTypes_t, JDAQTimeslice> in(inputFile, selector);
157 
158  counter_type counter = 0;
159 
160  for ( ; in.hasNext() && counter != inputFile.getLimit(); ++counter) {
161 
162  STATUS("event: " << setw(10) << counter << '\r'); DEBUG(endl);
163 
164  const JDAQTimeslice* timeslice = in.next();
165 
166  for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
167 
168  if (router.hasModule(frame->getModuleID())) {
169 
170  TH2D* h2 = manager[frame->getModuleID()];
171  const JModule& module = router.getModule(frame->getModuleID());
172 
173  JSuperFrame2D_t& buffer = JSuperFrame2D_t::demultiplex(*frame, module);
174 
175  if (option != none_t) {
176  for (JSuperFrame2D_t::iterator i = buffer.begin(); i != buffer.end(); ++i) {
177  if (option == join_t)
178  i->join(match);
179  else if (option == filter_t)
180  i->filter(match);
181  }
182  }
183 
184  JSuperFrame1D_t& data = JSuperFrame1D_t::multiplex(buffer);
185 
186  vector<double> t0(NUMBER_OF_PMTS, -numeric_limits<double>::max());
187 
188  for (JSuperFrame1D_t::const_iterator p = data.begin(), __end = data.end() - 1; p != __end; ) {
189 
190  JSuperFrame1D_t::const_iterator q = p;
191 
192  double ct_max = -1.0;
193  double dt_min = numeric_limits<double>::max();
194 
195  while (++q != __end && q->getT() - p->getT() < T_ns.getUpperLimit()) {
196 
197  const double ct = getDot(module.getPMT(p->getPMT()), module.getPMT(q->getPMT()));
198  const double dt = q->getT() - p->getT();
199 
200  if (ct > ct_max) {
201  ct_max = ct;
202  }
203 
204  if (dt < dt_min) {
205  dt_min = dt;
206  }
207  }
208 
209  if (multiplicity(distance(p,q)) && ct_max < ctMax && dt_min > T_ns.getLowerLimit()) {
210 
211  for (JSuperFrame1D_t::const_iterator __p = p; __p != q; ++__p) {
212 
213  if (__p->getT() > t0[__p->getPMT()] + deadTime_ns) {
214  h2->Fill(__p->getPMT(), __p->getToT());
215  }
216 
217  t0[__p->getPMT()] = __p->getT();
218  }
219  }
220 
221  p = q;
222  }
223  }
224  }
225  }
226 
227  STATUS(endl);
228 
229  //---------------------------------------------
230  // store histograms
231  //---------------------------------------------
232 
233  TFile out(outputFile.c_str(), "recreate");
234 
235  putObject(&out, JMeta(argc, argv));
236 
237  out << manager;
238 
239  out.Close();
240  out.Write();
241 }
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:71
Utility class to parse command line options.
Definition: JParser.hh:1410
Data structure for all trigger parameters.
double getDot(const JNeutrinoDirection &first, const JNeutrinoDirection &second)
Dot product.
Definition: JAstronomy.hh:409
Auxiliary class to manage set of histograms.
#define STATUS(A)
Definition: JMessage.hh:61
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:611
Long64_t counter_type
Type definition for counter.
Definition: JCounter.hh:24
JSuperFrame2D< hit_type > JSuperFrame2D_t
Definition: JDataFilter.cc:93
string outputFile
JLimit JLimit_t
Type definition of limit.
Definition: JLimit.hh:214
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1836
Data time slice.
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
int debug
debug level
Definition: JSirene.cc:59
#define FATAL(A)
Definition: JMessage.hh:65
static const char *const _2SToT
Histogram naming.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:72
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
JTriggerParameters getTriggerParameters(const JMultipleFileScanner_t &file_list)
Get trigger parameters.
double hit_type
Definition: JDataFilter.cc:89
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:60
JSuperFrame1D< hit_type > JSuperFrame1D_t
Definition: JDataFilter.cc:92
bool putObject(TDirectory *dir, const T &object)
Write object to ROOT directory.