xref: /petsc/src/dm/dt/tests/ex14.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1d6685f55SMatthew G. Knepley const char help[] = "Tests properties of probability distributions";
2d6685f55SMatthew G. Knepley 
3d6685f55SMatthew G. Knepley #include <petscdt.h>
4d6685f55SMatthew G. Knepley 
5d6685f55SMatthew G. Knepley /* Checks that
6d6685f55SMatthew G. Knepley    - the PDF integrates to 1
7d6685f55SMatthew G. Knepley    - the incomplete integral of the PDF is the CDF at many points
8d6685f55SMatthew G. Knepley */
9d6685f55SMatthew G. Knepley static PetscErrorCode VerifyDistribution(const char name[], PetscBool pos, PetscProbFunc pdf, PetscProbFunc cdf)
10d6685f55SMatthew G. Knepley {
11d6685f55SMatthew G. Knepley   const PetscInt digits = 14;
12d6685f55SMatthew G. Knepley   PetscReal      lower = pos ? 0. : -10., upper = 10.;
13d6685f55SMatthew G. Knepley   PetscReal      integral, integral2;
14d6685f55SMatthew G. Knepley   PetscInt       i;
15d6685f55SMatthew G. Knepley 
16d6685f55SMatthew G. Knepley   PetscFunctionBeginUser;
175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscDTTanhSinhIntegrate((void (*)(const PetscReal[], void *, PetscReal *)) pdf, lower, upper, digits, NULL, &integral));
18d6685f55SMatthew G. Knepley   PetscCheck(PetscAbsReal(integral - 1.0) < 100*PETSC_MACHINE_EPSILON, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PDF %s must integrate to 1, not %g", name, integral);
19d6685f55SMatthew G. Knepley   for (i = 0; i <= 10; ++i) {
20d6685f55SMatthew G. Knepley     const PetscReal x = i;
21d6685f55SMatthew G. Knepley 
225f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscDTTanhSinhIntegrate((void (*)(const PetscReal[], void *, PetscReal *)) pdf, lower, x, digits, NULL, &integral));
235f80ce2aSJacob Faibussowitsch     CHKERRQ(cdf(&x, NULL, &integral2));
24d6685f55SMatthew G. Knepley     PetscCheck(PetscAbsReal(integral - integral2) < PETSC_SQRT_MACHINE_EPSILON, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Integral of PDF %s %g != %g CDF at x = %g", name, integral, integral2, x);
25d6685f55SMatthew G. Knepley   }
26d6685f55SMatthew G. Knepley   PetscFunctionReturn(0);
27d6685f55SMatthew G. Knepley }
28d6685f55SMatthew G. Knepley 
29d6685f55SMatthew G. Knepley static PetscErrorCode TestDistributions()
30d6685f55SMatthew G. Knepley {
31d6685f55SMatthew G. Knepley   PetscProbFunc  pdf[]  = {PetscPDFMaxwellBoltzmann1D, PetscPDFMaxwellBoltzmann2D, PetscPDFMaxwellBoltzmann3D, PetscPDFGaussian1D};
32d6685f55SMatthew G. Knepley   PetscProbFunc  cdf[]  = {PetscCDFMaxwellBoltzmann1D, PetscCDFMaxwellBoltzmann2D, PetscCDFMaxwellBoltzmann3D, PetscCDFGaussian1D};
33d6685f55SMatthew G. Knepley   PetscBool      pos[]  = {PETSC_TRUE,                 PETSC_TRUE,                 PETSC_TRUE,                 PETSC_FALSE};
34d6685f55SMatthew G. Knepley   const char    *name[] = {"Maxwell-Boltzmann 1D",     "Maxwell-Boltzmann 2D",     "Maxwell-Boltzmann 3D",     "Gaussian"};
35d6685f55SMatthew G. Knepley 
36d6685f55SMatthew G. Knepley   PetscFunctionBeginUser;
37d6685f55SMatthew G. Knepley   for (PetscInt i = 0; i < (PetscInt) (sizeof(pdf)/sizeof(PetscProbFunc)); ++i) {
385f80ce2aSJacob Faibussowitsch     CHKERRQ(VerifyDistribution(name[i], pos[i], pdf[i], cdf[i]));
39d6685f55SMatthew G. Knepley   }
40d6685f55SMatthew G. Knepley   PetscFunctionReturn(0);
41d6685f55SMatthew G. Knepley }
42d6685f55SMatthew G. Knepley 
43d6685f55SMatthew G. Knepley static PetscErrorCode TestSampling()
44d6685f55SMatthew G. Knepley {
45d6685f55SMatthew G. Knepley   PetscProbFunc  cdf[2]     = {PetscCDFMaxwellBoltzmann1D, PetscCDFMaxwellBoltzmann2D};
46d6685f55SMatthew G. Knepley   PetscProbFunc  sampler[2] = {PetscPDFSampleGaussian1D,   PetscPDFSampleGaussian2D};
47d6685f55SMatthew G. Knepley   PetscInt       dim[2]     = {1, 2};
48d6685f55SMatthew G. Knepley   PetscRandom    rnd;
49d6685f55SMatthew G. Knepley   Vec            v;
50d6685f55SMatthew G. Knepley   PetscScalar   *a;
51d6685f55SMatthew G. Knepley   PetscReal      alpha, confidenceLevel = 0.05;
52d6685f55SMatthew G. Knepley   PetscInt       n = 1000, s, i, d;
53d6685f55SMatthew G. Knepley 
54d6685f55SMatthew G. Knepley   PetscFunctionBeginUser;
555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF, &rnd));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetInterval(rnd, 0, 1.));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(rnd));
58d6685f55SMatthew G. Knepley 
59d6685f55SMatthew G. Knepley   for (s = 0; s < 2; ++s) {
605f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCreateSeq(PETSC_COMM_SELF, n*dim[s], &v));
615f80ce2aSJacob Faibussowitsch     CHKERRQ(VecSetBlockSize(v, dim[s]));
625f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(v, &a));
63d6685f55SMatthew G. Knepley     for (i = 0; i < n; ++i) {
64d6685f55SMatthew G. Knepley       PetscReal r[3], o[3];
65d6685f55SMatthew G. Knepley 
665f80ce2aSJacob Faibussowitsch       for (d = 0; d < dim[s]; ++d) CHKERRQ(PetscRandomGetValueReal(rnd, &r[d]));
675f80ce2aSJacob Faibussowitsch       CHKERRQ(sampler[s](r, NULL, o));
68d6685f55SMatthew G. Knepley       for (d = 0; d < dim[s]; ++d) a[i*dim[s]+d] = o[d];
69d6685f55SMatthew G. Knepley     }
705f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(v, &a));
715f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscProbComputeKSStatistic(v, cdf[s], &alpha));
72d6685f55SMatthew G. Knepley     PetscCheck(alpha < confidenceLevel, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "KS finds sampling does not match the distribution at confidence level %.2g", confidenceLevel);
735f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&v));
74d6685f55SMatthew G. Knepley   }
755f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&rnd));
76d6685f55SMatthew G. Knepley   PetscFunctionReturn(0);
77d6685f55SMatthew G. Knepley }
78d6685f55SMatthew G. Knepley 
79d6685f55SMatthew G. Knepley int main(int argc, char **argv)
80d6685f55SMatthew G. Knepley {
81d6685f55SMatthew G. Knepley 
82*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL, help));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(TestDistributions());
845f80ce2aSJacob Faibussowitsch   CHKERRQ(TestSampling());
85*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
86*b122ec5aSJacob Faibussowitsch   return 0;
87d6685f55SMatthew G. Knepley }
88d6685f55SMatthew G. Knepley 
89d6685f55SMatthew G. Knepley /*TEST
90d6685f55SMatthew G. Knepley 
91d6685f55SMatthew G. Knepley   test:
92d6685f55SMatthew G. Knepley     suffix: 0
93d6685f55SMatthew G. Knepley     requires: ks
94d6685f55SMatthew G. Knepley     args:
95d6685f55SMatthew G. Knepley 
96d6685f55SMatthew G. Knepley TEST*/
97