Jpp  16.0.3
the software that should make you happy
 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 
15 
16 #include "JLang/JLangToolkit.hh"
17 #include "JPhysics/JConstants.hh"
18 
19 #include "JDetector/JDetector.hh"
23 
24 #include "JROOT/JRootToolkit.hh"
25 
26 #include "JTools/JRange.hh"
27 
28 #include "JSupport/JMeta.hh"
29 #include "JSupport/JSupport.hh"
32 
34 #include "JCalibrate/JFitToT.hh"
35 
36 #include "Jeep/JPrint.hh"
37 #include "Jeep/JParser.hh"
38 #include "Jeep/JMessage.hh"
39 
40 namespace {
41 
42  /**
43  * Wild card for histogram naming.\n
44  * The wild card will replaced by the module identifier.
45  */
46  static const char* const WILDCARD = "%";
47 }
48 
49 
50 /**
51  * \file
52  *
53  * Auxiliary program to fit time-over-threshold distributions.
54  * The fitted parameters are the PMT gain and gain spread.
55  * \author mkarel
56  */
57 int main(int argc, char **argv)
58 {
59  using namespace std;
60  using namespace JPP;
61  using namespace KM3NETDAQ;
62 
63  typedef JRange<double> JRange_t;
65 
66  string inputFile;
68  string detectorFile;
69  string pmtFile;
70  bool writeFits;
71 
72  double Wmin = 5e3; //!< Minimal histogram weight in fit range
73  double LeftMargin = 10.0; //!< Fit range left of time-over-threshold peak
74  double RightMargin = 10.0; //!< Fit range right of time-over-threshold peak
75  double gradientThreshold = -0.005; //!< Minimal right-to-left gradient value in time-over-threshold peak-search,\n
76  //!< normalized with respect to the total number of histogram entries
77  JRange_t gainLimits (FITTOT_GAIN_MIN,
78  FITTOT_GAIN_MAX); //!< gain limits
79  JRange_t gainSpreadLimits(FITTOT_GAINSPREAD_MIN,
80  FITTOT_GAINSPREAD_MAX); //!< gainspread limits
81 
82  JRange_t ToTfitRange;
83  string option;
84 
85  string regexp;
86  int debug;
87 
88  try {
89 
90  JProperties properties;
91 
92  properties.insert(gmake_property(Wmin));
93  properties.insert(gmake_property(LeftMargin));
94  properties.insert(gmake_property(RightMargin));
95  properties.insert(gmake_property(gradientThreshold));
96  properties.insert(gmake_property(gainLimits));
97  properties.insert(gmake_property(gainSpreadLimits));
98 
99  JParser<> zap("Auxiliary program to fit time-over-threshold distributions.");
100 
101  zap['f'] = make_field(inputFile, "input file (output from JCalibrateToT).");
102  zap['o'] = make_field(outputFile, "output file.") = "fit.root";
103  zap['a'] = make_field(detectorFile, "detector file.");
104  zap['P'] = make_field(pmtFile, "specify PMT file name that can be given to JTriggerEfficiency.") = "";
105  zap['w'] = make_field(writeFits, "write fit results.");
106  zap['O'] = make_field(option, "ROOT fit options, see TH1::Fit") = "";
107  zap['@'] = make_field(properties, "fit properties") = JPARSER::initialised();
108  zap['x'] = make_field(ToTfitRange, "ROOT fit range (time-over-threshold).") = JRange_t(0.0, 35.0);
109  zap['R'] = make_field(regexp, "regular expression for histogram name.") = MAKE_CSTRING(WILDCARD << _2SToT);
110  zap['d'] = make_field(debug, "debug.") = 1;
111 
112  zap(argc, argv);
113  }
114  catch(const exception &error) {
115  FATAL(error.what() << endl);
116  }
117 
118 
119  //----------------------------------------------------------
120  // load detector file and PMT parameters
121  //----------------------------------------------------------
122 
124 
125  try {
126  load(detectorFile, detector);
127  }
128  catch(const JException& error) {
129  FATAL(error);
130  }
131 
132 
133  JPMTParametersMap parametersMap;
134 
135  if (!pmtFile.empty()) {
136 
137  try {
138  parametersMap.load(pmtFile.c_str());
139  }
140  catch(const JException& error) {}
141  }
142 
143  parametersMap.comment.add(JMeta(argc, argv));
144 
145  // Set ROOT fit options
146  if (option.find('R') == string::npos) { option += 'R'; }
147  if (option.find('S') == string::npos) { option += 'S'; }
148  if (option.find('Q') == string::npos && debug < JEEP::debug_t) { option += 'Q'; }
149 
150 
151  TFile* in = TFile::Open(inputFile.c_str(), "exist");
152 
153  if (in == NULL || !in->IsOpen()) {
154  FATAL("File: " << inputFile << " not opened." << endl);
155  }
156 
157 
158  TH2D gain("gain", NULL,
159  100, 0.0, 2.0,
160  100, 0.0, 1.0);
161  TH1D chi2("chi2", NULL, 100, 0.0, 5.0);
162 
163  outputFile.open();
164 
165  if (!outputFile.is_open()) {
166  FATAL("Error opening file " << outputFile << endl);
167  }
168 
169  outputFile.put(JMeta(argc, argv));
170 
171  //----------------------------------------------------------
172  // loop over modules in detector file to fit histogram
173  //----------------------------------------------------------
174 
175  for (JDetector::iterator module = detector.begin(); module != detector.end(); ++module) {
176 
177  NOTICE("Module " << setw(10) << module->getID() << ' ' << module->getLocation() << endl);
178 
179  if (module->empty()) {
180  continue;
181  }
182 
183  //----------------------------------------------------------
184  // get histogram for this module
185  //----------------------------------------------------------
186 
187  TH2D* h2s = (TH2D*) in->Get(replace(regexp, WILDCARD, MAKE_STRING(module->getID())).c_str());
188 
189  if (h2s == NULL) {
190 
191  WARNING("No histogram " << module->getID() << "; skip fit." << endl);
192 
193  continue;
194  }
195 
196  const double ny = h2s->GetYaxis()->GetNbins();
197  const double ymin = h2s->GetYaxis()->GetXmin();
198  const double ymax = h2s->GetYaxis()->GetXmax();
199 
200  TH1D chi2_1d (MAKE_CSTRING(module->getID() << ".1chi2"), NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS-0.5);
201  TH1D gain_1d (MAKE_CSTRING(module->getID() << ".1gain"), NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS-0.5);
202  TH1D gainspread_1d(MAKE_CSTRING(module->getID() << ".1gainspread"), NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS-0.5);
203 
204  for (int ix = 1; ix <= h2s->GetXaxis()->GetNbins(); ++ix) {
205 
206  const JPMTIdentifier id(module->getID(), ix-1);
207 
208  const string name = MAKE_STRING(module->getID() << '.' << FILL(2,'0') <<
209  id.getPMTAddress() << FITTOT_SUFFIX);
210 
211  DEBUG("Processing histogram " << name << endl);
212 
213  TH1D h1(name.c_str(), NULL, ny, ymin, ymax);
214 
215  h1.Sumw2();
216 
217  for (int iy = h2s->GetNbinsY(); iy >= 1; --iy) {
218 
219  const double w = h1.GetBinWidth (iy);
220  const double y = h2s->GetBinContent(ix, iy);
221 
222  h1.SetBinContent(iy, y/w);
223  h1.SetBinError (iy, sqrt(y)/w);
224  }
225 
226  JRange_t range(ToTfitRange);
227 
228  const double weight = h1.Integral(h1.FindBin(range.getLowerLimit()),
229  h1.FindBin(range.getUpperLimit()));
230 
231  if (weight <= Wmin) {
232 
233  WARNING("Insufficient histogram entries for PMT " << id <<
234  "(" << weight << " < " << Wmin << "); skip fit." << endl);
235 
236  if (writeFits) {
237  outputFile.put(h1);
238  }
239 
240  continue;
241  }
242 
243  //----------------------------------------------------------
244  // find the mode of the time-over-threshold distribution
245  //----------------------------------------------------------
246 
247  double max = 0.0;
248  double ToTmode = JDETECTOR::TIME_OVER_THRESHOLD_NS;
249 
250  JPMTParameters& parameters = parametersMap[id];
251 
252  DEBUG(RIGHT(10) << "ToT" << RIGHT(10) << "Counts" << RIGHT(10) << "gradient" <<
253  RIGHT(10) << "threshold" << RIGHT(10) << "Mode" << RIGHT(10) << "Max" << endl);
254 
255  for (int i = h1.GetBinCenter(ny-1);
256  i > h1.GetBinCenter(h1.FindBin(parameters.mean_ns + parameters.sigma_ns));
257  --i) {
258 
259  const double x = h1.GetBinCenter (i);
260  const double y = h1.GetBinContent(i);
261 
262  const double gradient = ( (h1.GetBinContent(i-1) - h1.GetBinContent(i+1)) /
263  (h1.GetBinCenter (i+1) - h1.GetBinCenter (i-1)) );
264 
265  DEBUG(FIXED(10,1) << x << FIXED(10,1) << y << FIXED(10,1) << gradient <<
266  FIXED(10,1) << -fabs(gradientThreshold) * h1.Integral() << FIXED(10,1) << ToTmode << FIXED(10,1) << max << endl);
267 
268 
269  if (y > max) {
270 
271  ToTmode = x;
272  max = y;
273 
274  } else if (gradient < -fabs(gradientThreshold) * h1.Integral()) {
275 
276  break;
277  }
278  }
279 
280  //----------------------------------------------------------
281  // set fit range
282  //----------------------------------------------------------
283 
284  if (!range.is_valid()) {
285  range.setRange(ToTmode - LeftMargin,
286  ToTmode + RightMargin);
287  }
288 
289  range.setRange(std::max(h1.GetXaxis()->GetXmin(), range.getLowerLimit()),
290  std::min(h1.GetXaxis()->GetXmax(), range.getUpperLimit()));
291 
292  DEBUG(LEFT(10) << "Fit-range: " << FIXED(20,1) << range << endl);
293 
294  //----------------------------------------------------------
295  // perform minimum chi-square fit
296  //----------------------------------------------------------
297 
298  JFitToT fit(parameters, range);
299 
301  gainLimits.getLowerLimit(), gainLimits.getUpperLimit());
302  setParLimits(fit, fit.getModelParameter(&JFitToT::JFitToTParameters::gainSpread),
303  gainSpreadLimits.getLowerLimit(), gainSpreadLimits.getUpperLimit());
304 
305  const double initGain = fit.getCPU().getNPE(ToTmode);
306  const double initGainSpread = initGain / 3.0;
307 
308  if (!gainLimits(initGain)) { // Invalid initial gain value
309 
310  WARNING("Invalid initial gain for PMT " << id << "; set default values." << endl);
311 
312  parameters.gain = gainLimits .constrain(initGain);
313  parameters.gainSpread = gainSpreadLimits.constrain(initGainSpread);
314 
315  if (writeFits) {
316  outputFile.put(h1);
317  }
318 
319  continue;
320  }
321 
322  fit.gain = initGain;
323  fit.gainSpread = initGainSpread;
324 
325  const TFitResultPtr result = fit(h1, option);
326 
327  if (result.Get() == NULL) {
328  FATAL("Invalid TFitResultPtr " << h1.GetName() << endl);
329  }
330 
331  if (int(result) < 0 || !result->IsValid()) { // Fit failed
332 
333  WARNING("Fit failure for PMT " << id << "; set initialization values." << endl);
334 
335  h1.GetListOfFunctions()->Clear();
336 
337  parameters.gain = initGain;
338  parameters.gainSpread = initGainSpread;
339 
340  if (writeFits) {
341  outputFile.put(h1);
342  }
343 
344  continue;
345  }
346 
347  if ((ymin < range.getLowerLimit() || ymax > range.getUpperLimit()) &&
348  (option.find("0") == string::npos && option.find("N") == string::npos)) {
349 
350  // add function with full range
351 
352  TF1* f1 = new TF1(MAKE_CSTRING(FITTOT_FNAME << "_fullRange"),
353  &fit,
355  ymin,
356  ymax,
357  JFitToT::getNumberOfModelParameters());
358  f1->SetParameters(fit.GetParameters());
359  f1->SetLineStyle(kDotted);
360  f1->SetNpx(1000);
361 
362  h1.GetListOfFunctions()->Add(f1);
363  }
364 
365  // store fit results
366 
367  const double reducedChi2 = result->Chi2() / result->Ndf();
368 
369  const JFitToTParameters values(fit.GetParameters());
370  const JFitToTParameters errors(fit.GetParErrors());
371 
372  chi2_1d. SetBinContent(ix, reducedChi2);
373  gain_1d. SetBinContent(ix, values.gain);
374  gain_1d. SetBinError (ix, errors.gain);
375  gainspread_1d.SetBinContent(ix, values.gainSpread);
376  gainspread_1d.SetBinError (ix, errors.gainSpread);
377 
378  gain.Fill(fit.gain, fit.gainSpread);
379  chi2.Fill(TMath::Log10(reducedChi2));
380 
381  NOTICE("PMT " << id << " chi2 / NDF " << result->Chi2() << ' ' << result->Ndf() << endl);
382 
383  parameters.gain = fit.gain;
384  parameters.gainSpread = fit.gainSpread;
385 
386  if (writeFits) {
387  outputFile.put(h1);
388  }
389  }
390 
391  if (writeFits) {
392  outputFile.put(chi2_1d);
393  outputFile.put(gain_1d);
394  outputFile.put(gainspread_1d);
395  }
396  }
397 
398  for (JSingleFileScanner<JMeta> in(inputFile); in.hasNext(); ) {
399 
400  const JMeta& meta = *in.next();
401 
402  parametersMap.comment.add(meta);
403 
404  outputFile.put(meta);
405  }
406 
407  for (JRootFileReader<JDAQHeader> in(inputFile.c_str()); in.hasNext(); ) {
408  outputFile.put(*in.next());
409  }
410 
411 
412  if (pmtFile != "") {
413  parametersMap.store(pmtFile.c_str());
414  }
415 
416  outputFile.put(gain);
417  outputFile.put(chi2);
418 
419  outputFile.close();
420 }
421 
422 
Auxiliary class for ROOT I/O of application specific meta data.
Definition: JMeta.hh:70
Object writing to file.
Utility class to parse command line options.
Definition: JParser.hh:1500
General exception.
Definition: JException.hh:23
double gain
gain [unit]
data_type w[N+1][M+1]
Definition: JPolint.hh:757
#define WARNING(A)
Definition: JMessage.hh:65
const Double_t getModelParameter(const int i) const
Get model parameter.
Definition: JFitToT.hh:146
double getValue(const JScale_t scale)
Get numerical value corresponding to scale.
Definition: JScale.hh:47
debug
Definition: JMessage.hh:29
T getLowerLimit() const
Get lower limit.
Definition: JRange.hh:202
int main(int argc, char *argv[])
Definition: Main.cc:15
ROOT TTree parameter settings of various packages.
TString replace(const TString &target, const TRegexp &regexp, const T &replacement)
Replace regular expression in input by given replacement.
Definition: JPrintResult.cc:63
#define gmake_property(A)
macro to convert (template) parameter to JPropertiesElement object
double gainSpread
gain spread [unit]
static const double FITTOT_GAIN_MAX
Default maximal gain.
Definition: JFitToT.hh:44
Detector data structure.
Definition: JDetector.hh:89
Double_t gainSpread
PMT gain spread.
Definition: JFitToT.hh:169
const double TIME_OVER_THRESHOLD_NS
Specification for time-over-threshold corresponding to a one photo-electron pulse.
static const double FITTOT_GAIN_MIN
Default minimal gain.
Definition: JFitToT.hh:43
Recording of objects on file according a format that follows from the file name extension.
then echo Enter input within $TIMEOUT_S seconds echo n User name
Definition: JCookie.sh:42
Utility class to parse parameter values.
Definition: JProperties.hh:496
*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
then for HISTOGRAM in h0 h1
Definition: JMatrixNZ.sh:71
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition: JParser.hh:66
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:446
string outputFile
Parametrisation of time-over-threshold distribution.
Definition: JFitToT.hh:179
Data structure for detector geometry and calibration.
Fit parameters for two-fold coincidence rate due to K40.
Definition: JFitToT.hh:53
#define MAKE_STRING(A)
Make string.
Definition: JPrint.hh:142
const JPMTAnalogueSignalProcessor & getCPU() const
Access method for the analogue signal processor.
Definition: JFitToT.hh:210
static const std::string FITTOT_SUFFIX
Definition: JFitToT.hh:36
Type list.
Definition: JTypeList.hh:22
Scanning of objects from a single file according a format that follows from the extension of each fil...
Type definition of range.
Definition: JHead.hh:39
I/O formatting auxiliaries.
double mean_ns
mean time-over-threshold of threshold-band hits [ns]
Detector file.
Definition: JHead.hh:224
double sigma_ns
time-over-threshold standard deviation of threshold-band hits [ns]
static const double FITTOT_GAINSPREAD_MAX
Default maximal gain spread.
Definition: JFitToT.hh:47
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
void setRange(const range_type &range)
Set range.
Definition: JRange.hh:146
return result
Definition: JPolint.hh:743
ROOT I/O of application specific meta data.
T getUpperLimit() const
Get upper limit.
Definition: JRange.hh:213
#define NOTICE(A)
Definition: JMessage.hh:64
virtual double getNPE(const double tot_ns) const override
Get number of photo-electrons.
Physics constants.
Auxiliary class for map of PMT parameters.
T constrain(argument_type x) const
Constrain value to range.
Definition: JRange.hh:350
Double_t gain
PMT gain.
Definition: JFitToT.hh:168
int debug
debug level
Definition: JSirene.cc:63
bool setParLimits(TF1 &f1, const Int_t index, Double_t xmin, Double_t xmax)
Set fit parameter limits.
General purpose messaging.
Auxiliary data structure for sequence of same character.
Definition: JManip.hh:328
#define FATAL(A)
Definition: JMessage.hh:67
void load(const char *file_name)
Load from input file.
z range($ZMAX-$ZMIN)< $MINIMAL_DZ." fi fi typeset -Z 4 STRING typeset -Z 2 FLOOR JPlot1D -f $
then $JPP_DIR examples JDetector JToT o $OUTPUT_FILE n N $NPE P gain
Definition: JToT.sh:47
bool is_valid() const
Check validity of range.
Definition: JRange.hh:311
void load(const std::string &file_name, JDetector &detector)
Load detector from input file.
Auxiliary class to define a range between two values.
Utility class to parse command line options.
then usage $script< input_file >< detector_file > fi set_variable OUTPUT_DIR set_variable SELECTOR JDAQTimesliceL1 set_variable DEBUG case set_variable DEBUG
static const char *const _2SToT
Histogram naming.
PMT analogue signal processor.
static const double FITTOT_GAINSPREAD_MIN
Default minimal gain spread.
Definition: JFitToT.hh:46
Object reading from a list of files.
possible values
Data structure for PMT parameters.
do set_variable DETECTOR_TXT $WORKDIR detector
KM3NeT DAQ constants, bit handling, etc.
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:42
ROOT file reader.
static const std::string FITTOT_FNAME
Definition: JFitToT.hh:37
std::vector< double > weight
Definition: JAlgorithm.hh:417