Jpp 19.3.0-rc.5
the software that should make you happy
Loading...
Searching...
No Matches
JLine3ZRegressor.hh
Go to the documentation of this file.
1#ifndef __JFIT__JLINE3ZREGRESSOR__
2#define __JFIT__JLINE3ZREGRESSOR__
3
4#include <string>
5#include <array>
6#include <memory>
7#include <limits>
8
10#include "JPhysics/JPDFTable.hh"
11#include "JPhysics/JNPETable.hh"
13#include "JPhysics/JGeane.hh"
17#include "JTools/JRange.hh"
20#include "JMath/JZero.hh"
21#include "JFit/JTimeRange.hh"
22#include "JFit/JLine3Z.hh"
23#include "JFit/JPMTW0.hh"
24#include "JFit/JSimplex.hh"
25#include "JFit/JGandalf.hh"
26#include "JFit/JMEstimator.hh"
27#include "JFit/JRegressor.hh"
28#include "Jeep/JMessage.hh"
29
30
31/**
32 * \file
33 * Data regression method for JFIT::JLine3Z.
34 * \author mdejong
35 */
36namespace JFIT {}
37namespace JPP { using namespace JFIT; }
38
39namespace JFIT {
40
48
49
50 /**
51 * Regressor function object for JLine3Z fit using JSimplex minimiser.
52 */
53 template<>
54 struct JRegressor<JLine3Z, JSimplex> :
55 public JAbstractRegressor<JLine3Z, JSimplex>
56 {
57 using JAbstractRegressor<JLine3Z, JSimplex>::operator();
58
59 /**
60 * Constructor.
61 *
62 * \param sigma time resolution of hit [ns]
63 */
64 JRegressor(double sigma) :
65 estimator(new JMEstimatorNormal())
66 {
67 this->sigma = sigma;
68 }
69
70
71 /**
72 * Fit function.
73 * This method is used to determine the chi2 of given hit with respect to trajectory of muon.
74 *
75 * The template argument <tt>JHit_t</tt> refers to a data structure which should have the following member methods:
76 * - double getX(); // [m]
77 * - double getY(); // [m]
78 * - double getZ(); // [m]
79 * - double getT(); // [ns]
80 *
81 * \param track track
82 * \param hit hit
83 * \return chi2
84 */
85 template<class JHit_t>
86 double operator()(const JLine3Z& track, const JHit_t& hit) const
87 {
88 using namespace JPP;
89
90 JPosition3D D(hit.getX(), hit.getY(), hit.getZ());
91
92 D.sub(track.getPosition());
93
94 const double z = D.getDot(track.getDirection());
95 const double R = sqrt(D.getLengthSquared() - z*z);
96
97 const double t1 = track.getT() + (z + R * getKappaC()) * getInverseSpeedOfLight();
98
99 const double u = (t1 - hit.getT()) / sigma;
100
101 return estimator->getRho(u) * hit.getW();
102 }
103
104 std::shared_ptr<JMEstimator> estimator; //!< M-Estimator function
105 double sigma; //!< Time resolution [ns]
106 };
107
108
109 /**
110 * Template specialisation for storage of PDF tables.
111 */
112 template<>
114 {
120
125
127
128 static const int NUMBER_OF_PDFS = 6; //!< Number of PDFs
129
132
133
134 /**
135 * Default constructor.
136 */
139
140
141 /**
142 * Constructor.
143 *
144 * The PDF file descriptor should contain the wild card character JPHYSICS::WILDCARD which
145 * will be replaced by the PDF types listed in JRegressorStorage<JLine3Z, JGandalf>::pdf_t.
146 *
147 * The <tt>TTS_ns</tt> corresponds to the additional time smearing applied to the PDFs.
148 *
149 * \param fileDescriptor PDF file descriptor
150 * \param T_ns time range [ns]
151 * \param TTS_ns TTS [ns]
152 * \param numberOfPoints number of points for Gauss-Hermite integration of TTS
153 * \param epsilon precision for Gauss-Hermite integration of TTS
154 */
155 JRegressorStorage(const std::string& fileDescriptor,
156 const JTimeRange& T_ns,
157 const double TTS_ns,
158 const int numberOfPoints = 25,
159 const double epsilon = 1.0e-10) :
160 T_ns(T_ns)
161 {
162 using namespace std;
163 using namespace JPP;
164
165 const JPDF_t::JSupervisor supervisor(new JPDF_t::JDefaultResult(JMATH::zero));
166
167 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
168
169 const string file_name = getFilename(fileDescriptor, pdf_t[i]);
170
171 _pdf[i].load(file_name.c_str());
172
173 _pdf[i].setExceptionHandler(supervisor);
174 }
175
176 // Add PDFs
177
178 for (int i = 1; i < NUMBER_OF_PDFS; i += 2) {
179
180 _pdf[ i ].add(_pdf[i-1]);
181 _pdf[ i ].compress(numeric_limits<double>::max(), T_ns);
182
183 JPDF_t buffer;
184
185 _pdf[i-1].swap(buffer);
186
187 if (TTS_ns > 0.0) {
188 _pdf[i].blur(TTS_ns, numberOfPoints, epsilon);
189 }
190
191 _npe[ i ] = JNPE_t(_pdf[i], T_ns);
192 }
193 }
194
195
196 /**
197 * Get PDFs.
198 *
199 * \return PDFs
200 */
201 const JPDFs_t& getPDF() const
202 {
203 return _pdf;
204 }
205
206
207 /**
208 * Get NPEs.
209 *
210 * \return NPEs
211 */
212 const JNPEs_t& getNPE() const
213 {
214 return _npe;
215 }
216
217
218 /**
219 * Transform PDFs and NPEs.
220 *
221 * \param transformer transformer
222 */
223 void transform(const transformer_type& transformer)
224 {
225 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
226 _pdf[i].transform(transformer);
227 _npe[i].transform(transformer);
228 }
229 }
230
231
232 /**
233 * PDF types.
234 */
235 static const JPDFType_t pdf_t[NUMBER_OF_PDFS];
236
237 JTimeRange T_ns; //!< Time window with respect to Cherenkov hypothesis [ns]
238
239 private:
240 JPDFs_t _pdf; //!< PDFs
241 JNPEs_t _npe; //!< NPEs
242 };
243
244
245 /**
246 * PDF types.
247 */
248 const JPDFType_t JRegressorStorage<JLine3Z, JGandalf>::pdf_t[] = { DIRECT_LIGHT_FROM_MUON,
249 SCATTERED_LIGHT_FROM_MUON,
250 DIRECT_LIGHT_FROM_DELTARAYS,
251 SCATTERED_LIGHT_FROM_DELTARAYS,
252 DIRECT_LIGHT_FROM_EMSHOWERS,
253 SCATTERED_LIGHT_FROM_EMSHOWERS };
254
255
256 /**
257 * Regressor function object for JLine3Z fit using JGandalf minimiser.
258 */
259 template<>
260 struct JRegressor<JLine3Z, JGandalf> :
261 public JAbstractRegressor<JLine3Z, JGandalf>,
262 public JRegressorStorage <JLine3Z, JGandalf>
263 {
264 using JAbstractRegressor<JLine3Z, JGandalf>::operator();
265
266 typedef JRegressorStorage<JLine3Z, JGandalf> storage_type;
267
268 /**
269 * Default constructor
270 */
271 JRegressor() :
272 storage_type(),
273 pdf(getPDF()),
274 npe(getNPE()),
275 estimator()
276 {}
277
278
279 /**
280 * Constructor.
281 *
282 * The PDF file descriptor should contain the wild card character JPHYSICS::WILDCARD which
283 * will be replaced by the PDF types listed in JRegressorStorage<JLine3Z, JGandalf>::pdf_t.
284 *
285 * The <tt>TTS_ns</tt> corresponds to the additional Gaussian time smearing applied to the PDFs.
286 *
287 * \param fileDescriptor PDF file descriptor
288 * \param T_ns time range [ns]
289 * \param TTS_ns TTS [ns]
290 * \param numberOfPoints number of points for Gauss-Hermite integration of TTS
291 * \param epsilon precision for Gauss-Hermite integration of TTS
292 */
293 JRegressor(const std::string& fileDescriptor,
294 const JTimeRange& T_ns,
295 const double TTS_ns,
296 const int numberOfPoints = 25,
297 const double epsilon = 1.0e-10) :
298 storage_type(fileDescriptor, T_ns, TTS_ns, numberOfPoints, epsilon),
299 pdf(getPDF()),
300 npe(getNPE()),
301 estimator(new JMEstimatorNormal())
302 {}
303
304
305 /**
306 * Constructor.
307 *
308 * \param storage PDF storage
309 */
310 JRegressor(const storage_type& storage) :
311 pdf(storage.getPDF()),
312 npe(storage.getNPE()),
313 estimator(new JMEstimatorNormal())
314 {
315 T_ns = storage.T_ns;
316 }
317
318
319 /**
320 * Fit function.
321 * This method is used to determine the chi2 and gradient of given hit with respect to trajectory of muon.
322 *
323 * The template argument <tt>JHit_t</tt> refers to a data structure which should have the following member methods:
324 * - double getX(); // [m]
325 * - double getY(); // [m]
326 * - double getZ(); // [m]
327 * - double getDX(); // [u]
328 * - double getDY(); // [u]
329 * - double getDZ(); // [u]
330 * - double getT(); // [ns]
331 *
332 * \param track track
333 * \param hit hit
334 * \return chi2 and gradient
335 */
336 template<class JHit_t>
337 result_type operator()(const JLine3Z& track, const JHit_t& hit) const
338 {
339 using namespace JPP;
340
341 JPosition3D D(hit.getX(), hit.getY(), hit.getZ());
342 JDirection3D U(hit.getDX(), hit.getDY(), hit.getDZ());
343
344 D.sub(track.getPosition());
345
346 const double z = D.getDot(track.getDirection());
347 const double x = D.getX() - z * track.getDX();
348 const double y = D.getY() - z * track.getDY();
349 const double R2 = D.getLengthSquared() - z*z;
350 const double R = sqrt(R2);
351
352 const double t1 = track.getT() + (z + R * getTanThetaC()) * getInverseSpeedOfLight();
353
354 U.rotate(JRotation3Z(-atan2(y,x))); // rotate PMT axis to x-z plane
355
356 const double theta = U.getTheta();
357 const double phi = fabs(U.getPhi());
358
359 const double E = gWater.getE(E_GeV, z);
360 const double dt = T_ns.constrain(hit.getT() - t1);
361
362 JPDF_t::result_type H0 = getH0(hit.getR(), dt);
363 JPDF_t::result_type H1 = getH1(E, R, theta, phi, dt);
364
365 if (H1.V >= Vmax_npe) {
366 H1 *= Vmax_npe / H1.V;
367 }
368
369 H1 += H0; // signal + background
370
371 result_type result;
372
373 result.chi2 = H1.getChi2() - H0.getChi2(); // Likelihood ratio
374
375 const double wc = 1.0 - getTanThetaC() * z / R; // d(ct1)/d(z)
376
377 result.gradient = JLine3Z(JLine1Z(JPosition3D(-getTanThetaC() * D.getX() / R, // d(ct1)/d(x)
378 -getTanThetaC() * D.getY() / R, // d(ct1)/d(y)
379 0.0),
380 getSpeedOfLight()), // d(ct1)/d(t)
381 JVersor3Z(wc * (D.getX() - D.getZ()*track.getDX()/track.getDZ()), // d(ct1)/d(dx)
382 wc * (D.getY() - D.getZ()*track.getDY()/track.getDZ()))); // d(ct1)/d(dy)
383
384 result.gradient.mul(getInverseSpeedOfLight() * (H1.getDerivativeOfChi2() -
385 H0.getDerivativeOfChi2())); // x d(chi2)/d(ct1)
386
387 return result;
388 }
389
390
391 /**
392 * Fit function.
393 * This method is used to determine the chi2 and gradient of given PMT with respect to trajectory of muon.
394 *
395 * \param track track
396 * \param pmt pmt
397 * \return chi2 and gradient
398 */
399 result_type operator()(const JLine3Z& track, const JPMTW0& pmt) const
400 {
401 using namespace JGEOMETRY3D;
402 using namespace JPHYSICS;
403
404 JPosition3D D(pmt.getPosition());
405 JDirection3D U(pmt.getDirection());
406
407 D.sub(track.getPosition());
408
409 const double z = D.getDot(track.getDirection());
410 const double x = D.getX() - z * track.getDX();
411 const double y = D.getY() - z * track.getDY();
412 const double R2 = D.getLengthSquared() - z*z;
413 const double R = sqrt(R2);
414
415 U.rotate(JRotation3Z(-atan2(y,x))); // rotate PMT axis to x-z plane
416
417 const double theta = U.getTheta();
418 const double phi = fabs(U.getPhi());
419
420 const double E = gWater.getE(E_GeV, z);
421
422 JNPE_t::result_type H0 = getH0(pmt.getR());
423 JNPE_t::result_type H1 = getH1(E, R, theta, phi);
424
425 if (H1.f >= Vmax_npe) {
426 H1 *= Vmax_npe / H1.f;
427 }
428
429 H1 += H0;
430
431 const bool hit = pmt.getN() != 0;
432 const double u = H1.getChi2(hit);
433
434 result_type result;
435
436 result.chi2 = estimator->getRho(u);
437
438 result.gradient = JLine3Z(JLine1Z(JPosition3D(-D.getX() / R, // d(R)/d(x)
439 -D.getY() / R, // d(R)/d(y)
440 0.0),
441 0.0),
442 JVersor3Z(-z * D.getX() / R, // d(R)/d(dx)
443 -z * D.getY() / R)); // d(R)/d(dy)
444
445 result.gradient.mul(estimator->getPsi(u));
446 result.gradient.mul(H1.getDerivativeOfChi2(hit)); // x d(chi2)/d(R)
447
448 return result;
449 }
450
451
452 /**
453 * Get background hypothesis value for time differentiated PDF.
454 *
455 * \param R_Hz rate [Hz]
456 * \param t1 time [ns]
457 * \return hypothesis value
458 */
459 JPDF_t::result_type getH0(const double R_Hz,
460 const double t1) const
461 {
462 return JPDF_t::result_type(R_Hz * 1e-9, t1, T_ns);
463 }
464
465
466 /**
467 * Get signal hypothesis value for time differentiated PDF.
468 *
469 * \param E muon energy at minimum distance of approach [GeV]
470 * \param R minimum distance of approach [m]
471 * \param theta PMT zenith angle [rad]
472 * \param phi PMT azimuth angle [rad]
473 * \param t1 arrival time relative to Cherenkov hypothesis [ns]
474 * \return hypothesis value
475 */
476 JPDF_t::result_type getH1(const double E,
477 const double R,
478 const double theta,
479 const double phi,
480 const double t1) const
481 {
482 using namespace std;
483 using namespace JPP;
484
485 JPDF_t::result_type h1 = zero;
486
487 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
488
489 if (!pdf[i].empty() && R <= pdf[i].getXmax()) {
490
491 JPDF_t::result_type y1 = pdf[i](max(R, pdf[i].getXmin()), theta, phi, t1);
492
493 // safety measures
494
495 if (y1.f <= 0.0) {
496 y1.f = 0.0;
497 y1.fp = 0.0;
498 }
499
500 if (y1.v <= 0.0) {
501 y1.v = 0.0;
502 }
503
504 // energy dependence
505
506 if (is_deltarays(pdf_t[i])) {
507 y1 *= getDeltaRaysFromMuon(E);
508 } else if (is_bremsstrahlung(pdf_t[i])) {
509 y1 *= E;
510 }
511
512 h1 += y1;
513 }
514 }
515
516 return h1;
517 }
518
519
520 /**
521 * Get background hypothesis value for time integrated PDF.
522 *
523 * \param R_Hz rate [Hz]
524 * \return hypothesis value
525 */
526 JNPE_t::result_type getH0(const double R_Hz) const
527 {
528 return JNPE_t::result_type(R_Hz * 1e-9 * T_ns.getLength(), 0.0);
529 }
530
531
532 /**
533 * Get signal hypothesis value for time integrated PDF.
534 *
535 * \param E muon energy at minimum distance of approach [GeV]
536 * \param R minimum distance of approach [m]
537 * \param theta PMT zenith angle [rad]
538 * \param phi PMT azimuth angle [rad]
539 * \return hypothesis value
540 */
541 JNPE_t::result_type getH1(const double E,
542 const double R,
543 const double theta,
544 const double phi) const
545 {
546 using namespace std;
547 using namespace JPP;
548
549 JNPE_t::result_type h1 = JMATH::zero;
550
551 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
552
553 if (!npe[i].empty() && R <= npe[i].getXmax()) {
554
555 try {
556
557 JNPE_t::result_type y1 = npe[i](max(R, npe[i].getXmin()), theta, phi);
558
559 // safety measures
560
561 if (y1.f > 0.0) {
562
563 // energy dependence
564
565 if (is_deltarays(pdf_t[i])) {
566 y1 *= getDeltaRaysFromMuon(E);
567 } else if (is_bremsstrahlung(pdf_t[i])) {
568 y1 *= E;
569 }
570
571 h1 += y1;
572 }
573 }
574 catch(JException& error) {
575 ERROR(error << endl);
576 }
577 }
578 }
579
580 return h1;
581 }
582
583
584 /**
585 * Get maximal road width of PDF.
586 *
587 * \return road width [m]
588 */
589 inline double getRmax() const
590 {
591 double xmax = 0.0;
592
593 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
594 if (!pdf[i].empty() && pdf[i].getXmax() > xmax) {
595 xmax = pdf[i].getXmax();
596 }
597 }
598
599 return xmax;
600 }
601
602
603 static double Vmax_npe; //!< Maximal integral of PDF [npe]
604
605 const JPDFs_t& pdf; //!< PDF
606 const JNPEs_t& npe; //!< PDF
607
608 double E_GeV = 0.0; //!< Energy of muon at vertex [GeV]
609
610 std::shared_ptr<JMEstimator> estimator; //!< M-Estimator function
611 };
612
613
614 /**
615 * Default values.
616 */
617 double JRegressor<JLine3Z, JGandalf>::Vmax_npe = std::numeric_limits<double>::max();
618}
619
620#endif
Various implementations of functional maps.
Energy loss of muon.
Maximum likelihood estimator (M-estimators).
General purpose messaging.
#define ERROR(A)
Definition JMessage.hh:66
Auxiliary methods for PDF calculations.
Numbering scheme for PDF types.
Physics constants.
Auxiliary class to define a range between two values.
General purpose data regression method.
int numberOfPoints
Definition JResultPDF.cc:22
Definition of zero value for any class.
Fit method based on the Levenberg-Marquardt method.
Definition JGandalf.hh:87
Data structure for fit of straight line paralel to z-axis.
Definition JLine1Z.hh:29
Data structure for fit of straight line in positive z-direction.
Definition JLine3Z.hh:40
JVersor3D getDirection(const JVector3D &pos) const
Get photon direction of Cherenkov light on PMT.
Definition JLine3Z.hh:255
double getT(const JVector3D &pos) const
Get arrival time of Cherenkov light at given position.
Definition JLine3Z.hh:234
Simple fit method based on Powell's algorithm, see reference: Numerical Recipes in C++,...
Definition JSimplex.hh:44
Data structure for direction in three dimensions.
const JDirection3D & getDirection() const
Get direction.
Data structure for position in three dimensions.
const JPosition3D & getPosition() const
Get position.
Rotation around Z-axis.
Data structure for normalised vector in positive z-direction.
Definition JVersor3Z.hh:41
double getDZ() const
Get z direction.
Definition JVersor3Z.hh:169
double getDY() const
Get y direction.
Definition JVersor3Z.hh:158
double getDX() const
Get x direction.
Definition JVersor3Z.hh:147
General exception.
Definition JException.hh:24
Custom class for integrated values of the PDF of the arrival time of Cherenkov light.
Definition JNPETable.hh:46
Multi-dimensional PDF table for arrival time of Cherenkov light.
Definition JPDFTable.hh:44
Interface for weight application and coordinate transformation of function.
const double sigma[]
const double xmax
const double epsilon
double getNPE(const Hit &hit)
Get true charge of hit.
Auxiliary classes and methods for linear and iterative data regression.
Definition JEnergy.hh:15
Auxiliary classes and methods for 3D geometrical objects and operations.
Definition JAngle3D.hh:19
static const JZero zero
Function object to assign zero value.
Definition JZero.hh:105
Auxiliary methods for light properties of deep-sea water.
double getDeltaRaysFromMuon(const double E, const JRange< double > T_GeV=JRange< double >(DELTARAY_TMIN, DELTARAY_TMAX))
Equivalent EM-shower energy due to delta-rays per unit muon track length.
double getKappaC()
Get average R-dependence of arrival time of Cherenkov light (a.k.a.
bool is_deltarays(const int pdf)
Test if given PDF type corresponds to Cherenkov light from delta-rays.
Definition JPDFTypes.hh:151
bool is_bremsstrahlung(const int pdf)
Test if given PDF type corresponds to Cherenkov light from Bremsstrahlung.
Definition JPDFTypes.hh:137
JPDFType_t
PDF types.
Definition JPDFTypes.hh:24
@ SCATTERED_LIGHT_FROM_DELTARAYS
scattered light from delta-rays
Definition JPDFTypes.hh:33
@ DIRECT_LIGHT_FROM_EMSHOWERS
direct light from EM showers
Definition JPDFTypes.hh:29
@ SCATTERED_LIGHT_FROM_EMSHOWERS
scattered light from EM showers
Definition JPDFTypes.hh:30
@ SCATTERED_LIGHT_FROM_MUON
scattered light from muon
Definition JPDFTypes.hh:27
@ DIRECT_LIGHT_FROM_DELTARAYS
direct light from delta-rays
Definition JPDFTypes.hh:32
@ DIRECT_LIGHT_FROM_MUON
direct light from muon
Definition JPDFTypes.hh:26
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Abstract class for global fit method.
Definition JRegressor.hh:79
Normal M-estimator.
Auxiliary class for handling PMT geometry, rate and response.
Definition JPMTW0.hh:24
int getN() const
Get number of hits.
Definition JPMTW0.hh:67
double getR() const
Get rate.
Definition JPMTW0.hh:56
JPHYSICS::JNPETable< double, double, JNPEMaplist_t > JNPE_t
time integrated PDF
const JNPEs_t & getNPE() const
Get NPEs.
std::array< JPDF_t, NUMBER_OF_PDFS > JPDFs_t
PDFs.
JPHYSICS::JPDFTable< JFunction1D_t, JPDFMaplist_t > JPDF_t
time dependent PDF
JTimeRange T_ns
Time window with respect to Cherenkov hypothesis [ns].
JTOOLS::JMAPLIST< JTOOLS::JPolint1FunctionalMap, JTOOLS::JPolint0FunctionalGridMap, JTOOLS::JPolint0FunctionalGridMap >::maplist JPDFMaplist_t
std::array< JNPE_t, NUMBER_OF_PDFS > JNPEs_t
NPEs.
JTOOLS::JMAPLIST< JTOOLS::JPolint1FunctionalMapH, JTOOLS::JPolint1FunctionalGridMap, JTOOLS::JPolint1FunctionalGridMap >::maplist JNPEMaplist_t
const JPDFs_t & getPDF() const
Get PDFs.
JRegressorStorage(const std::string &fileDescriptor, const JTimeRange &T_ns, const double TTS_ns, const int numberOfPoints=25, const double epsilon=1.0e-10)
Constructor.
void transform(const transformer_type &transformer)
Transform PDFs and NPEs.
Template data structure for storage of internal data.
Template definition of a data regressor of given model.
Definition JRegressor.hh:70
Auxiliary class to set-up Hit.
Definition JSirene.hh:58
Auxiliary class for recursive map list generation.
Definition JMapList.hh:109
Type definition of a zero degree polynomial interpolation based on a JGridMap implementation.
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 1st degree polynomial interpolation based on a JMap implementation.
Type definition of a spline interpolation method based on a JCollection with JResultPDF result type.