1 #ifndef __JFIT__JLINE3ZREGRESSOR__
2 #define __JFIT__JLINE3ZREGRESSOR__
36 namespace JPP {
using namespace JFIT; }
75 template<
class JHit_t>
85 const double R = sqrt(
D.getLengthSquared() -
z*
z);
89 const double u = (t1 - hit.getT()) /
sigma;
91 return estimator->getRho(
u) * hit.getW();
137 const double epsilon = 1.0e-10) :
150 pdf(storage.getPDF()),
172 template<
class JHit_t>
178 JDirection3D U(hit.getDX(), hit.getDY(), hit.getDZ());
183 const double x =
D.getX() -
z * track.
getDX();
184 const double y =
D.getY() -
z * track.
getDY();
185 const double R2 =
D.getLengthSquared() -
z*
z;
186 const double R = (R2 > Rmin_m*Rmin_m ? sqrt(R2) : Rmin_m);
190 U.rotate(JRotation3Z(-atan2(
y,
x)));
192 const double theta = U.getTheta();
193 const double phi = fabs(U.getPhi());
196 const double dt = T_ns.constrain(hit.getT() - t1);
201 if (H1.V >= Vmax_npe) {
202 H1 *= Vmax_npe / H1.V;
209 result.chi2 = H1.getChi2() - H0.getChi2();
218 wc * (
D.getY() -
D.getZ()*track.
getDY()/track.
getDZ())));
221 H0.getDerivativeOfChi2()));
237 using namespace JGEOMETRY3D;
238 using namespace JPHYSICS;
246 const double x =
D.getX() -
z * track.
getDX();
247 const double y =
D.getY() -
z * track.
getDY();
248 const double R2 =
D.getLengthSquared() -
z*
z;
249 const double R = (R2 > Rmin_m*Rmin_m ? sqrt(R2) : Rmin_m);
251 U.rotate(JRotation3Z(-atan2(
y,
x)));
253 const double theta = U.getTheta();
254 const double phi = fabs(U.getPhi());
261 if (H1.f >= Vmax_npe) {
262 H1 *= Vmax_npe / H1.f;
267 const bool hit = pmt.
getN() != 0;
268 const double u = H1.getChi2(hit);
272 result.chi2 = estimator->getRho(
u);
281 result.gradient.mul(estimator->getPsi(
u));
282 result.gradient.mul(H1.getDerivativeOfChi2(hit));
296 const double t1)
const
316 const double t1)
const
323 for (
int i = 0;
i != NUMBER_OF_PDFS; ++
i) {
325 if (!pdf[
i].empty() && R <= pdf[
i].getXmax()) {
380 const double phi)
const
387 for (
int i = 0;
i != NUMBER_OF_PDFS; ++
i) {
389 if (!npe[
i].empty() && R <= npe[
i].getXmax()) {
411 ERROR(error << endl);
429 for (
int i = 0;
i != NUMBER_OF_PDFS; ++
i) {
430 if (!pdf[
i].empty() && pdf[
i].getXmax() > xmax) {
431 xmax = pdf[
i].getXmax();
Template definition of a data regressor of given model.
General purpose data regression method.
double getR() const
Get rate.
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.
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 $script< input file >[option[primary[working directory]]] nWhere option can be E
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.
double getT(const JVector3D &pos) const
Get arrival time of Cherenkov light at given position.
const JDirection3D & getDirection() const
Get direction.
JRegressor(double sigma)
Constructor.
Template specialisation for storage for PDF tables.
Data structure for fit of straight line in positive z-direction.
static const JGeaneWater gWater
Function object for energy loss of muon in sea water.
static const JZero zero
Function object to assign zero value.
JNPE_t::result_type getH0(const double R_Hz) const
Get background hypothesis value for time integrated PDF.
Template data structure for storage for PDF tables.
Definition of zero value for any class.
static double Rmin_m
Minimal distance of [m].
double getDY() const
Get y direction.
Various implementations of functional maps.
static JTimeRange T_ns
Time window with respect to Cherenkov hypothesis [ns].
double getDeltaRaysFromMuon(const double E, const JRange< double > T_GeV=JRange< double >(DELTARAY_TMIN, DELTARAY_TMAX))
Equivalent EM-shower energy due to delta-rays per unit muon track length.
JRegressor(const JRegressorStorage< JLine3Z, JGandalf > &storage)
Constructor.
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.
JVersor3D getDirection(const JVector3D &pos) const
Get photon direction of Cherenkov light on PMT.
minimiser_type::result_type result_type
The template JSharedPointer class can be used to share a pointer to an object.
multifunction_t::result_type result_type
JLANG::JSharedPointer< JMEstimator > estimator
M-Estimator function.
const JPosition3D & getPosition() const
Get position.
double operator()(const JLine3Z &track, const JHit_t &hit) const
Fit function.
result_type operator()(const JLine3Z &track, const JPMTW0 &pmt) const
Fit function.
General purpose messaging.
double getRmax() const
Get maximal road width of PDF.
Fit method based on the Levenberg-Marquardt method.
then JCookie sh JDataQuality D $DETECTOR_ID R
virtual double getE(const double E, const double dx) const override
Get energy of muon after specified distance.
then usage $script[energy[distance[z of PMT]]] fi case set_variable z
const double getSpeedOfLight()
Get speed of light.
double sigma
Time resolution [ns].
int getN() const
Get number of hits.
Auxiliary class to define a range between two values.
Simple fit method based on Powell's algorithm, see reference: Numerical Recipes in C++...
static double Vmax_npe
Maximal integral of PDF [npe].
const double getInverseSpeedOfLight()
Get inverse speed of light.
double getDX() const
Get x direction.
bool is_deltarays(const int pdf)
Test if given PDF type corresponds to Cherenkov light from delta-rays.
Data structure for fit of straight line paralel to z-axis.
double getTanThetaC()
Get average tangent of Cherenkov angle of water corresponding to group velocity.
double getNPE(const Hit &hit)
Get true charge of hit.
JRegressorStorage< JLine3Z, JGandalf > JRegressorStorage_t
Data structure for position in three dimensions.
transformablemultifunction_type::result_type result_type
double getDZ() const
Get z direction.
Data structure for normalised vector in positive z-direction.
JRegressor()
Default constructor.
Abstract class for global fit method.
do echo Generating $dir eval D
Maximum likelihood estimator (M-estimators).
JTOOLS::JRange< double > JTimeRange
Type definition for time range (unit [s]).