Jpp 20.0.0-195-g190c9e876
the software that should make you happy
Loading...
Searching...
No Matches
JDrawPDF.cc File Reference

Auxiliary program to draw PDF of Cherenkov light from muon. More...

#include <string>
#include <iostream>
#include <iomanip>
#include "TROOT.h"
#include "TFile.h"
#include "TH1D.h"
#include "JPhysics/JPDF.hh"
#include "JPhysics/Antares.hh"
#include "JPhysics/KM3NeT.hh"
#include "JPhysics/JPDFSupportkit.hh"
#include "JPhysics/JGeane.hh"
#include "JTools/JFunction1D_t.hh"
#include "JTools/JQuantiles.hh"
#include "JTools/JAbstractHistogram.hh"
#include "JGeometry3D/JAngle3D.hh"
#include "Jeep/JProperties.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 draw PDF of Cherenkov light from muon.

Author
mdejong, bofearraigh

Definition in file JDrawPDF.cc.

Function Documentation

◆ main()

int main ( int argc,
char ** argv )

Definition at line 31 of file JDrawPDF.cc.

32{
33 using namespace std;
34 using namespace JPP;
35
37
38 string outputFile;
40 double epsilon;
41 JAbsorptionLength absorptionLength;
42 JScatteringLength scatteringLength;
43 JScatteringProbability scatteringProbability;
44 double E;
45 double R;
46 double z;
47 JAngle3D dir;
48 vector<int> function;
49 JHistogram_t histogram;
50 int debug;
51
52 try {
53
54 JProperties properties;
55
56 properties.insert(gmake_property(MODULE_RADIUS_M));
57 properties.insert(gmake_property(absorptionLength));
58 properties.insert(gmake_property(scatteringLength));
59 properties.insert(gmake_property(scatteringProbability));
60
61 JParser<> zap("Auxiliary program to draw PDF of Cherenkov light from muon.");
62
63 zap['@'] = make_field(properties, endl
64 << "possible options absorptionLength: " << get_keys(absorptionLength) << endl
65 << "possible options scatteringLength: " << get_keys(scatteringLength) << endl
66 << "possible options scatteringProbability: " << get_keys(scatteringProbability) << endl) = JPARSER::initialised();
67 zap['o'] = make_field(outputFile) = "";
68 zap['n'] = make_field(numberOfPoints, "points for integration") = 25;
69 zap['e'] = make_field(epsilon, "precision for integration") = 1.0e-10;
70 zap['E'] = make_field(E, "muon energy at vertex [GeV]") = 1.0;
71 zap['R'] = make_field(R, "distance of approach [m]");
72 zap['z'] = make_field(z, "PMT z-position [m]");
73 zap['D'] = make_field(dir, "(theta, phi) of PMT [rad]");
74 zap['F'] = make_field(function, "PDF type");
75 zap['H'] = make_field(histogram, "histogram binning") = JHistogram_t();
76 zap['d'] = make_field(debug) = 0;
77
78 zap(argc, argv);
79 }
80 catch(const exception &error) {
81 FATAL(error.what() << endl);
82 }
83
84
85 const JPDF_C
86 pdf(NAMESPACE::getPhotocathodeArea(),
87 NAMESPACE::getQE,
88 NAMESPACE::getAngularAcceptance,
89 JAbsorptionLength::getAbsorptionLength,
90 JScatteringLength::getScatteringLength,
91 JScatteringProbability::getScatteringProbability,
92 NAMESPACE::getAmbientPressure(),
96 epsilon);
97
98 const double z_0 = z - R/getTanThetaC(); // emission point z-position
99 const double E_0 = gWater.getE(E, z_0); // energy of muon at emission point
100
101
102 if (outputFile == "") {
103
104 cout << "enter time (^C to exit) > " << flush;
105
106 for (double dt; cin >> dt; ) {
107
108 for (vector<int>::const_iterator F = function.begin(); F != function.end(); ++F) {
109
110 cout << setw(2) << *F << ' '
111 << SCIENTIFIC(7,1) << E << ' '
112 << FIXED(5,1) << R << ' '
113 << FIXED(5,1) << z << ' '
114 << FIXED(5,2) << dir.getTheta() << ' '
115 << FIXED(5,2) << dir.getPhi() << ' '
116 << FIXED(5,1) << dt << ' '
117 << SCIENTIFIC(9,3) << pdf.getLightFromMuon(*F, E_0, R, dir.getTheta(), dir.getPhi(), dt) << endl;
118 }
119 }
120
121 return 0;
122 }
123
124
125 TFile out(outputFile.c_str(), "recreate");
126
127 const double t0 = 0.0; // time [ns]
128
129 if (!histogram.is_valid()) {
130
131 if (function.size() == 1 && function[0] == DIRECT_LIGHT_FROM_MUON) {
132
133 histogram = JHistogram_t(t0 - 20.0, t0 + 50.0);
134
135 histogram.setBinWidth(0.1);
136
137 } else {
138
139 histogram = JHistogram_t(t0 - 20.0, t0 + 500.0);
140
141 histogram.setBinWidth(0.5);
142 }
143 }
144
145 TH1D h0("h0", NULL, histogram.getNumberOfBins(), histogram.getLowerLimit(), histogram.getUpperLimit());
146
148
149 if (E_0 >= MASS_MUON* (1.0/SIN_THETA_C_WATER)) { // muon emission energy has to be above Cherenkov threshold
150
151 if (z_0 >= 0 && z_0 <= gWater(E)) { // emission point is between start and end point of muon
152
153 for (int i = 1; i <= h0.GetNbinsX(); ++i) {
154
155 const double dt = h0.GetBinCenter(i) - t0;
156
157 double value = 0.0;
158
159 for (vector<int>::const_iterator F = function.begin(); F != function.end(); ++F) {
160 value += pdf.getLightFromMuon(*F, E_0, R, dir.getTheta(), dir.getPhi(), dt);
161 }
162
163 h0.SetBinContent(i, value);
164
165 f1[dt] = value;
166 }
167 }
168 }
169
170 f1.setExceptionHandler(new JSplineFunction1S_t::JDefaultResult(JMATH::zero));
171 f1.compile();
172
173 try {
174
175 const double T_ns = 5; // [ns]
176
177 JQuantiles quantiles(f1);
178
179 const double t1 = quantiles.getX();
180 const double y = f1(t1 + T_ns).v - f1(t1 - T_ns).v;
181
182 DEBUG("E " << E << endl);
183 DEBUG("E_0 " << E_0 << endl);
184 DEBUG("R " << R << endl);
185 DEBUG("z " << z << endl);
186 DEBUG("theta " << dir.getTheta() << endl);
187 DEBUG("phi " << dir.getPhi() << endl);
188 DEBUG("int " << quantiles.getIntegral() << endl);
189 DEBUG("t1 " << t1 << endl);
190 DEBUG("max " << quantiles.getY() << endl);
191 DEBUG("FWHM " << quantiles.getFWHM() << endl);
192 DEBUG("int[] " << y << endl);
193 }
194 catch(const exception&) {}
195
196 out.Write();
197 out.Close();
198}
string outputFile
#define DEBUG(A)
Message macros.
Definition JMessage.hh:62
#define FATAL(A)
Definition JMessage.hh:67
int debug
debug level
Definition JSirene.cc:74
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Definition JParser.hh:2142
#define gmake_property(A)
macros to convert (template) parameter to JPropertiesElement object
int numberOfPoints
Definition JResultPDF.cc:22
Utility class to parse parameter values.
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
Utility class to parse command line options.
Definition JParser.hh:1698
Probability Density Functions of the time response of a PMT with an implementation of the JAbstractPM...
Definition JPDF.hh:2186
Quantile calculator for a given interpolating function.
Definition JQuantiles.hh:34
T getLowerLimit() const
Get lower limit.
Definition JRange.hh:202
T getUpperLimit() const
Get upper limit.
Definition JRange.hh:213
const array_type< JKey_t > & get_keys(const std::map< JKey_t, JValue_t, JComparator_t, JAllocator_t > &data)
Method to create array of keys of map.
static const JZero zero
Function object to assign zero value.
Definition JZero.hh:105
double getMinimalWavelength()
Get minimal wavelength for PDF evaluations.
double getMaximalWavelength()
Get maximal wavelength for PDF evaluations.
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
Empty structure for specification of parser element that is initialised (i.e. does not require input)...
Definition JParser.hh:68
Auxiliary data structure to customize absorption length.
Auxiliary data structure to customize scattering length.
Auxiliary data structure to customize scattering probability.
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.
void setExceptionHandler(const JSupervisor &supervisor)
Set the supervisor for handling of exceptions.
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