Jpp 19.3.0-rc.2
the software that should make you happy
Loading...
Searching...
No Matches
JRootFitToGauss.cc File Reference

Program to test ROOT fit. More...

#include <string>
#include <iostream>
#include <iomanip>
#include <vector>
#include <cmath>
#include "TROOT.h"
#include "TFile.h"
#include "TF1.h"
#include "TH1D.h"
#include "JROOT/JMinimizer.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 ROOT fit.

Author
mdejong

Definition in file JRootFitToGauss.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 30 of file JRootFitToGauss.cc.

31{
32 using namespace std;
33 using namespace JPP;
34
35 string outputFile;
36 int numberOfEvents;
37 JGauss gauss;
38 JGauss precision;
39 int debug;
40
41 try {
42
43 JParser<> zap("Program to test ROOT fit.");
44
45 zap['o'] = make_field(outputFile) = "";
46 zap['n'] = make_field(numberOfEvents) = 1000;
47 zap['@'] = make_field(gauss) = JGauss(0.0, 1.0, 1000.0, 100.0);
48 zap['e'] = make_field(precision) = JGauss(0.05, 0.05, 25.0, 25.0);
49 zap['d'] = make_field(debug) = 3;
50
51 zap(argc, argv);
52 }
53 catch(const exception& error) {
54 FATAL(error.what() << endl);
55 }
56
57
58 ASSERT(numberOfEvents > 0);
59
60 TF1 fs("fs", "exp(-0.5 * (x-[0])*(x-[0]) / ([1]*[1]))");
61 TF1 fb("fb", "1.0");
62
63 fs.FixParameter(0, gauss.mean);
64 fs.FixParameter(1, gauss.sigma);
65
66 const Int_t nx = 21;
67 const Double_t xmin = -5.0;
68 const Double_t xmax = +5.0;
69
70 JQuantile Q[] = { JQuantile("mean "),
71 JQuantile("sigma "),
72 JQuantile("signal "),
73 JQuantile("background") };
74
75 TH1D H[] = { TH1D("ha", "", 101, -0.1, +0.1),
76 TH1D("hb", "", 101, -0.1, +0.1),
77 TH1D("hc", "", 101, -100.0, +100.0),
78 TH1D("hd", "", 101, -100.0, +100.0) };
79
80 JTimer timer;
81
82 for (int i = 0; i != numberOfEvents; ++i) {
83
84 STATUS("event: " << setw(10) << i << '\r'); DEBUG(endl);
85
86 TH1D h0("h0", NULL, nx, xmin, xmax);
87
88 h0.Sumw2();
89
90 h0.FillRandom("fs", (Int_t) gauss.signal);
91 h0.FillRandom("fb", (Int_t) gauss.background);
92
93 TF1 f1("f1", "[2]*exp(-0.5 * (x-[0])*(x-[0]) / ([1]*[1])) / (TMath::Sqrt(2.0*TMath::Pi())*[1]) + [3]");
94
95 f1.SetParameter(0, h0.GetMean());
96 f1.SetParameter(1, h0.GetRMS());
97 f1.SetParameter(2, h0.GetEntries() - h0.GetMinimum() * h0.GetNbinsX());
98 f1.SetParameter(3, h0.GetMinimum());
99
100 f1.SetParError(0, 1.0e-3);
101 f1.SetParError(1, 1.0e-3);
102 f1.SetParError(2, 0.5);
103 f1.SetParError(3, 0.5);
104
105 string option = "NWL";
106
107 if (debug < debug_t && option.find('Q') == string::npos) {
108 option += "Q";
109 }
110
111 timer.start();
112
113 h0.Fit(&f1, option.c_str());
114
115 timer.stop();
116
117 const double Y[] = { f1.GetParameter(0) - gauss.mean,
118 f1.GetParameter(1) - gauss.sigma,
119 f1.GetParameter(2) * nx / (xmax - xmin) - gauss.signal,
120 f1.GetParameter(3) * nx - gauss.background };
121
122 for (int i = 0; i != sizeof(Q)/sizeof(Q[0]); ++i) {
123 Q[i].put (Y[i]);
124 H[i].Fill(Y[i]);
125 }
126 }
127
128 for (int i = 0; i != sizeof(Q)/sizeof(Q[0]); ++i) {
129 NOTICE((i == 0 ? longprint : shortprint) << Q[i]);
130 }
131
132 if (debug >= notice_t) {
133 timer.print(cout, true, micro_t);
134 }
135
136 if (outputFile != "") {
137
138 TFile out(outputFile.c_str(), "recreate");
139
140 for (int i = 0; i != sizeof(H)/sizeof(H[0]); ++i) {
141 out << H[i];
142 }
143
144 out.Write();
145 out.Close();
146 }
147
148 for (int i = 0; i != sizeof(Q)/sizeof(Q[0]); ++i) {
149 ASSERT(fabs(Q[i].getMean()) < Q[i].getSTDev());
150 }
151
152 ASSERT(Q[0].getSTDev() < precision.mean);
153 ASSERT(Q[1].getSTDev() < precision.sigma);
154 ASSERT(Q[2].getSTDev() < precision.signal);
155 ASSERT(Q[3].getSTDev() < precision.background);
156
157 return 0;
158}
double getMean(vector< double > &v)
get mean of vector content
string outputFile
std::ostream & longprint(std::ostream &out)
Set long printing.
Definition JManip.hh:172
std::ostream & shortprint(std::ostream &out)
Set short printing.
Definition JManip.hh:144
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#define STATUS(A)
Definition JMessage.hh:63
#define ASSERT(A,...)
Assert macro.
Definition JMessage.hh:90
#define NOTICE(A)
Definition JMessage.hh:64
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:72
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
Auxiliary class for CPU timing and usage.
Definition JTimer.hh:33
void print(std::ostream &out, const JScale_t scale=milli_t) const
Print timer data.
Definition JTimer.hh:172
void stop()
Stop timer.
Definition JTimer.hh:127
void start()
Start timer.
Definition JTimer.hh:106
Utility class to parse command line options.
Definition JParser.hh:1698
const JPolynome f1(1.0, 2.0, 3.0)
Function.
const double xmax
const double xmin
double gauss(const double x, const double sigma)
Gauss function (normalised to 1 at x = 0).
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Auxiliary data structure for running average, standard deviation and quantiles.
Definition JQuantile.hh:46
void put(const double x, const double w=1.0)
Put value.
Definition JQuantile.hh:133