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

Example program to test JTOOLS::JHistogram1D. More...

#include <string>
#include <iostream>
#include <iomanip>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TMath.h"
#include "TRandom.h"
#include "JTools/JHistogram1D.hh"
#include "JTools/JCollection.hh"
#include "JTools/JGrid.hh"
#include "JTools/JSet.hh"
#include "JTools/JQuadrature.hh"
#include "JTools/JFunction1D_t.hh"
#include "JTools/JToolsToolkit.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

Example program to test JTOOLS::JHistogram1D.

Author
mdejong

Definition in file JTools/JHistogram1D.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 57 of file JTools/JHistogram1D.cc.

58{
59 using namespace std;
60 using namespace JPP;
61
62 string outputFile;
63 int numberOfEvents;
64 int numberOfBins;
65 bool quadrature;
66 bool subtract;
67 int debug;
68
69 try {
70
71 JParser<> zap("Example program to test 1D histogram.");
72
73 zap['o'] = make_field(outputFile) = "histogram.root";
74 zap['n'] = make_field(numberOfEvents);
75 zap['N'] = make_field(numberOfBins);
76 zap['Q'] = make_field(quadrature);
77 zap['D'] = make_field(subtract);
78 zap['d'] = make_field(debug) = 3;
79
80 zap(argc, argv);
81 }
82 catch(const exception &error) {
83 FATAL(error.what() << endl);
84 }
85
86
87 if (numberOfEvents <= 0) {
88 FATAL("No events." << endl);
89 }
90
91
92 const Int_t nx = 1000;
93 const Double_t xmin = -5.0;
94 const Double_t xmax = +5.0;
95
96
97 typedef JElement2D<double, double> JElement_t;
98
100
101
102 // abscissa
103
104 if (!quadrature)
105 histogram.configure(make_grid(numberOfBins, xmin, xmax));
106 else
107 histogram.configure(JGaussHermite(numberOfBins));
108
109
110 // fill
111
112 for (int i = 0; i != numberOfEvents; ++i) {
113 histogram.fill(gRandom->Gaus(0.0, 1.0), 1.0);
114 }
115
116 histogram.div((double) numberOfEvents);
117
118
119 // interpolation
120
121 typedef JSplineFunction1S_t JFunction1D_t;
122 //typedef JPolintFunction1D<3, JElement_t, JCollection, JResultHesse<double> > JFunction1D_t;
123
124
125 JFunction1D_t F1;
126
127 integrate(histogram, F1);
128
129 F1.compile();
130
131
132 JFunction1D_t f1;
133
134 makePDF(histogram, f1);
135
136 f1.compile();
137
138 f1.setExceptionHandler(new JFunction1D_t::JDefaultResult(JMATH::zero));
139 F1.setExceptionHandler(new JFunction1D_t::JDefaultResult(JMATH::zero));
140
141
142 TFile out(outputFile.c_str(), "recreate");
143
144 TH1D h0("h0", NULL, nx, xmin, xmax);
145 TH1D h1("h1", NULL, nx, xmin, xmax);
146 TH1D h2("h2", NULL, nx, xmin, xmax);
147
148 TH1D i0("i0", NULL, nx, xmin, xmax);
149 TH1D i1("i1", NULL, nx, xmin, xmax);
150
151 for (int i = 1; i <= h0.GetNbinsX(); ++i) {
152
153 const Double_t x = h0.GetBinCenter(i);
154
155 h0.SetBinContent(i, g1(x));
156 i0.SetBinContent(i, G1(x));
157
158 try {
159 h1.SetBinContent(i, get_value(F1(x).fp));
160 h2.SetBinContent(i, get_value(f1(x)));
161 i1.SetBinContent(i, get_value(F1(x)));
162 }
163 catch(const exception& error) {
164 ERROR(error.what() << endl);
165 }
166 }
167
168 if (subtract) {
169 h1.Add(&h0, -1.0);
170 h2.Add(&h0, -1.0);
171 i1.Add(&i0, -1.0);
172 }
173
174 out.Write();
175 out.Close();
176}
string outputFile
#define ERROR(A)
Definition JMessage.hh:66
#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
Double_t g1(const Double_t x)
Function.
Definition JQuantiles.cc:25
Double_t G1(const Double_t x)
Integral of method g1.
Definition JQuantiles.cc:37
bool quadrature
Definition JResultPDF.cc:23
int numberOfBins
number of bins for average CDF integral of optical module
Definition JSirene.cc:73
Utility class to parse command line options.
Definition JParser.hh:1698
Numerical integrator for .
Histogram in 1D.
void fill(typename JClass< abscissa_type >::argument_type x, typename JClass< contents_type >::argument_type w)
Fill histogram.
JHistogram1D & div(const double value)
Scale contents.
static const JZero zero
Function object to assign zero value.
Definition JZero.hh:105
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
JElement_t::ordinate_type integrate(const JCollection< JElement_t, JDistance_t > &input, typename JMappable< JElement_t >::map_type &output)
Conversion of data points to integral values.
void makePDF(const JHistogram1D< JElement_t, JContainer_t, JDistance_t > &input, typename JMappable< JElement_t >::map_type &output)
Conversion of histogram to probability density function (PDF).
JGrid< JAbscissa_t > make_grid(const int nx, const JAbscissa_t Xmin, const JAbscissa_t Xmax)
Helper method for JGrid.
Definition JGrid.hh:209
2D Element.
Definition JPolint.hh:1131
Type definition of a spline interpolation method based on a JCollection with JResultPDF result type.