Jpp - the software that should make you happy
 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 "km3net-dataformat/online/JDAQ.hh"
#include "JDAQ/JDAQTimesliceIO.hh"
#include "JDetector/JDetector.hh"
#include "JDetector/JDetectorToolkit.hh"
#include "JDetector/JModuleRouter.hh"
#include "JTrigger/JHitR0.hh"
#include "JTrigger/JMatchL0.hh"
#include "JTrigger/JSuperFrame1D.hh"
#include "JTrigger/JSuperFrame2D.hh"
#include "JTrigger/JPreprocessor.hh"
#include "JSupport/JMultipleFileScanner.hh"
#include "JSupport/JSupport.hh"
#include "JSupport/JMeta.hh"
#include "JROOT/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 44 of file JCalibrateToT.cc.

45 {
46  using namespace std;
47  using namespace JPP;
48  using namespace KM3NETDAQ;
49 
51  JLimit_t& numberOfEvents = inputFile.getLimit();
52  string outputFile;
53  string detectorFile;
54  JRange<double> T_ns;
55  double ctMax;
56  JRange<int> multiplicity;
57  double deadTime_us;
58  JROOTClassSelector selector;
59  JPreprocessor option;
60  int debug;
61 
62  try {
63 
64  JParser<> zap("Monitoring of PMT time-over-threshold distributions.");
65 
66  zap['f'] = make_field(inputFile, "input file.");
67  zap['o'] = make_field(outputFile, "output file.") = "calibrate_tot.root";
68  zap['n'] = make_field(numberOfEvents) = JLimit::max();
69  zap['a'] = make_field(detectorFile, "detector file.");
70  zap['T'] = make_field(T_ns, "time window [ns].") = JRange<double>(10.0, 25.0);
71  zap['c'] = make_field(ctMax, "maximal cosine space angle between PMT axes.") = 0.0;
72  zap['M'] = make_field(multiplicity, "multiplicity range of hits on DOM.") = JRange<int>(1, 2);
73  zap['D'] = make_field(deadTime_us, "L0 dead time (us)") = 10;
74  zap['C'] = make_field(selector, "timeslice selector, e.g. JDAQTimesliceL1.") = getROOTClassSelection<JDAQTimesliceTypes_t>();
75  zap['O'] = make_field(option, "hit pre-processing option.") = JPreprocessor::getOptions();
76  zap['d'] = make_field(debug, "debug flag.") = 1;
77 
78  zap(argc, argv);
79  }
80  catch(const exception &error) {
81  FATAL(error.what() << endl);
82  }
83 
84  cout.tie(&cerr);
85 
86  //-----------------------------------------------------------
87  // check the input parameters
88  //-----------------------------------------------------------
89 
90  if (!T_ns.is_valid()) {
91  FATAL("Invalid time window [ns] " << T_ns << endl);
92  }
93 
94  if (selector == JDAQTimeslice ::Class_Name() ||
95  selector == JDAQTimesliceL1::Class_Name()) {
96 
98 
99  try {
100  parameters = getTriggerParameters(inputFile);
101  }
102  catch(const JException& error) {
103  FATAL("No trigger parameters from input." << endl);
104  }
105 
106  if ((selector == JDAQTimeslice ::Class_Name() && parameters.writeL1.prescale > 0) ||
107  (selector == JDAQTimesliceL1::Class_Name())) {
108 
109  if (parameters.TMaxLocal_ns < T_ns.getUpperLimit()) {
110  FATAL("Option -T <T_ns> = " << T_ns.getUpperLimit() << " is larger than in the trigger " << parameters.TMaxLocal_ns << endl);
111  }
112  }
113  }
114 
115  if (!multiplicity.is_valid()) { FATAL("Invalid multiplicity " << multiplicity << endl); }
116  if ( multiplicity.getLowerLimit() < 1) { FATAL("Invalid multiplicity " << multiplicity << endl); }
117 
118  //-----------------------------------------------------------
119  // load detector file
120  //-----------------------------------------------------------
121 
123 
124  try {
125  load(detectorFile, detector);
126  }
127  catch(const JException& error) {
128  FATAL(error);
129  }
130 
131  if (detector.empty()) {
132  FATAL("Empty detector." << endl);
133  }
134 
135  const JModuleRouter router(detector);
136 
137  //-----------------------------------------------------------
138  // initialise histograms
139  //-----------------------------------------------------------
140 
141  const double zmin = -0.5; // [ns]
142  const double zmax = 256.5; // [ns]
143  const int nz = (int) ((zmax-zmin) / 1.0);
144 
145  JManager<int, TH2D> manager(new TH2D(MAKE_CSTRING("%" << _2SToT), NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS-0.5, nz, zmin, zmax));
146 
147  typedef JHitR0 hit_type;
148  typedef JSuperFrame2D<hit_type> JSuperFrame2D_t;
149  typedef JSuperFrame1D<hit_type> JSuperFrame1D_t;
150 
151  const JMatchL0<hit_type> match(T_ns.getUpperLimit()); // time window self-coincidences [ns]
152 
153  const double deadTime_ns = deadTime_us * 1e3;
154 
156 
157  counter_type counter = 0;
158 
159  for ( ; in.hasNext() && counter != inputFile.getLimit(); ++counter) {
160 
161  STATUS("event: " << setw(10) << counter << '\r'); DEBUG(endl);
162 
163  const JDAQTimeslice* timeslice = in.next();
164 
165  for (JDAQTimeslice::const_iterator frame = timeslice->begin(); frame != timeslice->end(); ++frame) {
166 
167  if (router.hasModule(frame->getModuleID())) {
168 
169  TH2D* h2 = manager[frame->getModuleID()];
170  const JModule& module = router.getModule(frame->getModuleID());
171 
172  JSuperFrame2D_t& buffer = JSuperFrame2D_t::demultiplex(*frame, module);
173 
174  buffer.preprocess(option, match);
175 
176  JSuperFrame1D_t& data = JSuperFrame1D_t::multiplex(buffer);
177 
178  vector<double> t0(NUMBER_OF_PMTS, numeric_limits<double>::lowest());
179 
180  for (JSuperFrame1D_t::const_iterator p = data.begin(), __end = data.end() - 1; p != __end; ) {
181 
182  JSuperFrame1D_t::const_iterator q = p;
183 
184  double ct_max = -1.0;
185  double dt_min = numeric_limits<double>::max();
186 
187  while (++q != __end && q->getT() - p->getT() < T_ns.getUpperLimit()) {
188 
189  const double ct = getDot(module.getPMT(p->getPMT()), module.getPMT(q->getPMT()));
190  const double dt = q->getT() - p->getT();
191 
192  if (ct > ct_max) {
193  ct_max = ct;
194  }
195 
196  if (dt < dt_min) {
197  dt_min = dt;
198  }
199  }
200 
201  if (multiplicity(distance(p,q)) && ct_max < ctMax && dt_min > T_ns.getLowerLimit()) {
202 
203  for (JSuperFrame1D_t::const_iterator __p = p; __p != q; ++__p) {
204 
205  if (__p->getT() > t0[__p->getPMT()] + deadTime_ns) {
206  h2->Fill(__p->getPMT(), __p->getToT());
207  }
208 
209  t0[__p->getPMT()] = __p->getT();
210  }
211  }
212 
213  p = q;
214  }
215  }
216  }
217  }
218 
219  STATUS(endl);
220 
221  //---------------------------------------------
222  // store histograms
223  //---------------------------------------------
224 
225  TFile out(outputFile.c_str(), "recreate");
226 
227  putObject(&out, JMeta(argc, argv));
228 
229  out << manager;
230 
231  out.Close();
232  out.Write();
233 }
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:70
Utility class to parse command line options.
Definition: JParser.hh:1500
General exception.
Definition: JException.hh:23
Data structure for a composite optical module.
Definition: JModule.hh:57
std::vector< T >::difference_type distance(typename std::vector< T >::const_iterator first, typename PhysicsEvent::const_iterator< T > second)
Specialisation of STL distance.
double getDot(const JNeutrinoDirection &first, const JNeutrinoDirection &second)
Dot product.
Definition: JAstronomy.hh:409
L0 match criterion.
Definition: JMatchL0.hh:27
#define STATUS(A)
Definition: JMessage.hh:63
Detector data structure.
Definition: JDetector.hh:80
Auxiliary class to select ROOT class based on class name.
Router for direct addressing of module data in detector data structure.
*fatal Wrong number of arguments esac JCookie sh typeset Z DETECTOR typeset Z SOURCE_RUN typeset Z TARGET_RUN set_variable PARAMETERS_FILE $WORKDIR parameters
Definition: diff-Tuna.sh:38
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:151
Long64_t counter_type
Type definition for counter.
Auxiliary class for multiplexing object iterators.
string outputFile
1-dimensional frame with time calibrated data from one optical module.
Auxiliary class for defining the range of iterations of objects.
Definition: JLimit.hh:41
Auxiliary class to manage set of compatible ROOT objects (e.g. histograms) using unique keys...
Definition: JManager.hh:43
Detector file.
Definition: JHead.hh:196
Reduced data structure for L0 hit.
Definition: JHitR0.hh:25
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
Data time slice.
int debug
debug level
Definition: JSirene.cc:63
const JPMT & getPMT(const int index) const
Get PMT.
Definition: JModule.hh:211
#define FATAL(A)
Definition: JMessage.hh:67
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
General purpose class for object reading from a list of file names.
static const char *const _2SToT
Histogram naming.
2-dimensional frame with time calibrated data from one optical module.
const JLimit & getLimit() const
Get limit.
Definition: JLimit.hh:73
do set_variable DETECTOR_TXT $WORKDIR detector
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable INPUT_FILE $argv[2] eval JPrintDetector a $DETECTOR O IDENTIFIER eval JPrintDetector a $DETECTOR O SUMMARY source JAcoustics sh $DETECTOR_ID CHECK_EXIT_CODE typeset A TRIPODS get_tripods $WORKDIR tripod txt TRIPODS for EMITTER in
Definition: JCanberra.sh:38
Auxiliary class for specifying the way of pre-processing of hits.
JTriggerParameters getTriggerParameters(const JMultipleFileScanner_t &file_list)
Get trigger parameters.
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
bool putObject(TDirectory *dir, const T &object)
Write object to ROOT directory.