Jpp  18.4.0
the software that should make you happy
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Classes | Typedefs | Enumerations | Functions
JFIT Namespace Reference

Auxiliary classes and methods for linear and iterative data regression. More...

Classes

class  JEnergy
 Data structure for fit of energy. More...
 
struct  JRegressor< JEnergy >
 Regressor function object for fit of muon energy. More...
 
class  JEstimator
 Template definition of linear fit. More...
 
class  JGandalf
 Fit method based on the Levenberg-Marquardt method. More...
 
struct  JParameter_t
 Auxiliary data structure for fit parameter. More...
 
struct  JModifier_t
 Auxiliary data structure for editable parameter. More...
 
struct  JGradient
 Conjugate gradient fit. More...
 
struct  JK40
 Auxiliary class for handling light yields due to K40 decays. More...
 
struct  JK40Hit
 Auxiliary class for simultaneously handling light yields and response of module. More...
 
class  JEstimator< JLegendre< JOrdinate_t, N > >
 Linear fit of Legendre polynomial. More...
 
class  JLine1Z
 Data structure for fit of straight line paralel to z-axis. More...
 
class  JEstimator< JLine1Z >
 Linear fit of straight line parallel to z-axis to set of hits (objects with position and time). More...
 
class  JLine3EZ
 Data structure for fit of straight line in positive z-direction with energy. More...
 
class  JLine3Z
 Data structure for fit of straight line in positive z-direction. More...
 
struct  JRegressor< JLine3Z, JSimplex >
 Regressor function object for JLine3Z fit using JSimplex minimiser. More...
 
struct  JRegressor< JLine3Z, JGandalf >
 Regressor function object for JLine3Z fit using JGandalf minimiser. More...
 
class  JMatrixNZ
 Determination of the co-variance matrix of hits for a track along z-axis (JFIT::JLine1Z). More...
 
struct  JMEstimator
 Interface for maximum likelihood estimator (M-estimator). More...
 
struct  JMEstimatorNormal
 Normal M-estimator. More...
 
struct  JMEstimatorLorentzian
 Lorentzian M-estimator. More...
 
struct  JMEstimatorLinear
 Linear M-estimator. More...
 
struct  JMEstimatorNull
 Null M-estimator. More...
 
struct  JMEstimatorTukey
 Tukey's biweight M-estimator. More...
 
struct  JMEstimatorNormalWithBackground
 Normal M-estimator with background. More...
 
struct  JModel
 Auxiliary class to match data points with given model. More...
 
struct  JModel< JLine1Z >
 Template specialisation of class JModel to match hit with muon trajectory along z-axis. More...
 
struct  JModel< JEnergy >
 Template specialisation of class JModel to match hit with muon energy. More...
 
struct  JModel< JPoint4D >
 Template specialisation of class JModel to match hit with bright point. More...
 
struct  JNPE
 Auxiliary class for handling various light yields. More...
 
struct  JNPEHit
 Auxiliary class for simultaneously handling light yields and response of PMT. More...
 
struct  JPMTW0
 Auxiliary class for handling PMT geometry, rate and response. More...
 
class  JPoint3D
 Data structure for position fit. More...
 
class  JEstimator< JPoint3D >
 Linear fit of crossing point (position) between axes (objects with position and direction). More...
 
class  JPoint4D
 Data structure for vertex fit. More...
 
class  JEstimator< JPoint4D >
 Linear fit of bright point (position and time) between hits (objects with position and time). More...
 
struct  JRegressor< JPoint4D, JSimplex >
 Regressor function object for JPoint4D fit using JSimplex minimiser. More...
 
class  JAbstractMinimiser
 Abstract minimiser. More...
 
struct  JRegressor
 Template definition of a data regressor of given model. More...
 
struct  JAbstractRegressor
 Abstract class for global fit method. More...
 
class  JSimplex
 Simple fit method based on Powell's algorithm, see reference: Numerical Recipes in C++, W.H. More...
 
struct  JRegressorHelper
 Template data structure for PDF tables. More...
 
struct  JRegressorHelper< JLine3Z, JGandalf >
 Template specialisation for PDF tables. More...
 
struct  JRegressorStorage
 Template data structure for storage for PDF tables. More...
 
struct  JRegressorStorage< JLine3Z, JGandalf >
 Template specialisation for storage for PDF tables. More...
 
