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 */
VerifyDistribution(const char name[],PetscBool pos,PetscProbFn * pdf,PetscProbFn * cdf)9f8662bd6SBarry Smith static PetscErrorCode VerifyDistribution(const char name[], PetscBool pos, PetscProbFn *pdf, PetscProbFn *cdf)
10d71ae5a4SJacob Faibussowitsch {
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;
179566063dSJacob Faibussowitsch PetscCall(PetscDTTanhSinhIntegrate((void (*)(const PetscReal[], void *, PetscReal *))pdf, lower, upper, digits, NULL, &integral));
186497c311SBarry Smith PetscCheck(PetscAbsReal(integral - 1) < 100 * PETSC_MACHINE_EPSILON, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "PDF %s must integrate to 1, not %g", name, (double)integral);
19d6685f55SMatthew G. Knepley for (i = 0; i <= 10; ++i) {
20d6685f55SMatthew G. Knepley const PetscReal x = i;
21d6685f55SMatthew G. Knepley
229566063dSJacob Faibussowitsch PetscCall(PetscDTTanhSinhIntegrate((void (*)(const PetscReal[], void *, PetscReal *))pdf, lower, x, digits, NULL, &integral));
239566063dSJacob Faibussowitsch PetscCall(cdf(&x, NULL, &integral2));
2463a3b9bcSJacob Faibussowitsch 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, (double)integral, (double)integral2, (double)x);
25d6685f55SMatthew G. Knepley }
263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
27d6685f55SMatthew G. Knepley }
28d6685f55SMatthew G. Knepley
TestDistributions()29d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestDistributions()
30d71ae5a4SJacob Faibussowitsch {
31f8662bd6SBarry Smith PetscProbFn *pdf[] = {PetscPDFMaxwellBoltzmann1D, PetscPDFMaxwellBoltzmann2D, PetscPDFMaxwellBoltzmann3D, PetscPDFGaussian1D};
32f8662bd6SBarry Smith PetscProbFn *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;
37f8662bd6SBarry Smith for (PetscInt i = 0; i < (PetscInt)PETSC_STATIC_ARRAY_LENGTH(pdf); ++i) PetscCall(VerifyDistribution(name[i], pos[i], pdf[i], cdf[i]));
383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
39d6685f55SMatthew G. Knepley }
40d6685f55SMatthew G. Knepley
TestSampling()41d71ae5a4SJacob Faibussowitsch static PetscErrorCode TestSampling()
42d71ae5a4SJacob Faibussowitsch {
43f8662bd6SBarry Smith PetscProbFn *cdf[2] = {PetscCDFMaxwellBoltzmann1D, PetscCDFMaxwellBoltzmann2D};
44f8662bd6SBarry Smith PetscProbFn *sampler[2] = {PetscPDFSampleGaussian1D, PetscPDFSampleGaussian2D};
45d6685f55SMatthew G. Knepley PetscInt dim[2] = {1, 2};
46d6685f55SMatthew G. Knepley PetscRandom rnd;
47d6685f55SMatthew G. Knepley Vec v;
48d6685f55SMatthew G. Knepley PetscScalar *a;
49d6685f55SMatthew G. Knepley PetscReal alpha, confidenceLevel = 0.05;
50d6685f55SMatthew G. Knepley PetscInt n = 1000, s, i, d;
51d6685f55SMatthew G. Knepley
52d6685f55SMatthew G. Knepley PetscFunctionBeginUser;
539566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF, &rnd));
549566063dSJacob Faibussowitsch PetscCall(PetscRandomSetInterval(rnd, 0, 1.));
559566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(rnd));
56d6685f55SMatthew G. Knepley
57d6685f55SMatthew G. Knepley for (s = 0; s < 2; ++s) {
589566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, n * dim[s], &v));
599566063dSJacob Faibussowitsch PetscCall(VecSetBlockSize(v, dim[s]));
609566063dSJacob Faibussowitsch PetscCall(VecGetArray(v, &a));
61d6685f55SMatthew G. Knepley for (i = 0; i < n; ++i) {
62d6685f55SMatthew G. Knepley PetscReal r[3], o[3];
63d6685f55SMatthew G. Knepley
649566063dSJacob Faibussowitsch for (d = 0; d < dim[s]; ++d) PetscCall(PetscRandomGetValueReal(rnd, &r[d]));
659566063dSJacob Faibussowitsch PetscCall(sampler[s](r, NULL, o));
66d6685f55SMatthew G. Knepley for (d = 0; d < dim[s]; ++d) a[i * dim[s] + d] = o[d];
67d6685f55SMatthew G. Knepley }
689566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(v, &a));
699566063dSJacob Faibussowitsch PetscCall(PetscProbComputeKSStatistic(v, cdf[s], &alpha));
7063a3b9bcSJacob Faibussowitsch PetscCheck(alpha < confidenceLevel, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "KS finds sampling does not match the distribution at confidence level %.2g", (double)confidenceLevel);
719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v));
72d6685f55SMatthew G. Knepley }
739566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&rnd));
743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
75d6685f55SMatthew G. Knepley }
76d6685f55SMatthew G. Knepley
main(int argc,char ** argv)77d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
78d71ae5a4SJacob Faibussowitsch {
79327415f7SBarry Smith PetscFunctionBeginUser;
809566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
819566063dSJacob Faibussowitsch PetscCall(TestDistributions());
829566063dSJacob Faibussowitsch PetscCall(TestSampling());
839566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
84b122ec5aSJacob Faibussowitsch return 0;
85d6685f55SMatthew G. Knepley }
86d6685f55SMatthew G. Knepley
87d6685f55SMatthew G. Knepley /*TEST
88d6685f55SMatthew G. Knepley
89d6685f55SMatthew G. Knepley test:
90d6685f55SMatthew G. Knepley suffix: 0
91d6685f55SMatthew G. Knepley requires: ks
92d6685f55SMatthew G. Knepley args:
93*3886731fSPierre Jolivet output_file: output/empty.out
94d6685f55SMatthew G. Knepley
95d6685f55SMatthew G. Knepley TEST*/
96