xref: /petsc/src/dm/dt/tests/ex14.c (revision d6685f554fbda8d96c6a5d73ab0e7a4e21a05c51)
1*d6685f55SMatthew G. Knepley const char help[] = "Tests properties of probability distributions";
2*d6685f55SMatthew G. Knepley 
3*d6685f55SMatthew G. Knepley #include <petscdt.h>
4*d6685f55SMatthew G. Knepley 
5*d6685f55SMatthew G. Knepley /* Checks that
6*d6685f55SMatthew G. Knepley    - the PDF integrates to 1
7*d6685f55SMatthew G. Knepley    - the incomplete integral of the PDF is the CDF at many points
8*d6685f55SMatthew G. Knepley */
9*d6685f55SMatthew G. Knepley static PetscErrorCode VerifyDistribution(const char name[], PetscBool pos, PetscProbFunc pdf, PetscProbFunc cdf)
10*d6685f55SMatthew G. Knepley {
11*d6685f55SMatthew G. Knepley   const PetscInt digits = 14;
12*d6685f55SMatthew G. Knepley   PetscReal      lower = pos ? 0. : -10., upper = 10.;
13*d6685f55SMatthew G. Knepley   PetscReal      integral, integral2;
14*d6685f55SMatthew G. Knepley   PetscInt       i;
15*d6685f55SMatthew G. Knepley   PetscErrorCode ierr;
16*d6685f55SMatthew G. Knepley 
17*d6685f55SMatthew G. Knepley   PetscFunctionBeginUser;
18*d6685f55SMatthew G. Knepley   ierr = PetscDTTanhSinhIntegrate((void (*)(const PetscReal[], void *, PetscReal *)) pdf, lower, upper, digits, NULL, &integral);CHKERRQ(ierr);
19*d6685f55SMatthew 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);
20*d6685f55SMatthew G. Knepley   for (i = 0; i <= 10; ++i) {
21*d6685f55SMatthew G. Knepley     const PetscReal x = i;
22*d6685f55SMatthew G. Knepley 
23*d6685f55SMatthew G. Knepley     ierr = PetscDTTanhSinhIntegrate((void (*)(const PetscReal[], void *, PetscReal *)) pdf, lower, x, digits, NULL, &integral);CHKERRQ(ierr);
24*d6685f55SMatthew G. Knepley     ierr = cdf(&x, NULL, &integral2);CHKERRQ(ierr);
25*d6685f55SMatthew 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);
26*d6685f55SMatthew G. Knepley   }
27*d6685f55SMatthew G. Knepley   PetscFunctionReturn(0);
28*d6685f55SMatthew G. Knepley }
29*d6685f55SMatthew G. Knepley 
30*d6685f55SMatthew G. Knepley static PetscErrorCode TestDistributions()
31*d6685f55SMatthew G. Knepley {
32*d6685f55SMatthew G. Knepley   PetscProbFunc  pdf[]  = {PetscPDFMaxwellBoltzmann1D, PetscPDFMaxwellBoltzmann2D, PetscPDFMaxwellBoltzmann3D, PetscPDFGaussian1D};
33*d6685f55SMatthew G. Knepley   PetscProbFunc  cdf[]  = {PetscCDFMaxwellBoltzmann1D, PetscCDFMaxwellBoltzmann2D, PetscCDFMaxwellBoltzmann3D, PetscCDFGaussian1D};
34*d6685f55SMatthew G. Knepley   PetscBool      pos[]  = {PETSC_TRUE,                 PETSC_TRUE,                 PETSC_TRUE,                 PETSC_FALSE};
35*d6685f55SMatthew G. Knepley   const char    *name[] = {"Maxwell-Boltzmann 1D",     "Maxwell-Boltzmann 2D",     "Maxwell-Boltzmann 3D",     "Gaussian"};
36*d6685f55SMatthew G. Knepley   PetscErrorCode ierr;
37*d6685f55SMatthew G. Knepley 
38*d6685f55SMatthew G. Knepley   PetscFunctionBeginUser;
39*d6685f55SMatthew G. Knepley   for (PetscInt i = 0; i < (PetscInt) (sizeof(pdf)/sizeof(PetscProbFunc)); ++i) {
40*d6685f55SMatthew G. Knepley     ierr = VerifyDistribution(name[i], pos[i], pdf[i], cdf[i]);CHKERRQ(ierr);
41*d6685f55SMatthew G. Knepley   }
42*d6685f55SMatthew G. Knepley   PetscFunctionReturn(0);
43*d6685f55SMatthew G. Knepley }
44*d6685f55SMatthew G. Knepley 
45*d6685f55SMatthew G. Knepley static PetscErrorCode TestSampling()
46*d6685f55SMatthew G. Knepley {
47*d6685f55SMatthew G. Knepley   PetscProbFunc  cdf[2]     = {PetscCDFMaxwellBoltzmann1D, PetscCDFMaxwellBoltzmann2D};
48*d6685f55SMatthew G. Knepley   PetscProbFunc  sampler[2] = {PetscPDFSampleGaussian1D,   PetscPDFSampleGaussian2D};
49*d6685f55SMatthew G. Knepley   PetscInt       dim[2]     = {1, 2};
50*d6685f55SMatthew G. Knepley   PetscRandom    rnd;
51*d6685f55SMatthew G. Knepley   Vec            v;
52*d6685f55SMatthew G. Knepley   PetscScalar   *a;
53*d6685f55SMatthew G. Knepley   PetscReal      alpha, confidenceLevel = 0.05;
54*d6685f55SMatthew G. Knepley   PetscInt       n = 1000, s, i, d;
55*d6685f55SMatthew G. Knepley   PetscErrorCode ierr;
56*d6685f55SMatthew G. Knepley 
57*d6685f55SMatthew G. Knepley   PetscFunctionBeginUser;
58*d6685f55SMatthew G. Knepley   ierr = PetscRandomCreate(PETSC_COMM_SELF, &rnd);CHKERRQ(ierr);
59*d6685f55SMatthew G. Knepley   ierr = PetscRandomSetInterval(rnd, 0, 1.);CHKERRQ(ierr);
60*d6685f55SMatthew G. Knepley   ierr = PetscRandomSetFromOptions(rnd);CHKERRQ(ierr);
61*d6685f55SMatthew G. Knepley 
62*d6685f55SMatthew G. Knepley   for (s = 0; s < 2; ++s) {
63*d6685f55SMatthew G. Knepley     ierr = VecCreateSeq(PETSC_COMM_SELF, n*dim[s], &v);CHKERRQ(ierr);
64*d6685f55SMatthew G. Knepley     ierr = VecSetBlockSize(v, dim[s]);CHKERRQ(ierr);
65*d6685f55SMatthew G. Knepley     ierr = VecGetArray(v, &a);CHKERRQ(ierr);
66*d6685f55SMatthew G. Knepley     for (i = 0; i < n; ++i) {
67*d6685f55SMatthew G. Knepley       PetscReal r[3], o[3];
68*d6685f55SMatthew G. Knepley 
69*d6685f55SMatthew G. Knepley       for (d = 0; d < dim[s]; ++d) {ierr = PetscRandomGetValueReal(rnd, &r[d]);CHKERRQ(ierr);}
70*d6685f55SMatthew G. Knepley       ierr = sampler[s](r, NULL, o);CHKERRQ(ierr);
71*d6685f55SMatthew G. Knepley       for (d = 0; d < dim[s]; ++d) a[i*dim[s]+d] = o[d];
72*d6685f55SMatthew G. Knepley     }
73*d6685f55SMatthew G. Knepley     ierr = VecRestoreArray(v, &a);CHKERRQ(ierr);
74*d6685f55SMatthew G. Knepley     ierr = PetscProbComputeKSStatistic(v, cdf[s], &alpha);CHKERRQ(ierr);
75*d6685f55SMatthew 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);
76*d6685f55SMatthew G. Knepley     ierr = VecDestroy(&v);CHKERRQ(ierr);
77*d6685f55SMatthew G. Knepley   }
78*d6685f55SMatthew G. Knepley   ierr = PetscRandomDestroy(&rnd);CHKERRQ(ierr);
79*d6685f55SMatthew G. Knepley   PetscFunctionReturn(0);
80*d6685f55SMatthew G. Knepley }
81*d6685f55SMatthew G. Knepley 
82*d6685f55SMatthew G. Knepley int main(int argc, char **argv)
83*d6685f55SMatthew G. Knepley {
84*d6685f55SMatthew G. Knepley   PetscErrorCode ierr;
85*d6685f55SMatthew G. Knepley 
86*d6685f55SMatthew G. Knepley   ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr;
87*d6685f55SMatthew G. Knepley   ierr = TestDistributions();CHKERRQ(ierr);
88*d6685f55SMatthew G. Knepley   ierr = TestSampling();CHKERRQ(ierr);
89*d6685f55SMatthew G. Knepley   ierr = PetscFinalize();
90*d6685f55SMatthew G. Knepley   return ierr;
91*d6685f55SMatthew G. Knepley }
92*d6685f55SMatthew G. Knepley 
93*d6685f55SMatthew G. Knepley /*TEST
94*d6685f55SMatthew G. Knepley 
95*d6685f55SMatthew G. Knepley   test:
96*d6685f55SMatthew G. Knepley     suffix: 0
97*d6685f55SMatthew G. Knepley     requires: ks
98*d6685f55SMatthew G. Knepley     args:
99*d6685f55SMatthew G. Knepley 
100*d6685f55SMatthew G. Knepley TEST*/
101