class  JShower3EZ
 Data structure for fit of straight line in positive z-direction with energy. More...
 
struct  JRegressor< JShower3EZ, JSimplex >
 Regressor function object for JShower3EZ fit using JSimplex minimiser. More...
 
struct  JRegressor< JShower3EZ, JGandalf >
 Regressor function object for JShower3EZ fit using JGandalf minimiser. More...
 
class  JShower3Z
 Data structure for cascade in positive z-direction. More...
 
struct  JRegressor< JShowerEH, JSimplex >
 Regressor function object for JShowerEH fit using JSimplex minimiser. More...
 
struct  JRegressor< JPoint4D, JGandalf >
 Regressor function object for JPoint4D fit using JGandalf minimiser. More...
 
class  JShowerEH
 Data structure for fit of straight line in positive z-direction with energy. More...
 
struct  JRegressor< JEnergy, JSimplex >
 Regressor function object for JShower3EZ fit using JSimplex minimiser. More...
 
struct  JShowerNPE
 Auxiliary class for handling EM shower light yield. More...
 
struct  JShowerNPEHit
 Auxiliary class for simultaneously handling light yields and response of PMT. More...
 
class  JVectorNZ
 Determination of the time residual vector of hits for a track along z-axis (JFIT::JLine1Z). More...
 
class  JEnergyCorrection
 Auxiliary class for correction of energy determined by JEnergy.cc. More...
 
class  JFit
 Data structure for track fit results with history and optional associated values. More...
 
class  JEvt
 Data structure for set of track fit results. More...
 
struct  JEvent
 Auxiliary class for historical event. More...
 
struct  JHistory
 Container for historical events. More...
 
class  JShowerEnergyCorrection
 Auxiliary class for correction of energy determined by JShowerEnergy.cc. More...
 

Typedefs

typedef JTOOLS::JRange< double > JZRange
 

Enumerations

enum  JMEstimator_t {
  EM_NORMAL = 0, EM_LORENTZIAN = 1, EM_LINEAR = 2, EM_NULL = 3,
  EM_TUKEY = 4, EM_NORMALWITHBACKGROUND = 5
}
 Definition of the various M-Estimators available to use. More...
 

Functions

double getP (const double expval, bool hit)
 Get Poisson probability to observe a hit or not for given expectation value for the number of hits. More...
 
double getChi2 (const double P)
 Get chi2 corresponding to given probability. More...
 
double getChi2 (const double expval, bool hit)
 Get chi2 to observe a hit or not for given expectation value for the number of hits. More...
 
template<class JModel_t , class JHit_t >
double getChi2 (const JModel_t &model, const JHit_t &hit, const double sigma)
 Determine chi2 of a hit for a given model and time resolution. More...
 
template<class JModel_t , class T >
double getChi2 (const JModel_t &model, T __begin, T __end, const double sigma)
 Determine chi2 of data for given model and time resolution. More...
 
double getChi2 (const JVectorNZ &Y, const JMatrixNZ &V)
 Determine chi2 using full covariance matrix. More...
 
template<class T >
double getChi2 (const JLine1Z &track, T __begin, T __end, const JMatrixNZ &V)
 Determine chi2 of data for given track using full covariance matrix. More...
 
template<class T >
double getChi2 (const JLine1Z &track, T __begin, T __end, const double alpha, const double sigma)
 Determine chi2 of data for given track and angular and time resolution. More...
 
double getChi2 (const JVectorNZ &Y, const JMatrixNZ &V, const int i)
 Determine difference between chi2 with and without hit using full covariance matrix. More...
 
template<class JModel_t , class JFit_t , class T >
double getChi2 (const JModel_t &model, const JFit_t &fit, T __begin, T __end)
 Get chi2 of data for given model and fit function. More...
 
double getProbability (const size_t N, const size_t M, const JK40Rates &R_Hz, const double T_ns)
 Get probability due to random background. More...
 
JMEstimatorgetMEstimator (const int type)
 Get M-Estimator. More...
 

Detailed Description

Auxiliary classes and methods for linear and iterative data regression.

Author
mdejong
adomi
lquinn
mdejong, adomi

Typedef Documentation

typedef JTOOLS::JRange<double> JFIT::JZRange

Definition at line 21 of file JFit/JModel.hh.

Enumeration Type Documentation

