1#ifndef __JFIT__JLINE3ZREGRESSOR__
2#define __JFIT__JLINE3ZREGRESSOR__
37namespace JPP {
using namespace JFIT; }
85 template<
class JHit_t>
97 const double t1 =
track.getT() + (z + R * getKappaC()) * getInverseSpeedOfLight();
99 const double u = (t1 -
hit.getT()) / sigma;
101 return estimator->getRho(
u) *
hit.getW();
128 static const int NUMBER_OF_PDFS = 6;
156 const JTimeRange& T_ns,
159 const double epsilon = 1.0e-10) :
167 for (
int i = 0; i != NUMBER_OF_PDFS; ++i) {
169 const string file_name = getFilename(fileDescriptor, pdf_t[i]);
171 _pdf[i].load(file_name.c_str());
178 for (
int i = 1; i < NUMBER_OF_PDFS; i += 2) {
180 _pdf[ i ].add(_pdf[i-1]);
185 _pdf[i-1].swap(buffer);
191 _npe[ i ] =
JNPE_t(_pdf[i], T_ns);
225 for (
int i = 0; i != NUMBER_OF_PDFS; ++i) {
226 _pdf[i].transform(transformer);
227 _npe[i].transform(transformer);
235 static const JPDFType_t pdf_t[NUMBER_OF_PDFS];
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
296 const JTimeRange& T_ns,
299 const double epsilon = 1.0e-10) :
313 pdf(storage.getPDF()),
314 npe(storage.getNPE()),
338 template<
class JHit_t>
349 const double x = D.
getX() - z *
track.getDX();
350 const double y = D.
getY() - z *
track.getDY();
352 const double R = sqrt(R2);
354 const double t1 =
track.getT() + (z + R * getTanThetaC()) * getInverseSpeedOfLight();
359 const double phi = fabs(U.
getPhi());
361 const double E = gWater.getE(E_GeV, z);
362 const double dt = T_ns.constrain(
hit.getT() - t1);
367 if (
H1.V >= Vmax_npe) {
368 H1 *= Vmax_npe /
H1.V;
375 result.chi2 =
H1.getChi2() -
H0.getChi2();
377 const double wc = 1.0 - getTanThetaC() * z / R;
380 -getTanThetaC() * D.
getY() / R,
386 result.gradient.mul(getInverseSpeedOfLight() * (
H1.getDerivativeOfChi2() -
387 H0.getDerivativeOfChi2()));
412 const double x = D.
getX() - z *
track.getDX();
413 const double y = D.
getY() - z *
track.getDY();
415 const double R = sqrt(R2);
420 const double phi = fabs(U.
getPhi());
422 const double E = gWater.getE(E_GeV, z);
427 if (
H1.f >= Vmax_npe) {
428 H1 *= Vmax_npe /
H1.f;
433 const bool hit = pmt.
getN() != 0;
434 const double u =
H1.getChi2(
hit);
438 result.chi2 = estimator->getRho(
u);
447 result.gradient.mul(estimator->getPsi(
u));
448 result.gradient.mul(
H1.getDerivativeOfChi2(
hit));
462 const double t1)
const
482 const double t1)
const
489 for (
int i = 0; i != NUMBER_OF_PDFS; ++i) {
491 if (!pdf[i].empty() && R <= pdf[i].getXmax()) {
508 if (is_deltarays(pdf_t[i])) {
509 y1 *= JDeltaRays::getEnergyLossFromMuon(E);
510 }
else if (is_bremsstrahlung(pdf_t[i])) {
546 const double phi)
const
553 for (
int i = 0; i != NUMBER_OF_PDFS; ++i) {
555 if (!npe[i].empty() && R <= npe[i].getXmax()) {
567 if (is_deltarays(pdf_t[i])) {
568 y1 *= JDeltaRays::getEnergyLossFromMuon(E);
569 }
else if (is_bremsstrahlung(pdf_t[i])) {
595 for (
int i = 0; i != NUMBER_OF_PDFS; ++i) {
596 if (!pdf[i].empty() && pdf[i].getXmax() > xmax) {
597 xmax = pdf[i].getXmax();
Various implementations of functional maps.
Maximum likelihood estimator (M-estimators).
General purpose messaging.
Numbering scheme for PDF types.
Auxiliary class to define a range between two values.
General purpose data regression method.
Definition of zero value for any class.
Fit method based on the Levenberg-Marquardt method.
Data structure for fit of straight line paralel to z-axis.
Data structure for fit of straight line in positive z-direction.
Simple fit method based on Powell's algorithm, see reference: Numerical Recipes in C++,...
Data structure for direction in three dimensions.
JDirection3D & rotate(const JRotation3D &R)
Rotate.
const JDirection3D & getDirection() const
Get direction.
Data structure for position in three dimensions.
double getDot(const JAngle3D &angle) const
Get dot product.
const JPosition3D & getPosition() const
Get position.
double getY() const
Get y position.
double getLengthSquared() const
Get length squared.
double getZ() const
Get z position.
JVector3D & sub(const JVector3D &vector)
Subtract vector.
double getX() const
Get x position.
double getTheta() const
Get theta angle.
double getPhi() const
Get phi angle.
Data structure for normalised vector in positive z-direction.
Template definition of a multi-dimensional oscillation probability interpolation table.
Custom class for integrated values of the PDF of the arrival time of Cherenkov light.
multifunction_t::result_type result_type
Multi-dimensional PDF table for arrival time of Cherenkov light.
transformablemultifunction_type::result_type result_type
transformablemultifunction_type::transformer_type transformer_type
void blur(const double TTS, const int numberOfPoints=25, const double epsilon=1.0e-10, const double quantile=0.99)
Blur PDF.
Auxiliary classes and methods for linear and iterative data regression.
Auxiliary classes and methods for 3D geometrical objects and operations.
static const JZero zero
Function object to assign zero value.
Auxiliary methods for light properties of deep-sea water.
@ SCATTERED_LIGHT_FROM_DELTARAYS
scattered light from delta-rays
@ DIRECT_LIGHT_FROM_EMSHOWERS
direct light from EM showers
@ SCATTERED_LIGHT_FROM_EMSHOWERS
scattered light from EM showers
@ SCATTERED_LIGHT_FROM_MUON
scattered light from muon
@ DIRECT_LIGHT_FROM_DELTARAYS
direct light from delta-rays
@ DIRECT_LIGHT_FROM_MUON
direct light from muon
This name space includes all other name spaces (except KM3NETDAQ, KM3NET and ANTARES).
Abstract class for global fit method.
Data structure for return value of fit function.
Auxiliary class for handling PMT geometry, rate and response.
int getN() const
Get number of hits.
double getR() const
Get rate.
Template specialisation for storage of PDF tables.
JTOOLS::JMAPLIST< JTOOLS::JPolint1FunctionalMapH, JTOOLS::JPolint1FunctionalGridMap, JTOOLS::JPolint1FunctionalGridMap >::maplist JNPEMaplist_t
const JNPEs_t & getNPE() const
Get NPEs.
JTOOLS::JMAPLIST< JTOOLS::JPolint1FunctionalMap, JTOOLS::JPolint0FunctionalGridMap, JTOOLS::JPolint0FunctionalGridMap >::maplist JPDFMaplist_t
JRegressorStorage()
Default constructor.
JPHYSICS::JPDFTable< JFunction1D_t, JPDFMaplist_t > JPDF_t
time dependent PDF
JPDF_t::transformer_type transformer_type
JTimeRange T_ns
Time window with respect to Cherenkov hypothesis [ns].
std::array< JNPE_t, NUMBER_OF_PDFS > JNPEs_t
NPEs.
std::array< JPDF_t, NUMBER_OF_PDFS > JPDFs_t
PDFs.
const JPDFs_t & getPDF() const
Get PDFs.
JPHYSICS::JNPETable< double, double, JNPEMaplist_t > JNPE_t
time integrated PDF
JRegressorStorage(const std::string &fileDescriptor, const JTimeRange &T_ns, const double TTS_ns, const int numberOfPoints=25, const double epsilon=1.0e-10)
Constructor.
JTOOLS::JSplineFunction1S_t JFunction1D_t
void transform(const transformer_type &transformer)
Transform PDFs and NPEs.
Template data structure for storage of internal data.
result_type operator()(const JLine3Z &track, const JHit_t &hit) const
Fit function.
JRegressor(const std::string &fileDescriptor, const JTimeRange &T_ns, const double TTS_ns, const int numberOfPoints=25, const double epsilon=1.0e-10)
Constructor.
JRegressorStorage< JLine3Z, JGandalf > storage_type
JNPE_t::result_type getH0(const double R_Hz) const
Get background hypothesis value for time integrated PDF.
std::shared_ptr< JMEstimator > estimator
M-Estimator function.
result_type operator()(const JLine3Z &track, const JPMTW0 &pmt) const
Fit function.
static double Vmax_npe
Maximal integral of PDF [npe].
JRegressor(const storage_type &storage)
Constructor.
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.
JPDF_t::result_type getH0(const double R_Hz, const double t1) const
Get background hypothesis value for time differentiated PDF.
JRegressor()
Default constructor.
double getRmax() const
Get maximal road width of 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.
std::shared_ptr< JMEstimator > estimator
M-Estimator function.
double operator()(const JLine3Z &track, const JHit_t &hit) const
Fit function.
double sigma
Time resolution [ns].
JRegressor(double sigma)
Constructor.
Template definition of a data regressor of given model.
Auxiliary class to set-up Hit.