Jpp  16.0.3
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
JLine3ZRegressor.hh
Go to the documentation of this file.
1 #ifndef __JFIT__JLINE3ZREGRESSOR__
2 #define __JFIT__JLINE3ZREGRESSOR__
3 
4 #include <string>
5 
6 #include "JPhysics/JPDFTypes.hh"
7 #include "JPhysics/JPDFTable.hh"
8 #include "JPhysics/JNPETable.hh"
10 #include "JPhysics/JGeane.hh"
11 #include "JTools/JFunction1D_t.hh"
13 #include "JPhysics/JConstants.hh"
14 #include "JTools/JRange.hh"
17 #include "JMath/JZero.hh"
18 #include "JLang/JSharedPointer.hh"
19 #include "JFit/JTimeRange.hh"
20 #include "JFit/JLine3Z.hh"
21 #include "JFit/JPMTW0.hh"
22 #include "JFit/JSimplex.hh"
23 #include "JFit/JGandalf.hh"
24 #include "JFit/JMEstimator.hh"
25 #include "JFit/JRegressor.hh"
26 #include "Jeep/JMessage.hh"
27 
28 
29 /**
30  * \file
31  * Data regression method for JFIT::JLine3Z.
32  * \author mdejong
33  */
34 namespace JFIT {}
35 namespace JPP { using namespace JFIT; }
36 
37 namespace JFIT {
38 
46 
47 
48  /**
49  * Regressor function object for JLine3Z fit using JSimplex minimiser.
50  */
51  template<>
53  public JAbstractRegressor<JLine3Z, JSimplex>
54  {
56 
57  /**
58  * Constructor.
59  *
60  * \param sigma time resolution of hit [ns]
61  */
62  JRegressor(double sigma) :
63  estimator(new JMEstimatorNormal())
64  {
65  this->sigma = sigma;
66  }
67 
68 
69  /**
70  * Fit function.
71  * This method is used to determine the chi2 of given hit with respect to trajectory of muon.
72  *
73  * The template argument <tt>JHit_t</tt> refers to a data structure which should have the following member methods:
74  * - double getX(); // [m]
75  * - double getY(); // [m]
76  * - double getZ(); // [m]
77  * - double getT(); // [ns]
78  *
79  * \param track track
80  * \param hit hit
81  * \return chi2
82  */
83  template<class JHit_t>
84  double operator()(const JLine3Z& track, const JHit_t& hit) const
85  {
86  using namespace JPP;
87 
88  JPosition3D D(hit.getX(), hit.getY(), hit.getZ());
89 
90  D.sub(track.getPosition());
91 
92  const double z = D.getDot(track.getDirection());
93  const double R = sqrt(D.getLengthSquared() - z*z);
94 
95  const double t1 = track.getT() + (z + R * getKappaC()) * getInverseSpeedOfLight();
96 
97  const double u = (t1 - hit.getT()) / sigma;
98 
99  return estimator->getRho(u) * hit.getW();
100  }
101 
103  double sigma; //!< Time resolution [ns]
104  };
105 
106 
107  /**
108  * Regressor function object for JLine3Z fit using JGandalf minimiser.
109  */
110  template<>
112  public JAbstractRegressor<JLine3Z, JGandalf>
113  {
115 
121 
126 
127  /**
128  * Default constructor
129  */
131  {};
132 
133 
134  /**
135  * Constructor.
136  *
137  * The PDF file descriptor should contain the wild card character JPHYSICS::WILD_CARD which
138  * will be replaced by the corresponding PDF types listed in JRegressor<JLine3Z, JGandalf>::pdf_t.
139  *
140  * The <tt>TTS</tt> corresponds to the additional time smearing applied to the PDFs.
141  *
142  * \param fileDescriptor PDF file descriptor
143  * \param TTS TTS [ns]
144  * \param numberOfPoints number of points for Gauss-Hermite integration of TTS
145  * \param epsilon precision for Gauss-Hermite integration of TTS
146  */
147  JRegressor(const std::string& fileDescriptor,
148  const double TTS,
149  const int numberOfPoints = 25,
150  const double epsilon = 1.0e-10) :
151  E_GeV(0.0),
152  estimator(new JMEstimatorNormal())
153  {
154  using namespace std;
155  using namespace JPP;
156 
157  const JPDF_t::JSupervisor supervisor(new JPDF_t::JDefaultResult(JMATH::zero));
158 
159  for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
160 
161  try {
162 
163  const string file_name = getFilename(fileDescriptor, pdf_t[i]);
164 
165  NOTICE("loading PDF from file " << file_name << "... " << flush);
166 
167  pdf[i].load(file_name.c_str());
168 
169  NOTICE("OK" << endl);
170 
171  pdf[i].setExceptionHandler(supervisor);
172  }
173  catch(const JException& error) {
174  FATAL(error.what() << endl);
175  }
176  }
177 
178  // Add PDFs
179 
180  for (int i = 1; i < NUMBER_OF_PDFS; i += 2) {
181 
182  pdf[ i ].add(pdf[i-1]);
183 
184  JPDF_t buffer;
185 
186  pdf[i-1].swap(buffer);
187 
188  if (TTS > 0.0) {
189  pdf[i].blur(TTS, numberOfPoints, epsilon);
190  } else if (TTS < 0.0) {
191  ERROR("Illegal value of TTS [ns]: " << TTS << endl);
192  }
193 
194  npe[ i ] = JNPE_t(pdf[i]);
195  }
196  }
197 
198 
199  /**
200  * Fit function.
201  * This method is used to determine the chi2 and gradient of given hit with respect to trajectory of muon.
202  *
203  * The template argument <tt>JHit_t</tt> refers to a data structure which should have the following member methods:
204  * - double getX(); // [m]
205  * - double getY(); // [m]
206  * - double getZ(); // [m]
207  * - double getDX(); // [u]
208  * - double getDY(); // [u]
209  * - double getDZ(); // [u]
210  * - double getT(); // [ns]
211  *
212  * \param track track
213  * \param hit hit
214  * \return chi2 and gradient
215  */
216  template<class JHit_t>
217  result_type operator()(const JLine3Z& track, const JHit_t& hit) const
218  {
219  using namespace JPP;
220 
221  JPosition3D D(hit.getX(), hit.getY(), hit.getZ());
222  JDirection3D U(hit.getDX(), hit.getDY(), hit.getDZ());
223 
224  D.sub(track.getPosition());
225 
226  const double z = D.getDot(track.getDirection());
227  const double x = D.getX() - z * track.getDX();
228  const double y = D.getY() - z * track.getDY();
229  const double R2 = D.getLengthSquared() - z*z;
230  const double R = (R2 > Rmin_m*Rmin_m ? sqrt(R2) : Rmin_m);
231 
232  const double t1 = track.getT() + (z + R * getTanThetaC()) * getInverseSpeedOfLight();
233 
234  U.rotate(JRotation3Z(-atan2(y,x))); // rotate PMT axis to x-z plane
235 
236  const double theta = U.getTheta();
237  const double phi = fabs(U.getPhi());
238 
239  const double E = gWater.getE(E_GeV, z);
240  const double dt = T_ns.constrain(hit.getT() - t1);
241 
242  JPDF_t::result_type H0 = getH0(hit.getR(), dt);
243  JPDF_t::result_type H1 = getH1(E, R, theta, phi, dt);
244 
245  if (H1.V >= Vmax_npe) {
246  H1 *= Vmax_npe / H1.V;
247  }
248 
249  H1 += H0; // signal + background
250 
252 
253  result.chi2 = H1.getChi2() - H0.getChi2(); // Likelihood ratio
254 
255  const double wc = 1.0 - getTanThetaC() * z / R; // d(ct1)/d(z)
256 
257  result.gradient = JLine3Z(JLine1Z(JPosition3D(-getTanThetaC() * D.getX() / R, // d(ct1)/d(x)
258  -getTanThetaC() * D.getY() / R, // d(ct1)/d(y)
259  0.0),
260  getSpeedOfLight()), // d(ct1)/d(t)
261  JVersor3Z(wc * (D.getX() - D.getZ()*track.getDX()/track.getDZ()), // d(ct1)/d(dx)
262  wc * (D.getY() - D.getZ()*track.getDY()/track.getDZ()))); // d(ct1)/d(dy)
263 
264  result.gradient.mul(getInverseSpeedOfLight() * (H1.getDerivativeOfChi2() -
265  H0.getDerivativeOfChi2())); // x d(chi2)/d(ct1)
266 
267  return result;
268  }
269 
270 
271  /**
272  * Fit function.
273  * This method is used to determine the chi2 and gradient of given PMT with respect to trajectory of muon.
274  *
275  * \param track track
276  * \param pmt pmt
277  * \return chi2 and gradient
278  */
279  result_type operator()(const JLine3Z& track, const JPMTW0& pmt) const
280  {
281  using namespace JGEOMETRY3D;
282  using namespace JPHYSICS;
283 
284  JPosition3D D(pmt.getPosition());
285  JDirection3D U(pmt.getDirection());
286 
287  D.sub(track.getPosition());
288 
289  const double z = D.getDot(track.getDirection());
290  const double x = D.getX() - z * track.getDX();
291  const double y = D.getY() - z * track.getDY();
292  const double R2 = D.getLengthSquared() - z*z;
293  const double R = (R2 > Rmin_m*Rmin_m ? sqrt(R2) : Rmin_m);
294 
295  U.rotate(JRotation3Z(-atan2(y,x))); // rotate PMT axis to x-z plane
296 
297  const double theta = U.getTheta();
298  const double phi = fabs(U.getPhi());
299 
300  const double E = gWater.getE(E_GeV, z);
301 
302  JNPE_t::result_type H0 = getH0(pmt.getR());
303  JNPE_t::result_type H1 = getH1(E, R, theta, phi);
304 
305  if (H1.f >= Vmax_npe) {
306  H1 *= Vmax_npe / H1.f;
307  }
308 
309  H1 += H0;
310 
311  const bool hit = pmt.getN() != 0;
312  const double u = H1.getChi2(hit);
313 
315 
316  result.chi2 = estimator->getRho(u);
317 
318  result.gradient = JLine3Z(JLine1Z(JPosition3D(-D.getX() / R, // d(R)/d(x)
319  -D.getY() / R, // d(R)/d(y)
320  0.0),
321  0.0),
322  JVersor3Z(-z * D.getX() / R, // d(R)/d(dx)
323  -z * D.getY() / R)); // d(R)/d(dy)
324 
325  result.gradient.mul(estimator->getPsi(u));
326  result.gradient.mul(H1.getDerivativeOfChi2(hit)); // x d(chi2)/d(R)
327 
328  return result;
329  }
330 
331 
332  /**
333  * Get background hypothesis value for time differentiated PDF.
334  *
335  * \param R_Hz rate [Hz]
336  * \param t1 time [ns]
337  * \return hypothesis value
338  */
339  JPDF_t::result_type getH0(const double R_Hz,
340  const double t1) const
341  {
342  return JPDF_t::result_type(R_Hz * 1e-9, t1, T_ns);
343  }
344 
345 
346  /**
347  * Get signal hypothesis value for time differentiated PDF.
348  *
349  * \param E muon energy at minimum distance of approach [GeV]
350  * \param R minimum distance of approach [m]
351  * \param theta PMT zenith angle [rad]
352  * \param phi PMT azimuth angle [rad]
353  * \param t1 arrival time relative to Cherenkov hypothesis [ns]
354  * \return hypothesis value
355  */
357  const double R,
358  const double theta,
359  const double phi,
360  const double t1) const
361  {
362  using namespace std;
363  using namespace JPP;
364 
366 
367  for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
368 
369  if (!pdf[i].empty() && R <= pdf[i].getXmax()) {
370 
371  JPDF_t::result_type y1 = pdf[i](max(R, pdf[i].getXmin()), theta, phi, t1);
372 
373  // safety measures
374 
375  if (y1.f <= 0.0) {
376  y1.f = 0.0;
377  y1.fp = 0.0;
378  }
379 
380  if (y1.v <= 0.0) {
381  y1.v = 0.0;
382  }
383 
384  // energy dependence
385 
386  if (is_deltarays(pdf_t[i])) {
387  y1 *= getDeltaRaysFromMuon(E);
388  } else if (is_bremsstrahlung(pdf_t[i])) {
389  y1 *= E;
390  }
391 
392  h1 += y1;
393  }
394  }
395 
396  return h1;
397  }
398 
399 
400  /**
401  * Get background hypothesis value for time integrated PDF.
402  *
403  * \param R_Hz rate [Hz]
404  * \return hypothesis value
405  */
406  JNPE_t::result_type getH0(const double R_Hz) const
407  {
408  return JNPE_t::result_type(R_Hz * 1e-9 * T_ns.getLength(), 0.0);
409  }
410 
411 
412  /**
413  * Get signal hypothesis value for time integrated PDF.
414  *
415  * \param E muon energy at minimum distance of approach [GeV]
416  * \param R minimum distance of approach [m]
417  * \param theta PMT zenith angle [rad]
418  * \param phi PMT azimuth angle [rad]
419  * \return hypothesis value
420  */
422  const double R,
423  const double theta,
424  const double phi) const
425  {
426  using namespace std;
427  using namespace JPP;
428 
430 
431  for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
432 
433  if (!npe[i].empty() && R <= npe[i].getXmax()) {
434 
435  try {
436 
437  JNPE_t::result_type y1 = npe[i](max(R, npe[i].getXmin()), theta, phi);
438 
439  // safety measures
440 
441  if (y1.f > 0.0) {
442 
443  // energy dependence
444 
445  if (is_deltarays(pdf_t[i])) {
446  y1 *= getDeltaRaysFromMuon(E);
447  } else if (is_bremsstrahlung(pdf_t[i])) {
448  y1 *= E;
449  }
450 
451  h1 += y1;
452  }
453  }
454  catch(JException& error) {
455  ERROR(error << endl);
456  }
457  }
458  }
459 
460  return h1;
461  }
462 
463 
464  /**
465  * Get maximal road width of PDF.
466  *
467  * \return road width [m]
468  */
469  inline double getRmax() const
470  {
471  double xmax = 0.0;
472 
473  for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
474  if (!pdf[i].empty() && pdf[i].getXmax() > xmax) {
475  xmax = pdf[i].getXmax();
476  }
477  }
478 
479  return xmax;
480  }
481 
482 
483  /**
484  * Set maximal road width of PDF.
485  *
486  * \param Rmax road width [m]
487  */
488  inline void setRmax(const double Rmax)
489  {
490  for (int i = 0; i != NUMBER_OF_PDFS; ++i) {
491  if (!pdf[i].empty()) {
492  pdf[i].compress(JTOOLS::JRange<typename JFunction1D_t::abscissa_type>(0.0, Rmax));
493  }
494  }
495  }
496 
497 
498  static JTimeRange T_ns; //!< Time window with respect to Cherenkov hypothesis [ns]
499  static double Vmax_npe; //!< Maximal integral of PDF [npe]
500  static double Rmin_m; //!< Minimal distance of [m]
501 
502  static const int NUMBER_OF_PDFS = 6;
503 
504  static const JPDFType_t pdf_t[NUMBER_OF_PDFS];
505 
506  JPDF_t pdf[NUMBER_OF_PDFS]; //!< PDF
507  JNPE_t npe[NUMBER_OF_PDFS]; //!< PDF
508  double E_GeV; //!< Energy of muon at vertex [GeV]
509 
511  };
512 
513 
514  /**
515  * PDF types.
516  */
523 
524 
525  /**
526  * Default values.
527  */
529  double JRegressor<JLine3Z, JGandalf>::Vmax_npe = std::numeric_limits<double>::max();
531 }
532 
533 #endif
Template definition of a data regressor of given model.
Definition: JRegressor.hh:68
Type definition of a spline interpolation method based on a JCollection with JResultPDF result type...
General exception.
Definition: JException.hh:23
General purpose data regression method.
double getR() const
Get rate.
Definition: JPMTW0.hh:56
JTOOLS::JSplineFunction1S_t JFunction1D_t
JPDF_t::result_type getH0(const double R_Hz, const double t1) const
Get background hypothesis value for time differentiated PDF.
JPDF_t::result_type getH1(const double E, const double R, const double theta, const double phi, const double t1) const
Get signal hypothesis value for time differentiated PDF.
Auxiliary methods for PDF calculations.
JNPE_t::result_type getH1(const double E, const double R, const double theta, const double phi) const
Get signal hypothesis value for time integrated PDF.
then usage E
Definition: JMuonPostfit.sh:35
JRegressor(const std::string &fileDescriptor, const double TTS, const int numberOfPoints=25, const double epsilon=1.0e-10)
Constructor.
double getKappaC()
Get average R-dependence of arrival time of Cherenkov light (a.k.a.
Auxiliary class for handling PMT geometry, rate and response.
Definition: JPMTW0.hh:22
Energy loss of muon.
double getT(const JVector3D &pos) const
Get arrival time of Cherenkov light at given position.
Definition: JLine3Z.hh:234
const JDirection3D & getDirection() const
Get direction.
JTOOLS::JRange< double > JTimeRange
Type definition for time range (unit [ns]).
JRegressor(double sigma)
Constructor.
direct light from EM showers
Definition: JPDFTypes.hh:32
void setRmax(const double Rmax)
Set maximal road width of PDF.
void blur(const double TTS, const int numberOfPoints=25, const double epsilon=1.0e-10, const double quantile=0.99)
Blur PDF.
Definition: JPDFTable.hh:110
then for HISTOGRAM in h0 h1
Definition: JMatrixNZ.sh:71
Data structure for fit of straight line in positive z-direction.
Definition: JLine3Z.hh:36
static const JGeaneWater gWater
Function object for energy loss of muon in sea water.
Definition: JGeane.hh:323
static const JZero zero
Function object to assign zero value.
Definition: JZero.hh:105
JNPE_t::result_type getH0(const double R_Hz) const
Get background hypothesis value for time integrated PDF.
Type definition of a 1st degree polynomial interpolation based on a JMap implementation.
Definition of zero value for any class.
static double Rmin_m
Minimal distance of [m].
direct light from muon
Definition: JPDFTypes.hh:29
double getDY() const
Get y direction.
Definition: JVersor3Z.hh:158
Various implementations of functional maps.
static JTimeRange T_ns
Time window with respect to Cherenkov hypothesis [ns].
Numbering scheme for PDF types.
result_type operator()(const JLine3Z &track, const JHit_t &hit) const
Fit function.
bool is_bremsstrahlung(const int pdf)
Test if given PDF type corresponds to Cherenkov light from Bremsstrahlung.
Definition: JPDFTypes.hh:163
JVersor3D getDirection(const JVector3D &pos) const
Get photon direction of Cherenkov light on PMT.
Definition: JLine3Z.hh:255
scattered light from muon
Definition: JPDFTypes.hh:30
Type definition of a 1st degree polynomial interpolation based on a JGridMap implementation.
minimiser_type::result_type result_type
Definition: JRegressor.hh:80
The template JSharedPointer class can be used to share a pointer to an object.
JTOOLS::JMAPLIST< JTOOLS::JPolint1FunctionalMapH, JTOOLS::JPolint1FunctionalGridMap, JTOOLS::JPolint1FunctionalGridMap >::maplist JNPEMaplist_t
Type definition of a 1st degree polynomial interpolation based on a JMap implementation.
return result
Definition: JPolint.hh:743
double E_GeV
Energy of muon at vertex [GeV].
scattered light from delta-rays
Definition: JPDFTypes.hh:36
#define NOTICE(A)
Definition: JMessage.hh:64
#define ERROR(A)
Definition: JMessage.hh:66
Physics constants.
scattered light from EM showers
Definition: JPDFTypes.hh:33
Normal M-estimator.
Definition: JMEstimator.hh:66
JLANG::JSharedPointer< JMEstimator > estimator
M-Estimator function.
const JPosition3D & getPosition() const
Get position.
Definition: JPosition3D.hh:130
double operator()(const JLine3Z &track, const JHit_t &hit) const
Fit function.
JPDFType_t
PDF types.
Definition: JPDFTypes.hh:27
result_type operator()(const JLine3Z &track, const JPMTW0 &pmt) const
Fit function.
JPHYSICS::JNPETable< double, double, JNPEMaplist_t > JNPE_t
time integrated PDF
Auxiliary class for recursive map list generation.
Definition: JMapList.hh:108
then usage $script[distance] fi case set_variable R
Definition: JDrawLED.sh:43
Range of values.
Definition: JRange.hh:38
General purpose messaging.
double getRmax() const
Get maximal road width of PDF.
#define FATAL(A)
Definition: JMessage.hh:67
Fit method based on the Levenberg-Marquardt method.
Definition: JGandalf.hh:52
direct light from delta-rays
Definition: JPDFTypes.hh:35
virtual double getE(const double E, const double dx) const override
Get energy of muon after specified distance.
Definition: JGeane.hh:255
virtual const char * what() const override
Get error message.
Definition: JException.hh:48
const double getSpeedOfLight()
Get speed of light.
double sigma
Time resolution [ns].
int getN() const
Get number of hits.
Definition: JPMTW0.hh:67
Auxiliary class to define a range between two values.
Simple fit method based on Powell&#39;s algorithm, see reference: Numerical Recipes in C++...
Definition: JSimplex.hh:42
static double Vmax_npe
Maximal integral of PDF [npe].
const double getInverseSpeedOfLight()
Get inverse speed of light.
double getDX() const
Get x direction.
Definition: JVersor3Z.hh:147
Type definition of a zero degree polynomial interpolation based on a JGridMap implementation.
int numberOfPoints
Definition: JResultPDF.cc:22
bool is_deltarays(const int pdf)
Test if given PDF type corresponds to Cherenkov light from delta-rays.
Definition: JPDFTypes.hh:177
Data structure for fit of straight line paralel to z-axis.
Definition: JLine1Z.hh:27
double getDeltaRaysFromMuon(const double E)
Equivalent EM-shower energy due to delta-rays per unit muon track length.
Definition: JPDFToolkit.hh:67
std::string getFilename(const std::string &file_name)
Get file name part, i.e. part after last JEEP::PATHNAME_SEPARATOR if any.
Definition: JeepToolkit.hh:88
double getTanThetaC()
Get average tangent of Cherenkov angle of water corresponding to group velocity.
JPHYSICS::JPDFTable< JFunction1D_t, JPDFMaplist_t > JPDF_t
time dependent PDF
Data structure for position in three dimensions.
Definition: JPosition3D.hh:36
transformablemultifunction_type::result_type result_type
Definition: JPDFTable.hh:49
double getDZ() const
Get z direction.
Definition: JVersor3Z.hh:169
JLANG::JSharedPointer< JMEstimator > estimator
M-Estimator function.
double u[N+1]
Definition: JPolint.hh:755
Data structure for normalised vector in positive z-direction.
Definition: JVersor3Z.hh:39
JTOOLS::JMAPLIST< JTOOLS::JPolint1FunctionalMap, JTOOLS::JPolint0FunctionalGridMap, JTOOLS::JPolint0FunctionalGridMap >::maplist JPDFMaplist_t
Abstract class for global fit method.
Definition: JRegressor.hh:75
do echo Generating $dir eval D
Definition: JDrawLED.sh:53
Maximum likelihood estimator (M-estimators).