Definition of the various M-Estimators available to use.

Enumerator
EM_NORMAL 
EM_LORENTZIAN 
EM_LINEAR 
EM_NULL 
EM_TUKEY 
EM_NORMALWITHBACKGROUND 

Definition at line 185 of file JMEstimator.hh.

Function Documentation

double JFIT::getP ( const double  expval,
bool  hit 
)
inline

Get Poisson probability to observe a hit or not for given expectation value for the number of hits.

Parameters
expvalexpectation value
hithit
Returns
probability

Definition at line 41 of file JFitToolkit.hh.

42  {
43  if (hit)
44  return -expm1(-expval); // 1 - P(0)
45  else
46  return exp(-expval); // P(0)
47  }
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable STRING $argv[2] set_array QUANTILES set_variable FORMULA *[0] exp(-0.5 *(x-[1])*(x-[1])/([2]*[2]))" set_variable MODULE `getModule -a $DETECTOR -L "$STRING 0"` source JAcousticsToolkit.sh typeset -A TRIPODS get_tripods $WORKDIR/tripod.txt TRIPODS XMEAN
double JFIT::getChi2 ( const double  P)
inline

Get chi2 corresponding to given probability.

Parameters
Pprobability
Returns
chi2

Definition at line 56 of file JFitToolkit.hh.

57  {
58  return -log(P);
59  }
then set_variable DIR else fatal Wrong number of arguments fi for INPUT_FILE in ls rt $DIR stage * log
then $DIR JPlotNPE PDG P
Definition: JPlotNPE-PDG.sh:62
double JFIT::getChi2 ( const double  expval,
bool  hit 
)
inline

Get chi2 to observe a hit or not for given expectation value for the number of hits.

Parameters
expvalexpectation value
hithit
Returns
probability

Definition at line 69 of file JFitToolkit.hh.

70  {
71  if (hit)
72  return -log1p(-exp(-expval));
73  else
74  return expval;
75  }
then fatal Wrong number of arguments fi set_variable DETECTOR $argv[1] set_variable STRING $argv[2] set_array QUANTILES set_variable FORMULA *[0] exp(-0.5 *(x-[1])*(x-[1])/([2]*[2]))" set_variable MODULE `getModule -a $DETECTOR -L "$STRING 0"` source JAcousticsToolkit.sh typeset -A TRIPODS get_tripods $WORKDIR/tripod.txt TRIPODS XMEAN
template<class JModel_t , class JHit_t >
double JFIT::getChi2 ( const JModel_t &  model,
const JHit_t &  hit,
const double  sigma 
)
inline

Determine chi2 of a hit for a given model and time resolution.

The template argument JModel_t refers to a data structure which should provide for the following method:

  double getT(const JHit_t& hit) const;    // get expected time of hit
Parameters
modelmodel
hithit
sigmasigma
Returns
chi2

Definition at line 93 of file JFitToolkit.hh.

94  {
95  const double ds = (hit.getT() - model.getT(hit)) / sigma;
96 
97  return ds * ds;
98  }
const double sigma[]
Definition: JQuadrature.cc:74
template<class JModel_t , class T >
double JFIT::getChi2 ( const JModel_t &  model,
T  __begin,
T  __end,
const double  sigma 
)
inline

Determine chi2 of data for given model and time resolution.

The template argument JModel_t refers to a data structure which should provide for the following method:

  double getT(const value_type& hit) const;    // expected time of hit

where value_type corresponds to the value type if the input data.

Parameters
modelmodel
__beginbegin of data
__endend of data
sigmatime resolution [ns]
Returns
chi2

Definition at line 118 of file JFitToolkit.hh.

