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