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