Jpp  16.0.0-rc.1
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JResultPDF.cc
Go to the documentation of this file.
1 
2 #include <string>
3 #include <iostream>
4 #include <iomanip>
5 
6 #include "TROOT.h"
7 #include "TFile.h"
8 #include "TH1D.h"
9 #include "TMath.h"
10 #include "TString.h"
11 
12 #include "JTools/JFunction1D_t.hh"
13 #include "JTools/JQuadrature.hh"
14 #include "JTools/JGrid.hh"
15 
16 #include "Jeep/JParser.hh"
17 #include "Jeep/JMessage.hh"
18 
19 
20 // command line options
21 
24 int debug;
25 
26 
27 namespace {
28 
29  /**
30  * Test function f1().
31  *
32  * \param x abscissa value
33  * \return ordinate value
34  */
35  inline Double_t f1(const Double_t x)
36  {
37  return TMath::Gaus(x, 0.0, 1.0, kTRUE);
38  }
39 
40 
41  /**
42  * Derivative of f1().
43  *
44  * \param x abscissa value
45  * \return derivative value
46  */
47  inline Double_t g1(const Double_t x)
48  {
49  return f1(x) * -x;
50  }
51 
52 
53  /**
54  * Integral of f1().
55  *
56  * \param x abscissa value
57  * \return integral value
58  */
59  inline Double_t G1(const Double_t x)
60  {
61  return 0.5 * (1.0 + TMath::Erf(sqrt(0.5)*x));
62  }
63 
64 
65  /**
66  * Method to test interpolations for PDF.
67  */
68  template<class JFunction1D_t>
69  void test(const char* title)
70  {
71  using namespace JPP;
72 
73  JFunction1D_t buffer;
74 
75  const int nx = 1000;
76  const double xmin = -3.5;
77  const double xmax = +3.5;
78 
79  if (!quadrature)
80  buffer.configure(make_grid(numberOfPoints, xmin, xmax), f1);
81  else
82  buffer.configure(JGaussHermite(numberOfPoints), f1);
83 
84  buffer.compile();
85 
86  TH1D f0(TString("f0[") + title + "]", NULL, nx, xmin, xmax);
87  TH1D g0(TString("g0[") + title + "]", NULL, nx, xmin, xmax);
88  TH1D G0(TString("G0[") + title + "]", NULL, nx, xmin, xmax);
89 
90  for (int i = 1; i <= f0.GetNbinsX(); ++i) {
91 
92  const Double_t x = f0.GetBinCenter(i);
93 
94  try {
95 
96  const typename JFunction1D_t::result_type result = buffer(x);
97 
98  f0.SetBinContent(i, result.f);
99  g0.SetBinContent(i, result.fp);
100  G0.SetBinContent(i, result.v);
101  }
102  catch(const std::exception& error) {
103  //ERROR(error.what() << endl);
104  }
105  }
106 
107  f0.Write();
108  g0.Write();
109  G0.Write();
110  }
111 }
112 
113 /**
114  * \file
115  *
116  * Example program to test interpolations for PDF.
117  * \author mdejong
118  */
119 int main(int argc, char **argv)
120 {
121  using namespace std;
122  using namespace JPP;
123 
124  string outputFile;
125 
126  try {
127 
128  JParser<> zap("Example program to test interpolations for PDF.");
129 
130  zap['o'] = make_field(outputFile);
131  zap['N'] = make_field(numberOfPoints) = 15;
132  zap['Q'] = make_field(quadrature);
133  zap['d'] = make_field(debug) = 3;
134 
135  zap(argc, argv);
136  }
137  catch(const exception &error) {
138  FATAL(error.what() << endl);
139  }
140 
141 
142  TFile out(outputFile.c_str(), "recreate");
143 
144  const int nx = 1000;
145  const double xmin = -3.5;
146  const double xmax = +3.5;
147 
148  TH1D f0(TString("f0[") + "true" + "]", NULL, nx, xmin, xmax);
149  TH1D g0(TString("g0[") + "true" + "]", NULL, nx, xmin, xmax);
150  TH1D G0(TString("G0[") + "true" + "]", NULL, nx, xmin, xmax);
151 
152  for (int i = 1; i <= f0.GetNbinsX(); ++i) {
153 
154  const Double_t x = f0.GetBinCenter(i);
155 
156  f0.SetBinContent(i, f1(x));
157  g0.SetBinContent(i, g1(x));
158  G0.SetBinContent(i, G1(x));
159  }
160 
161  {
162  const int N = 4;
163 
164  typedef JPolintFunction1D<N, JPolintElement2S<double, double>, JCollection, JResultPDF<double> > JFunction1D_t;
165 
166  test<JFunction1D_t>("Polint");
167  }
168 
169  {
170  typedef JSplineFunction1S_t JFunction1D_t;
171 
172  test<JFunction1D_t>("Spline");
173  }
174 
175  out.Write();
176  out.Close();
177 }
Utility class to parse command line options.
Definition: JParser.hh:1500
int main(int argc, char *argv[])
Definition: Main.cc:15
then JShowerPostfit f $INPUT_FILE o $OUTPUT_FILE N
string outputFile
bool quadrature
Definition: JResultPDF.cc:23
Double_t G1(const Double_t x)
Integral of method g1.
Definition: JQuantiles.cc:37
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
return result
Definition: JPolint.hh:743
int debug
debug level
Definition: JSirene.cc:63
General purpose messaging.
#define FATAL(A)
Definition: JMessage.hh:67
Auxiliary classes for numerical integration.
Utility class to parse command line options.
int numberOfPoints
Definition: JResultPDF.cc:22
JGrid< JAbscissa_t > make_grid(const int nx, const JAbscissa_t Xmin, const JAbscissa_t Xmax)
Helper method for JGrid.
Definition: JGrid.hh:177
Double_t g1(const Double_t x)
Function.
Definition: JQuantiles.cc:25