Jpp
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JGandalfFitToGauss.cc File Reference

Program to test JFIT::JGandalf algorithm. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include <cmath>
#include "TROOT.h"
#include "TFile.h"
#include "TF1.h"
#include "TH1D.h"
#include "TMath.h"
#include "JFit/JGandalf.hh"
#include "JTools/JElement.hh"
#include "JTools/JConstants.hh"
#include "JTools/JQuantile.hh"
#include "JROOT/JRootToolkit.hh"
#include "JMath/JGauss.hh"
#include "Jeep/JTimer.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

Program to test JFIT::JGandalf algorithm.

Author
mdejong

Definition in file JGandalfFitToGauss.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 72 of file JGandalfFitToGauss.cc.

73 {
74  using namespace std;
75  using namespace JPP;
76 
77  string outputFile;
78  int numberOfEvents;
79  JGauss gauss;
80  JGauss precision;
81  int debug;
82 
83  try {
84 
85  JParser<> zap("Program to test JGandalf algorithm.");
86 
87  zap['o'] = make_field(outputFile) = "";
88  zap['n'] = make_field(numberOfEvents) = 100;
89  zap['@'] = make_field(gauss) = JGauss(0.0, 1.0, 1000.0, 100.0);
90  zap['e'] = make_field(precision) = JGauss(0.05, 0.05, 25.0, 25.0);
91  zap['d'] = make_field(debug) = 3;
92 
93  zap(argc, argv);
94  }
95  catch(const exception& error) {
96  FATAL(error.what() << endl);
97  }
98 
99  using namespace JFIT;
100 
101  ASSERT(numberOfEvents > 0);
102 
104 
105  TF1 fs("fs", "exp(-0.5 * (x-[0])*(x-[0]) / ([1]*[1]))");
106  TF1 fb("fb", "1.0");
107 
108  fs.FixParameter(0, gauss.mean);
109  fs.FixParameter(1, gauss.sigma);
110 
111  const Int_t nx = 21;
112  const Double_t xmin = -5.0;
113  const Double_t xmax = +5.0;
114 
115  JQuantile Q[] = { JQuantile("mean "),
116  JQuantile("sigma "),
117  JQuantile("signal "),
118  JQuantile("background") };
119 
120  TH1D H[] = { TH1D("ha", "", 101, -0.1, +0.1),
121  TH1D("hb", "", 101, -0.1, +0.1),
122  TH1D("hc", "", 101, -100.0, +100.0),
123  TH1D("hd", "", 101, -100.0, +100.0) };
124 
125  JTimer timer;
126 
127  for (int i = 0; i != numberOfEvents; ++i) {
128 
129  STATUS("event: " << setw(10) << i << '\r'); DEBUG(endl);
130 
131  TH1D h0("h0", NULL, nx, xmin, xmax);
132 
133  h0.Sumw2();
134 
135  h0.FillRandom("fs", (Int_t) gauss.signal);
136  h0.FillRandom("fb", (Int_t) gauss.background);
137 
138  vector<JElement_t> data;
139 
140  for (Int_t i = 1; i <= h0.GetNbinsX(); ++i) {
141  data.push_back(JElement_t(h0.GetBinCenter (i),
142  h0.GetBinContent(i)));
143  }
144 
145  JGandalf<JGauss> fit;
146 
147  fit.parameters.resize(4);
148 
149  fit.parameters[0] = &JGauss::mean;
150  fit.parameters[1] = &JGauss::sigma;
151  fit.parameters[2] = &JGauss::signal;
152  fit.parameters[3] = &JGauss::background;
153 
154  fit.value = JGauss(h0.GetMean(),
155  h0.GetRMS(),
156  h0.GetEntries() * (xmax - xmin) / nx - h0.GetMinimum(),
157  h0.GetMinimum());
158 
159  DEBUG("Start value " << fit.value << endl);
160 
161  timer.start();
162 
163  const double chi2 = fit(g1, data.begin(), data.end());
164 
165  timer.stop();
166 
167  DEBUG("Final value " << fit.value << endl);
168  DEBUG("Chi2 " << chi2 << endl);
169 
170  const double Y[] = { fit.value.mean - gauss.mean,
171  fit.value.sigma - gauss.sigma,
172  fit.value.signal * nx / (xmax - xmin) - gauss.signal,
173  fit.value.background * nx - gauss.background };
174 
175  for (int i = 0; i != sizeof(Q)/sizeof(Q[0]); ++i) {
176  Q[i].put (Y[i]);
177  H[i].Fill(Y[i]);
178  }
179  }
180 
181  for (int i = 0; i != sizeof(Q)/sizeof(Q[0]); ++i) {
182  NOTICE((i == 0 ? longprint : shortprint) << Q[i]);
183  }
184 
185  if (debug >= notice_t) {
186  timer.print(cout, true, micro_t);
187  }
188 
189  if (outputFile != "") {
190 
191  TFile out(outputFile.c_str(), "recreate");
192 
193  for (int i = 0; i != sizeof(H)/sizeof(H[0]); ++i) {
194  out << H[i];
195  }
196 
197  out.Write();
198  out.Close();
199  }
200 
201  for (int i = 0; i != sizeof(Q)/sizeof(Q[0]); ++i) {
202  ASSERT(fabs(Q[i].getMean()) < Q[i].getSTDev());
203  }
204 
205  ASSERT(Q[0].getSTDev() < precision.mean);
206  ASSERT(Q[1].getSTDev() < precision.sigma);
207  ASSERT(Q[2].getSTDev() < precision.signal);
208  ASSERT(Q[3].getSTDev() < precision.background);
209 
210  return 0;
211 }
static const double H
Planck constant [eV s].
Definition: JConstants.hh:25
Utility class to parse command line options.
Definition: JParser.hh:1410
micro
Definition: JScale.hh:29
JModel_t value
Definition: JGandalf.hh:265
#define STATUS(A)
Definition: JMessage.hh:61
notice
Definition: JMessage.hh:30
std::vector< parameter_type > parameters
Definition: JGandalf.hh:267
double mean
Definition: JGauss.hh:169
string outputFile
double getMean(vector< double > &v)
get mean of vector content
Quantile calculator.
Definition: JQuantile.hh:34
#define ASSERT(A)
Assert macro.
Definition: JMessage.hh:72
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1836
Auxiliary class for CPU timing and usage.
Definition: JTimer.hh:32
#define NOTICE(A)
Definition: JMessage.hh:62
int debug
debug level
Definition: JSirene.cc:59
#define FATAL(A)
Definition: JMessage.hh:65
Fit method based on the Levenberg-Marquardt method.
Definition: JGandalf.hh:45
std::ostream & longprint(std::ostream &out)
Set long printing.
Definition: JPrint.hh:189
double sigma
Definition: JGauss.hh:170
Data structure for Gaussian function on top of a flat background.
Definition: JGauss.hh:30
double signal
Definition: JGauss.hh:171
std::ostream & shortprint(std::ostream &out)
Set short printing.
Definition: JPrint.hh:161
double background
Definition: JGauss.hh:172
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:60
Double_t g1(const Double_t x)
Function.
Definition: JQuantiles.cc:25