Jpp  17.3.0
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JPlotCDF.cc
Go to the documentation of this file.
1 
2 #include <string>
3 #include <iostream>
4 #include <fstream>
5 #include <iomanip>
6 #include <vector>
7 #include <map>
8 
9 #include "TROOT.h"
10 #include "TFile.h"
11 #include "TH1D.h"
12 #include "TProfile.h"
13 #include "TRandom3.h"
14 
15 #include "JTools/JFunction1D_t.hh"
18 #include "JPhysics/JCDFTable.hh"
19 #include "JPhysics/JPDFTypes.hh"
20 #include "JPhysics/JPDFToolkit.hh"
21 #include "JGeometry3D/JAngle3D.hh"
22 #include "JROOT/JManager.hh"
23 #include "JROOT/JRootToolkit.hh"
24 #include "JGizmo/JGizmoToolkit.hh"
25 #include "Jeep/JTimer.hh"
26 #include "Jeep/JParser.hh"
27 #include "Jeep/JMessage.hh"
28 
29 
30 /**
31  * \file
32  *
33  * Program to verify generation of arrival times of Cherenkov photons from a muon using tabulated CDF.
34  * \author mdejong
35  */
36 int main(int argc, char **argv)
37 {
38  using namespace std;
39  using namespace JPP;
40 
41  typedef JAbstractHistogram<double> JHistogram_t;
42 
43  vector<string> inputFile;
44  string outputFile;
45  double E;
46  double R;
47  JAngle3D dir;
48  int numberOfEvents;
49  JHistogram_t histogram;
50  int debug;
51 
52  try {
53 
54  JParser<> zap("Program to verify generation of arrival times of Cherenkov photons from a muon using tabulated CDF.");
55 
56  zap['f'] = make_field(inputFile);
57  zap['o'] = make_field(outputFile) = "";
58  zap['E'] = make_field(E, "muon energy [GeV]") = 1.0;
59  zap['R'] = make_field(R, "distance of approach [m]");
60  zap['D'] = make_field(dir, "(theta, phi) of PMT [rad]");
61  zap['n'] = make_field(numberOfEvents) = 0;
62  zap['H'] = make_field(histogram, "histogram binning") = JHistogram_t();
63  zap['d'] = make_field(debug) = 2;
64 
65  zap(argc, argv);
66  }
67  catch(const exception& error) {
68  FATAL(error.what() << endl);
69  }
70 
71 
72  typedef JHermiteSplineFunction1D_t JFunction1D_t;
73  //typedef JSplineFunction1D_t JFunction1D_t;
74  //typedef JPolint1Function1D_t JFunction1D_t;
75  typedef JMAPLIST<JPolint1FunctionalMap,
76  JPolint1FunctionalGridMap,
77  JPolint1FunctionalGridMap>::maplist JMapList_t;
78  typedef JCDFTable<JFunction1D_t, JMapList_t> JCDF_t;
79 
80  const int N = inputFile.size();
81 
82  int type[N];
83  JCDF_t cdf [N];
84 
85  try {
86 
87  for (int i = 0; i != N; ++i) {
88 
89  NOTICE("loading input from file " << inputFile[i] << "... " << flush);
90 
91  type[i] = getPDFType(inputFile[i]);
92 
93  cdf [i].load(inputFile[i].c_str());
94 
95  NOTICE("OK" << endl);
96  }
97  }
98  catch(const JException& error) {
99  FATAL(error.what() << endl);
100  }
101 
102 
103  if (outputFile == "") {
104 
105  for ( ; ; ) {
106 
107  double x;
108 
109  cout << "> " << flush;
110  cin >> x;
111 
112  cout << " --> ";
113 
114  for (int i = 0; i != N; ++i) {
115 
116  try {
117 
118  double npe = cdf[i].getNPE (R, dir.getTheta(), dir.getPhi());
119  double t = cdf[i].getTime(R, dir.getTheta(), dir.getPhi(), x);
120 
121  if (is_bremsstrahlung(type[i])) {
122  npe *= E;
123  } else if (is_deltarays(type[i])) {
124  npe *= getDeltaRaysFromMuon(E);
125  }
126 
127  cout << ' ' << FIXED(6,2) << t << ' ' << FIXED(5,2) << npe;
128  }
129  catch(const exception& error) {
130  ERROR(error.what() << endl);
131  }
132  }
133 
134  cout << endl;
135  }
136 
137  return 0;
138  }
139 
140 
141  int function = 0;
142 
143  if (inputFile.size() == 1 && getPDFType(inputFile[0]) == DIRECT_LIGHT_FROM_MUON) {
144  function = 1;
145  }
146 
147  TH1D* h0 = NULL;
148 
149  //const double z0 = -100.0; // [m]
150  //const double t0 = (-z0 + R * getTanThetaC()) / C; // time [ns]
151  const double t0 = 0.0; // time [ns]
152 
153  if (!histogram.is_valid()) {
154 
156 
157  double x = -100.0;
158 
159  {
160  for ( ; x < -10.0; x += 5.0) { X.push_back(t0 + x); }
161  for ( ; x < +20.0; x += 1.0) { X.push_back(t0 + x); }
162  for ( ; x < +50.0; x += 2.0) { X.push_back(t0 + x); }
163  }
164  if (function != 1) {
165  for ( ; x < +100.0; x += 5.0) { X.push_back(t0 + x); }
166  for ( ; x < +250.0; x += 10.0) { X.push_back(t0 + x); }
167  for ( ; x < +500.0; x += 25.0) { X.push_back(t0 + x); }
168  for ( ; x < +900.0; x += 50.0) { X.push_back(t0 + x); }
169  }
170 
171  h0 = new TH1D("h0", NULL, X.size() - 1, X.data());
172 
173  } else {
174 
175  h0 = new TH1D("h0", NULL, histogram.getNumberOfBins(), histogram.getLowerLimit(), histogram.getUpperLimit());
176  }
177 
178  h0->Sumw2();
179 
180  JManager<int, TH1D> H1(new TH1D("h1[%]", NULL, 100000, 0.0, +1.0));
181 
182  for (int i = 0; i != N; ++i) {
183 
184  for (int j = 1; j <= H1->GetNbinsX(); ++j) {
185 
186  try {
187 
188  const double x = H1->GetBinCenter(j);
189  const double t = cdf[i].getTime(R, dir.getTheta(), dir.getPhi(), x);
190 
191  H1[i]->SetBinContent(j, t);
192  }
193  catch(const exception& error) {
194  ERROR(error.what() << endl);
195  }
196  }
197  }
198 
199  if (numberOfEvents > 0) {
200 
201  JTimer timer;
202 
203  timer.reset();
204 
205  for (int counter = 0; counter != numberOfEvents; ++counter) {
206 
207  if (counter%1000== 0) {
208  STATUS(setw(10) << counter << '\r'); DEBUG(endl);
209  }
210 
211  timer.start();
212 
213  for (int i = 0; i != N; ++i) {
214 
215  try {
216 
217  double npe = cdf[i].getNPE(R, dir.getTheta(), dir.getPhi());
218 
219  if (is_bremsstrahlung(type[i])) {
220  npe *= E;
221  } else if (is_deltarays(type[i])) {
222  npe *= getDeltaRaysFromMuon(E);
223  }
224 
225  for (int j = gRandom->Poisson(npe); j != 0; --j) {
226 
227  const double x = gRandom->Rndm();
228  const double t = cdf[i].getTime(R, dir.getTheta(), dir.getPhi(), x);
229 
230  h0->Fill(t + t0);
231  }
232  }
233  catch(const exception& error) {
234  NOTICE(error.what() << endl);
235  }
236  }
237 
238  timer.stop();
239  }
240  STATUS(endl);
241 
242  const double W = 1.0 / (double) numberOfEvents;
243 
244  if (debug > JEEP::status_t) {
245  timer.print(cout, W);
246  }
247 
248  convertToPDF(*h0, "WE", W);
249  }
250 
251  TFile out(outputFile.c_str(), "recreate");
252 
253  out << *h0 << H1;
254 
255  out.Write();
256  out.Close();
257 }
Utility class to parse command line options.
Definition: JParser.hh:1517
int main(int argc, char *argv[])
Definition: Main.cc:15
Auxiliary methods for PDF calculations.
then usage E
Definition: JMuonPostfit.sh:35
#define STATUS(A)
Definition: JMessage.hh:63
then JShowerPostfit f $INPUT_FILE o $OUTPUT_FILE N
Dynamic ROOT object management.
Auxiliary data structure for floating point format specification.
Definition: JManip.hh:446
string outputFile
direct light from muon
Definition: JPDFTypes.hh:26
Various implementations of functional maps.
Numbering scheme for PDF types.
bool is_bremsstrahlung(const int pdf)
Test if given PDF type corresponds to Cherenkov light from Bremsstrahlung.
Definition: JPDFTypes.hh:137
status
Definition: JMessage.hh:30
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition: JParser.hh:1993
void convertToPDF(TH1 &h1, const std::string &option="NW", const double factor=1.0)
Convert 1D histogram to PDF.
#define NOTICE(A)
Definition: JMessage.hh:64
#define ERROR(A)
Definition: JMessage.hh:66
General purpose messaging.
#define FATAL(A)
Definition: JMessage.hh:67
then for FUNCTION in pdf npe cdf
Definition: JPlotNPE-PDG.sh:73
int getPDFType(const std::string &file_name)
Get PDF type.
Definition: JPDFTypes.hh:77
then JCookie sh JDataQuality D $DETECTOR_ID R
Definition: JDataQuality.sh:41
Utility class to parse command line options.
bool is_deltarays(const int pdf)
Test if given PDF type corresponds to Cherenkov light from delta-rays.
Definition: JPDFTypes.hh:151
double getDeltaRaysFromMuon(const double E)
Equivalent EM-shower energy due to delta-rays per unit muon track length.
Definition: JPDFToolkit.hh:67
no fit printf nominal n $STRING awk v X
int j
Definition: JPolint.hh:703
int debug
debug level
#define DEBUG(A)
Message macros.
Definition: JMessage.hh:62