Jpp test-rotations-old-533-g2bdbdb559
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/JPDFTypes.hh"
11#include "JPhysics/JPDFTable.hh"
12#include "JPhysics/JNPETable.hh"
14#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 */
249 DIRECT_LIGHT_FROM_MUON,
250 SCATTERED_LIGHT_FROM_MUON,
251 DIRECT_LIGHT_FROM_DELTARAYS,
252 SCATTERED_LIGHT_FROM_DELTARAYS,
253 DIRECT_LIGHT_FROM_EMSHOWERS,
254 SCATTERED_LIGHT_FROM_EMSHOWERS
255 };
256
257
258 /**
259 * Regressor function object for JLine3Z fit using JGandalf minimiser.
260 */
261 template<>
262 struct JRegressor<JLine3Z, JGandalf> :
263 public JAbstractRegressor<JLine3Z, JGandalf>,
264 public JRegressorStorage <JLine3Z, JGandalf>
265 {
266 using JAbstractRegressor<JLine3Z, JGandalf>::operator();
267
268 typedef JRegressorStorage<JLine3Z, JGandalf> storage_type;
269
270 /**
271 * Default constructor
272 */
273 JRegressor() :
274 storage_type(),
275 pdf(getPDF()),
276 npe(getNPE()),
277 estimator()
278 {}
279
280
281 /**
282 * Constructor.
283 *
284 * The PDF file descriptor should contain the wild card character JPHYSICS::WILDCARD which
285 * will be replaced by the PDF types listed in JRegressorStorage<JLine3Z, JGandalf>::pdf_t.
286 *
287 * The <tt>TTS_ns</tt> corresponds to the additional Gaussian time smearing applied to the PDFs.
288 *
289 * \param fileDescriptor PDF file descriptor
290 * \param T_ns time range [ns]
291 * \param TTS_ns TTS [ns]
292 * \param numberOfPoints number of points for Gauss-Hermite integration of TTS
293 * \param epsilon precision for Gauss-Hermite integration of TTS
294 */
295 JRegressor(const std::string& fileDescriptor,
296 const JTimeRange& T_ns,
297 const double TTS_ns,
298 const int numberOfPoints = 25,
299 const double epsilon = 1.0e-10) :
300 storage_type(fileDescriptor, T_ns, TTS_ns, numberOfPoints, epsilon),
301 pdf(getPDF()),
302 npe(getNPE()),
303 estimator(new JMEstimatorNormal())
304 {}
305
306
307 /**
308 * Constructor.
309 *
310 * \param storage PDF storage
311 */
312 JRegressor(const storage_type& storage) :
313 pdf(storage.getPDF()),
314 npe(storage.getNPE()),
315 estimator(new JMEstimatorNormal())
316 {
317 T_ns = storage.T_ns;
318 }
319
320
321 /**
322 * Fit function.
323 * This method is used to determine the chi2 and gradient of given hit with respect to trajectory of muon.
324 *
325 * The template argument <tt>JHit_t</tt> refers to a data structure which should have the following member methods:
326 * - double getX(); // [m]
327 * - double getY(); // [m]
328 * - double getZ(); // [m]
329 * - double getDX(); // [u]
330 * - double getDY(); // [u]
331 * - double getDZ(); // [u]
332 * - double getT(); // [ns]
333 *
334 * \param track track
335 * \param hit hit
336 * \return chi2 and gradient
337 */
338 template<class JHit_t>
339 result_type operator()(const JLine3Z& track, const JHit_t& hit) const
340 {
341 using namespace JPP;
342
343 JPosition3D D(hit.getX(), hit.getY(), hit.getZ());
344 JDirection3D U(hit.getDX(), hit.getDY(), hit.getDZ());
345
346 D.sub(track.getPosition());
347
348 const double z = D.getDot(track.getDirection());
349 const double x = D.getX() - z * track.getDX();
350 const double y = D.getY() - z * track.getDY();
351 const double R2 = D.getLengthSquared() - z*z;
352 const double R = sqrt(R2);
353
354 const double t1 = track.getT() + (z + R * getTanThetaC()) * getInverseSpeedOfLight();
355
356 U.rotate(JRotation3Z(-atan2(y,x))); // rotate PMT axis to x-z plane
357
358 const double theta = U.getTheta();
359 const double phi = fabs(U.getPhi());
360
361 const double E = gWater.getE(E_GeV, z);
362 const double dt = T_ns.constrain(hit.getT() - t1);
363
364 JPDF_t::result_type H0 = getH0(hit.getR(), dt);
365 JPDF_t::result_type H1 = getH1(E, R, theta, phi, dt);
366
367 if (H1.V >= Vmax_npe) {
368 H1 *= Vmax_npe / H1.V;
369 }
370
371 H1 += H0; // signal + background
372
373 result_type result;
374
375 result.chi2 = H1.getChi2() - H0.getChi2(); // Likelihood ratio
376
377 const double wc = 1.0 - getTanThetaC() * z / R; // d(ct1)/d(z)
378
379 result.gradient = JLine3Z(JLine1Z(JPosition3D(-getTanThetaC() * D.getX() / R, // d(ct1)/d(x)
380 -getTanThetaC() * D.getY() / R, // d(ct1)/d(y)
381 0.0),
382 getSpeedOfLight()), // d(ct1)/d(t)
383 JVersor3Z(wc * (D.getX() - D.getZ()*track.getDX()/track.getDZ()), // d(ct1)/d(dx)
384 wc * (D.getY() - D.getZ()*track.getDY()/track.getDZ()))); // d(ct1)/d(dy)
385
386 result.gradient.mul(getInverseSpeedOfLight() * (H1.getDerivativeOfChi2() -
387 H0.getDerivativeOfChi2())); // x d(chi2)/d(ct1)
388
389 return result;
390 }
391
392
393 /**
394 * Fit function.
395 * This method is used to determine the chi2 and gradient of given PMT with respect to trajectory of muon.
396 *
397 * \param track track
398 * \param pmt pmt
399 * \return chi2 and gradient
400 */
401 result_type operator()(const JLine3Z& track, const JPMTW0& pmt) const
402 {
403 using namespace JGEOMETRY3D;
404 using namespace JPHYSICS;
405
406 JPosition3D D(pmt.getPosition());
407 JDirection3D U(pmt.getDirection());
408
409 D.sub(track.getPosition());
410
411 const double z = D.getDot(track.getDirection());
412 const double x = D.getX() - z * track.getDX();
413 const double y = D.getY() - z * track.getDY();
414 const double R2 = D.getLengthSquared() - z*z;
415 const double R = sqrt(R2);
416
417 U.rotate(JRotation3Z(-atan2(y,x))); // rotate PMT axis to x-z plane
418
419 const double theta = U.getTheta();
420 const double phi = fabs(U.getPhi());
421
422 const double E = gWater.getE(E_GeV, z);
423
424 JNPE_t::result_type H0 = getH0(pmt.getR());
425 JNPE_t::result_type H1 = getH1(E, R, theta, phi);
426
427 if (H1.f >= Vmax_npe) {
428 H1 *= Vmax_npe / H1.f;
429 }
430
431 H1 += H0;
432
433 const bool hit = pmt.getN() != 0;
434 const double u = H1.getChi2(hit);
435
436 result_type result;
437
438 result.chi2 = estimator->getRho(u);
439
440 result.gradient = JLine3Z(JLine1Z(JPosition3D(-D.getX() / R, // d(R)/d(x)
441 -D.getY() / R, // d(R)/d(y)
442 0.0),
443 0.0),
444 JVersor3Z(-z * D.getX() / R, // d(R)/d(dx)
445 -z * D.getY() / R)); // d(R)/d(dy)
446
447 result.gradient.mul(estimator->getPsi(u));
448 result.gradient.mul(H1.getDerivativeOfChi2(hit)); // x d(chi2)/d(R)
449
450 return result;
451 }
452
453
454 /**
455 * Get background hypothesis value for time differentiated PDF.
456 *
457 * \param R_Hz rate [Hz]
458 * \param t1 time [ns]
459 * \return hypothesis value
460 */
461 JPDF_t::result_type getH0(const double R_Hz,
462 const double t1) const
463 {
464 return JPDF_t::result_type(R_Hz * 1e-9, t1, T_ns);
465 }
466
467
468 /**
469 * Get signal hypothesis value for time differentiated PDF.
470 *
471 * \param E muon energy at minimum distance of approach [GeV]
472 * \param R minimum distance of approach [m]
473 * \param theta PMT zenith angle [rad]
474 * \param phi PMT azimuth angle [rad]
475 * \param t1 arrival time relative to Cherenkov hypothesis [ns]
476 * \return hypothesis value
477 */
478 JPDF_t::result_type getH1(const double E,
479 const double R,
480 const double theta,
481 const double phi,
482 const double t1) const
483 {
484 using namespace std;
485 using namespace JPP;
486
487 JPDF_t::result_type h1 = zero;
488
489 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
490
491 if (!pdf[i].empty() && R <= pdf[i].getXmax()) {
492
493 JPDF_t::result_type y1 = pdf[i](max(R, pdf[i].getXmin()), theta, phi, t1);
494
495 // safety measures
496
497 if (y1.f <= 0.0) {
498 y1.f = 0.0;
499 y1.fp = 0.0;
500 }
501
502 if (y1.v <= 0.0) {
503 y1.v = 0.0;
504 }
505
506 // energy dependence
507
508 if (is_deltarays(pdf_t[i])) {
509 y1 *= JDeltaRays::getEnergyLossFromMuon(E);
510 } else if (is_bremsstrahlung(pdf_t[i])) {
511 y1 *= E;
512 }
513
514 h1 += y1;
515 }
516 }
517
518 return h1;
519 }
520
521
522 /**
523 * Get background hypothesis value for time integrated PDF.
524 *
525 * \param R_Hz rate [Hz]
526 * \return hypothesis value
527 */
528 JNPE_t::result_type getH0(const double R_Hz) const
529 {
530 return JNPE_t::result_type(R_Hz * 1e-9 * T_ns.getLength(), 0.0);
531 }
532
533
534 /**
535 * Get signal hypothesis value for time integrated PDF.
536 *
537 * \param E muon energy at minimum distance of approach [GeV]
538 * \param R minimum distance of approach [m]
539 * \param theta PMT zenith angle [rad]
540 * \param phi PMT azimuth angle [rad]
541 * \return hypothesis value
542 */
543 JNPE_t::result_type getH1(const double E,
544 const double R,
545 const double theta,
546 const double phi) const
547 {
548 using namespace std;
549 using namespace JPP;
550
551 JNPE_t::result_type h1 = JMATH::zero;
552
553 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
554
555 if (!npe[i].empty() && R <= npe[i].getXmax()) {
556
557 try {
558
559 JNPE_t::result_type y1 = npe[i](max(R, npe[i].getXmin()), theta, phi);
560
561 // safety measures
562
563 if (y1.f > 0.0) {
564
565 // energy dependence
566
567 if (is_deltarays(pdf_t[i])) {
568 y1 *= JDeltaRays::getEnergyLossFromMuon(E);
569 } else if (is_bremsstrahlung(pdf_t[i])) {
570 y1 *= E;
571 }
572
573 h1 += y1;
574 }
575 }
576 catch(JException& error) {
577 ERROR(error << endl);
578 }
579 }
580 }
581
582 return h1;
583 }
584
585
586 /**
587 * Get maximal road width of PDF.
588 *
589 * \return road width [m]
590 */
591 inline double getRmax() const
592 {
593 double xmax = 0.0;
594
595 for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
596 if (!pdf[i].empty() && pdf[i].getXmax() > xmax) {
597 xmax = pdf[i].getXmax();
598 }
599 }
600
601 return xmax;
602 }
603
604
605 static double Vmax_npe; //!< Maximal integral of PDF [npe]
606
607 const JPDFs_t& pdf; //!< PDF
608 const JNPEs_t& npe; //!< PDF
609
610 double E_GeV = 0.0; //!< Energy of muon at vertex [GeV]
611
612 std::shared_ptr<JMEstimator> estimator; //!< M-Estimator function
613 };
614
615
616 /**
617 * Default values.
618 */
619 double JRegressor<JLine3Z, JGandalf>::Vmax_npe = std::numeric_limits<double>::max();
620}
621
622#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
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.
double getNPE(const Hit &hit)
Get true charge of hit.
Auxiliary classes and methods for linear and iterative data regression.
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 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:60
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.