npstat is hosted by Hepforge, IPPP Durham
NPStat  5.10.0
Distributions1D.hh
Go to the documentation of this file.
1 #ifndef NPSTAT_DISTRIBUTIONS1D_HH_
2 #define NPSTAT_DISTRIBUTIONS1D_HH_
3 
4 /*!
5 // \file Distributions1D.hh
6 //
7 // \brief A number of useful 1-d continuous statistical distributions
8 //
9 // Author: I. Volobouev
10 //
11 // November 2009
12 */
13 
15 
16 namespace npstat {
17  /**
18  // The uniform distribution is defined here by a constant density
19  // equal to 1 between 0 and 1 and equal to 0 everywhere else
20  */
22  {
23  public:
24  inline Uniform1D(const double location, const double scale)
26  inline virtual Uniform1D* clone() const {return new Uniform1D(*this);}
27 
28  inline virtual ~Uniform1D() {}
29 
30  // Methods needed for I/O
31  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
32  virtual bool write(std::ostream& os) const;
33 
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);
37 
38  protected:
39  virtual bool isEqual(const AbsDistribution1D& r) const
41 
42  private:
44 
45  inline Uniform1D(const double location, const double scale,
46  const std::vector<double>& /* params */)
48  inline static int nParameters() {return 0;}
49 
50  double unscaledDensity(double x) const;
51  double unscaledCdf(double x) const;
52  double unscaledQuantile(double x) const;
53 
54  inline double unscaledExceedance(const double x) const
55  {return 1.0 - unscaledCdf(x);}
56  };
57 
58  /** Isosceles triangle distribution: 1 - |x| supported on [-1, 1] */
60  {
61  public:
62  inline IsoscelesTriangle1D(const double location, const double scale)
64  inline virtual IsoscelesTriangle1D* clone() const
65  {return new IsoscelesTriangle1D(*this);}
66 
67  inline virtual ~IsoscelesTriangle1D() {}
68 
69  // Methods needed for I/O
70  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
71  virtual bool write(std::ostream& os) const;
72 
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);
77 
78  protected:
79  virtual bool isEqual(const AbsDistribution1D& r) const
81 
82  private:
84 
85  inline IsoscelesTriangle1D(const double location, const double scale,
86  const std::vector<double>& /* params */)
88  inline static int nParameters() {return 0;}
89 
90  double unscaledDensity(double x) const;
91  double unscaledCdf(double x) const;
92  double unscaledQuantile(double x) const;
93 
94  inline double unscaledExceedance(const double x) const
95  {return 1.0 - unscaledCdf(x);}
96  };
97 
98  /** Exponential distribution. "scale" is the decay time. */
100  {
101  public:
102  inline Exponential1D(const double location, const double scale)
104  inline virtual Exponential1D* clone() const
105  {return new Exponential1D(*this);}
106 
107  inline virtual ~Exponential1D() {}
108 
109  // Methods needed for I/O
110  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
111  virtual bool write(std::ostream& os) const;
112 
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);
116 
117  protected:
118  virtual bool isEqual(const AbsDistribution1D& r) const
120 
121  private:
123 
124  inline Exponential1D(const double location, const double scale,
125  const std::vector<double>& /* params */)
127  inline static int nParameters() {return 0;}
128 
129  double unscaledDensity(double x) const;
130  double unscaledCdf(double x) const;
131  double unscaledQuantile(double x) const;
132  double unscaledExceedance(double x) const;
133  };
134 
135  /** Laplace distribution */
137  {
138  public:
139  inline Laplace1D(const double location, const double scale)
141  inline virtual Laplace1D* clone() const
142  {return new Laplace1D(*this);}
143 
144  inline virtual ~Laplace1D() {}
145 
146  // Methods needed for I/O
147  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
148  virtual bool write(std::ostream& os) const;
149 
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);
153 
154  protected:
155  virtual bool isEqual(const AbsDistribution1D& r) const
157 
158  private:
160 
161  inline Laplace1D(const double location, const double scale,
162  const std::vector<double>& /* params */)
164  inline static int nParameters() {return 0;}
165 
166  double unscaledDensity(double x) const;
167  double unscaledCdf(double x) const;
168  double unscaledQuantile(double x) const;
169  double unscaledExceedance(double x) const;
170  };
171 
172  /** Logistic distribution */
174  {
175  public:
176  inline Logistic1D(const double location, const double scale)
178  inline virtual Logistic1D* clone() const
179  {return new Logistic1D(*this);}
180 
181  inline virtual ~Logistic1D() {}
182 
183  // Methods needed for I/O
184  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
185  virtual bool write(std::ostream& os) const;
186 
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);
190 
191  protected:
192  virtual bool isEqual(const AbsDistribution1D& r) const
194 
195  private:
197 
198  inline Logistic1D(const double location, const double scale,
199  const std::vector<double>& /* params */)
201  inline static int nParameters() {return 0;}
202 
203  double unscaledDensity(double x) const;
204  double unscaledCdf(double x) const;
205  double unscaledQuantile(double x) const;
206  double unscaledExceedance(double x) const;
207  };
208 
209  /**
210  // A distribution whose density has a simple quadratic shape.
211  // The support is from 0 to 1, and the coefficients "a" and "b"
212  // are the coefficients for the Legendre polynomials of 1st
213  // and 2nd degree translated to the support region. Note that
214  // only those values of "a" and "b" that guarantee non-negativity
215  // of the density are allowed, otherwise the code will generate
216  // a run-time error.
217  */
219  {
220  public:
221  Quadratic1D(double location, double scale, double a, double b);
222  inline virtual Quadratic1D* clone() const
223  {return new Quadratic1D(*this);}
224 
225  inline virtual ~Quadratic1D() {}
226 
227  inline double a() const {return a_;}
228  inline double b() const {return b_;}
229 
230  // Methods needed for I/O
231  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
232  virtual bool write(std::ostream& os) const;
233 
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);
237 
238  protected:
239  virtual bool isEqual(const AbsDistribution1D&) const;
240 
241  private:
243 
244  Quadratic1D(double location, double scale,
245  const std::vector<double>& params);
246  inline static int nParameters() {return 2;}
247 
248  void verifyNonNegative();
249  double unscaledDensity(double x) const;
250  double unscaledCdf(double x) const;
251  double unscaledQuantile(double x) const;
252 
253  inline double unscaledExceedance(const double x) const
254  {return 1.0 - unscaledCdf(x);}
255 
256  double a_;
257  double b_;
258  };
259 
260  /**
261  // A distribution whose density logarithm has a simple quadratic
262  // shape. The support is from 0 to 1, and the coefficients "a" and "b"
263  // are the coefficients for the Legendre polynomials of 1st and 2nd
264  // degree translated to the support region.
265  */
267  {
268  public:
269  LogQuadratic1D(double location, double scale, double a, double b);
270  inline virtual LogQuadratic1D* clone() const
271  {return new LogQuadratic1D(*this);}
272 
273  inline virtual ~LogQuadratic1D() {}
274 
275  inline double a() const {return a_;}
276  inline double b() const {return b_;}
277 
278  // Methods needed for I/O
279  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
280  virtual bool write(std::ostream& os) const;
281 
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);
285 
286  protected:
287  virtual bool isEqual(const AbsDistribution1D&) const;
288 
289  private:
291 
292  LogQuadratic1D(double location, double scale,
293  const std::vector<double>& params);
294  inline static int nParameters() {return 2;}
295 
296  void normalize();
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;
301 
302  inline double unscaledExceedance(const double x) const
303  {return 1.0 - unscaledCdf(x);}
304 
305  long double ref_;
306  long double range_;
307  double a_;
308  double b_;
309  double k_;
310  double s_;
311  double norm_;
312  };
313 
314  /** The Gaussian (or Normal) distribution */
316  {
317  public:
318  Gauss1D(double location, double scale);
319  inline virtual Gauss1D* clone() const {return new Gauss1D(*this);}
320 
321  inline virtual ~Gauss1D() {}
322 
323  // Methods needed for I/O
324  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
325  virtual bool write(std::ostream& os) const;
326 
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);
330 
331  protected:
332  virtual bool isEqual(const AbsDistribution1D& r) const
334 
335  private:
337 
338  Gauss1D(const double location, const double scale,
339  const std::vector<double>& params);
340  inline static int nParameters() {return 0;}
341 
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;
348  double xmin_;
349  double xmax_;
350  };
351 
352  /** Gaussian distribution truncated at some number of sigmas */
354  {
355  public:
356  TruncatedGauss1D(double location, double scale, double nsigma);
357  inline virtual TruncatedGauss1D* clone() const
358  {return new TruncatedGauss1D(*this);}
359 
360  inline virtual ~TruncatedGauss1D() {}
361 
362  inline double nsigma() const {return nsigma_;}
363 
364  // Methods needed for I/O
365  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
366  virtual bool write(std::ostream& os) const;
367 
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);
371 
372  protected:
373  virtual bool isEqual(const AbsDistribution1D&) const;
374 
375  private:
377 
378  TruncatedGauss1D(double location, double scale,
379  const std::vector<double>& params);
380  inline static int nParameters() {return 1;}
381 
382  void initialize();
383  double unscaledDensity(double x) const;
384  double unscaledCdf(double x) const;
385  double unscaledQuantile(double x) const;
386  unsigned unscaledRandom(AbsRandomGenerator& g,
387  double* generatedRandom) const;
388 
389  inline double unscaledExceedance(const double x) const
390  {return unscaledCdf(-x);}
391 
392  double nsigma_;
393  double norm_;
394  double cdf0_;
395  };
396 
397  /**
398  // Gaussian distribution on the [0, 1] interval, mirrored at the boundaries.
399  // This is the Green's function of the diffusion equation on [0, 1]. The
400  // interval can be shifted and scaled as for the uniform distribution.
401  */
403  {
404  public:
405  MirroredGauss1D(double location, double scale,
406  double meanOn0_1, double sigmaOn0_1);
407 
408  inline virtual MirroredGauss1D* clone() const
409  {return new MirroredGauss1D(*this);}
410 
411  inline virtual ~MirroredGauss1D() {}
412 
413  inline double meanOn0_1() const {return mu0_;}
414  inline double sigmaOn0_1() const {return sigma0_;}
415 
416  // Methods needed for I/O
417  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
418  virtual bool write(std::ostream& os) const;
419 
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);
423 
424  protected:
425  virtual bool isEqual(const AbsDistribution1D&) const;
426 
427  private:
429 
430  MirroredGauss1D(double location, double scale,
431  const std::vector<double>& params);
432  inline static int nParameters() {return 2;}
433 
434  void validate();
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;
440 
441  double mu0_;
442  double sigma0_;
443  };
444 
445  /**
446  // Bifurcated Gaussian distribution. Different sigmas
447  // can be used on the left and on the right, with constructor
448  // parameter "leftSigmaFraction" specifying the ratio of
449  // the left sigma to the sum of sigmas (this ratio must be
450  // between 0 and 1). Different truncations in terms of the
451  // number of sigmas can be used as well.
452  */
454  {
455  public:
456  BifurcatedGauss1D(double location, double scale,
457  double leftSigmaFraction,
458  double nSigmasLeft, double nSigmasRight);
459  inline virtual BifurcatedGauss1D* clone() const
460  {return new BifurcatedGauss1D(*this);}
461 
462  inline virtual ~BifurcatedGauss1D() {}
463 
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_;}
470 
471  // Methods needed for I/O
472  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
473  virtual bool write(std::ostream& os) const;
474 
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);
478 
479  protected:
480  virtual bool isEqual(const AbsDistribution1D&) const;
481 
482  private:
483  BifurcatedGauss1D(double location, double scale);
484 
486 
487  BifurcatedGauss1D(double location, double scale,
488  const std::vector<double>& params);
489  inline static int nParameters() {return 3;}
490 
491  void initialize();
492  double unscaledDensity(double x) const;
493  double unscaledCdf(double x) const;
494  double unscaledQuantile(double x) const;
495  double unscaledExceedance(double x) const;
496 
497  double leftSigma_;
498  double rightSigma_;
499  double nSigmasLeft_;
500  double nSigmasRight_;
501  double norm_;
502  double leftCdfFrac_;
503  double cdf0Left_;
504  double cdf0Right_;
505  };
506 
507  /** Symmetric beta distribution */
509  {
510  public:
511  SymmetricBeta1D(double location, double scale, double power);
512 
513  inline virtual SymmetricBeta1D* clone() const
514  {return new SymmetricBeta1D(*this);}
515 
516  inline virtual ~SymmetricBeta1D() {}
517 
518  inline double power() const {return n_;}
519 
520  // Methods needed for I/O
521  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
522  virtual bool write(std::ostream& os) const;
523 
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);
527 
528  protected:
529  virtual bool isEqual(const AbsDistribution1D&) const;
530 
531  private:
533 
534  SymmetricBeta1D(double location, double scale,
535  const std::vector<double>& params);
536  inline static int nParameters() {return 1;}
537 
538  double unscaledDensity(double x) const;
539  double unscaledCdf(double x) const;
540  double unscaledQuantile(double x) const;
541 
542  inline double unscaledExceedance(const double x) const
543  {return unscaledCdf(-x);}
544 
545  double calculateNorm();
546 
547  double n_;
548  double norm_;
549  int intpow_;
550  };
551 
552  /** Beta1D density is proportional to x^(apha-1) * (1-x)^(beta-1) */
554  {
555  public:
556  Beta1D(double location, double scale, double alpha, double beta);
557  inline virtual Beta1D* clone() const
558  {return new Beta1D(*this);}
559 
560  inline virtual ~Beta1D() {}
561 
562  inline double alpha() const {return alpha_;}
563  inline double beta() const {return beta_;}
564 
565  // Methods needed for I/O
566  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
567  virtual bool write(std::ostream& os) const;
568 
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);
572 
573  protected:
574  virtual bool isEqual(const AbsDistribution1D&) const;
575 
576  private:
578 
579  Beta1D(double location, double scale,
580  const std::vector<double>& params);
581  inline static int nParameters() {return 2;}
582 
583  double unscaledDensity(double x) const;
584  double unscaledCdf(double x) const;
585  double unscaledExceedance(double x) const;
586  double unscaledQuantile(double x) const;
587 
588  double calculateNorm() const;
589 
590  double alpha_;
591  double beta_;
592  double norm_;
593  };
594 
595  /** Shiftable gamma distribution */
597  {
598  public:
599  Gamma1D(double location, double scale, double alpha);
600  inline virtual Gamma1D* clone() const
601  {return new Gamma1D(*this);}
602 
603  inline virtual ~Gamma1D() {}
604 
605  inline double alpha() const {return alpha_;}
606 
607  // Methods needed for I/O
608  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
609  virtual bool write(std::ostream& os) const;
610 
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);
614 
615  protected:
616  virtual bool isEqual(const AbsDistribution1D&) const;
617 
618  private:
620 
621  Gamma1D(double location, double scale,
622  const std::vector<double>& params);
623  inline static int nParameters() {return 1;}
624 
625  static unsigned unscaledRng(double alpha, AbsRandomGenerator& g,
626  double* generatedRandom);
627 
628  double unscaledDensity(double x) const;
629  double unscaledCdf(double x) const;
630  double unscaledExceedance(double x) const;
631  double unscaledQuantile(double x) const;
632  unsigned unscaledRandom(AbsRandomGenerator& g,
633  double* generatedRandom) const;
634  void initialize();
635 
636  double alpha_;
637  double norm_;
638  double uplim_;
639  };
640 
641  /**
642  // Pareto distribution. Location parameter is location of 0, scale
643  // parameter is the distance between 0 and the start of the density
644  // (like the normal Pareto distribution location parameter).
645  */
647  {
648  public:
649  Pareto1D(double location, double scale, double powerParameter);
650  inline virtual Pareto1D* clone() const {return new Pareto1D(*this);}
651 
652  inline virtual ~Pareto1D() {}
653 
654  inline double powerParameter() const {return c_;}
655 
656  // Methods needed for I/O
657  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
658  virtual bool write(std::ostream& os) const;
659 
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);
663 
664  protected:
665  virtual bool isEqual(const AbsDistribution1D&) const;
666 
667  private:
669 
670  Pareto1D(double location, double scale,
671  const std::vector<double>& params);
672  inline static int nParameters() {return 1;}
673 
674  void initialize();
675  double unscaledDensity(double x) const;
676  double unscaledCdf(double x) const;
677  double unscaledQuantile(double x) const;
678  double unscaledExceedance(double x) const;
679 
680  double c_;
681  double support_;
682  };
683 
684  /**
685  // Uniform distribution with Pareto tail attached to the right, where
686  // the support of the uniform would normally end. Location parameter
687  // is location of 0, scale parameter is the width of the uniform part
688  // (like the normal Pareto distribution location parameter).
689  */
691  {
692  public:
693  UniPareto1D(double location, double scale, double powerParameter);
694  inline virtual UniPareto1D* clone() const {return new UniPareto1D(*this);}
695 
696  inline virtual ~UniPareto1D() {}
697 
698  inline double powerParameter() const {return c_;}
699 
700  // Methods needed for I/O
701  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
702  virtual bool write(std::ostream& os) const;
703 
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);
707 
708  protected:
709  virtual bool isEqual(const AbsDistribution1D&) const;
710 
711  private:
713 
714  UniPareto1D(double location, double scale,
715  const std::vector<double>& params);
716  inline static int nParameters() {return 1;}
717 
718  void initialize();
719  double unscaledDensity(double x) const;
720  double unscaledCdf(double x) const;
721  double unscaledQuantile(double x) const;
722  double unscaledExceedance(double x) const;
723 
724  double c_;
725  double support_;
726  double amplitude_;
727  };
728 
729  /** "Huber" distribution */
731  {
732  public:
733  Huber1D(double location, double scale, double tailWeight);
734  inline virtual Huber1D* clone() const {return new Huber1D(*this);}
735 
736  inline virtual ~Huber1D() {}
737 
738  inline double tailWeight() const {return tailWeight_;}
739  inline double tailStart() const {return a_;}
740 
741  // Methods needed for I/O
742  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
743  virtual bool write(std::ostream& os) const;
744 
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);
748 
749  protected:
750  virtual bool isEqual(const AbsDistribution1D&) const;
751 
752  private:
754 
755  Huber1D(double location, double scale,
756  const std::vector<double>& params);
757  inline static int nParameters() {return 1;}
758 
759  void initialize();
760  double unscaledDensity(double x) const;
761  double unscaledCdf(double x) const;
762  double unscaledQuantile(double x) const;
763 
764  inline double unscaledExceedance(const double x) const
765  {return unscaledCdf(-x);}
766 
767  double weight(double a) const;
768 
769  double tailWeight_;
770  double a_;
771  double normfactor_;
772  double support_;
773  double cdf0_;
774  };
775 
776  /** Cauchy (or Breit-Wigner) distribution */
778  {
779  public:
780  Cauchy1D(double location, double scale);
781  inline virtual Cauchy1D* clone() const {return new Cauchy1D(*this);}
782 
783  inline virtual ~Cauchy1D() {}
784 
785  // Methods needed for I/O
786  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
787  virtual bool write(std::ostream& os) const;
788 
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);
792 
793  protected:
794  virtual bool isEqual(const AbsDistribution1D& r) const
796 
797  private:
799 
800  Cauchy1D(const double location, const double scale,
801  const std::vector<double>& params);
802  inline static int nParameters() {return 0;}
803 
804  double unscaledDensity(double x) const;
805  double unscaledCdf(double x) const;
806  double unscaledQuantile(double x) const;
807 
808  inline double unscaledExceedance(const double x) const
809  {return unscaledCdf(-x);}
810 
811  double support_;
812  };
813 
814  /**
815  // Log-normal distribution represented by its mean, standard
816  // deviation, and skewness. This representation is more useful
817  // than other representations encountered in statistical literature.
818  */
820  {
821  public:
822  LogNormal(double mean, double stdev, double skewness);
823  inline virtual LogNormal* clone() const {return new LogNormal(*this);}
824 
825  inline virtual ~LogNormal() {}
826 
827  inline double skewness() const {return skew_;}
828 
829  bool isGaussian() const;
830  double getDelta() const;
831  double getGamma() const;
832  double getXi() const;
833 
834  // Methods needed for I/O
835  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
836  virtual bool write(std::ostream& os) const;
837 
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);
841 
842  protected:
843  virtual bool isEqual(const AbsDistribution1D&) const;
844 
845  private:
847 
848  LogNormal(double location, double scale,
849  const std::vector<double>& params);
850  inline static int nParameters() {return 1;}
851 
852  void initialize();
853  double unscaledDensity(double x) const;
854  double unscaledCdf(double x) const;
855  double unscaledExceedance(double x) const;
856  double unscaledQuantile(double x) const;
857 
858  double skew_;
859 
860  double w_;
861  double logw_;
862  double s_;
863  double xi_;
864  double emgamovd_;
865  };
866 
867  /**
868  // Moyal Distribution (originally derived by Moyal as an approximation
869  // to the Landau distribution)
870  */
872  {
873  public:
874  Moyal1D(double location, double scale);
875  inline virtual Moyal1D* clone() const {return new Moyal1D(*this);}
876 
877  inline virtual ~Moyal1D() {}
878 
879  // Methods needed for I/O
880  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
881  virtual bool write(std::ostream& os) const;
882 
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);
886 
887  protected:
888  virtual bool isEqual(const AbsDistribution1D& r) const
890 
891  private:
893 
894  Moyal1D(double location, double scale,
895  const std::vector<double>& params);
896 
897  inline static int nParameters() {return 0;}
898 
899  double unscaledDensity(double x) const;
900  double unscaledCdf(double x) const;
901  double unscaledQuantile(double x) const;
902  double unscaledExceedance(double x) const;
903 
904  double xmax_;
905  double xmin_;
906  };
907 
908  /** Student's t-distribution */
910  {
911  public:
912  StudentsT1D(double location, double scale, double nDegreesOfFreedom);
913  inline virtual StudentsT1D* clone() const
914  {return new StudentsT1D(*this);}
915 
916  inline virtual ~StudentsT1D() {}
917 
918  inline double nDegreesOfFreedom() const {return nDoF_;}
919 
920  // Methods needed for I/O
921  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
922  virtual bool write(std::ostream& os) const;
923 
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);
927 
928  protected:
929  virtual bool isEqual(const AbsDistribution1D&) const;
930 
931  private:
933 
934  StudentsT1D(double location, double scale,
935  const std::vector<double>& params);
936  inline static int nParameters() {return 1;}
937 
938  void initialize();
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;
944 
945  double nDoF_;
946  double normfactor_;
947  double power_;
948  double bignum_;
949  };
950 
951  /** Distribution defined by an interpolation table */
953  {
954  public:
955  // The "data" array gives (unnormalized) density values at
956  // equidistant intervals. data[0] is density at 0.0, and
957  // data[dataLen-1] is density at 1.0. If "dataLen" is less
958  // than 2, uniform distribution will be created. Internally,
959  // the data is kept in double precision.
960  //
961  // "interpolationDegree" must be less than 4 and less than "dataLen".
962  //
963  template <typename Real>
964  Tabulated1D(double location, double scale,
965  const Real* data, unsigned dataLen,
966  unsigned interpolationDegree);
967 
968  inline Tabulated1D(const double location, const double scale,
969  const std::vector<double>& table,
970  const unsigned interpolationDegree)
972  {
973  const unsigned long sz = table.size();
974  initialize(sz ? &table[0] : (double*)0, sz, interpolationDegree);
975  }
976 
977  inline virtual Tabulated1D* clone() const
978  {return new Tabulated1D(*this);}
979 
980  inline virtual ~Tabulated1D() {}
981 
982  inline unsigned interpolationDegree() const {return deg_;}
983  inline unsigned tableLength() const {return len_;}
984  inline const double* tableData() const {return &table_[0];}
985 
986  // Methods needed for I/O
987  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
988  virtual bool write(std::ostream& os) const;
989 
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);
993 
994  protected:
995  virtual bool isEqual(const AbsDistribution1D&) const;
996 
997  private:
999 
1000  // The following constructor creates interpolator
1001  // of maximum degree possible
1002  Tabulated1D(double location, double scale,
1003  const std::vector<double>& params);
1004  inline static int nParameters() {return -1;}
1005 
1006  double unscaledDensity(double x) const;
1007  double unscaledCdf(double x) const;
1008  double unscaledQuantile(double x) const;
1009  double unscaledExceedance(double x) const;
1010 
1011  template <typename Real> void initialize(
1012  const Real* data, unsigned dataLen, unsigned interpolationDegree);
1013 
1014  void normalize();
1015  double interpolate(double x) const;
1016  double intervalInteg(unsigned intervalNumber) const;
1017  double interpIntegral(double x0, double x1) const;
1018 
1019  std::vector<double> table_;
1020  std::vector<double> cdf_;
1021  std::vector<double> exceed_;
1022  double step_;
1023  unsigned len_;
1024  unsigned deg_;
1025  };
1026 
1027  /**
1028  // Another interpolated distribution. For this one, we will assume
1029  // that the coordinates correspond to 1-d histogram bin centers.
1030  */
1032  {
1033  public:
1034  // The "data" array gives density values at equidistant intervals.
1035  // data[0] is density at 0.5/dataLen, and data[dataLen-1] is density
1036  // at 1.0 - 0.5/dataLen.
1037  //
1038  // "interpolationDegree" must be less than 2.
1039  //
1040  template <typename Real>
1041  BinnedDensity1D(double location, double scale,
1042  const Real* data, unsigned dataLen,
1043  unsigned interpolationDegree);
1044 
1045  inline BinnedDensity1D(const double location, const double scale,
1046  const std::vector<double>& table,
1047  const unsigned interpolationDegree)
1049  {
1050  const unsigned long sz = table.size();
1051  initialize(sz ? &table[0] : (double*)0, sz, interpolationDegree);
1052  }
1053 
1054  inline virtual BinnedDensity1D* clone() const
1055  {return new BinnedDensity1D(*this);}
1056 
1057  inline virtual ~BinnedDensity1D() {}
1058 
1059  inline unsigned interpolationDegree() const {return deg_;}
1060  inline unsigned tableLength() const {return len_;}
1061  inline const double* tableData() const {return &table_[0];}
1062 
1063  // Methods needed for I/O
1064  virtual gs::ClassId classId() const {return gs::ClassId(*this);}
1065  virtual bool write(std::ostream& os) const;
1066 
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);
1070 
1071  protected:
1072  virtual bool isEqual(const AbsDistribution1D&) const;
1073 
1074  private:
1076 
1077  // The following constructor creates interpolator
1078  // of maximum degree possible
1079  BinnedDensity1D(double location, double scale,
1080  const std::vector<double>& params);
1081  inline static int nParameters() {return -1;}
1082 
1083  double unscaledDensity(double x) const;
1084  double unscaledCdf(double x) const;
1085  double unscaledQuantile(double x) const;
1086 
1087  inline double unscaledExceedance(const double x) const
1088  {return 1.0 - unscaledCdf(x);}
1089 
1090  template <typename Real> void initialize(
1091  const Real* data, unsigned dataLen, unsigned interpolationDegree);
1092 
1093  void normalize();
1094  double interpolate(double x) const;
1095 
1096  std::vector<double> table_;
1097  std::vector<double> cdf_;
1098  double step_;
1099  unsigned len_;
1100  unsigned deg_;
1101  unsigned firstNonZeroBin_;
1102  unsigned lastNonZeroBin_;
1103  };
1104 }
1105 
1106 #include "npstat/stat/Distributions1D.icc"
1107 
1108 #endif // NPSTAT_DISTRIBUTIONS1D_HH_
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: Distributions1D.hh:22
virtual Uniform1D * clone() const
Definition: Distributions1D.hh:26
virtual gs::ClassId classId() const
Definition: Distributions1D.hh:31
virtual bool isEqual(const AbsDistribution1D &r) const
Definition: Distributions1D.hh:39
Definition: AbsArrayProjector.hh:14
Definition: AbsDistribution1D.hh:31
Definition: AbsRandomGenerator.hh:27