Jpp test-rotations-old-533-g2bdbdb559
the software that should make you happy
Loading...
Searching...
No Matches
JPlotPDG.cc File Reference

Auxiliary program to plot PDF of Cherenkov light from EM-shower using interpolation tables. More...

#include <string>
#include <iostream>
#include <fstream>
#include <iomanip>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "JTools/JFunction1D_t.hh"
#include "JTools/JFunctionalMap_t.hh"
#include "JTools/JAbstractHistogram.hh"
#include "JPhysics/JPDFTransformer.hh"
#include "JPhysics/JPDFTable.hh"
#include "JPhysics/JPDFTypes.hh"
#include "JPhysics/JGeanz.hh"
#include "JGeometry3D/JAngle3D.hh"
#include "Jeep/JPrint.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

Auxiliary program to plot PDF of Cherenkov light from EM-shower using interpolation tables.

Author
mdejong

Definition in file JPlotPDG.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 30 of file JPlotPDG.cc.

31{
32 using namespace std;
33 using namespace JPP;
34
36
37 vector<string> inputFile;
38 string outputFile;
39 double E;
40 double D;
41 double cd;
42 JAngle3D dir;
43 double TTS_ns;
44 size_t profile;
45 bool cms;
46 JHistogram_t histogram;
47 int debug;
48
49 try {
50
51 JParser<> zap("Auxiliary program to plot PDF of Cherenkov light from EM-shower using interpolation tables.");
52
53 zap['f'] = make_field(inputFile);
54 zap['o'] = make_field(outputFile) = "";
55 zap['E'] = make_field(E, "shower energy [GeV]") = 1.0;
56 zap['R'] = make_field(D, "distance [m]");
57 zap['c'] = make_field(cd, "cosine emission angle");
58 zap['D'] = make_field(dir, "(theta, phi) of PMT [rad]");
59 zap['T'] = make_field(TTS_ns, "PMT time smearing [ns]") = 0.0; // [ns]
60 zap['N'] = make_field(profile, "steps shower profile") = 0;
61 zap['C'] = make_field(cms, "center shower profile");
62 zap['H'] = make_field(histogram, "histogram binning") = JHistogram_t();
63 zap['d'] = make_field(debug) = 0;
64
65 zap(argc, argv);
66 }
67 catch(const exception &error) {
68 FATAL(error.what() << endl);
69 }
70
71
72 typedef JSplineFunction1S_t JFunction1D_t;
76 JPolint1FunctionalGridMap>::maplist JMapList_t;
78
79 const int N = inputFile.size();
80
81 int type[N];
82 JPDF_t pdf [N];
83
84 try {
85
86 for (int i = 0; i != N; ++i) {
87
88 NOTICE("loading input from file " << inputFile[i] << "... " << flush);
89
90 type[i] = getPDFType(inputFile[i]);
91
92 pdf [i].load(inputFile[i].c_str());
93
94 pdf [i].setExceptionHandler(new JFunction1D_t::JDefaultResult(JMATH::zero));
95
96 if (TTS_ns > 0.0) {
97 pdf[i].blur(TTS_ns);
98 }
99
100 NOTICE("OK" << endl);
101 }
102 }
103 catch(const JException& error) {
104 FATAL(error.what() << endl);
105 }
106
107
108 if (outputFile == "") {
109
110 cout << "enter time (^C to exit) > " << flush;
111
112 for (double dt; cin >> dt; ) {
113
114 for (int i = 0; i != N; ++i) {
115
116 JFunction1D_t::result_type y = pdf[i](D, cd, dir.getTheta(), dir.getPhi(), dt);
117
118 cout << setw(2) << type[i] << ' '
119 << SCIENTIFIC(7,1) << E << ' '
120 << FIXED(5,1) << D << ' '
121 << FIXED(5,2) << cd << ' '
122 << FIXED(5,2) << dir.getTheta() << ' '
123 << FIXED(5,2) << dir.getPhi() << ' '
124 << FIXED(5,1) << dt << ' '
125 << SCIENTIFIC(9,3) << get_value(y) << ' '
126 << SCIENTIFIC(9,3) << get_derivative(y) << ' '
127 << SCIENTIFIC(9,3) << get_integral(y) << ' '
128 << SCIENTIFIC(9,3) << get_total_integral(y) << endl;
129 }
130 }
131
132 return 0;
133 }
134
135
136 TFile out(outputFile.c_str(), "recreate");
137
138 int function = 0;
139
140 if (inputFile.size() == 1 &&
141 inputFile.begin()->find(getLabel(SCATTERED_LIGHT_FROM_EMSHOWER)) == string::npos) {
142 function = 1;
143 }
144
145 //const double t0 = D * getIndexOfRefraction() / C; // time [ns]
146 const double t0 = 0.0; // time [ns]
147
148 if (!histogram.is_valid()) {
149
150 if (function == 1) {
151
152 histogram = JHistogram_t(t0 - 20.0, t0 + 50.0);
153
154 histogram.setBinWidth(0.1);
155
156 } else {
157
158 histogram = JHistogram_t(t0 - 20.0, t0 + 500.0);
159
160 histogram.setBinWidth(0.5);
161 }
162 }
163
164 TH1D h0("h0", NULL, histogram.getNumberOfBins(), histogram.getLowerLimit(), histogram.getUpperLimit());
165
166 const double D0 = D;
167 const double z0 = D0 * cd;
168 const double R = D0 * sqrt((1.0 + cd) * (1.0 - cd));
169
170 double W = 1.0;
171
173
174 if (profile >= 2) {
175
176 const double z1 = (cms ? geanz.getMaximum(E) :0.0);
177
178 W = 1.0 / (double) profile;
179
180 for (double x = 0.5*W; x < 1.0; x += W) {
181 Z.push_back(geanz.getLength(E, x) - z1);
182 }
183
184 } else {
185
186 W = 1.0;
187 Z = { 0.0 };
188 }
189
190 for (int i = 1; i <= h0.GetNbinsX(); ++i) {
191
192 const double dt = h0.GetBinCenter(i) - t0;
193
194 JFunction1D_t::result_type Y = JMATH::zero;
195
196 for (const double z1 : Z) {
197
198 const double z = z0 - z1;
199
200 D = sqrt(R*R + z*z);
201 cd = z / D;
202
203 const double t1 = dt - (z1 + (D - D0) * getIndexOfRefraction()) / C;
204
205 for (int j = 0; j != N; ++j) {
206 Y += pdf[j](D, cd, dir.getTheta(), dir.getPhi(), t1) * E * W;
207 }
208
209 DEBUG("sum: "
210 << FIXED(7,2) << dt << ' '
211 << FIXED(7,2) << z1 << ' '
212 << FIXED(7,2) << t1 << ' '
213 << SCIENTIFIC(12,3) << get_value(Y) << endl);
214 }
215 DEBUG(endl);
216
217 h0.SetBinContent(i, get_value(Y));
218 }
219
220 out.Write();
221 out.Close();
222}
string outputFile
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#define NOTICE(A)
Definition JMessage.hh:64
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:72
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
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
virtual const char * what() const override
Get error message.
Definition JException.hh:64
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
std::string getLabel(const JLocation &location)
Get module label for monitoring and other applications.
Definition JLocation.hh:247
static const JZero zero
Function object to assign zero value.
Definition JZero.hh:105
int getPDFType(const std::string &file_name)
Get PDF type.
Definition JPDFTypes.hh:77
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
JResultEvaluator< JResult_t >::result_type get_derivative(const JResult_t &value)
Helper method to convert function value to derivative value.
Definition JResult.hh:1011
int j
Definition JPolint.hh:801
JResultEvaluator< JResult_t >::result_type get_total_integral(const JResult_t &value)
Helper method to convert function value to total integral value.
Definition JResult.hh:1037
JResultEvaluator< JResult_t >::result_type get_integral(const JResult_t &value)
Helper method to convert function value to partial integral value.
Definition JResult.hh:1024
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