119  {
120  double chi2 = 0.0;
121 
122  for (T hit = __begin; hit != __end; ++hit) {
123  chi2 += getChi2(model, *hit, sigma);
124  }
125 
126  return chi2;
127  }
const double sigma[]
Definition: JQuadrature.cc:74
do set_variable OUTPUT_DIRECTORY $WORKDIR T
then if[[!-f $DETECTOR]] then JDetector sh $DETECTOR fi cat $WORKDIR trigger_parameters txt<< EOFtrigger3DMuon.enabled=1;trigger3DMuon.numberOfHits=5;trigger3DMuon.gridAngle_deg=1;ctMin=0.0;TMaxLocal_ns=15.0;EOF set_variable TRIGGEREFFICIENCY_TRIGGERED_EVENTS_ONLY INPUT_FILES=() for((i=1;$i<=$NUMBER_OF_RUNS;++i));do JSirene.sh $DETECTOR $JPP_DATA/genhen.km3net_wpd_V2_0.evt.gz $WORKDIR/sirene_ ${i}.root JTriggerEfficiency.sh $DETECTOR $DETECTOR $WORKDIR/sirene_ ${i}.root $WORKDIR/trigger_efficiency_ ${i}.root $WORKDIR/trigger_parameters.txt $JPP_DATA/PMT_parameters.txt INPUT_FILES+=($WORKDIR/trigger_efficiency_ ${i}.root) done for ANGLE_DEG in $ANGLES_DEG[*];do set_variable SIGMA_NS 3.0 set_variable OUTLIERS 3 set_variable OUTPUT_FILE $WORKDIR/matrix\[${ANGLE_DEG}\deg\].root $JPP_DIR/examples/JReconstruction-f"$INPUT_FILES[*]"-o $OUTPUT_FILE-S ${SIGMA_NS}-A ${ANGLE_DEG}-O ${OUTLIERS}-d ${DEBUG}--!fiif[[$OPTION=="plot"]];then if((0));then for H1 in h0 h1;do JPlot1D-f"$WORKDIR/matrix["${^ANGLES_DEG}" deg].root:${H1}"-y"1 2e3"-Y-L TR-T""-\^"number of events [a.u.]"-> o chi2
Definition: JMatrixNZ.sh:106
double getChi2(const double P)
Get chi2 corresponding to given probability.
Definition: JFitToolkit.hh:56
double JFIT::getChi2 ( const JVectorNZ Y,
const JMatrixNZ V 
)
inline

Determine chi2 using full covariance matrix.

Parameters
Yresiduals
Vcovariance matrix
Returns
chi2

Definition at line 137 of file JFitToolkit.hh.

139  {
140  return V.getDot(Y);
141  }
static double getDot(const variance &first, const variance &second)
Get dot product.
Definition: JMatrixNZ.hh:196
template<class T >
double JFIT::getChi2 ( const JLine1Z track,
T  __begin,
T  __end,
const JMatrixNZ V 
)
inline

Determine chi2 of data for given track using full covariance matrix.

Parameters
tracktrack
__beginbegin of data
__endend of data
Vcovariance matrix
Returns
chi2

Definition at line 154 of file JFitToolkit.hh.

158  {
159  const JVectorNZ Y(track, __begin, __end);
160 
161  return getChi2(Y, V);
162  }
then fatal Wrong number of arguments fi set_variable STRING $argv[1] set_variable DETECTORXY_TXT $WORKDIR $DETECTORXY_TXT tail read X Y CHI2 RMS printf optimum n $X $Y $CHI2 $RMS awk v Y
Determination of the time residual vector of hits for a track along z-axis (JFIT::JLine1Z).
Definition: JVectorNZ.hh:21
double getChi2(const double P)
Get chi2 corresponding to given probability.
Definition: JFitToolkit.hh:56
template<class T >
double JFIT::getChi2 ( const JLine1Z track,
T  __begin,
T  __end,
const double  alpha,
const double  sigma 
)
inline

Determine chi2 of data for given track and angular and time resolution.

Parameters
tracktrack
__beginbegin of data
__endend of data
alphaangular resolution [deg]
sigmatime resolution [ns]
Returns
chi2

Definition at line 176 of file JFitToolkit.hh.

177  {
178  JMatrixNZ V(track, __begin, __end, alpha, sigma);
179 
180  V.invert();
181 
182  return getChi2(track, __begin, __end, V);
183  }
V(JDAQEvent-JTriggerReprocessor)*1.0/(JDAQEvent+1.0e-10)
const double sigma[]
Definition: JQuadrature.cc:74
Determination of the co-variance matrix of hits for a track along z-axis (JFIT::JLine1Z).
Definition: JMatrixNZ.hh:28
double getChi2(const double P)
Get chi2 corresponding to given probability.
Definition: JFitToolkit.hh:56
double JFIT::getChi2 ( const JVectorNZ Y,
const JMatrixNZ V,
const int  i 
)
inline

Determine difference between chi2 with and without hit using full covariance matrix.

