PDF types.
PDF types.
35{
38
39 string fileDescriptor;
41 string option;
42 double E;
43 double cd;
45 double TMax_ns;
47
48 try {
49
51
54 zap[
'O'] =
make_field(option) =
"KM3NeT",
"Antares";
55 zap[
'E'] =
make_field(E,
"muon/shower energy [GeV]");
56 zap[
'c'] =
make_field(cd,
"cosine emission angle");
57 zap[
'D'] =
make_field(dir,
"(theta, phi) of PMT [rad]");
58 zap[
'T'] =
make_field(TMax_ns,
"L1 coincidence time window [ns]") = 20.0;
60
62
63 zap(argc, argv);
64 }
65 catch(const exception &error) {
66 FATAL(error.what() << endl);
67 }
68
69
70 const double epsilon = 1.0e-6;
72
73
74
75
77
78 if (option == "KM3NeT") {
79
111
112 } else if (option == "Antares") {
113
117
118 } else {
119
120 FATAL(
"Fatal error at detector option " << option << endl);
121 };
122
123
124
125
127
130 }
131
132
134
135
136 TH1D h0m("L0m", NULL, 300, 1.0, 151.0);
137 TH1D h1m("L1m", NULL, 300, 1.0, 151.0);
138
139 TH1D h0s("L0s", NULL, 300, 1.0, 151.0);
140 TH1D h1s("L1s", NULL, 300, 1.0, 151.0);
141
142
143 {
149
150
151
152
154 DIRECT_LIGHT_FROM_MUON,
155 SCATTERED_LIGHT_FROM_MUON,
156 DIRECT_LIGHT_FROM_EMSHOWERS,
157 SCATTERED_LIGHT_FROM_EMSHOWERS,
158 DIRECT_LIGHT_FROM_DELTARAYS,
159 SCATTERED_LIGHT_FROM_DELTARAYS
160 };
161
162 const double TTS = 2.0;
164 const double epsilon = 1.0e-10;
165
166 const int NUMBER_OF_PDFS = sizeof(type)/sizeof(type[0]);
167
168 JPDF_t pdf[NUMBER_OF_PDFS];
169
170 const JPDF_t::JSupervisor supervisor(
new JPDF_t::JDefaultResult(
JMATH::zero));
171
172 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
173
174 try {
175
176 const string file_name = getFilename(fileDescriptor, type[i]);
177
178 NOTICE(
"loading PDF from file " << file_name <<
"... " << flush);
179
180 pdf[i].load(file_name.c_str());
181
183
184 pdf[i].setExceptionHandler(supervisor);
186 }
189 }
190 }
191
192 const double Tmin = -15.0;
193 const double Tmax = +250.0;
194 const double dt = 1.0;
195
196 for (int ix = 1; ix <= h0m.GetNbinsX(); ++ix) {
197
198 const double R = h0m.GetBinCenter(ix);
199
200 double V = 0.0;
201
202 for (vector<JAngle3D>::const_iterator pmt =
PMT.begin(); pmt !=
PMT.end(); ++pmt) {
203
204 const double theta = pi.constrain(pmt->getTheta());
205 const double phi = pi.constrain(fabs(pmt->getPhi()));
206
207 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
208
209 double yi = pdf[i](R, theta, phi, Tmax).V;
210
212 yi *= E;
214 yi *= JDeltaRays::getEnergyLossFromMuon(E);
215 }
216
217 V += yi;
218 }
219
220 h0m.SetBinContent(ix, 1.0 - TMath::PoissonI(0,V));
221 }
222 }
223
224 const int NUMBER_OF_PMTS =
PMT.size();
225
226 double Vi[NUMBER_OF_PMTS];
227
228 for (int ix = 1; ix <= h1m.GetNbinsX(); ++ix) {
229
230 const double R = h1m.GetBinCenter(ix);
231
232 double Y = 0.0;
233
234 for (
double x = Tmin;
x <= Tmax;
x += dt) {
235
236 double V = 0.0;
237
238 for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
239
240 Vi[pmt] = 0.0;
241
242 const double theta = pi.constrain(
PMT[pmt].getTheta());
243 const double phi = pi.constrain(fabs(
PMT[pmt].getPhi()));
244
245 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
246
247 double yi[] = {
248 pdf[i](R, theta, phi,
x).v,
249 pdf[i](R, theta, phi, x + TMax_ns).v
250 };
251
253 yi[0] *= E;
254 yi[1] *= E;
256 yi[0] *= JDeltaRays::getEnergyLossFromMuon(E);
257 yi[1] *= JDeltaRays::getEnergyLossFromMuon(E);
258 }
259
260 Vi[pmt] += yi[1] - yi[0];
261 }
262
263 V += Vi[pmt];
264 }
265
266 for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
267
268 const double theta = pi.constrain(
PMT[pmt].getTheta());
269 const double phi = pi.constrain(fabs(
PMT[pmt].getPhi()));
270
272
273 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
274
275 double yi = pdf[i](R, theta, phi,
x).f;
276
278 yi *= E;
280 yi *= JDeltaRays::getEnergyLossFromMuon(E);
281 }
282
284 }
285
286 if (y > 0.0) {
287 Y +=
y * (1.0 - TMath::PoissonI(0, V - Vi[pmt])) * dt;
288 }
289 }
290 }
291
292 h1m.SetBinContent(ix, 1.0 - TMath::PoissonI(0,Y));
293 }
294 }
295
296
297 {
304
305
306
307
309 DIRECT_LIGHT_FROM_EMSHOWER,
310 SCATTERED_LIGHT_FROM_EMSHOWER
311 };
312
313 const double TTS = 2.0;
315 const double epsilon = 1.0e-10;
316
317 const int NUMBER_OF_PDFS = sizeof(type)/sizeof(type[0]);
318
319 JPDF_t pdf[NUMBER_OF_PDFS];
320
321 const JPDF_t::JSupervisor supervisor(
new JPDF_t::JDefaultResult(
JMATH::zero));
322
323 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
324
325 try {
326
327 const string file_name = getFilename(fileDescriptor, type[i]);
328
329 NOTICE(
"loading PDF from file " << file_name <<
"... " << flush);
330
331 pdf[i].load(file_name.c_str());
332
334
335 pdf[i].setExceptionHandler(supervisor);
337 }
340 }
341 }
342
343 const double Tmin = -15.0;
344 const double Tmax = +250.0;
345 const double dt = 1.0;
346
347 for (int ix = 1; ix <= h0s.GetNbinsX(); ++ix) {
348
349 const double D = h0s.GetBinCenter(ix);
350
351 double V = 0.0;
352
353 for (vector<JAngle3D>::const_iterator pmt =
PMT.begin(); pmt !=
PMT.end(); ++pmt) {
354
355 const double theta = pi.constrain(pmt->getTheta());
356 const double phi = pi.constrain(fabs(pmt->getPhi()));
357
358 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
359
360 double yi = pdf[i](D, cd, theta, phi, Tmax).V;
361
362 yi *= E;
363
364 V += yi;
365 }
366 }
367
368 h0s.SetBinContent(ix, 1.0 - TMath::PoissonI(0,V));
369 }
370
371
372 const int NUMBER_OF_PMTS =
PMT.size();
373
374 double Vi[NUMBER_OF_PMTS];
375
376 for (int ix = 1; ix <= h1s.GetNbinsX(); ++ix) {
377
378 const double D = h1s.GetBinCenter(ix);
379
380 double Y = 0.0;
381
382 for (
double x = Tmin;
x <= Tmax;
x += dt) {
383
384 double V = 0.0;
385
386 for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
387
388 Vi[pmt] = 0.0;
389
390 const double theta = pi.constrain(
PMT[pmt].getTheta());
391 const double phi = pi.constrain(fabs(
PMT[pmt].getPhi()));
392
393 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
394
395 double yi[] = {
396 pdf[i](D, cd, theta, phi,
x).v,
397 pdf[i](D, cd, theta, phi, x + TMax_ns).v
398 };
399
400 yi[0] *= E;
401 yi[1] *= E;
402
403 Vi[pmt] += yi[1] - yi[0];
404 }
405
406 V += Vi[pmt];
407 }
408
409 for (int pmt = 0; pmt != NUMBER_OF_PMTS; ++pmt) {
410
411 const double theta = pi.constrain(
PMT[pmt].getTheta());
412 const double phi = pi.constrain(fabs(
PMT[pmt].getPhi()));
413
415
416 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
417
418 double yi = pdf[i](D, cd, theta, phi,
x).f;
419
420 yi *= E;
421
423 }
424
425 if (y > 0.0) {
426 Y +=
y * (1.0 - TMath::PoissonI(0, V - Vi[pmt])) * dt;
427 }
428 }
429 }
430
431 h1s.SetBinContent(ix, 1.0 - TMath::PoissonI(0,Y));
432 }
433 }
434
435
436 out.Write();
437 out.Close();
438}
JDAQPMTIdentifier PMT
Command line options.
#define make_field(A,...)
macro to convert parameter to JParserTemplateElement object
Data structure for angles in three dimensions.
Data structure for direction in three dimensions.
JDirection3D & rotate(const JRotation3D &R)
Rotate.
virtual const char * what() const override
Get error message.
Utility class to parse command line options.
Multi-dimensional PDF table for arrival time of Cherenkov light.
static const JZero zero
Function object to assign zero value.
bool is_deltarays(const int pdf)
Test if given PDF type corresponds to Cherenkov light from delta-rays.
bool is_bremsstrahlung(const int pdf)
Test if given PDF type corresponds to Cherenkov light from Bremsstrahlung.
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Empty structure for specification of parser element that is not initialised (i.e. does require input)...