Jpp 19.3.0-rc.1
the software that should make you happy
Loading...
Searching...
No Matches
JPlotPDG.cc
Go to the documentation of this file.
1
2#include <string>
3#include <iostream>
4#include <fstream>
5#include <iomanip>
6
7#include "TROOT.h"
8#include "TFile.h"
9#include "TH1D.h"
10
15#include "JPhysics/JPDFTable.hh"
16#include "JPhysics/JPDFTypes.hh"
18#include "Jeep/JPrint.hh"
19#include "Jeep/JParser.hh"
20#include "Jeep/JMessage.hh"
21
22
23/**
24 * \file
25 *
26 * Program to plot PDF of Cherenkov light from EM-shower or scattered light from muon using interpolation tables.
27 * \author mdejong
28 */
29int main(int argc, char **argv)
30{
31 using namespace std;
32 using namespace JPP;
33
35
36 vector<string> inputFile;
37 string outputFile;
38 double E;
39 double D;
40 double cd;
41 JAngle3D dir;
42 double TTS_ns;
43 JHistogram_t histogram;
44 int debug;
45
46 try {
47
48 JParser<> zap("Program to plot PDF of Cherenkov light from EM-shower or scattered light from muon using interpolation tables.");
49
50 zap['f'] = make_field(inputFile);
51 zap['o'] = make_field(outputFile) = "";
52 zap['E'] = make_field(E, "shower energy [GeV]") = 1.0;
53 zap['R'] = make_field(D, "distance [m]");
54 zap['c'] = make_field(cd, "cosine emission angle");
55 zap['D'] = make_field(dir, "(theta, phi) of PMT [rad]");
56 zap['T'] = make_field(TTS_ns, "PMT time smearing [ns]") = 0.0; // [ns]
57 zap['H'] = make_field(histogram, "histogram binning") = JHistogram_t();
58 zap['d'] = make_field(debug) = 0;
59
60 zap(argc, argv);
61 }
62 catch(const exception &error) {
63 FATAL(error.what() << endl);
64 }
65
66
67 typedef JSplineFunction1S_t JFunction1D_t;
71 JPolint1FunctionalGridMap>::maplist JMapList_t;
73
74 const int N = inputFile.size();
75
76 int type[N];
77 JPDF_t pdf [N];
78
79 try {
80
81 for (int i = 0; i != N; ++i) {
82
83 NOTICE("loading input from file " << inputFile[i] << "... " << flush);
84
85 type[i] = getPDFType(inputFile[i]);
86
87 pdf [i].load(inputFile[i].c_str());
88
89 pdf [i].setExceptionHandler(new JFunction1D_t::JDefaultResult(JMATH::zero));
90
91 if (TTS_ns > 0.0) {
92 pdf[i].blur(TTS_ns);
93 }
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 (double dt; cin >> dt; ) {
106
107 for (int i = 0; i != N; ++i) {
108
109 JFunction1D_t::result_type y = pdf[i](D, cd, dir.getTheta(), dir.getPhi(), dt);
110
111 cout << setw(2) << type[i] << ' '
112 << SCIENTIFIC(7,1) << E << ' '
113 << FIXED(5,1) << D << ' '
114 << FIXED(5,2) << cd << ' '
115 << FIXED(5,2) << dir.getTheta() << ' '
116 << FIXED(5,2) << dir.getPhi() << ' '
117 << FIXED(5,1) << dt << ' '
118 << SCIENTIFIC(9,3) << get_value(y) << ' '
119 << SCIENTIFIC(9,3) << get_derivative(y) << ' '
120 << SCIENTIFIC(9,3) << get_integral(y) << ' '
121 << SCIENTIFIC(9,3) << get_total_integral(y) << endl;
122 }
123 }
124
125 return 0;
126 }
127
128
129 TFile out(outputFile.c_str(), "recreate");
130
131 int function = 0;
132
133 if (inputFile.size() == 1 &&
134 inputFile.begin()->find(getLabel(SCATTERED_LIGHT_FROM_EMSHOWER)) == string::npos) {
135 function = 1;
136 }
137
138 //const double t0 = D * getIndexOfRefraction() / C; // time [ns]
139 const double t0 = 0.0; // time [ns]
140
141 if (!histogram.is_valid()) {
142
143 if (function == 1) {
144
145 histogram = JHistogram_t(t0 - 20.0, t0 + 50.0);
146
147 histogram.setBinWidth(0.1);
148
149 } else {
150
151 histogram = JHistogram_t(t0 - 20.0, t0 + 500.0);
152
153 histogram.setBinWidth(0.5);
154 }
155 }
156
157 TH1D h0("h0", NULL, histogram.getNumberOfBins(), histogram.getLowerLimit(), histogram.getUpperLimit());
158
159 for (int i = 1; i <= h0.GetNbinsX(); ++i) {
160
161 const double dt = h0.GetBinCenter(i) - t0;
162
163 JFunction1D_t::result_type Y = JMATH::zero;
164
165 for (int j = 0; j != N; ++j) {
166 Y += pdf[j](D, cd, dir.getTheta(), dir.getPhi(), dt) * E;
167 }
168
169 h0.SetBinContent(i, get_value(Y));
170 }
171
172 out.Write();
173 out.Close();
174}
string outputFile
Various implementations of functional maps.
General purpose messaging.
#define NOTICE(A)
Definition JMessage.hh:64
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:72
Numbering scheme for PDF types.
Utility class to parse command line options.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
int main(int argc, char **argv)
Definition JPlotPDG.cc:29
I/O formatting auxiliaries.
Data structure for angles in three dimensions.
Definition JAngle3D.hh:35
double getTheta() const
Get theta angle.
Definition JAngle3D.hh:86
double getPhi() const
Get phi angle.
Definition JAngle3D.hh:97
General exception.
Definition JException.hh:24
Utility class to parse command line options.
Definition JParser.hh:1698
Multi-dimensional PDF table for arrival time of Cherenkov light.
Definition JPDFTable.hh:44
T getLowerLimit() const
Get lower limit.
Definition JRange.hh:202
T getUpperLimit() const
Get upper limit.
Definition JRange.hh:213
static const JZero zero
Function object to assign zero value.
Definition JZero.hh:105
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Auxiliary data structure for floating point format specification.
Definition JManip.hh:448
Simple data structure for histogram binning.
void setBinWidth(const abscissa_type dx, int option=0)
Set bin width.
bool is_valid() const
Check validity of histogram binning.
int getNumberOfBins() const
Get number of bins.
Auxiliary class for recursive map list generation.
Definition JMapList.hh:109
Type definition of a 1st degree polynomial interpolation based on a JGridMap implementation.
Type definition of a 1st degree polynomial interpolation based on a JMap implementation.
Type definition of a spline interpolation method based on a JCollection with JResultPDF result type.
Auxiliary data structure for floating point format specification.
Definition JManip.hh:488