Jpp  16.0.2
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Functions
JPlotCDF.cc File Reference

Program to verify generation of arrival times of Cherenkov photons from a muon using tabulated CDF. More...

#include <string>
#include <iostream>
#include <fstream>
#include <iomanip>
#include <vector>
#include <map>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "TProfile.h"
#include "TRandom3.h"
#include "JTools/JFunction1D_t.hh"
#include "JTools/JFunctionalMap_t.hh"
#include "JTools/JAbstractHistogram.hh"
#include "JPhysics/JPDFTransformer.hh"
#include "JPhysics/JCDFTable.hh"
#include "JPhysics/JPDFTypes.hh"
#include "JGeometry3D/JAngle3D.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

Program to verify generation of arrival times of Cherenkov photons from a muon using tabulated CDF.

Author
mdejong

Definition in file JPlotCDF.cc.

Function Documentation

int main ( int  argc,
char **  argv 
)

Definition at line 33 of file JPlotCDF.cc.

34 {
35  using namespace std;
36  using namespace JPP;
37 
38  typedef JAbstractHistogram<double> JHistogram_t;
39 
40  string inputFile;
41  string outputFile;
42  double E;
43  double R;
44  JAngle3D dir;
45  int numberOfEvents;
46  JHistogram_t histogram;
47  int debug;
48 
49  try {
50 
51  JParser<> zap("Program to verify generation of arrival times of Cherenkov photons from a muon using tabulated CDF.");
52 
53  zap['f'] = make_field(inputFile);
54  zap['o'] = make_field(outputFile) = "";
55  zap['E'] = make_field(E, "muon energy [GeV]") = 1.0;
56  zap['R'] = make_field(R, "distance of approach [m]");
57  zap['D'] = make_field(dir, "(theta, phi) of PMT [rad]");
58  zap['N'] = make_field(numberOfEvents) = 0;
59  zap['H'] = make_field(histogram, "histogram binning") = JHistogram_t();
60  zap['d'] = make_field(debug) = 0;
61 
62  zap(argc, argv);
63  }
64  catch(const exception& error) {
65  FATAL(error.what() << endl);
66  }
67 
68 
69  typedef JSplineFunction1D_t JFunction1D_t;
70  //typedef JPolint1Function1D_t JFunction1D_t;
71  typedef JMAPLIST<JPolint1FunctionalMap,
72  JPolint1FunctionalGridMap,
73  JPolint1FunctionalGridMap>::maplist JMapList_t;
74  typedef JCDFTable<JFunction1D_t, JMapList_t> JCDF_t;
75 
76  JCDF_t cdf;
77 
78 
79  try {
80 
81  NOTICE("loading input from file " << inputFile << "... " << flush);
82 
83  cdf.load(inputFile.c_str());
84 
85  NOTICE("OK" << endl);
86  }
87  catch(const JException& error) {
88  FATAL(error.what() << endl);
89  }
90 
91 
92  if (outputFile == "") {
93 
94  for ( ; ; ) {
95 
96  double x;
97 
98  cout << "> " << flush;
99  cin >> x;
100 
101  try {
102 
103  const double npe = cdf.getNPE (R, dir.getTheta(), dir.getPhi()) * E;
104  const double t = cdf.getTime(R, dir.getTheta(), dir.getPhi(), x);
105 
106  cout << " --> " << t << ' ' << npe;
107  }
108  catch(const exception& error) {
109  ERROR(error.what() << endl);
110  }
111 
112  cout << endl;
113  }
114 
115  } else {
116 
117 
118  TFile out(outputFile.c_str(), "recreate");
119 
120  int function = 1;
121 
122  if (inputFile.find(getLabel(DIRECT_LIGHT_FROM_MUON)) == string::npos) {
123  function = 0;
124  }
125 
126  //const double z0 = -100.0; // [m]
127  //const double t0 = (-z0 + R * getTanThetaC()) / C; // time [ns]
128  const double t0 = 0.0; // time [ns]
129 
130  if (!histogram.is_valid()) {
131 
132  if (function == 0) {
133 
134  histogram = JHistogram_t(t0 - 20.0, t0 + 50.0);
135 
136  histogram.setBinWidth(0.1);
137 
138  } else {
139 
140  histogram = JHistogram_t(t0 - 20.0, t0 + 500.0);
141 
142  histogram.setBinWidth(0.5);
143  }
144  }
145 
146  TH1D h0("h0", NULL, histogram.getNumberOfBins(), histogram.getLowerLimit(), histogram.getUpperLimit());
147  TH1D h1("h1", NULL, 100000, 0.0, +1.0);
148 
149  h0.Sumw2();
150 
151  for (int j = 1; j <= h1.GetNbinsX(); ++j) {
152 
153  try {
154 
155  const double x = h1.GetBinCenter(j);
156  const double t = cdf.getTime(R, dir.getTheta(), dir.getPhi(), x);
157 
158  h1.SetBinContent(j, t);
159  }
160  catch(const exception& error) {
161  ERROR(error.what() << endl);
162  }
163  }
164 
165  if (numberOfEvents > 0) {
166 
167 
168  JTimer timer;
169 
170  timer.reset();
171 
172  for (int i = 0; i != numberOfEvents; ++i) {
173 
174  if (i%1000== 0) { NOTICE(i << '\r'); }
175 
176  timer.start();
177 
178  try {
179 
180  const double npe = cdf.getNPE(R, dir.getTheta(), dir.getPhi()) * E;
181 
182  for (int j = gRandom->Poisson(npe); j != 0; --j) {
183 
184  const double x = gRandom->Rndm();
185  const double t = cdf.getTime(R, dir.getTheta(), dir.getPhi(), x);
186 
187  h0.Fill(t + t0);
188  }
189  }
190  catch(const exception& error) {
191  NOTICE(error.what() << endl);
192  }
193 
194  timer.stop();
195  }
196  NOTICE(endl);
197 
198  const double w = 1.0 / (double) numberOfEvents;
199 
200  if (debug > JEEP::notice_t) {
201  timer.print(cout, w);
202  }
203 
204  // normalise histogram (dP/dt)
205 
206  for (int j = 1; j <= h0.GetNbinsX(); ++j) {
207 
208  const double y = h0.GetBinContent(j);
209  const double z = h0.GetBinError (j);
210  const double dx = h0.GetBinWidth (j);
211 
212  h0.SetBinContent(j, y / (dx*numberOfEvents*E));
213  h0.SetBinError (j, z / (dx*numberOfEvents*E));
214  }
215  }
216 
217  out.Write();
218  out.Close();
219  }
220 }
Utility class to parse command line options.
Definition: JParser.hh:1500
data_type w[N+1][M+1]
Definition: JPolint.hh:757
then usage E
Definition: JMuonPostfit.sh:35
std::string getLabel(const JLocation &location)
Get module label for monitoring and other applications.
Definition: JLocation.hh:246
notice
Definition: JMessage.hh:32
then for HISTOGRAM in h0 h1
Definition: JMatrixNZ.sh:71
string outputFile
direct light from muon
Definition: JPDFTypes.hh:29
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1961
#define NOTICE(A)
Definition: JMessage.hh:64
#define ERROR(A)
Definition: JMessage.hh:66
int debug
debug level
Definition: JSirene.cc:63
then usage $script[distance] fi case set_variable R
Definition: JDrawLED.sh:43
#define FATAL(A)
Definition: JMessage.hh:67
then for FUNCTION in pdf npe cdf
Definition: JPlotNPE-PDG.sh:73
int j
Definition: JPolint.hh:682