xref: /petsc/src/dm/dt/interface/dtprob.c (revision ffa8c5705e8ab2cf85ee1d14dbe507a6e2eb5283)
1 #include <petscdt.h> /*I "petscdt.h" I*/
2 
3 #include <petscvec.h>
4 #include <petscdraw.h>
5 #include <petsc/private/petscimpl.h>
6 
7 const char * const DTProbDensityTypes[] = {"constant", "gaussian", "maxwell_boltzmann", "DTProb Density", "DTPROB_DENSITY_", NULL};
8 
9 /*@
10   PetscPDFMaxwellBoltzmann1D - PDF for the Maxwell-Boltzmann distribution in 1D
11 
12   Not collective
13 
14   Input Parameter:
15 . x - Speed in [0, \infty]
16 
17   Output Parameter:
18 . p - The probability density at x
19 
20   Level: beginner
21 
22 .seealso: PetscCDFMaxwellBoltzmann1D
23 @*/
24 PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
25 {
26   p[0] = PetscSqrtReal(2. / PETSC_PI) * PetscExpReal(-0.5 * PetscSqr(x[0]));
27   return 0;
28 }
29 
30 /*@
31   PetscCDFMaxwellBoltzmann1D - CDF for the Maxwell-Boltzmann distribution in 1D
32 
33   Not collective
34 
35   Input Parameter:
36 . x - Speed in [0, \infty]
37 
38   Output Parameter:
39 . p - The probability density at x
40 
41   Level: beginner
42 
43 .seealso: PetscPDFMaxwellBoltzmann1D
44 @*/
45 PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
46 {
47   p[0] = PetscErfReal(x[0] / PETSC_SQRT2);
48   return 0;
49 }
50 
51 /*@
52   PetscPDFMaxwellBoltzmann2D - PDF for the Maxwell-Boltzmann distribution in 2D
53 
54   Not collective
55 
56   Input Parameter:
57 . x - Speed in [0, \infty]
58 
59   Output Parameter:
60 . p - The probability density at x
61 
62   Level: beginner
63 
64 .seealso: PetscCDFMaxwellBoltzmann2D
65 @*/
66 PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
67 {
68   p[0] = x[0] * PetscExpReal(-0.5 * PetscSqr(x[0]));
69   return 0;
70 }
71 
72 /*@
73   PetscCDFMaxwellBoltzmann2D - CDF for the Maxwell-Boltzmann distribution in 2D
74 
75   Not collective
76 
77   Input Parameter:
78 . x - Speed in [0, \infty]
79 
80   Output Parameter:
81 . p - The probability density at x
82 
83   Level: beginner
84 
85 .seealso: PetscPDFMaxwellBoltzmann2D
86 @*/
87 PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
88 {
89   p[0] = 1. - PetscExpReal(-0.5 * PetscSqr(x[0]));
90   return 0;
91 }
92 
93 /*@
94   PetscPDFMaxwellBoltzmann3D - PDF for the Maxwell-Boltzmann distribution in 3D
95 
96   Not collective
97 
98   Input Parameter:
99 . x - Speed in [0, \infty]
100 
101   Output Parameter:
102 . p - The probability density at x
103 
104   Level: beginner
105 
106 .seealso: PetscCDFMaxwellBoltzmann3D
107 @*/
108 PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
109 {
110   p[0] = PetscSqrtReal(2. / PETSC_PI) * PetscSqr(x[0]) * PetscExpReal(-0.5 * PetscSqr(x[0]));
111   return 0;
112 }
113 
114 /*@
115   PetscCDFMaxwellBoltzmann3D - CDF for the Maxwell-Boltzmann distribution in 3D
116 
117   Not collective
118 
119   Input Parameter:
120 . x - Speed in [0, \infty]
121 
122   Output Parameter:
123 . p - The probability density at x
124 
125   Level: beginner
126 
127 .seealso: PetscPDFMaxwellBoltzmann3D
128 @*/
129 PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
130 {
131   p[0] = PetscErfReal(x[0] / PETSC_SQRT2) - PetscSqrtReal(2. / PETSC_PI) * x[0] * PetscExpReal(-0.5 * PetscSqr(x[0]));
132   return 0;
133 }
134 
135 /*@
136   PetscPDFGaussian1D - PDF for the Gaussian distribution in 1D
137 
138   Not collective
139 
140   Input Parameter:
141 . x - Coordinate in [-\infty, \infty]
142 
143   Output Parameter:
144 . p - The probability density at x
145 
146   Level: beginner
147 
148 .seealso: PetscPDFMaxwellBoltzmann3D
149 @*/
150 PetscErrorCode PetscPDFGaussian1D(const PetscReal x[], const PetscReal scale[], PetscReal p[])
151 {
152   const PetscReal sigma = scale ? scale[0] : 1.;
153   p[0] = PetscSqrtReal(1. / (2.*PETSC_PI)) * PetscExpReal(-0.5 * PetscSqr(x[0] / sigma)) / sigma;
154   return 0;
155 }
156 
157 PetscErrorCode PetscCDFGaussian1D(const PetscReal x[], const PetscReal scale[], PetscReal p[])
158 {
159   const PetscReal sigma = scale ? scale[0] : 1.;
160   p[0] = 0.5 * (1. + PetscErfReal(x[0] / PETSC_SQRT2 / sigma));
161   return 0;
162 }
163 
164 /*@
165   PetscPDFSampleGaussian1D - Sample uniformly from a Gaussian distribution in 1D
166 
167   Not collective
168 
169   Input Parameters:
170 + p - A uniform variable on [0, 1]
171 - dummy - ignored
172 
173   Output Parameter:
174 . x - Coordinate in [-\infty, \infty]
175 
176   Level: beginner
177 
178   Note: http://www.mimirgames.com/articles/programming/approximations-of-the-inverse-error-function/
179   https://stackoverflow.com/questions/27229371/inverse-error-function-in-c
180 @*/
181 PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
182 {
183   const PetscReal q = 2*p[0] - 1.;
184   const PetscInt  maxIter = 100;
185   PetscReal       ck[100], r = 0.;
186   PetscInt        k, m;
187 
188   PetscFunctionBeginHot;
189   /* Transform input to [-1, 1] since the code below computes the inverse error function */
190   for (k = 0; k < maxIter; ++k) ck[k] = 0.;
191   ck[0] = 1;
192   r = ck[0]* (PetscSqrtReal(PETSC_PI) / 2.) * q;
193   for (k = 1; k < maxIter; ++k) {
194     const PetscReal temp = 2.*k + 1.;
195 
196     for (m = 0; m <= k-1; ++m) {
197       PetscReal denom = (m + 1.)*(2.*m + 1.);
198 
199       ck[k] += (ck[m]*ck[k-1-m])/denom;
200     }
201     r += (ck[k]/temp) * PetscPowReal((PetscSqrtReal(PETSC_PI)/2.) * q, 2.*k+1.);
202   }
203   /* Scale erfinv() by \sqrt{\pi/2} */
204   x[0] = PetscSqrtReal(PETSC_PI * 0.5) * r;
205   PetscFunctionReturn(0);
206 }
207 
208 /*@
209   PetscPDFGaussian2D - PDF for the Gaussian distribution in 2D
210 
211   Not collective
212 
213   Input Parameters:
214 + x - Coordinate in [-\infty, \infty]^2
215 - dummy - ignored
216 
217   Output Parameter:
218 . p - The probability density at x
219 
220   Level: beginner
221 
222 .seealso: PetscPDFSampleGaussian2D(), PetscPDFMaxwellBoltzmann3D()
223 @*/
224 PetscErrorCode PetscPDFGaussian2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
225 {
226   p[0] = (1. / PETSC_PI) * PetscExpReal(-0.5 * (PetscSqr(x[0]) + PetscSqr(x[1])));
227   return 0;
228 }
229 
230 /*@
231   PetscPDFSampleGaussian2D - Sample uniformly from a Gaussian distribution in 2D
232 
233   Not collective
234 
235   Input Parameters:
236 + p - A uniform variable on [0, 1]^2
237 - dummy - ignored
238 
239   Output Parameter:
240 . x - Coordinate in [-\infty, \infty]^2
241 
242   Level: beginner
243 
244   Note: https://en.wikipedia.org/wiki/Box%E2%80%93Muller_transform
245 
246 .seealso: PetscPDFGaussian2D(), PetscPDFMaxwellBoltzmann3D()
247 @*/
248 PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
249 {
250   const PetscReal mag = PetscSqrtReal(-2.0 * PetscLogReal(p[0]));
251   x[0] = mag * PetscCosReal(2.0 * PETSC_PI * p[1]);
252   x[1] = mag * PetscSinReal(2.0 * PETSC_PI * p[1]);
253   return 0;
254 }
255 
256 /*@
257   PetscPDFConstant1D - PDF for the uniform distribution in 1D
258 
259   Not collective
260 
261   Input Parameter:
262 . x - Coordinate in [-1, 1]
263 
264   Output Parameter:
265 . p - The probability density at x
266 
267   Level: beginner
268 
269 .seealso: PetscCDFConstant1D()
270 @*/
271 PetscErrorCode PetscPDFConstant1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
272 {
273   p[0] = x[0] > -1. && x[0] <= 1. ? 0.5 : 0.;
274   return 0;
275 }
276 
277 PetscErrorCode PetscCDFConstant1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
278 {
279   p[0] = x[0] <= -1. ? 0. : (x[0] > 1. ? 1. : 0.5*(x[0] + 1.));
280   return 0;
281 }
282 
283 /*@
284   PetscPDFSampleConstant1D - Sample uniformly from a uniform distribution on [-1, 1] in 1D
285 
286   Not collective
287 
288   Input Parameter:
289 . p - A uniform variable on [0, 1]
290 
291   Output Parameter:
292 . x - Coordinate in [-1, 1]
293 
294   Level: beginner
295 
296 .seealso: PetscPDFConstant1D(), PetscCDFConstant1D()
297 @*/
298 PetscErrorCode PetscPDFSampleConstant1D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
299 {
300   x[0] = 2.*p[1] - 1.;
301   return 0;
302 }
303 
304 /*@C
305   PetscProbCreateFromOptions - Return the probability distribution specified by the argumetns and options
306 
307   Not collective
308 
309   Input Parameters:
310 + dim    - The dimension of sample points
311 . prefix - The options prefix, or NULL
312 - name   - The option name for the probility distribution type
313 
314   Output Parameters:
315 + pdf - The PDF of this type
316 . cdf - The CDF of this type
317 - sampler - The PDF sampler of this type
318 
319   Level: intermediate
320 
321 .seealso: PetscPDFMaxwellBoltzmann1D(), PetscPDFGaussian1D(), PetscPDFConstant1D()
322 @*/
323 PetscErrorCode PetscProbCreateFromOptions(PetscInt dim, const char prefix[], const char name[], PetscProbFunc *pdf, PetscProbFunc *cdf, PetscProbFunc *sampler)
324 {
325   DTProbDensityType den = DTPROB_DENSITY_GAUSSIAN;
326   PetscErrorCode    ierr;
327 
328   PetscFunctionBegin;
329   ierr = PetscOptionsBegin(PETSC_COMM_SELF, prefix, "PetscProb Options", "DT");PetscCall(ierr);
330   PetscCall(PetscOptionsEnum(name, "Method to compute PDF <constant, gaussian>", "", DTProbDensityTypes, (PetscEnum) den, (PetscEnum *) &den, NULL));
331   ierr = PetscOptionsEnd();PetscCall(ierr);
332 
333   if (pdf) {PetscValidPointer(pdf, 4); *pdf = NULL;}
334   if (cdf) {PetscValidPointer(cdf, 5); *cdf = NULL;}
335   if (sampler) {PetscValidPointer(sampler, 6); *sampler = NULL;}
336   switch (den) {
337     case DTPROB_DENSITY_CONSTANT:
338       switch (dim) {
339         case 1:
340           if (pdf) *pdf = PetscPDFConstant1D;
341           if (cdf) *cdf = PetscCDFConstant1D;
342           if (sampler) *sampler = PetscPDFSampleConstant1D;
343           break;
344         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]);
345       }
346       break;
347     case DTPROB_DENSITY_GAUSSIAN:
348       switch (dim) {
349         case 1:
350           if (pdf) *pdf = PetscPDFGaussian1D;
351           if (cdf) *cdf = PetscCDFGaussian1D;
352           if (sampler) *sampler = PetscPDFSampleGaussian1D;
353           break;
354         case 2:
355           if (pdf) *pdf = PetscPDFGaussian2D;
356           if (sampler) *sampler = PetscPDFSampleGaussian2D;
357           break;
358         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]);
359       }
360       break;
361     case DTPROB_DENSITY_MAXWELL_BOLTZMANN:
362       switch (dim) {
363         case 1:
364           if (pdf) *pdf = PetscPDFMaxwellBoltzmann1D;
365           if (cdf) *cdf = PetscCDFMaxwellBoltzmann1D;
366           break;
367         case 2:
368           if (pdf) *pdf = PetscPDFMaxwellBoltzmann2D;
369           if (cdf) *cdf = PetscCDFMaxwellBoltzmann2D;
370           break;
371         case 3:
372           if (pdf) *pdf = PetscPDFMaxwellBoltzmann3D;
373           if (cdf) *cdf = PetscCDFMaxwellBoltzmann3D;
374           break;
375         default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]);
376       }
377       break;
378     default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Density type %s is not supported", DTProbDensityTypes[PetscMax(0, PetscMin(den, DTPROB_NUM_DENSITY))]);
379   }
380   PetscFunctionReturn(0);
381 }
382 
383 #ifdef PETSC_HAVE_KS
384 EXTERN_C_BEGIN
385 #include <KolmogorovSmirnovDist.h>
386 EXTERN_C_END
387 #endif
388 
389 /*@C
390   PetscProbComputeKSStatistic - Compute the Kolmogorov-Smirnov statistic for the empirical distribution for an input vector, compared to an analytic CDF.
391 
392   Collective on v
393 
394   Input Parameters:
395 + v   - The data vector, blocksize is the sample dimension
396 - cdf - The analytic CDF
397 
398   Output Parameter:
399 . alpha - The KS statisic
400 
401   Level: advanced
402 
403   Note: The Kolmogorov-Smirnov statistic for a given cumulative distribution function $F(x)$ is
404 $
405 $    D_n = \sup_x \left| F_n(x) - F(x) \right|
406 $
407 where $\sup_x$ is the supremum of the set of distances, and the empirical distribution
408 function $F_n(x)$ is discrete, and given by
409 $
410 $    F_n =  # of samples <= x / n
411 $
412 The empirical distribution function $F_n(x)$ is discrete, and thus had a ``stairstep''
413 cumulative distribution, making $n$ the number of stairs. Intuitively, the statistic takes
414 the largest absolute difference between the two distribution functions across all $x$ values.
415 
416 .seealso: PetscProbFunc
417 @*/
418 PetscErrorCode PetscProbComputeKSStatistic(Vec v, PetscProbFunc cdf, PetscReal *alpha)
419 {
420 #if !defined(PETSC_HAVE_KS)
421   SETERRQ(PetscObjectComm((PetscObject) v), PETSC_ERR_SUP, "No support for Kolmogorov-Smirnov test.\nReconfigure using --download-ks");
422 #else
423   PetscViewer        viewer = NULL;
424   PetscViewerFormat  format;
425   PetscDraw          draw;
426   PetscDrawLG        lg;
427   PetscDrawAxis      axis;
428   const PetscScalar *a;
429   PetscReal         *speed, Dn = PETSC_MIN_REAL;
430   PetscBool          isascii = PETSC_FALSE, isdraw = PETSC_FALSE, flg;
431   PetscInt           n, p, dim, d;
432   PetscMPIInt        size;
433   const char        *names[2] = {"Analytic", "Empirical"};
434   char               title[PETSC_MAX_PATH_LEN];
435   PetscOptions       options;
436   const char        *prefix;
437   MPI_Comm           comm;
438 
439   PetscFunctionBegin;
440   PetscCall(PetscObjectGetComm((PetscObject) v, &comm));
441   PetscCall(PetscObjectGetOptionsPrefix((PetscObject) v, &prefix));
442   PetscCall(PetscObjectGetOptions((PetscObject) v, &options));
443   PetscCall(PetscOptionsGetViewer(comm, options, prefix, "-ks_monitor", &viewer, &format, &flg));
444   if (flg) {
445     PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii));
446     PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERDRAW,  &isdraw));
447   }
448   if (isascii) {
449     PetscCall(PetscViewerPushFormat(viewer, format));
450   } else if (isdraw) {
451     PetscCall(PetscViewerPushFormat(viewer, format));
452     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
453     PetscCall(PetscDrawLGCreate(draw, 2, &lg));
454     PetscCall(PetscDrawLGSetLegend(lg, names));
455   }
456 
457   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) v), &size));
458   PetscCall(VecGetLocalSize(v, &n));
459   PetscCall(VecGetBlockSize(v, &dim));
460   n   /= dim;
461   /* TODO Parallel algorithm is harder */
462   if (size == 1) {
463     PetscCall(PetscMalloc1(n, &speed));
464     PetscCall(VecGetArrayRead(v, &a));
465     for (p = 0; p < n; ++p) {
466       PetscReal mag = 0.;
467 
468       for (d = 0; d < dim; ++d) mag += PetscSqr(PetscRealPart(a[p*dim+d]));
469       speed[p] = PetscSqrtReal(mag);
470     }
471     PetscCall(PetscSortReal(n, speed));
472     PetscCall(VecRestoreArrayRead(v, &a));
473     for (p = 0; p < n; ++p) {
474       const PetscReal x = speed[p], Fn = ((PetscReal) p) / n;
475       PetscReal       F, vals[2];
476 
477       PetscCall(cdf(&x, NULL, &F));
478       Dn = PetscMax(PetscAbsReal(Fn - F), Dn);
479       if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "x: %g F: %g Fn: %g Dn: %.2g\n", x, F, Fn, Dn));
480       if (isdraw)  {vals[0] = F; vals[1] = Fn; PetscCall(PetscDrawLGAddCommonPoint(lg, x, vals));}
481     }
482     if (speed[n-1] < 6.) {
483       const PetscReal k = (PetscInt) (6. - speed[n-1]) + 1, dx = (6. - speed[n-1])/k;
484       for (p = 0; p <= k; ++p) {
485         const PetscReal x = speed[n-1] + p*dx, Fn = 1.0;
486         PetscReal       F, vals[2];
487 
488         PetscCall(cdf(&x, NULL, &F));
489         Dn = PetscMax(PetscAbsReal(Fn - F), Dn);
490         if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "x: %g F: %g Fn: %g Dn: %.2g\n", x, F, Fn, Dn));
491         if (isdraw)  {vals[0] = F; vals[1] = Fn; PetscCall(PetscDrawLGAddCommonPoint(lg, x, vals));}
492       }
493     }
494     PetscCall(PetscFree(speed));
495   }
496   if (isdraw) {
497     PetscCall(PetscDrawLGGetAxis(lg, &axis));
498     PetscCall(PetscSNPrintf(title, PETSC_MAX_PATH_LEN, "Kolmogorov-Smirnov Test (Dn %.2g)", Dn));
499     PetscCall(PetscDrawAxisSetLabels(axis, title, "x", "CDF(x)"));
500     PetscCall(PetscDrawLGDraw(lg));
501     PetscCall(PetscDrawLGDestroy(&lg));
502   }
503   if (viewer) {
504     PetscCall(PetscViewerFlush(viewer));
505     PetscCall(PetscViewerPopFormat(viewer));
506     PetscCall(PetscViewerDestroy(&viewer));
507   }
508   *alpha = KSfbar((int) n, (double) Dn);
509   PetscFunctionReturn(0);
510 #endif
511 }
512