1 #ifndef NPSTAT_DISTRIBUTIONS1D_HH_
2 #define NPSTAT_DISTRIBUTIONS1D_HH_
31 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
32 virtual bool write(std::ostream& os)
const;
34 static inline const char* classname() {
return "npstat::Uniform1D";}
35 static inline unsigned version() {
return 1;}
36 static Uniform1D* read(
const gs::ClassId&
id, std::istream& in);
46 const std::vector<double>& )
48 inline static int nParameters() {
return 0;}
50 double unscaledDensity(
double x)
const;
51 double unscaledCdf(
double x)
const;
52 double unscaledQuantile(
double x)
const;
54 inline double unscaledExceedance(
const double x)
const
55 {
return 1.0 - unscaledCdf(x);}
70 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
71 virtual bool write(std::ostream& os)
const;
73 static inline const char* classname()
74 {
return "npstat::IsoscelesTriangle1D";}
75 static inline unsigned version() {
return 1;}
76 static IsoscelesTriangle1D* read(
const gs::ClassId&
id, std::istream& in);
86 const std::vector<double>& )
88 inline static int nParameters() {
return 0;}
90 double unscaledDensity(
double x)
const;
91 double unscaledCdf(
double x)
const;
92 double unscaledQuantile(
double x)
const;
94 inline double unscaledExceedance(
const double x)
const
95 {
return 1.0 - unscaledCdf(x);}
110 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
111 virtual bool write(std::ostream& os)
const;
113 static inline const char* classname() {
return "npstat::Exponential1D";}
114 static inline unsigned version() {
return 1;}
115 static Exponential1D* read(
const gs::ClassId&
id, std::istream& in);
125 const std::vector<double>& )
127 inline static int nParameters() {
return 0;}
129 double unscaledDensity(
double x)
const;
130 double unscaledCdf(
double x)
const;
131 double unscaledQuantile(
double x)
const;
132 double unscaledExceedance(
double x)
const;
147 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
148 virtual bool write(std::ostream& os)
const;
150 static inline const char* classname() {
return "npstat::Laplace1D";}
151 static inline unsigned version() {
return 1;}
152 static Laplace1D* read(
const gs::ClassId&
id, std::istream& in);
162 const std::vector<double>& )
164 inline static int nParameters() {
return 0;}
166 double unscaledDensity(
double x)
const;
167 double unscaledCdf(
double x)
const;
168 double unscaledQuantile(
double x)
const;
169 double unscaledExceedance(
double x)
const;
184 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
185 virtual bool write(std::ostream& os)
const;
187 static inline const char* classname() {
return "npstat::Logistic1D";}
188 static inline unsigned version() {
return 1;}
189 static Logistic1D* read(
const gs::ClassId&
id, std::istream& in);
199 const std::vector<double>& )
201 inline static int nParameters() {
return 0;}
203 double unscaledDensity(
double x)
const;
204 double unscaledCdf(
double x)
const;
205 double unscaledQuantile(
double x)
const;
206 double unscaledExceedance(
double x)
const;
227 inline double a()
const {
return a_;}
228 inline double b()
const {
return b_;}
231 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
232 virtual bool write(std::ostream& os)
const;
234 static inline const char* classname() {
return "npstat::Quadratic1D";}
235 static inline unsigned version() {
return 3;}
236 static Quadratic1D* read(
const gs::ClassId&
id, std::istream& in);
245 const std::vector<double>& params);
246 inline static int nParameters() {
return 2;}
248 void verifyNonNegative();
249 double unscaledDensity(
double x)
const;
250 double unscaledCdf(
double x)
const;
251 double unscaledQuantile(
double x)
const;
253 inline double unscaledExceedance(
const double x)
const
254 {
return 1.0 - unscaledCdf(x);}
275 inline double a()
const {
return a_;}
276 inline double b()
const {
return b_;}
279 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
280 virtual bool write(std::ostream& os)
const;
282 static inline const char* classname() {
return "npstat::LogQuadratic1D";}
283 static inline unsigned version() {
return 3;}
284 static LogQuadratic1D* read(
const gs::ClassId&
id, std::istream& in);
293 const std::vector<double>& params);
294 inline static int nParameters() {
return 2;}
297 long double quadInteg(
long double x)
const;
298 double unscaledDensity(
double x)
const;
299 double unscaledCdf(
double x)
const;
300 double unscaledQuantile(
double x)
const;
302 inline double unscaledExceedance(
const double x)
const
303 {
return 1.0 - unscaledCdf(x);}
324 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
325 virtual bool write(std::ostream& os)
const;
327 static inline const char* classname() {
return "npstat::Gauss1D";}
328 static inline unsigned version() {
return 1;}
329 static Gauss1D* read(
const gs::ClassId&
id, std::istream& in);
339 const std::vector<double>& params);
340 inline static int nParameters() {
return 0;}
342 double unscaledDensity(
double x)
const;
343 double unscaledCdf(
double x)
const;
344 double unscaledQuantile(
double x)
const;
345 double unscaledExceedance(
double x)
const;
346 unsigned unscaledRandom(AbsRandomGenerator& g,
347 double* generatedRandom)
const;
362 inline double nsigma()
const {
return nsigma_;}
365 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
366 virtual bool write(std::ostream& os)
const;
368 static inline const char* classname() {
return "npstat::TruncatedGauss1D";}
369 static inline unsigned version() {
return 1;}
370 static TruncatedGauss1D* read(
const gs::ClassId&
id, std::istream& in);
379 const std::vector<double>& params);
380 inline static int nParameters() {
return 1;}
383 double unscaledDensity(
double x)
const;
384 double unscaledCdf(
double x)
const;
385 double unscaledQuantile(
double x)
const;
387 double* generatedRandom)
const;
389 inline double unscaledExceedance(
const double x)
const
390 {
return unscaledCdf(-x);}
406 double meanOn0_1,
double sigmaOn0_1);
413 inline double meanOn0_1()
const {
return mu0_;}
414 inline double sigmaOn0_1()
const {
return sigma0_;}
417 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
418 virtual bool write(std::ostream& os)
const;
420 static inline const char* classname() {
return "npstat::MirroredGauss1D";}
421 static inline unsigned version() {
return 1;}
422 static MirroredGauss1D* read(
const gs::ClassId&
id, std::istream& in);
431 const std::vector<double>& params);
432 inline static int nParameters() {
return 2;}
435 double unscaledDensity(
double x)
const;
436 double unscaledCdf(
double x)
const;
437 double unscaledQuantile(
double x)
const;
438 double unscaledExceedance(
double x)
const;
439 long double ldCdf(
double x)
const;
457 double leftSigmaFraction,
458 double nSigmasLeft,
double nSigmasRight);
464 inline double leftSigmaFraction()
const
465 {
return leftSigma_/(leftSigma_ + rightSigma_);}
466 inline double nSigmasLeft()
const
467 {
return nSigmasLeft_;}
468 inline double nSigmasRight()
const
469 {
return nSigmasRight_;}
472 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
473 virtual bool write(std::ostream& os)
const;
475 static inline const char* classname() {
return "npstat::BifurcatedGauss1D";}
476 static inline unsigned version() {
return 1;}
477 static BifurcatedGauss1D* read(
const gs::ClassId&
id, std::istream& in);
488 const std::vector<double>& params);
489 inline static int nParameters() {
return 3;}
492 double unscaledDensity(
double x)
const;
493 double unscaledCdf(
double x)
const;
494 double unscaledQuantile(
double x)
const;
495 double unscaledExceedance(
double x)
const;
500 double nSigmasRight_;
518 inline double power()
const {
return n_;}
521 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
522 virtual bool write(std::ostream& os)
const;
524 static inline const char* classname() {
return "npstat::SymmetricBeta1D";}
525 static inline unsigned version() {
return 1;}
526 static SymmetricBeta1D* read(
const gs::ClassId&
id, std::istream& in);
535 const std::vector<double>& params);
536 inline static int nParameters() {
return 1;}
538 double unscaledDensity(
double x)
const;
539 double unscaledCdf(
double x)
const;
540 double unscaledQuantile(
double x)
const;
542 inline double unscaledExceedance(
const double x)
const
543 {
return unscaledCdf(-x);}
545 double calculateNorm();
558 {
return new Beta1D(*
this);}
560 inline virtual ~
Beta1D() {}
562 inline double alpha()
const {
return alpha_;}
563 inline double beta()
const {
return beta_;}
566 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
567 virtual bool write(std::ostream& os)
const;
569 static inline const char* classname() {
return "npstat::Beta1D";}
570 static inline unsigned version() {
return 1;}
571 static Beta1D* read(
const gs::ClassId&
id, std::istream& in);
580 const std::vector<double>& params);
581 inline static int nParameters() {
return 2;}
583 double unscaledDensity(
double x)
const;
584 double unscaledCdf(
double x)
const;
585 double unscaledExceedance(
double x)
const;
586 double unscaledQuantile(
double x)
const;
588 double calculateNorm()
const;
605 inline double alpha()
const {
return alpha_;}
608 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
609 virtual bool write(std::ostream& os)
const;
611 static inline const char* classname() {
return "npstat::Gamma1D";}
612 static inline unsigned version() {
return 1;}
613 static Gamma1D* read(
const gs::ClassId&
id, std::istream& in);
622 const std::vector<double>& params);
623 inline static int nParameters() {
return 1;}
626 double* generatedRandom);
628 double unscaledDensity(
double x)
const;
629 double unscaledCdf(
double x)
const;
630 double unscaledExceedance(
double x)
const;
631 double unscaledQuantile(
double x)
const;
633 double* generatedRandom)
const;
654 inline double powerParameter()
const {
return c_;}
657 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
658 virtual bool write(std::ostream& os)
const;
660 static inline const char* classname() {
return "npstat::Pareto1D";}
661 static inline unsigned version() {
return 1;}
662 static Pareto1D* read(
const gs::ClassId&
id, std::istream& in);
671 const std::vector<double>& params);
672 inline static int nParameters() {
return 1;}
675 double unscaledDensity(
double x)
const;
676 double unscaledCdf(
double x)
const;
677 double unscaledQuantile(
double x)
const;
678 double unscaledExceedance(
double x)
const;
698 inline double powerParameter()
const {
return c_;}
701 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
702 virtual bool write(std::ostream& os)
const;
704 static inline const char* classname() {
return "npstat::UniPareto1D";}
705 static inline unsigned version() {
return 1;}
706 static UniPareto1D* read(
const gs::ClassId&
id, std::istream& in);
715 const std::vector<double>& params);
716 inline static int nParameters() {
return 1;}
719 double unscaledDensity(
double x)
const;
720 double unscaledCdf(
double x)
const;
721 double unscaledQuantile(
double x)
const;
722 double unscaledExceedance(
double x)
const;
738 inline double tailWeight()
const {
return tailWeight_;}
739 inline double tailStart()
const {
return a_;}
742 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
743 virtual bool write(std::ostream& os)
const;
745 static inline const char* classname() {
return "npstat::Huber1D";}
746 static inline unsigned version() {
return 1;}
747 static Huber1D* read(
const gs::ClassId&
id, std::istream& in);
756 const std::vector<double>& params);
757 inline static int nParameters() {
return 1;}
760 double unscaledDensity(
double x)
const;
761 double unscaledCdf(
double x)
const;
762 double unscaledQuantile(
double x)
const;
764 inline double unscaledExceedance(
const double x)
const
765 {
return unscaledCdf(-x);}
767 double weight(
double a)
const;
786 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
787 virtual bool write(std::ostream& os)
const;
789 static inline const char* classname() {
return "npstat::Cauchy1D";}
790 static inline unsigned version() {
return 1;}
791 static Cauchy1D* read(
const gs::ClassId&
id, std::istream& in);
801 const std::vector<double>& params);
802 inline static int nParameters() {
return 0;}
804 double unscaledDensity(
double x)
const;
805 double unscaledCdf(
double x)
const;
806 double unscaledQuantile(
double x)
const;
808 inline double unscaledExceedance(
const double x)
const
809 {
return unscaledCdf(-x);}
822 LogNormal(
double mean,
double stdev,
double skewness);
827 inline double skewness()
const {
return skew_;}
829 bool isGaussian()
const;
830 double getDelta()
const;
831 double getGamma()
const;
832 double getXi()
const;
835 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
836 virtual bool write(std::ostream& os)
const;
838 static inline const char* classname() {
return "npstat::LogNormal";}
839 static inline unsigned version() {
return 1;}
840 static LogNormal* read(
const gs::ClassId&
id, std::istream& in);
849 const std::vector<double>& params);
850 inline static int nParameters() {
return 1;}
853 double unscaledDensity(
double x)
const;
854 double unscaledCdf(
double x)
const;
855 double unscaledExceedance(
double x)
const;
856 double unscaledQuantile(
double x)
const;
880 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
881 virtual bool write(std::ostream& os)
const;
883 static inline const char* classname() {
return "npstat::Moyal1D";}
884 static inline unsigned version() {
return 1;}
885 static Moyal1D* read(
const gs::ClassId&
id, std::istream& in);
895 const std::vector<double>& params);
897 inline static int nParameters() {
return 0;}
899 double unscaledDensity(
double x)
const;
900 double unscaledCdf(
double x)
const;
901 double unscaledQuantile(
double x)
const;
902 double unscaledExceedance(
double x)
const;
918 inline double nDegreesOfFreedom()
const {
return nDoF_;}
921 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
922 virtual bool write(std::ostream& os)
const;
924 static inline const char* classname() {
return "npstat::StudentsT1D";}
925 static inline unsigned version() {
return 1;}
926 static StudentsT1D* read(
const gs::ClassId&
id, std::istream& in);
935 const std::vector<double>& params);
936 inline static int nParameters() {
return 1;}
939 double effectiveSupport()
const;
940 double unscaledDensity(
double x)
const;
941 double unscaledCdf(
double x)
const;
942 double unscaledExceedance(
double x)
const;
943 double unscaledQuantile(
double x)
const;
963 template <
typename Real>
965 const Real* data,
unsigned dataLen,
966 unsigned interpolationDegree);
969 const std::vector<double>& table,
970 const unsigned interpolationDegree)
973 const unsigned long sz = table.size();
974 initialize(sz ? &table[0] : (
double*)0, sz, interpolationDegree);
982 inline unsigned interpolationDegree()
const {
return deg_;}
983 inline unsigned tableLength()
const {
return len_;}
984 inline const double* tableData()
const {
return &table_[0];}
987 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
988 virtual bool write(std::ostream& os)
const;
990 static inline const char* classname() {
return "npstat::Tabulated1D";}
991 static inline unsigned version() {
return 1;}
992 static Tabulated1D* read(
const gs::ClassId&
id, std::istream& in);
1003 const std::vector<double>& params);
1004 inline static int nParameters() {
return -1;}
1006 double unscaledDensity(
double x)
const;
1007 double unscaledCdf(
double x)
const;
1008 double unscaledQuantile(
double x)
const;
1009 double unscaledExceedance(
double x)
const;
1011 template <
typename Real>
void initialize(
1012 const Real* data,
unsigned dataLen,
unsigned interpolationDegree);
1015 double interpolate(
double x)
const;
1016 double intervalInteg(
unsigned intervalNumber)
const;
1017 double interpIntegral(
double x0,
double x1)
const;
1019 std::vector<double> table_;
1020 std::vector<double> cdf_;
1021 std::vector<double> exceed_;
1040 template <
typename Real>
1042 const Real* data,
unsigned dataLen,
1043 unsigned interpolationDegree);
1046 const std::vector<double>& table,
1047 const unsigned interpolationDegree)
1050 const unsigned long sz = table.size();
1051 initialize(sz ? &table[0] : (
double*)0, sz, interpolationDegree);
1059 inline unsigned interpolationDegree()
const {
return deg_;}
1060 inline unsigned tableLength()
const {
return len_;}
1061 inline const double* tableData()
const {
return &table_[0];}
1064 virtual gs::ClassId
classId()
const {
return gs::ClassId(*
this);}
1065 virtual bool write(std::ostream& os)
const;
1067 static inline const char* classname() {
return "npstat::BinnedDensity1D";}
1068 static inline unsigned version() {
return 1;}
1069 static BinnedDensity1D* read(
const gs::ClassId&
id, std::istream& in);
1080 const std::vector<double>& params);
1081 inline static int nParameters() {
return -1;}
1083 double unscaledDensity(
double x)
const;
1084 double unscaledCdf(
double x)
const;
1085 double unscaledQuantile(
double x)
const;
1087 inline double unscaledExceedance(
const double x)
const
1088 {
return 1.0 - unscaledCdf(x);}
1090 template <
typename Real>
void initialize(
1091 const Real* data,
unsigned dataLen,
unsigned interpolationDegree);
1094 double interpolate(
double x)
const;
1096 std::vector<double> table_;
1097 std::vector<double> cdf_;
1101 unsigned firstNonZeroBin_;
1102 unsigned lastNonZeroBin_;
1106 #include "npstat/stat/Distributions1D.icc"
Factories for 1-d distributions for use in interpretive language environments.
Definition: AbsDistribution1D.hh:165
double scale() const
Definition: AbsDistribution1D.hh:183
virtual bool isEqual(const AbsDistribution1D &other) const
Definition: AbsDistribution1D.hh:242
double location() const
Definition: AbsDistribution1D.hh:180
Definition: Distributions1D.hh:554
virtual bool isEqual(const AbsDistribution1D &) const
virtual Beta1D * clone() const
Definition: Distributions1D.hh:557
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:566
Definition: Distributions1D.hh:454
virtual bool isEqual(const AbsDistribution1D &) const
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:472
virtual BifurcatedGauss1D * clone() const
Definition: Distributions1D.hh:459
Definition: Distributions1D.hh:1032
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:1064
virtual bool isEqual(const AbsDistribution1D &) const
virtual BinnedDensity1D * clone() const
Definition: Distributions1D.hh:1054
Definition: Distributions1D.hh:778
virtual bool isEqual(const AbsDistribution1D &r) const
Definition: Distributions1D.hh:794
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:786
virtual Cauchy1D * clone() const
Definition: Distributions1D.hh:781
Definition: Distributions1D.hh:100
virtual bool isEqual(const AbsDistribution1D &r) const
Definition: Distributions1D.hh:118
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:110
virtual Exponential1D * clone() const
Definition: Distributions1D.hh:104
Definition: Distributions1D.hh:597
virtual bool isEqual(const AbsDistribution1D &) const
virtual Gamma1D * clone() const
Definition: Distributions1D.hh:600
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:608
Definition: Distributions1D.hh:316
virtual Gauss1D * clone() const
Definition: Distributions1D.hh:319
virtual bool isEqual(const AbsDistribution1D &r) const
Definition: Distributions1D.hh:332
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:324
Definition: Distributions1D.hh:731
virtual bool isEqual(const AbsDistribution1D &) const
virtual Huber1D * clone() const
Definition: Distributions1D.hh:734
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:742
Definition: Distributions1D.hh:60
virtual bool isEqual(const AbsDistribution1D &r) const
Definition: Distributions1D.hh:79
virtual IsoscelesTriangle1D * clone() const
Definition: Distributions1D.hh:64
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:70
Definition: Distributions1D.hh:137
virtual bool isEqual(const AbsDistribution1D &r) const
Definition: Distributions1D.hh:155
virtual Laplace1D * clone() const
Definition: Distributions1D.hh:141
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:147
Definition: Distributions1D.hh:820
virtual LogNormal * clone() const
Definition: Distributions1D.hh:823
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:835
virtual bool isEqual(const AbsDistribution1D &) const
Definition: Distributions1D.hh:267
virtual bool isEqual(const AbsDistribution1D &) const
virtual LogQuadratic1D * clone() const
Definition: Distributions1D.hh:270
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:279
Definition: Distributions1D.hh:174
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:184
virtual bool isEqual(const AbsDistribution1D &r) const
Definition: Distributions1D.hh:192
virtual Logistic1D * clone() const
Definition: Distributions1D.hh:178
Definition: Distributions1D.hh:403
virtual bool isEqual(const AbsDistribution1D &) const
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:417
virtual MirroredGauss1D * clone() const
Definition: Distributions1D.hh:408
Definition: Distributions1D.hh:872
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:880
virtual bool isEqual(const AbsDistribution1D &r) const
Definition: Distributions1D.hh:888
virtual Moyal1D * clone() const
Definition: Distributions1D.hh:875
Definition: Distributions1D.hh:647
virtual bool isEqual(const AbsDistribution1D &) const
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:657
virtual Pareto1D * clone() const
Definition: Distributions1D.hh:650
Definition: Distributions1D.hh:219
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:231
virtual Quadratic1D * clone() const
Definition: Distributions1D.hh:222
virtual bool isEqual(const AbsDistribution1D &) const
Definition: Distribution1DFactory.hh:35
Definition: Distributions1D.hh:910
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:921
virtual bool isEqual(const AbsDistribution1D &) const
virtual StudentsT1D * clone() const
Definition: Distributions1D.hh:913
Definition: Distributions1D.hh:509
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:521
virtual bool isEqual(const AbsDistribution1D &) const
virtual SymmetricBeta1D * clone() const
Definition: Distributions1D.hh:513
Definition: Distributions1D.hh:953
virtual Tabulated1D * clone() const
Definition: Distributions1D.hh:977
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:987
virtual bool isEqual(const AbsDistribution1D &) const
Definition: Distributions1D.hh:354
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:365
virtual TruncatedGauss1D * clone() const
Definition: Distributions1D.hh:357
virtual bool isEqual(const AbsDistribution1D &) const
Definition: Distributions1D.hh:691
virtual UniPareto1D * clone() const
Definition: Distributions1D.hh:694
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:701
virtual bool isEqual(const AbsDistribution1D &) const
Definition: AbsArrayProjector.hh:14
Definition: AbsDistribution1D.hh:31
Definition: AbsRandomGenerator.hh:27