1#ifndef __JFIT__JLINE3ZREGRESSOR__
2#define __JFIT__JLINE3ZREGRESSOR__
9#include "JPhysics/JConstants.hh"
10#include "JPhysics/JPDFTypes.hh"
11#include "JPhysics/JPDFTable.hh"
12#include "JPhysics/JNPETable.hh"
13#include "JPhysics/JDeltaRays.hh"
14#include "JPhysics/JGeane.hh"
15#include "JTools/JFunction1D_t.hh"
16#include "JTools/JFunctionalMap_t.hh"
17#include "JTools/JRange.hh"
18#include "JGeometry3D/JPosition3D.hh"
19#include "JGeometry3D/JDirection3D.hh"
20#include "JMath/JZero.hh"
28#include "Jeep/JMessage.hh"
37namespace JPP {
using namespace JFIT; }
41 using JPHYSICS::JPDFType_t;
42 using JPHYSICS::DIRECT_LIGHT_FROM_MUON;
43 using JPHYSICS::SCATTERED_LIGHT_FROM_MUON;
44 using JPHYSICS::DIRECT_LIGHT_FROM_DELTARAYS;
45 using JPHYSICS::SCATTERED_LIGHT_FROM_DELTARAYS;
46 using JPHYSICS::DIRECT_LIGHT_FROM_EMSHOWERS;
47 using JPHYSICS::SCATTERED_LIGHT_FROM_EMSHOWERS;
85 template<
class JHit_t>
90 JPosition3D D(hit.getX(), hit.getY(), hit.getZ());
92 D.sub(track.getPosition());
95 const double R = sqrt(D.getLengthSquared() - z*z);
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();
116 typedef JTOOLS::JMAPLIST<JTOOLS::JPolint1FunctionalMap,
117 JTOOLS::JPolint0FunctionalGridMap,
119 typedef JPHYSICS::JPDFTable<JFunction1D_t, JPDFMaplist_t>
JPDF_t;
121 typedef JTOOLS::JMAPLIST<JTOOLS::JPolint1FunctionalMapH,
122 JTOOLS::JPolint1FunctionalGridMap,
124 typedef JPHYSICS::JNPETable<double, double, JNPEMaplist_t>
JNPE_t;
128 static const int NUMBER_OF_PDFS = 6;
130 typedef std::array<JPDF_t, NUMBER_OF_PDFS>
JPDFs_t;
131 typedef std::array<JNPE_t, NUMBER_OF_PDFS>
JNPEs_t;
156 const JTimeRange& T_ns,
158 const int numberOfPoints = 25,
159 const double epsilon = 1.0e-10) :
165 const JPDF_t::JSupervisor supervisor(
new JPDF_t::JDefaultResult(JMATH::zero));
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());
173 _pdf[i].setExceptionHandler(supervisor);
178 for (
int i = 1; i < NUMBER_OF_PDFS; i += 2) {
180 _pdf[ i ].add(_pdf[i-1]);
181 _pdf[ i ].compress(numeric_limits<double>::max(), T_ns);
185 _pdf[i-1].swap(buffer);
188 _pdf[i].blur(TTS_ns, numberOfPoints, epsilon);
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,
298 const int numberOfPoints = 25,
299 const double epsilon = 1.0e-10) :
300 storage_type(fileDescriptor, T_ns, TTS_ns, numberOfPoints, epsilon),
313 pdf(storage.getPDF()),
314 npe(storage.getNPE()),
338 template<
class JHit_t>
343 JPosition3D D(hit.getX(), hit.getY(), hit.getZ());
344 JDirection3D U(hit.getDX(), hit.getDY(), hit.getDZ());
346 D.sub(track.getPosition());
349 const double x = D.getX() - z * track.getDX();
350 const double y = D.getY() - z * track.getDY();
351 const double R2 = D.getLengthSquared() - z*z;
352 const double R = sqrt(R2);
354 const double t1 = track.
getT() + (z + R * getTanThetaC()) * getInverseSpeedOfLight();
356 U.rotate(JRotation3Z(-atan2(y,x)));
358 const double theta = U.getTheta();
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);
364 JPDF_t::result_type H0 = getH0(hit.getR(), dt);
365 JPDF_t::result_type H1 = getH1(E, R, theta, phi, dt);
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,
383 JVersor3Z(wc * (D.getX() - D.getZ()*track.getDX()/track.getDZ()),
384 wc * (D.getY() - D.getZ()*track.getDY()/track.getDZ())));
386 result.
gradient.mul(getInverseSpeedOfLight() * (H1.getDerivativeOfChi2() -
387 H0.getDerivativeOfChi2()));
403 using namespace JGEOMETRY3D;
404 using namespace JPHYSICS;
406 JPosition3D D(pmt.getPosition());
407 JDirection3D U(pmt.getDirection());
409 D.sub(track.getPosition());
412 const double x = D.getX() - z * track.getDX();
413 const double y = D.getY() - z * track.getDY();
414 const double R2 = D.getLengthSquared() - z*z;
415 const double R = sqrt(R2);
417 U.rotate(JRotation3Z(-atan2(y,x)));
419 const double theta = U.getTheta();
420 const double phi = fabs(U.getPhi());
422 const double E = gWater.getE(E_GeV, z);
424 JNPE_t::result_type H0 = getH0(pmt.
getR());
425 JNPE_t::result_type H1 = getH1(E, R, theta, phi);
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);
444 JVersor3Z(-z * D.getX() / R,
447 result.
gradient.mul(estimator->getPsi(u));
448 result.
gradient.mul(H1.getDerivativeOfChi2(hit));
461 JPDF_t::result_type
getH0(
const double R_Hz,
462 const double t1)
const
464 return JPDF_t::result_type(R_Hz * 1e-9, t1, T_ns);
478 JPDF_t::result_type
getH1(
const double E,
482 const double t1)
const
487 JPDF_t::result_type h1 = zero;
489 for (
int i = 0; i != NUMBER_OF_PDFS; ++i) {
491 if (!pdf[i].empty() && R <= pdf[i].getXmax()) {
493 JPDF_t::result_type y1 = pdf[i](max(R, pdf[i].getXmin()), theta, phi, t1);
508 if (is_deltarays(pdf_t[i])) {
509 y1 *= JDeltaRays::getEnergyLossFromMuon(E);
510 }
else if (is_bremsstrahlung(pdf_t[i])) {
528 JNPE_t::result_type
getH0(
const double R_Hz)
const
530 return JNPE_t::result_type(R_Hz * 1e-9 * T_ns.getLength(), 0.0);
543 JNPE_t::result_type
getH1(
const double E,
546 const double phi)
const
551 JNPE_t::result_type h1 = JMATH::zero;
553 for (
int i = 0; i != NUMBER_OF_PDFS; ++i) {
555 if (!npe[i].empty() && R <= npe[i].getXmax()) {
559 JNPE_t::result_type y1 = npe[i](max(R, npe[i].getXmin()), theta, phi);
567 if (is_deltarays(pdf_t[i])) {
568 y1 *= JDeltaRays::getEnergyLossFromMuon(E);
569 }
else if (is_bremsstrahlung(pdf_t[i])) {
576 catch(JException& error) {
577 ERROR(error << endl);
595 for (
int i = 0; i != NUMBER_OF_PDFS; ++i) {
596 if (!pdf[i].empty() && pdf[i].getXmax() > xmax) {
597 xmax = pdf[i].getXmax();
Maximum likelihood estimator (M-estimators).
General purpose data regression method.
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.
JVersor3D getDirection(const JVector3D &pos) const
Get photon direction of Cherenkov light on PMT.
double getT(const JVector3D &pos) const
Get arrival time of Cherenkov light at given position.
Simple fit method based on Powell's algorithm, see reference: Numerical Recipes in C++,...
Auxiliary classes and methods for linear and iterative data regression.
Abstract class for global fit method.
Data structure for return value of fit function.
JModel_t gradient
partial derivatives of chi2
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.