Jpp  15.0.1
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JTools/JHistogram2D.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 "TH2D.h"
9 #include "TMath.h"
10 #include "TRandom.h"
11 
12 #include "JTools/JHistogram1D_t.hh"
15 #include "JTools/JCollection.hh"
16 #include "JTools/JMultiGrid.hh"
17 #include "JTools/JMultiSet.hh"
18 #include "JTools/JQuadrature.hh"
19 #include "JTools/JFunction1D_t.hh"
21 #include "JTools/JToolsToolkit.hh"
22 #include "JTools/JMultiPDF.hh"
23 
24 #include "Jeep/JParser.hh"
25 #include "Jeep/JMessage.hh"
26 
27 
28 namespace {
29  /**
30  * Function g1.
31  *
32  * \param x x value
33  * \return function value
34  */
35  inline Double_t g1(const Double_t x)
36  {
37  return TMath::Gaus(x, 0.0, 1.0, kTRUE);
38  }
39 
40 
41  /**
42  * Integral of g1.
43  *
44  * \param x x value
45  * \return integral value
46  */
47  inline Double_t G1(const Double_t x)
48  {
49  return 0.5 * (1.0 + TMath::Erf(sqrt(0.5)*x));
50  }
51 
52  /**
53  * Function g2.
54  *
55  * \param x x value
56  * \param y y value
57  * \return function value
58  */
59  inline Double_t g2(const Double_t x, const Double_t y)
60  {
61  return g1(x) * g1(y);
62  }
63 
64 
65  /**
66  * Integral of g2.
67  *
68  * \param x x value
69  * \param y y value
70  * \return integral value
71  */
72  inline Double_t G2(const Double_t x, const Double_t y)
73  {
74  return G1(x) * G1(y);
75  }
76 }
77 
78 
79 /**
80  * \file
81  *
82  * Example program to test JTOOLS::JMultiHistogram.
83  * \author mdejong
84  */
85 int main(int argc, char **argv)
86 {
87  using namespace std;
88 
89  string outputFile;
90  int numberOfEvents;
91  int numberOfBins;
92  bool quadrature;
93  bool subtract;
94  int debug;
95 
96  try {
97 
98  JParser<> zap("Example program to test 2D histogram.");
99 
100  zap['o'] = make_field(outputFile) = "histogram.root";
101  zap['n'] = make_field(numberOfEvents);
102  zap['N'] = make_field(numberOfBins);
103  zap['Q'] = make_field(quadrature);
104  zap['D'] = make_field(subtract);
105  zap['d'] = make_field(debug) = 3;
106 
107  zap(argc, argv);
108  }
109  catch(const exception &error) {
110  FATAL(error.what() << endl);
111  }
112 
113 
114  if (numberOfEvents <= 0) {
115  FATAL("No events." << endl);
116  }
117 
118  using namespace JPP;
119 
120  const Int_t nx = 100;
121  const Double_t xmin = -3.0;
122  const Double_t xmax = +3.0;
123 
124  typedef JMultiHistogram<JHistogram1D_t,
125  JMapList<JHistogramMap_t> > JHistogram2D_t;
126 
127  JHistogram2D_t histogram;
128 
129 
130  // abscissa
131 
132  if (!quadrature)
133  configure(histogram, make_grid(numberOfBins + 1, xmin, xmax));
134  else
135  configure(histogram, make_set(JGaussHermite(numberOfBins + 1)));
136 
137 
138  // fill
139 
140  for (int i = 0; i != numberOfEvents; ++i) {
141 
142  const double x = gRandom->Gaus(0.0, 1.0);
143  const double y = gRandom->Gaus(0.0, 1.0);
144 
145  histogram.fill(x, y, 1.0);
146  }
147 
148 
149  histogram.div((double) numberOfEvents);
150 
151 
152  // interpolation based on function values
153 
154  JMultiPDF<JPolint2Function1D_t,
155  JMapList<JPolint2FunctionalMap> > f2(histogram);
156 
157 
158  // interpolation based on integral values
159 
160  accumulate(histogram);
161 
162  JMultiFunction<JPolint3Function1H_t,
163  JMapList<JPolint3FunctionalMapH> > F2;
164 
165  copy(histogram, F2);
166 
167  F2.compile();
168 
169 
170  TFile out(outputFile.c_str(), "recreate");
171 
172  TH2D h0("h0", NULL, nx, xmin, xmax, nx, xmin, xmax);
173  TH2D h1("h1", NULL, nx, xmin, xmax, nx, xmin, xmax);
174  TH2D h2("h2", NULL, nx, xmin, xmax, nx, xmin, xmax);
175 
176 
177  for (int i = 1; i <= h0.GetXaxis()->GetNbins(); ++i) {
178  for (int j = 1; j <= h0.GetYaxis()->GetNbins(); ++j) {
179 
180  const Double_t x = h0.GetXaxis()->GetBinCenter(i);
181  const Double_t y = h0.GetYaxis()->GetBinCenter(j);
182 
183  h0.SetBinContent(i, j, g2(x,y));
184 
185  try {
186  h1.SetBinContent(i, j, F2(x,y).fp.fp);
187  h2.SetBinContent(i, j, f2(x,y));
188  }
189  catch(const exception& error) {
190  //ERROR(error.what() << endl);
191  }
192  }
193  }
194 
195  if (subtract) {
196  h1.Add(&h0, -1.0);
197  h2.Add(&h0, -1.0);
198  }
199 
200  out.Write();
201  out.Close();
202 }
203 
Utility class to parse command line options.
Definition: JParser.hh:1500
int main(int argc, char *argv[])
Definition: Main.cc:15
then for HISTOGRAM in h0 h1
Definition: JMatrixNZ.sh:71
string outputFile
Various implementations of functional maps.
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
This include file contains various recursive methods to operate on multi-dimensional collections...
JHistogram1D< JElement2D< double, double >, JCollection > JHistogram1D_t
Type definition of a 1 dimensional histogram based on a JCollection.
General purpose class for a collection of sorted elements.
int debug
debug level
Definition: JSirene.cc:63
T make_set(T __begin, T __end, JResult_t std::iterator_traits< T >::value_type::*value, const JComparator_t &comparator)
Method to exclude outliers from already sorted data.
Definition: JVectorize.hh:169
General purpose messaging.
#define FATAL(A)
Definition: JMessage.hh:67
Auxiliary classes for numerical integration.
void configure(const T &value, const JAbstractCollection< JAbscissa_t > &bounds, JBool< false > option)
Configuration of value.
Utility class to parse command line options.
void copy(const Head &from, JHead &to)
Copy header from from to to.
Definition: JHead.cc:139
int j
Definition: JPolint.hh:666
JGrid< JAbscissa_t > make_grid(const int nx, const JAbscissa_t Xmin, const JAbscissa_t Xmax)
Helper method for JGrid.
Definition: JGrid.hh:177
int numberOfBins
number of bins for average CDF integral of optical module
Definition: JSirene.cc:64
void accumulate(T &value, JBool< false > option)
Accumulation of value.
Double_t g1(const Double_t x)
Function.
Definition: JQuantiles.cc:25