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