Parameters
Yresiduals
Vcovariance matrix
iindex of excluded hit
Returns
chi2

Definition at line 194 of file JFitToolkit.hh.

197  {
198  double chi2 = 0.0;
199 
200  for (size_t j = 0; j != V.size(); ++j) {
201  chi2 += V(i,j) * Y[j];
202  }
203 
204  return chi2*chi2 / V(i,i);
205  }
V(JDAQEvent-JTriggerReprocessor)*1.0/(JDAQEvent+1.0e-10)
size_t size() const
Get dimension of matrix.
Definition: JMatrixND.hh:199
then if[[!-f $DETECTOR]] then JDetector sh $DETECTOR fi cat $WORKDIR trigger_parameters txt<< EOFtrigger3DMuon.enabled=1;trigger3DMuon.numberOfHits=5;trigger3DMuon.gridAngle_deg=1;ctMin=0.0;TMaxLocal_ns=15.0;EOF set_variable TRIGGEREFFICIENCY_TRIGGERED_EVENTS_ONLY INPUT_FILES=() for((i=1;$i<=$NUMBER_OF_RUNS;++i));do JSirene.sh $DETECTOR $JPP_DATA/genhen.km3net_wpd_V2_0.evt.gz $WORKDIR/sirene_ ${i}.root JTriggerEfficiency.sh $DETECTOR $DETECTOR $WORKDIR/sirene_ ${i}.root $WORKDIR/trigger_efficiency_ ${i}.root $WORKDIR/trigger_parameters.txt $JPP_DATA/PMT_parameters.txt INPUT_FILES+=($WORKDIR/trigger_efficiency_ ${i}.root) done for ANGLE_DEG in $ANGLES_DEG[*];do set_variable SIGMA_NS 3.0 set_variable OUTLIERS 3 set_variable OUTPUT_FILE $WORKDIR/matrix\[${ANGLE_DEG}\deg\].root $JPP_DIR/examples/JReconstruction-f"$INPUT_FILES[*]"-o $OUTPUT_FILE-S ${SIGMA_NS}-A ${ANGLE_DEG}-O ${OUTLIERS}-d ${DEBUG}--!fiif[[$OPTION=="plot"]];then if((0));then for H1 in h0 h1;do JPlot1D-f"$WORKDIR/matrix["${^ANGLES_DEG}" deg].root:${H1}"-y"1 2e3"-Y-L TR-T""-\^"number of events [a.u.]"-> o chi2
Definition: JMatrixNZ.sh:106
int j
Definition: JPolint.hh:792
template<class JModel_t , class JFit_t , class T >
double JFIT::getChi2 ( const JModel_t &  model,
const JFit_t fit,
T  __begin,
T  __end 
)
inline

Get chi2 of data for given model and fit function.

The template argument JFit_t refers to a data structure which should provide for the function object operator:

  double operator()(const JModel_t& model, const value_type&) const;    // chi2

where JModel_t refers to the given model and value_type to the value type if the input data. The return value should correspond to the chi2 of the hit.

Parameters
modelmodel
fitfit function
__beginbegin of data
__endend of data

Definition at line 227 of file JFitToolkit.hh.

