Jpp  18.0.0-rc.3
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JMultiPDF.cc File Reference

Example program to test conditional probability density function $P(x_0,x_1\|x_2,x_3)$ using JTOOLS::JMultiPDF and JTOOLS::JMultiFunction. More...

#include <string>
#include <iostream>
#include <iomanip>
#include "TMath.h"
#include "TRandom.h"
#include "JTools/JHistogram1D_t.hh"
#include "JTools/JHistogramMap_t.hh"
#include "JTools/JMultiHistogram.hh"
#include "JTools/JFunction1D_t.hh"
#include "JTools/JMultiFunction.hh"
#include "JTools/JMultiPDF.hh"
#include "JTools/JToolsToolkit.hh"
#include "JTools/JQuadrature.hh"
#include "JTools/JQuantile.hh"
#include "JIO/JFileStreamIO.hh"
#include "JLang/JObjectIO.hh"
#include "Jeep/JPrint.hh"
#include "Jeep/JTimer.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 conditional probability density function $P(x_0,x_1\|x_2,x_3)$ using JTOOLS::JMultiPDF and JTOOLS::JMultiFunction.

Author
mdejong

Definition in file JMultiPDF.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 69 of file JMultiPDF.cc.

70 {
71  using namespace std;
72 
73  string inputFile;
74  string outputFile;
75  int numberOfEvents;
76  int numberOfBins;
77  int debug;
78 
79  try {
80 
81  JParser<> zap("Example program to test conditional probability density function.");
82 
83  zap['f'] = make_field(inputFile) = "";
84  zap['o'] = make_field(outputFile) = "";
85  zap['n'] = make_field(numberOfEvents);
86  zap['N'] = make_field(numberOfBins) = 15;
87  zap['d'] = make_field(debug) = 2;
88 
89  zap(argc, argv);
90  }
91  catch(const exception &error) {
92  FATAL(error.what() << endl);
93  }
94 
95 
96  if (numberOfEvents <= 0) {
97  FATAL("No events." << endl);
98  }
99 
100  using namespace JPP;
101 
102 
103  const double xmin = -10.0;
104  const double xmax = +10.0;
105  const double dx = (xmax - xmin) / (numberOfBins - 1);
106 
107 
108  if (outputFile != "") {
109 
110  typedef JMultiHistogram<JHistogram1D_t,
111  JMAPLIST<JHistogramGridMap_t,
112  JHistogramGridMap_t,
113  JHistogramMap_t>::maplist> JMultiHistogram_t;
114 
115  JMultiHistogram_t histogram;
116 
117  const JGaussHermite bounds(numberOfBins);
118 
119  for (double x0 = xmin; x0 <= xmax + 0.5*dx; x0 += dx) {
120  for (double x1 = xmin; x1 <= xmax + 0.5*dx; x1 += dx) {
121  for (JGaussHermite::const_iterator x2 = bounds.begin(); x2 != bounds.end(); ++x2) {
122  for (JGaussHermite::const_iterator x3 = bounds.begin(); x3 != bounds.end(); ++x3) {
123  histogram[x0][x1][x0 + x2->getX()][x1 + x3->getX()] = 0.0;
124  }
125  }
126  }
127  }
128 
129  // fill
130 
131  JTimer timer("histogram");
132 
133  for (int i = 0; i != numberOfEvents; ++i) {
134 
135  if (i%1000 == 0) {
136  STATUS("event: " << setw(10) << i << '\r'); DEBUG(endl);
137  }
138 
139  const double x0 = gRandom->Uniform(xmin, xmax);
140  const double x1 = gRandom->Uniform(xmin, xmax);
141  const double x2 = x0 + gRandom->Gaus(0.0, 1.0);
142  const double x3 = x1 + gRandom->Gaus(0.0, 1.0);
143  const double w = 1.0;
144 
145  timer.start();
146 
147  histogram.fill(x0, x1, x2, x3, w);
148 
149  timer.stop();
150  }
151  STATUS(endl);
152 
153  if (debug >= status_t) {
154  timer.print(cout, true, micro_t);
155  }
156 
157  try {
158 
159  NOTICE("storing output to file " << outputFile << "... " << flush);
160 
161  JLANG::store<JIO::JFileStreamWriter>(outputFile.c_str(), histogram);
162 
163  NOTICE("OK" << endl);
164  }
165  catch(const JException& error) {
166  FATAL(error.what() << endl);
167  }
168  }
169 
170 
171  if (inputFile != "") {
172 
173  typedef JMultiHistogram<JHistogram1D_t,
174  JMAPLIST<JHistogramMap_t>::maplist> JHistogram2D_t;
175 
176  typedef JMultiHistogram<JHistogram2D_t,
177  JMAPLIST<JHistogramGridMap_t,
178  JHistogramGridMap_t>::maplist> JMultiHistogram_t;
179 
180  JMultiHistogram_t histogram;
181 
182  try {
183 
184  NOTICE("loading input to file " << inputFile << "... " << flush);
185 
186  JLANG::load<JIO::JFileStreamReader>(inputFile.c_str(), histogram);
187 
188  NOTICE("OK" << endl);
189  }
190  catch(const JException& error) {
191  FATAL(error.what() << endl);
192  }
193 
194 
195  typedef JMultiPDF<JPolint1Function1D_t,
196  JMAPLIST<JPolint1FunctionalMap>::maplist> JFunction2D_t; // PDF
197 
198  typedef JMultiFunction<JFunction2D_t,
199  JMAPLIST<JPolint1FunctionalGridMap,
200  JPolint1FunctionalGridMap>::maplist> JMultiFunction_t; // interpolation
201 
202 
203  JMultiFunction_t pdf(histogram);
204 
205 
206  // test
207 
208  JQuantile Q;
209 
210  for (int i = 0; i != numberOfEvents; ++i) {
211 
212  const double x0 = gRandom->Uniform(xmin, xmax);
213  const double x1 = gRandom->Uniform(xmin, xmax);
214  const double x2 = x0 + gRandom->Gaus(0.0, 1.0);
215  const double x3 = x1 + gRandom->Gaus(0.0, 1.0);
216 
217  try {
218 
219  const double u = g4 (x0, x1, x2, x3);
220  const double v = pdf(x0, x1, x2, x3);
221 
222  Q.put(u - v);
223  }
224  catch(const std::exception& error) {}
225  }
226 
227 
228  JQuantile V;
229 
230  for (JMultiFunction_t::super_const_iterator i = pdf.super_begin(); i != pdf.super_end(); ++i) {
231  V.put(getIntegral((*i).getValue().getMultiFunction()));
232  }
233 
234  cout << "normalisation "
235  << FIXED(6,4) << V.getMean() << " +/- "
236  << FIXED(6,4) << V.getSTDev() << endl;
237  cout << "efficiency "
238  << FIXED(5,3) << (double) Q.getCount() / (double) numberOfEvents << endl;
239 
240  cout << "mean " << SCIENTIFIC(10,2) << Q.getMean() << endl;
241  cout << "RMS " << SCIENTIFIC(10,2) << Q.getSTDev() << endl;
242  }
243 }
const double xmax
Definition: JQuadrature.cc:24
Utility class to parse command line options.
Definition: JParser.hh:1514
data_type w[N+1][M+1]
Definition: JPolint.hh:778
Q(UTCMax_s-UTCMin_s)-livetime_s
micro
Definition: JScale.hh:29
#define STATUS(A)
Definition: JMessage.hh:63
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:446
V(JDAQEvent-JTriggerReprocessor)*1.0/(JDAQEvent+1.0e-10)
string outputFile
status
Definition: JMessage.hh:30
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1989
#define NOTICE(A)
Definition: JMessage.hh:64
JHistogram1D< JElement2D< double, double >, JCollection > JHistogram1D_t
Type definition of a 1 dimensional histogram based on a JCollection.
#define FATAL(A)
Definition: JMessage.hh:67
const double xmin
Definition: JQuadrature.cc:23
data_type v[N+1][M+1]
Definition: JPolint.hh:777
double u[N+1]
Definition: JPolint.hh:776
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:484
int numberOfBins
number of bins for average CDF integral of optical module
Definition: JSirene.cc:69
JContainer_t::ordinate_type getIntegral(const JContainer_t &input)
Get integral of input data points.
int debug
debug level
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62