Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JFitToT.cc
Go to the documentation of this file.
1 #include <string>
2 #include <iostream>
3 #include <iomanip>
4 
5 #include "TROOT.h"
6 #include "TFile.h"
7 #include "TH1D.h"
8 #include "TH2D.h"
9 #include "TF1.h"
10 #include "TFitResult.h"
11 #include "TMath.h"
12 
13 #include "JLang/JLangToolkit.hh"
14 #include "JTools/JConstants.hh"
16 
17 #include "JDetector/JDetector.hh"
21 
22 #include "JROOT/JRootToolkit.hh"
23 #include "JTools/JRange.hh"
24 #include "JSupport/JMeta.hh"
25 
27 #include "JCalibrate/JFitToT.hh"
28 
29 #include "Jeep/JPrint.hh"
30 #include "Jeep/JParser.hh"
31 #include "Jeep/JMessage.hh"
32 
33 namespace {
34 
35  /**
36  * Wild card for histogram naming.\n
37  * The wild card will replaced by the module identifier.
38  */
39  static const char* const WILDCARD = "%";
40 }
41 
42 
43 /**
44  * \file
45  *
46  * Auxiliary program to fit time-over-threshold distributions.
47  * The fitted parameters are the PMT gain and gain spread.
48  * \author mkarel
49  */
50 int main(int argc, char **argv)
51 {
52  using namespace std;
53  using namespace JPP;
54  using namespace KM3NETDAQ;
55 
56  typedef JRange<double> JRange_t;
57 
58  string inputFile;
59  string outputFile;
60  string detectorFile;
61  string pmtFile;
62  bool writeFits;
63  JRange_t fitRange;
64  string option;
65  double WMin;
66  string regexp;
67  int debug;
68 
69  try {
70 
71  JParser<> zap("Auxiliary program to fit time-over-threshold distributions.");
72 
73  zap['f'] = make_field(inputFile, "input file (output from JCalibrateToT).");
74  zap['o'] = make_field(outputFile, "output file.") = "fit.root";
75  zap['a'] = make_field(detectorFile, "detector file.");
76  zap['P'] = make_field(pmtFile, "specify PMT file name that can be given to JTriggerEfficiency.") = "";
77  zap['w'] = make_field(writeFits, "write fit results.");
78  zap['x'] = make_field(fitRange, "ROOT fit range (time-over-threshold).") = JRange_t(0.0, 35.0);
79  zap['O'] = make_field(option, "ROOT fit options, see TH1::Fit.") = "";
80  zap['W'] = make_field(WMin, "minimal histogram contents.") = 1000.0;
81  zap['R'] = make_field(regexp, "regular expression for histogram name.") = MAKE_CSTRING(WILDCARD << _2SToT);
82  zap['d'] = make_field(debug, "debug.") = 1;
83 
84  zap(argc, argv);
85  }
86  catch(const exception &error) {
87  FATAL(error.what() << endl);
88  }
89 
90 
91  //----------------------------------------------------------
92  // load detector file and PMT parameters
93  //----------------------------------------------------------
94 
96 
97  try {
98  load(detectorFile, detector);
99  }
100  catch(const JException& error) {
101  FATAL(error);
102  }
103 
104 
106 
107  if (pmtFile != "") {
108  try {
109  parameters.load(pmtFile.c_str());
110  }
111  catch(const JException& error) {}
112  }
113 
114  parameters.comment.add(JMeta(argc, argv));
115 
116  // Set ROOT fit options
117  if (option.find('R') == string::npos) { option += 'R'; }
118  if (option.find('S') == string::npos) { option += 'S'; }
119  if (option.find('Q') == string::npos && debug < JEEP::debug_t) { option += 'Q'; }
120 
121 
122  TFile* in = TFile::Open(inputFile.c_str(), "exist");
123 
124  if (in == NULL || !in->IsOpen()) {
125  FATAL("File: " << inputFile << " not opened." << endl);
126  }
127 
128 
129  TH2D gain("gain", NULL,
130  100, 0.0, 2.0,
131  100, 0.0, 1.0);
132  TH1D chi2("chi2", NULL, 100, 0.0, 5.0);
133 
134 
135  TFile out(outputFile.c_str(), "recreate");
136 
137  //----------------------------------------------------------
138  // loop over modules in detector file to fit histogram
139  //----------------------------------------------------------
140 
141  for (JDetector::iterator module = detector.begin(); module != detector.end(); ++module) {
142 
143  DEBUG("Module " << module->getID() << endl);
144 
145  //----------------------------------------------------------
146  // get histogram for this module
147  //----------------------------------------------------------
148 
149  TH2D* h2s = (TH2D*) in->Get(replace(regexp, WILDCARD, MAKE_STRING(module->getID())).c_str());
150 
151  if (h2s == NULL) {
152 
153  WARNING("No histogram " << module->getID() << "; skip fit." << endl);
154 
155  continue;
156  }
157 
158  const int ny = h2s->GetYaxis()->GetNbins();
159  const double ymin = h2s->GetYaxis()->GetXmin();
160  const double ymax = h2s->GetYaxis()->GetXmax();
161 
162  TH1D chi2_1d (MAKE_CSTRING(module->getID() << ".1chi2"), NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS-0.5);
163  TH1D gain_1d (MAKE_CSTRING(module->getID() << ".1gain"), NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS-0.5);
164  TH1D gainspread_1d(MAKE_CSTRING(module->getID() << ".1gainspread"), NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS-0.5);
165 
166  for (int ix = 1; ix <= h2s->GetXaxis()->GetNbins(); ++ix) {
167 
168  const JPMTIdentifier id(module->getID(), ix-1);
169 
170  parameters[id].setPMTParameters(parameters.getDefaultPMTParameters());
171 
172  TH1D h1(MAKE_CSTRING(module->getID() << "." << id.getPMTAddress() << ".1ToT"), NULL, ny, ymin, ymax);
173 
174  h1.Sumw2();
175 
176  for (int iy = h2s->GetNbinsY(); iy >= 1; --iy) {
177 
178  const double w = h1.GetBinWidth (iy);
179  const double y = h2s->GetBinContent(ix, iy);
180 
181  h1.SetBinContent(iy, y/w);
182  h1.SetBinError (iy, sqrt(y)/w);
183  }
184 
185  // Normalize histogram
186  const double Wtot = h1.Integral();
187  h1.Scale(1.0/Wtot);
188 
189 
190  //----------------------------------------------------------
191  // set-up fitting procedure
192  //----------------------------------------------------------
193 
194  JFitToT fit(parameters.getPMTParameters(id), fitRange);
195  const TFitResultPtr result = fit(h1, WMin/Wtot, option.c_str());
196 
197  if ((int(result) < 0) || (!result->IsValid())) { // Fit failed
198 
199  WARNING("Fit failure for PMT " << id << "; set default values." << endl);
200 
201  h1.GetListOfFunctions()->Clear();
202 
203  const double initG = fit.getInitGain();
204  const double initGS = fit.getInitGainSpread();
205 
206  parameters[id].gain = ( initG < FITTOT_GAIN_MIN ? FITTOT_GAIN_MIN :
207  (initG > FITTOT_GAIN_MAX ? FITTOT_GAIN_MAX : initG) );
208  parameters[id].gainSpread = ( initGS < FITTOT_GAINSPREAD_MIN ? FITTOT_GAINSPREAD_MIN :
209  (initGS > FITTOT_GAINSPREAD_MAX ? FITTOT_GAINSPREAD_MAX : initGS) );
210 
211  if (writeFits) {
212  out << h1;
213  }
214 
215  continue;
216  }
217 
218  if ((ymin < fitRange.getLowerLimit() || ymax > fitRange.getUpperLimit()) &&
219  (option.find("0") == string::npos && option.find("N") == string::npos)) {
220 
221  // add function with full range
222 
223  TF1* f1 = new TF1(MAKE_CSTRING(FITTOT_FNAME << "_fullRange"),
224  &fit,
226  ymin,
227  ymax,
228  JFitToT::getNumberOfModelParameters());
229  f1->SetParameters(fit.GetParameters());
230  f1->SetLineStyle(kDotted);
231  f1->SetNpx(1000);
232 
233  h1.GetListOfFunctions()->Add(f1);
234  }
235 
236  // store fit results
237 
238  const double reducedChi2 = result->Chi2() / result->Ndf();
239 
240  const JFitToTParameters values(fit.GetParameters());
241  const JFitToTParameters errors(fit.GetParErrors());
242 
243  chi2_1d. SetBinContent(ix, reducedChi2);
244  gain_1d. SetBinContent(ix, values.gain);
245  gain_1d. SetBinError (ix, errors.gain);
246  gainspread_1d.SetBinContent(ix, values.gainSpread);
247  gainspread_1d.SetBinError (ix, errors.gainSpread);
248 
249  gain.Fill(fit.gain, fit.gainSpread);
250  chi2.Fill(TMath::Log10(reducedChi2));
251 
252  NOTICE("PMT " << id << " chi2 / NDF " << result->Chi2() << ' ' << result->Ndf() << endl);
253 
254  parameters[id].gain = fit.gain;
255  parameters[id].gainSpread = fit.gainSpread;
256 
257  if (writeFits) {
258  out << h1;
259  }
260  }
261 
262  if (writeFits) {
263  out << chi2_1d << gain_1d << gainspread_1d;
264  }
265  }
266 
267  if (pmtFile != "") {
268  parameters.store(pmtFile.c_str());
269  }
270 
271  out << gain << chi2;
272 
273  out.Write();
274  out.Close();
275 }
276 
277 
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:1493
General exception.
Definition: JException.hh:23
data_type w[N+1][M+1]
Definition: JPolint.hh:708
#define WARNING(A)
Definition: JMessage.hh:65
double getValue(const JScale_t scale)
Get numerical value corresponding to scale.
Definition: JScale.hh:47
debug
Definition: JMessage.hh:29
TString replace(const TString &target, const TRegexp &regexp, const T &replacement)
Replace regular expression in input by given replacement.
Definition: JPrintResult.cc:63
static const double FITTOT_GAIN_MAX
Maximal gain.
Definition: JFitToT.hh:37
Detector data structure.
Definition: JDetector.hh:80
static const double FITTOT_GAIN_MIN
Minimal gain.
Definition: JFitToT.hh:36
esac print_variable DETECTOR INPUT_FILE OUTPUT_FILE CDF for TYPE in
Definition: JSirene.sh:45
*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:708
then for HISTOGRAM in h0 h1
Definition: JMatrixNZ.sh:69
string outputFile
Parametrisation of time-over-threshold distribution.
Definition: JFitToT.hh:173
Data structure for detector geometry and calibration.
Fit parameters for two-fold coincidence rate due to K40.
Definition: JFitToT.hh:50
#define MAKE_STRING(A)
Make string.
Definition: JPrint.hh:699
then echo Variable JPP_DIR undefined exit fi source $JPP_DIR setenv sh $JPP_DIR set_variable DEBUG set_variable NPE set_variable FIT_RANGE set_variable OUTPUT_DIR tmp set_variable OUTPUT_JSON $OUTPUT_DIR HVtuning json set_variable OUTPUT_ROOT $OUTPUT_DIR HVtuning root set_variable FIT_OPTIONS RME set_variable PMT_DEFAULT gain
Definition: JTuneHV.sh:30
Constants.
I/O formatting auxiliaries.
Detector file.
Definition: JHead.hh:130
static const double FITTOT_GAINSPREAD_MAX
Maximal gain spread.
Definition: JFitToT.hh:40
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1954
return result
Definition: JPolint.hh:695
ROOT I/O of application specific meta data.
#define NOTICE(A)
Definition: JMessage.hh:64
Auxiliary class for map of PMT parameters.
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
int debug
debug level
Definition: JSirene.cc:61
General purpose messaging.
#define FATAL(A)
Definition: JMessage.hh:67
JRange< Double_t > JRange_t
Definition: JFitToT.hh:34
Auxiliary class to define a range between two values.
Utility class to parse command line options.
static const char *const _2SToT
Histogram naming.
PMT analogue signal processor.
static const double FITTOT_GAINSPREAD_MIN
Minimal gain spread.
Definition: JFitToT.hh:39
KM3NeT DAQ constants, bit handling, etc.
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
static const std::string FITTOT_FNAME
Definition: JFitToT.hh:45
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62
int main(int argc, char *argv[])
Definition: Main.cpp:15