Jpp
JFitK40.cc
Go to the documentation of this file.
1 
2 #include <string>
3 #include <iostream>
4 #include <iomanip>
5 #include <map>
6 #include <cmath>
7 
8 #include "TROOT.h"
9 #include "TFile.h"
10 #include "TH1D.h"
11 #include "TH2D.h"
12 #include "TF2.h"
13 #include "TMath.h"
14 #include "TFitResult.h"
15 #include "Math/MinimizerOptions.h"
16 
17 #include "JTools/JConstants.hh"
18 #include "JDAQ/JDAQ.hh"
19 #include "JDetector/JDetector.hh"
22 
23 #include "JROOT/JRootToolkit.hh"
24 #include "JTools/JRange.hh"
25 #include "JTools/JQuantile.hh"
26 #include "JSupport/JMeta.hh"
27 
29 #include "JCalibrate/JFitK40.hh"
30 
31 #include "Jeep/JProperties.hh"
32 #include "Jeep/JPrint.hh"
33 #include "Jeep/JParser.hh"
34 #include "Jeep/JMessage.hh"
35 
36 namespace {
37 
38  double STDEV = 2.0; //!< number of standard deviations
39  double QE_MIN_VALUE = 0.01; //!< minimal QE value (below, set to zero)
40  double QE_MAX_ERROR = 0.5; //!< maximal QE error (below, set to zero)
41  double MINIMAL_RATE_HZ = 1.0e-2; //!< minimal coincidence rate (bin-by-bin)
42 
43  /**
44  * Check validity of QE measurement.
45  *
46  * The QE measurement is valid when:
47  * - fitted value is at least STDEV standard deviations larger than QE_MIN_VALUE; and
48  * - fitted error is smaller than QE_MAX_ERROR.
49  *
50  * \param qe_value QE value
51  * \param qe_error QE error
52  * \return true if valid; else false
53  */
54  inline bool is_valid(const double qe_value, const double qe_error)
55  {
56  return (qe_value > QE_MIN_VALUE + STDEV * qe_error &&
57  qe_error < QE_MAX_ERROR);
58  }
59 }
60 
61 
62 /**
63  * \file
64  *
65  * Auxiliary program to fit PMT parameters from JMergeCalibrateK40.cc output.
66  * \author mdejong
67  */
68 int main(int argc, char **argv)
69 {
70  using namespace std;
71  using namespace JPP;
72  using namespace KM3NETDAQ;
73 
74  typedef JRange<double> JRange_t;
75  typedef multimap<int, int> JTDC_t;
77 
78  const double epsilon = 1.0e-10; // minimal error for ROOT plotting
79  JFitK40Parameters fitk40; // setting of internal fit parameters
80 
81  string inputFile;
82  string outputFile;
83  string detectorFile;
84  string pmtFile;
85  JTDC_t TDC;
86  bool reverse;
87  bool overwriteDetector;
88  bool writeFits;
89  bool fitAngDist;
90  bool fitModel;
91  double mu;
92  JRange_t X;
93  JRange_t Y;
94  string option;
95  int debug;
96 
97  try {
98 
99  JProperties properties;
100 
101  properties.insert(gmake_property(STDEV));
102  properties.insert(gmake_property(QE_MIN_VALUE));
103  properties.insert(gmake_property(QE_MAX_ERROR));
104  properties.insert(gmake_property(MINIMAL_RATE_HZ));
105  properties.insert(gmake_property(fitk40.Rate_Hz));
106  properties.insert(gmake_property(fitk40.p1));
107  properties.insert(gmake_property(fitk40.p2));
108  properties.insert(gmake_property(fitk40.p3));
109  properties.insert(gmake_property(fitk40.p4));
110 
111  JParser<> zap("Auxiliary program to fit PMT parameters from JMergeCalibrateK40 output.");
112 
113  zap['@'] = make_field(properties) = JPARSER::initialised();
114  zap['f'] = make_field(inputFile, "input file (output from JMergeCalibrateK40).");
115  zap['o'] = make_field(outputFile, "output file.") = "fit.root";
116  zap['a'] = make_field(detectorFile, "detector file.");
117  zap['P'] = make_field(pmtFile, "specify PMT file name that can be given to JTriggerEfficiency.") = "";
118  zap['!'] = make_field(TDC,
119  "Fix time offset(s) of PMT(s) of certain module(s), e.g."
120  "\n-! \"808969848 0 808982077 23\" will fix offset of PMT 0 of module 808969848 and of PMT 23 of module 808982077."
121  "\nSame PMT offset can be fixed for all optical modules, e.g."
122  "\n-! \"-1 0 -1 23\" will fix offsets of PMTs 0 and 23 for all optical modules.")
124  zap['r'] = make_field(reverse, "reverse TDC constraints due to option -! <TDC>.");
125  zap['A'] = make_field(overwriteDetector, "overwrite detector file provided through '-a' with correct time offsets.");
126  zap['w'] = make_field(writeFits, "write fit results.");
127  zap['D'] = make_field(fitAngDist, "fit angular distribution; fix normalisation.");
128  zap['M'] = make_field(fitModel, "fit angular distribution; fix PMT QEs = 1.0.");
129  zap['E'] = make_field(mu, "expectation value for npe given two-fold coincidence") = 0.0;
130  zap['x'] = make_field(X, "ROOT fit range (PMT pairs).") = JRange_t();
131  zap['y'] = make_field(Y, "ROOT fit range (time residual).") = JRange_t();
132  zap['O'] = make_field(option, "ROOT fit option, see TH1::Fit.") = "";
133  zap['d'] = make_field(debug, "debug.") = 1;
134 
135  zap(argc, argv);
136  }
137  catch(const exception &error) {
138  FATAL(error.what() << endl);
139  }
140 
141 
142  ROOT::Math::MinimizerOptions::SetDefaultMaxFunctionCalls(1000000);
143 
144  //----------------------------------------------------------
145  // load detector file and PMT parameters
146  //----------------------------------------------------------
147 
148  JDetector detector;
149 
150  try {
151  load(detectorFile, detector);
152  }
153  catch(const JException& error) {
154  FATAL(error);
155  }
156 
157 
158  JPMTParametersMap parameters;
159 
160  if (pmtFile != "") {
161  try {
162  parameters.load(pmtFile.c_str());
163  }
164  catch(const JException& error) {
165  //ERROR(error.what());
166  }
167  }
168 
169  parameters.comment.add(JMeta(argc, argv));
170 
171  {
172  //----------------------------------------------------------
173  // module identifier -1 corresponds to any module
174  //----------------------------------------------------------
175 
176  const range_type range = TDC.equal_range(-1); // module wild card
177 
178  if (range.first != range.second) {
179 
180  const JTDC_t buffer(range.first, range.second); // temporary copy
181 
182  TDC.erase(range.first, range.second); // remove entries
183 
184  for (JDetector::iterator module = detector.begin(); module != detector.end(); ++module) {
185 
186  for (JTDC_t::const_iterator i = buffer.begin(); i != buffer.end(); ++i) {
187 
188  //----------------------------------------------------------
189  // PMT identifier -1 corresponds to any PMT in module
190  //----------------------------------------------------------
191 
192  if (i->second == -1) {
193 
194  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
195  TDC.insert(make_pair(module->getID(), pmt));
196  }
197 
198  } else {
199 
200  TDC.insert(make_pair(module->getID(), i->second));
201  }
202  }
203  }
204  }
205  }
206 
207  if (reverse) {
208 
209  //----------------------------------------------------------
210  // reverse TDC constraints
211  //----------------------------------------------------------
212 
213  JTDC_t buffer;
214 
215  for (JTDC_t::const_iterator __p = TDC.begin(); __p != TDC.end(); ) {
216 
217  JTDC_t::const_iterator __q = __p;
218 
219  for ( ; __q != TDC.end() && __q->first == __p->first; ++__q) {}
220 
221  for (int i = 0; i != NUMBER_OF_PMTS; ++i) {
222 
223  JTDC_t::const_iterator __i = __p;
224 
225  for ( ; __i != __q && __i->second != i; ++__i) {}
226 
227  if (__i == __q) {
228  buffer.insert(std::make_pair(__p->first,i));
229  }
230  }
231 
232  __p = __q;
233  }
234 
235  TDC.swap(buffer);
236  }
237 
238  for (JTDC_t::const_iterator tdc = TDC.begin(); tdc != TDC.end(); ++tdc) {
239  DEBUG("PMT " << setw(10) << tdc->first << ' ' << setw(2) << tdc->second << " constrain t0." << endl);
240  }
241 
242  for (JTDC_t::const_iterator tdc = TDC.begin(); tdc != TDC.end(); ++tdc) {
243  if (tdc->first < 0) {
244  FATAL("Illegal module identifier: " << tdc->first << endl);
245  }
246  if (tdc->second < 0 || tdc->second >= NUMBER_OF_PMTS) {
247  FATAL("Illegal TDC (= readout channel) identifier: " << tdc->second << " {0, .., " << NUMBER_OF_PMTS - 1 << "}" << endl);
248  }
249  }
250 
251 
252  if (option.find('O') == string::npos) { option += '0'; }
253  if (option.find('S') == string::npos) { option += 'S'; }
254  if (option.find('Q') == string::npos && debug < JEEP::debug_t) { option += 'Q'; }
255 
256 
257  TFile* in = TFile::Open(inputFile.c_str(), "exist");
258 
259  if (in == NULL || !in->IsOpen()) {
260  FATAL("File: " << inputFile << " not opened." << endl);
261  }
262 
263 
264  // Histograms with global fit results.
265  TH1D hc("chi2", NULL, 500, 0.0, 5.0);
266 
267  TH1D p1("p1", NULL, 501, -5.0, +5.0);
268  TH1D p2("p2", NULL, 501, -5.0, +5.0);
269  TH1D p3("p3", NULL, 501, -5.0, +5.0);
270  TH1D p4("p4", NULL, 501, -5.0, +5.0);
271  TH1D bg("bg", NULL, 501, -0.1, +0.1);
272  TH1D cc("cc", NULL, 501, -0.1, +0.1);
273 
274 
275  TFile out(outputFile.c_str(), "recreate");
276 
277  //----------------------------------------------------------
278  // loop over modules in detector file to fit histogram
279  //----------------------------------------------------------
280 
281  for (JDetector::iterator module = detector.begin(); module != detector.end(); ++module) {
282 
283  NOTICE("Module " << setw(10) << module->getID() << endl);
284 
285  //----------------------------------------------------------
286  // get histogram for this module
287  //----------------------------------------------------------
288 
289  TH2D* h2s = (TH2D*) in->Get(MAKE_CSTRING(module->getID() << _2R));
290 
291  if (h2s == NULL || h2s->GetEntries() == 0) {
292 
293  NOTICE("No data for module " << module->getID() << "; set corresponding QEs to 0." << endl);
294 
295  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
296  parameters[JPMTIdentifier(module->getID(), pmt)].QE = 0.0;
297  }
298 
299  continue;
300  }
301 
302  //----------------------------------------------------------
303  // get constraints
304  //----------------------------------------------------------
305 
306  const range_type range = TDC.equal_range(module->getID());
307 
308  //----------------------------------------------------------
309  // set-up fitting procedure
310  //----------------------------------------------------------
311 
312  JFitK40 fit(*module,
313  h2s->GetXaxis()->GetXmin(),
314  h2s->GetXaxis()->GetXmax(),
315  h2s->GetYaxis()->GetXmin(),
316  h2s->GetYaxis()->GetXmax(),
317  range.first == range.second);
318 
319  fit.setModelParameters(fitk40.getModelParameters());
320 
321  for (JTDC_t::const_iterator i = range.first; i != range.second; ++i) {
322  fixParameter(fit, fit.getModelParameter(i->second, &JPMTParameters_t::t0));
323  }
324 
325  if (fitModel) {
326  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
327  fixParameter(fit, fit.getModelParameter(pmt, &JPMTParameters_t::QE));
328  }
329 
330  } else {
331 
332  fixParameter(fit, fit.getModelParameter(&JFitK40::Rate_Hz));
333 
334  if (!fitAngDist) {
335  fixParameter(fit, fit.getModelParameter(&JFitK40::p1));
336  fixParameter(fit, fit.getModelParameter(&JFitK40::p2));
337  fixParameter(fit, fit.getModelParameter(&JFitK40::p3));
338  fixParameter(fit, fit.getModelParameter(&JFitK40::p4));
339  }
340  }
341 
342  //----------------------------------------------------------
343  // check validity of data
344  //----------------------------------------------------------
345 
346  vector<int> buffer(NUMBER_OF_PMTS, 1);
347 
348  for (int ix = 1; ix <= h2s->GetXaxis()->GetNbins(); ++ix) {
349 
350  Double_t value = 0.0;
351  Double_t error = 0.0;
352 
353  for (int iy = 1; iy <= h2s->GetYaxis()->GetNbins(); ++iy) {
354  value += h2s->GetBinContent(ix,iy);
355  error += h2s->GetBinError(ix,iy) * h2s->GetBinError(ix,iy);
356  }
357 
358  error = sqrt(error);
359 
360  if (value < MINIMAL_RATE_HZ + STDEV*error) {
361 
362  const JFitK40::pair_type pair = fit.getPair(ix-1);
363 
364  buffer[pair.first] += 1;
365  buffer[pair.second] += 1;
366  }
367  }
368 
369  try {
370 
371  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
372 
373  if (buffer[pmt] == NUMBER_OF_PMTS) {
374 
375  WARNING("PMT " << setw(10) << module->getID() << ' ' << setw(2) << pmt << " too low a rate; switch off." << endl);
376 
377  fit.disablePMT(pmt);
378  }
379  }
380  }
381  catch(const JException& error) {
382 
383  NOTICE("Skip fit for module " << module->getID() << "; set corresponding QEs to 0." << endl);
384 
385  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
386  parameters[JPMTIdentifier(module->getID(), pmt)].QE = 0.0;
387  }
388 
389  continue;
390  }
391 
392  // constrain fit to specified range(s)
393 
394  if (X != JRange_t()) {
395  h2s->GetXaxis()->SetRangeUser(max(X.getLowerLimit(), h2s->GetXaxis()->GetXmin()),
396  min(X.getUpperLimit(), h2s->GetXaxis()->GetXmax()));
397  }
398  if (Y != JRange_t()) {
399  h2s->GetYaxis()->SetRangeUser(max(Y.getLowerLimit(), h2s->GetYaxis()->GetXmin()),
400  min(Y.getUpperLimit(), h2s->GetYaxis()->GetXmax()));
401  }
402 
403 
404  TFitResultPtr result = fit(*h2s, option);
405 
406 
407  bool refit = false;
408 
409  try {
410 
411  const JFitK40Parameters values(fit.GetParameters());
412  const JFitK40Parameters errors(fit.GetParErrors());
413 
414  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
415 
416  if (!result->IsParameterFixed(fit.getModelParameter(pmt, &JPMTParameters_t::QE)) && !is_valid(values.getQE(pmt), errors.getQE(pmt))) {
417 
418  WARNING("PMT " << setw(10) << module->getID() << ' ' << setw(2) << pmt << ' '
419  << "measured QE "
420  << FIXED(5,2) << values.getQE(pmt) << " +/- "
421  << FIXED(5,2) << errors.getQE(pmt)
422  << " compatible with zero -> set QE = 0 and t0 = 0; -> refit." << endl);
423 
424 
425  fit.disablePMT(pmt);
426 
427  refit = true;
428  }
429  }
430  }
431  catch(const JException& error) {
432 
433  NOTICE("Skip fit for module " << module->getID() << "; set corresponding QEs to 0." << endl);
434 
435  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
436  parameters[JPMTIdentifier(module->getID(), pmt)].QE = 0.0;
437  }
438 
439  continue;
440  }
441 
442 
443  if (refit) {
444  result = fit(*h2s, option);
445  }
446 
447 
448  const JFitK40Parameters values(fit.GetParameters());
449  const JFitK40Parameters errors(fit.GetParErrors());
450 
451  hc.Fill(TMath::Log10(result->Chi2() / result->Ndf()));
452 
453  // To force drawing of fit functions bin-by-bin compatible to the corresponding histogram data,
454  // the function values will be stored as histograms, rather than as associated function objects.
455 
456  TH2D h2f(MAKE_CSTRING(module->getID() << _2F),
457  NULL,
458  h2s->GetXaxis()->GetNbins(), h2s->GetXaxis()->GetXmin (), h2s->GetXaxis()->GetXmax (),
459  h2s->GetYaxis()->GetNbins(), h2s->GetYaxis()->GetXmin (), h2s->GetYaxis()->GetXmax ());
460 
461  for (int ix = 1; ix <= h2f.GetXaxis()->GetNbins(); ++ix) {
462  for (int iy = 1; iy <= h2f.GetYaxis()->GetNbins(); ++iy) {
463 
464  const Double_t x = h2f.GetXaxis()->GetBinCenter(ix);
465  const Double_t y = h2f.GetYaxis()->GetBinCenter(iy);
466 
467  h2f.SetBinContent(ix, iy, fit.Eval(x,y));
468  h2f.SetBinError (ix, iy, 0.0);
469  }
470  }
471 
472 
473  if (debug >= JEEP::notice_t) {
474 
475  cout << "Fit result chi2 / NDF " << result->Chi2() << ' ' << result->Ndf() << ' ' << (result->IsValid() ? "" : "failed") << endl;
476 
477  cout << "Rate_Hz= " << FIXED(8,4) << values.Rate_Hz << (result->IsParameterFixed(fit.getModelParameter(&JFitK40::Rate_Hz)) ? " *** fixed *** " : "") << endl;
478 
479  cout << "p1= " << FIXED(8,4) << values.p1 << (result->IsParameterFixed(fit.getModelParameter(&JFitK40::p1)) ? " *** fixed *** " : "") << endl;
480  cout << "p2= " << FIXED(8,4) << values.p2 << (result->IsParameterFixed(fit.getModelParameter(&JFitK40::p2)) ? " *** fixed *** " : "") << endl;
481  cout << "p3= " << FIXED(8,4) << values.p3 << (result->IsParameterFixed(fit.getModelParameter(&JFitK40::p2)) ? " *** fixed *** " : "") << endl;
482  cout << "p4= " << FIXED(8,4) << values.p4 << (result->IsParameterFixed(fit.getModelParameter(&JFitK40::p3)) ? " *** fixed *** " : "") << endl;
483 
484  cout << "Background (correlated): " << FIXED(8,5) << values.bg << (result->IsParameterFixed(fit.getModelParameter(&JFitK40::bg)) ? " *** fixed *** " : "") << endl;
485  cout << "Background (uncorrelated): " << FIXED(8,5) << values.cc << (result->IsParameterFixed(fit.getModelParameter(&JFitK40::cc)) ? " *** fixed *** " : "") << endl;
486 
487  cout << " PMT t0 [ns] TTS [ns] R" << endl;
488 
489  JQuantile Q[3];
490 
491  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
492 
493  cout << ' '
494  << setw(2) << pmt << ' '
495  << FIXED(7,2) << fit.getT0 (pmt) << ' '
496  << FIXED(7,2) << fit.getTTS(pmt) << ' '
497  << FIXED(7,3) << fit.getQE (pmt);
498 
499  cout << (result->IsParameterFixed(fit.getModelParameter(pmt, &JPMTParameters_t::QE)) ? " *** fixed QE *** " : "");
500  cout << (result->IsParameterFixed(fit.getModelParameter(pmt, &JPMTParameters_t::TTS)) ? " *** fixed TTS *** " : "");
501  cout << (result->IsParameterFixed(fit.getModelParameter(pmt, &JPMTParameters_t::t0)) ? " *** fixed t0 *** " : "");
502  cout << (fit.is_disabled (pmt) ? " (disabled)" : "");
503  cout << (fit.is_average_t0(pmt) ? " (averaged)" : "");
504  cout << endl;
505 
506  Q[0].put(fit.getT0 (pmt));
507  Q[1].put(fit.getTTS(pmt));
508  Q[2].put(fit.getQE (pmt));
509  }
510 
511  cout << "<> "
512  << FIXED(7,2) << Q[0].getMean() << ' '
513  << FIXED(7,2) << Q[1].getMean() << ' '
514  << FIXED(7,3) << Q[2].getMean() << endl;
515  cout << "+/- "
516  << FIXED(7,2) << Q[0].getSTDev() << ' '
517  << FIXED(7,2) << Q[1].getSTDev() << ' '
518  << FIXED(7,3) << Q[2].getSTDev() << endl;
519  }
520 
521  h2s->Write();
522  h2f .Write();
523 
524  if (writeFits) {
525 
526  // Histograms with fit results.
527 
528  p1.Fill(values.p1);
529  p2.Fill(values.p2);
530  p3.Fill(values.p3);
531  p4.Fill(values.p4);
532  bg.Fill(values.bg);
533  cc.Fill(values.cc);
534 
535  TH1D h1t(MAKE_CSTRING(module->getID() << ".1t0"), NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS - 0.5);
536  TH1D h1s(MAKE_CSTRING(module->getID() << ".1TTS"), NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS - 0.5);
537  TH1D h1q(MAKE_CSTRING(module->getID() << ".1QE"), NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS - 0.5);
538 
539  h1t.Sumw2();
540  h1s.Sumw2();
541  h1q.Sumw2();
542 
543  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
544 
545  h1t.SetBinContent(pmt + 1, values.getT0 (pmt));
546  h1s.SetBinContent(pmt + 1, values.getTTS(pmt));
547  h1q.SetBinContent(pmt + 1, values.getQE (pmt));
548 
549  h1t.SetBinError (pmt + 1, max(errors.getT0 (pmt), epsilon));
550  h1s.SetBinError (pmt + 1, max(errors.getTTS(pmt), epsilon));
551  h1q.SetBinError (pmt + 1, max(errors.getQE (pmt), epsilon));
552  }
553 
554  out << h1t << h1s << h1q;
555  }
556 
557 
558  // save output
559 
560  for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
561 
562  parameters[JPMTIdentifier(module->getID(), pmt)].QE = values.getQE(pmt);
563 
564  if (overwriteDetector) {
565  module->getPMT(pmt).addT0(fit.getT0(pmt));
566  }
567  }
568  }
569 
570 
571  if (writeFits) {
572  out << hc << p1 << p2 << p3 << p4 << bg << cc;
573  }
574 
575  out.Close();
576 
577 
578  if (overwriteDetector) {
579 
580  NOTICE("Store calibration data on file " << detectorFile << endl);
581 
582  detector.comment.add(JMeta(argc, argv));
583 
584  store(detectorFile, detector);
585  }
586 
587  if (pmtFile != "") {
588 
589  //----------------------------------------------------------
590  // at this stage, fitted QE corresponds to hit probability
591  //----------------------------------------------------------
592 
593  try {
594  parameters.convertHitProbabilityToQE(mu);
595  }
596  catch(const JException& error) {
597  FATAL(error.what());
598  }
599 
600  ofstream out(pmtFile.c_str());
601 
602  out << parameters << endl;
603 
604  out.close();
605  }
606 
607  return 0;
608 }
JPMTParametersMap.hh
JMeta.hh
JDAQ.hh
FIXED
Auxiliary data structure for floating point format specification.
Definition: JPrint.hh:481
JMessage.hh
JPrint.hh
JCALIBRATE::_2F
static const char *const _2F
Name extension for 2F rate.
Definition: JCalibrateK40.hh:37
JEEP::notice_t
notice
Definition: JMessage.hh:32
JAANET::is_valid
bool is_valid(const T &value)
Check validity of given value.
Definition: JHead.hh:823
JPARSER::initialised
Empty structure for specification of parser element that is initialised (i.e.
Definition: JParser.hh:63
JDETECTOR::load
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
Definition: JDetectorToolkit.hh:456
JFitK40.hh
std::vector< int >
KM3NETDAQ::NUMBER_OF_PMTS
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
JPARSER::JParser
Utility class to parse command line options.
Definition: JParser.hh:1493
NOTICE
#define NOTICE(A)
Definition: JMessage.hh:64
JQuantile.hh
main
int main(int argc, char **argv)
Definition: JFitK40.cc:68
JPP
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Definition: JAAnetToolkit.hh:37
MAKE_CSTRING
#define MAKE_CSTRING(A)
Make C-string.
Definition: JPrint.hh:708
JRange.hh
WARNING
#define WARNING(A)
Definition: JMessage.hh:65
debug
int debug
debug level
Definition: JSirene.cc:59
JConstants.hh
JTOOLS::result
return result
Definition: JPolint.hh:695
JRootToolkit.hh
std::pair
Definition: JSTDTypes.hh:15
JParser.hh
JDetectorToolkit.hh
JDETECTOR::store
void store(const JString &file_name, const JDetector &detector)
Store detector to output file.
Definition: JDetectorToolkit.hh:532
std::multimap
Definition: JSTDTypes.hh:17
JCALIBRATE::pair_type
JCombinatorics::pair_type pair_type
Definition: JCalibrateK40.hh:39
make_field
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1954
gmake_property
#define gmake_property(A)
macro to convert (template) parameter to JPropertiesElement object
Definition: JProperties.hh:1319
JEEP::JProperties
Utility class to parse parameter values.
Definition: JProperties.hh:496
JROOT::fixParameter
bool fixParameter(TF1 &f1, const JFitParameter_t &parameter)
Fix fit parameter.
Definition: JRootToolkit.hh:249
DEBUG
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
std
Definition: jaanetDictionary.h:36
KM3NETDAQ
KM3NeT DAQ data structures and auxiliaries.
Definition: DataQueue.cc:39
JProperties.hh
p1
TPaveText * p1
Definition: JDrawModule3D.cc:35
JCalibrateK40.hh
JDetector.hh
JCALIBRATE::_2R
static const char *const _2R
Name extension for 2D rate.
Definition: JCalibrateK40.hh:36
FATAL
#define FATAL(A)
Definition: JMessage.hh:67
outputFile
string outputFile
Definition: JDAQTimesliceSelector.cc:37
JEEP::debug_t
debug
Definition: JMessage.hh:29