228  {
229  double chi2 = 0.0;
230 
231  for (T hit = __begin; hit != __end; ++hit) {
232  chi2 += fit(model, *hit);
233  }
234 
235  return chi2;
236  }
do set_variable OUTPUT_DIRECTORY $WORKDIR T
then if[[!-f $DETECTOR]] then JDetector sh $DETECTOR fi cat $WORKDIR trigger_parameters txt<< EOFtrigger3DMuon.enabled=1;trigger3DMuon.numberOfHits=5;trigger3DMuon.gridAngle_deg=1;ctMin=0.0;TMaxLocal_ns=15.0;EOF set_variable TRIGGEREFFICIENCY_TRIGGERED_EVENTS_ONLY INPUT_FILES=() for((i=1;$i<=$NUMBER_OF_RUNS;++i));do JSirene.sh $DETECTOR $JPP_DATA/genhen.km3net_wpd_V2_0.evt.gz $WORKDIR/sirene_ ${i}.root JTriggerEfficiency.sh $DETECTOR $DETECTOR $WORKDIR/sirene_ ${i}.root $WORKDIR/trigger_efficiency_ ${i}.root $WORKDIR/trigger_parameters.txt $JPP_DATA/PMT_parameters.txt INPUT_FILES+=($WORKDIR/trigger_efficiency_ ${i}.root) done for ANGLE_DEG in $ANGLES_DEG[*];do set_variable SIGMA_NS 3.0 set_variable OUTLIERS 3 set_variable OUTPUT_FILE $WORKDIR/matrix\[${ANGLE_DEG}\deg\].root $JPP_DIR/examples/JReconstruction-f"$INPUT_FILES[*]"-o $OUTPUT_FILE-S ${SIGMA_NS}-A ${ANGLE_DEG}-O ${OUTLIERS}-d ${DEBUG}--!fiif[[$OPTION=="plot"]];then if((0));then for H1 in h0 h1;do JPlot1D-f"$WORKDIR/matrix["${^ANGLES_DEG}" deg].root:${H1}"-y"1 2e3"-Y-L TR-T""-\^"number of events [a.u.]"-> o chi2
Definition: JMatrixNZ.sh:106
double JFIT::getProbability ( const size_t  N,
const size_t  M,
const JK40Rates &  R_Hz,
const double  T_ns 
)
inline

Get probability due to random background.

Parameters
Nnumber of active PMTs
Mmultiplicity of hits
R_Hzrates [Hz]
T_nstime window [ns]
Returns
probability

Definition at line 248 of file JFitToolkit.hh.

252  {
253  using namespace JPP;
254 
255  if (N == 0) {
256  return (M == 0 ? 1.0 : 0.0);
257  }
258 
259  const double T = T_ns * 1.0e-9;
260  const double p = -expm1(-R_Hz.getSinglesRate() * N * T); // 1 - P(0) due to singles rate
261 
262  double P = 0.0;
263 
264  for (multiplicity_type i = R_Hz.getLowerL1Multiplicity(); i <= R_Hz.getUpperL1Multiplicity(); ++i) {
265  P -= expm1(-R_Hz.getMultiplesRate(i) * T); // 1 - P(0) due to multiples rates
266  }
267 
268  P *= (1.0 + p); // 1 - P(0) due to multiples rates M-1 combined with singles rate
269 
270  return (poisson(M, R_Hz.getSinglesRate() * N * T) - // P(M) due to singles rate
271  expm1(-R_Hz.getMultiplesRate(M) * T) - // P(M) due to multiples rate
272  expm1(-R_Hz.getMultiplesRate(M-1) * T) * p) / (1.0 + P); // P(M-1) * P(1) due to multiples rate combined with singles rate
273  }
do set_variable OUTPUT_DIRECTORY $WORKDIR T
size_t multiplicity_type
Type definition of multiplicity.
Definition: JK40Rates.hh:33
then usage $script< input file >[option[primary[working directory]]] nWhere option can be N
Definition: JMuonPostfit.sh:40
double poisson(const size_t n, const double mu)
Poisson probability density distribition.
then $DIR JPlotNPE PDG P
Definition: JPlotNPE-PDG.sh:62
JMEstimator* JFIT::getMEstimator ( const int  type)
inline

Get M-Estimator.

Note that for M-estimators with an additional parameter, defaults are used.

Parameters
typetype
Returns
pointer to newly created M-Estimator (may be NULL)

Definition at line 203 of file JMEstimator.hh.

204  {
205  switch (type) {
206 
207  case EM_NORMAL:
208  return new JMEstimatorNormal();
209 
210  case EM_LORENTZIAN:
211  return new JMEstimatorLorentzian();
212 
213  case EM_LINEAR:
214  return new JMEstimatorLinear();
215 
216  case EM_NULL:
217  return new JMEstimatorNull();
218 
219  case EM_TUKEY:
220  return new JMEstimatorTukey(5.0);
221 
223  return new JMEstimatorNormalWithBackground(1.0e-5);
224 
225  default:
226  return NULL;
227  }
228  }
Null M-estimator.
Definition: JMEstimator.hh:92
Linear M-estimator.
Definition: JMEstimator.hh:78
Normal M-estimator.
Definition: JMEstimator.hh:52
Tukey&#39;s biweight M-estimator.
Definition: JMEstimator.hh:105
Normal M-estimator with background.
Definition: JMEstimator.hh:152
Lorentzian M-estimator.
Definition: JMEstimator.hh:65