Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JFitSinglesRates.cc
Go to the documentation of this file.
1 #include <string>
2 #include <iostream>
3 #include <iomanip>
4 #include <map>
5 
6 #include "TROOT.h"
7 #include "TFile.h"
8 #include "TH1D.h"
9 #include "TH2D.h"
10 #include "TF1.h"
11 #include "TF2.h"
12 #include "TMath.h"
13 
14 #include "JTools/JConstants.hh"
15 #include "JDAQ/JDAQ.hh"
16 #include "JDAQ/JDAQSummaryslice.hh"
17 #include "JDetector/JDetector.hh"
19 
20 #include "JTools/JRange.hh"
21 
22 #include "Jeep/JParser.hh"
23 #include "Jeep/JMessage.hh"
24 
25 /**
26  * \file
27  *
28  * Auxiliary program to fit singles rate distributions. The fitted parameters are the mean singles rate and the spread.
29  * \author mkarel
30  */
31 int main(int argc, char **argv)
32 {
33  using namespace std;
34  using namespace JPP;
35 
36  string inputFile;
37  string outputFile;
38  string detectorFile;
39  bool setDefaultLimits;
40  double peakFraction;
41  string option;
42  bool writeFits;
43  int debug;
44 
45  //Parser module handles input parameters
46  try {
47 
48  JParser<> zap("Auxiliary program to fit singles rate distributions.");
49 
50  zap['f'] = make_field(inputFile);
51  zap['o'] = make_field(outputFile) = "fit.root";
52  zap['a'] = make_field(detectorFile);
53  zap['L'] = make_field(setDefaultLimits);
54  zap['p'] = make_field(peakFraction) = 0.1;
55  zap['O'] = make_field(option) = "R W M";
56  zap['w'] = make_field(writeFits);
57  zap['d'] = make_field(debug) = 1;
58 
59  zap(argc, argv);
60  }
61  catch(const exception &error) {
62  FATAL(error.what() << endl);
63  }
64 
65 
66  using namespace KM3NETDAQ;
67 
68  JDetector detector;
69 
70  try {
71  load(detectorFile, detector);
72  }
73  catch(const JException& error) {
74  FATAL(error);
75  }
76 
77  TFile* in = TFile::Open(inputFile.c_str(), "exist");
78 
79  if (in == NULL || !in->IsOpen()) {
80  FATAL("File: " << inputFile << " not opened." << endl);
81  }
82 
83  TFile out(outputFile.c_str(), "recreate");
84 
85  for (JDetector::iterator module = detector.begin(); module != detector.end(); ++module) {
86 
87  DEBUG("Module " << module->getID() << endl);
88  const JString title = JString(module->getID());
89 
90  TH2D* h2s = (TH2D*) in->Get((title+".2S").c_str());
91 
92  if (h2s == NULL) {
93 
94  WARNING("No data in histogram " << title << endl);
95 
96  continue;
97  }
98 
99  const double factor = 1.0/1000 ;
100 
101  TH1D mean((title+".1mean").c_str(), NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS-0.5) ;
102  TH1D sigma((title+".1sigma").c_str(), NULL, NUMBER_OF_PMTS, -0.5, NUMBER_OF_PMTS-0.5) ;
103 
104  for (int i = 1; i <= h2s->GetXaxis()->GetNbins(); ++i) {
105 
106  TH1D slice((title+Form(".%i.1S", i-1)).c_str(), NULL, h2s->GetYaxis()->GetNbins(), JDAQRate::getData(factor)) ;
107 
108  slice.Sumw2() ;
109 
110  for (int j = 1 ; j <= h2s->GetNbinsY() ; ++j) {
111 
112  slice.SetBinContent(j, h2s->GetBinContent(i, j)) ;
113  slice.SetBinError (j, sqrt(h2s->GetBinContent(i, j))) ;
114 
115  }
116 
117  if (slice.Integral() <= 0.0) {
118 
119  WARNING("No data in PMT " << i-1 << " of module " << title << ", skipping fit" << endl) ;
120 
121  continue ;
122  }
123 
124  slice.Scale(1.0/slice.Integral()) ;
125 
126  // determine fit range
127  int binmax = slice.GetMaximumBin() ;
128  double ratemax = slice.GetBinContent(binmax) ;
129 
130  int bin_l = binmax ;
131  int bin_r = binmax ;
132  while ((slice.GetBinContent(bin_l)==0 || slice.GetBinContent(bin_l) >= peakFraction*ratemax) && bin_l > 0) { bin_l-- ; } ; bin_l++ ;
133  while ((slice.GetBinContent(bin_r)==0 || slice.GetBinContent(bin_r) >= peakFraction*ratemax) && bin_r < slice.GetNbinsX()) { bin_r++ ; } ; bin_r-- ;
134 
135  JRange<double> fitrange(slice.GetXaxis()->GetBinCenter(bin_l), slice.GetXaxis()->GetBinCenter(bin_r)) ;
136 
137  TF1 f1("f1", "[0]*exp(-0.5*(x-[1])*(x-[1])/([2]*[2]))/(sqrt(2*pi)*[2])");
138 
139  // default parameters for start of fit
140  f1.SetParameter(0, 0.2) ;
141  f1.SetParameter(1, 6.0) ; // [kHz]
142  f1.SetParameter(2, 0.3) ; // [kHz]
143 
144  if (setDefaultLimits) {
145  f1.SetParLimits(0, 0.0, 10.0);
146  f1.SetParLimits(1, 0.0, 20.0); // [kHz]
147  f1.SetParLimits(2, 0.0, 1.0); // [kHz]
148  }
149 
150  if (debug < JEEP::debug_t) {
151  option += " Q";
152  }
153 
154  DEBUG("Fit histogram " << slice.GetName() << " in range " << fitrange << " kHz " << endl ) ;
155 
156  slice.Fit(&f1, option.c_str(), "same", max(0.0, fitrange.first), fitrange.second) ;
157 
158  DEBUG( f1.GetParameter(0)<<" "<<f1.GetParameter(1)<<" "<<f1.GetParameter(2)<<endl ) ;
159 
160  mean.SetBinContent(i, f1.GetParameter(1)) ;
161  sigma.SetBinContent(i, f1.GetParameter(2)) ;
162 
163  if (writeFits) {
164  slice.Write() ;
165  }
166 
167  }
168 
169  mean.Write() ;
170  sigma.Write() ;
171 
172  }
173 
174  out.Close() ;
175 
176 }
177 
178 
Utility class to parse command line options.
Definition: JParser.hh:1410
#define WARNING(A)
Definition: JMessage.hh:63
debug
Definition: JMessage.hh:27
string outputFile
Data structure for detector geometry and calibration.
Constants.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1836
void load(const JString &file_name, JDetector &detector)
Load detector from input file.
int debug
debug level
Definition: JSirene.cc:59
General purpose messaging.
#define FATAL(A)
Definition: JMessage.hh:65
Auxiliary class to define a range between two values.
Utility class to parse command line options.
KM3NeT DAQ constants, bit handling, etc.
static const int NUMBER_OF_PMTS
Total number of PMTs in module.
Definition: JDAQ.hh:26
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:60
int main(int argc, char *argv[])
Definition: Main.cpp:15