xref: /petsc/src/dm/dt/interface/dtprob.c (revision dce8aeba1c9b69b19f651c53d8a6b674bd7e9cbd)
1d6685f55SMatthew G. Knepley #include <petscdt.h> /*I "petscdt.h" I*/
2d6685f55SMatthew G. Knepley 
3d6685f55SMatthew G. Knepley #include <petscvec.h>
4d6685f55SMatthew G. Knepley #include <petscdraw.h>
5d6685f55SMatthew G. Knepley #include <petsc/private/petscimpl.h>
6d6685f55SMatthew G. Knepley 
7d6685f55SMatthew G. Knepley const char *const DTProbDensityTypes[] = {"constant", "gaussian", "maxwell_boltzmann", "DTProb Density", "DTPROB_DENSITY_", NULL};
8d6685f55SMatthew G. Knepley 
9d6685f55SMatthew G. Knepley /*@
10d6685f55SMatthew G. Knepley   PetscPDFMaxwellBoltzmann1D - PDF for the Maxwell-Boltzmann distribution in 1D
11d6685f55SMatthew G. Knepley 
12d6685f55SMatthew G. Knepley   Not collective
13d6685f55SMatthew G. Knepley 
14d6685f55SMatthew G. Knepley   Input Parameter:
15d6685f55SMatthew G. Knepley . x - Speed in [0, \infty]
16d6685f55SMatthew G. Knepley 
17d6685f55SMatthew G. Knepley   Output Parameter:
18d6685f55SMatthew G. Knepley . p - The probability density at x
19d6685f55SMatthew G. Knepley 
20d6685f55SMatthew G. Knepley   Level: beginner
21d6685f55SMatthew G. Knepley 
22*dce8aebaSBarry Smith .seealso: `PetscCDFMaxwellBoltzmann1D()`
23d6685f55SMatthew G. Knepley @*/
24d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
25d71ae5a4SJacob Faibussowitsch {
26d6685f55SMatthew G. Knepley   p[0] = PetscSqrtReal(2. / PETSC_PI) * PetscExpReal(-0.5 * PetscSqr(x[0]));
27d6685f55SMatthew G. Knepley   return 0;
28d6685f55SMatthew G. Knepley }
29d6685f55SMatthew G. Knepley 
30d6685f55SMatthew G. Knepley /*@
31d6685f55SMatthew G. Knepley   PetscCDFMaxwellBoltzmann1D - CDF for the Maxwell-Boltzmann distribution in 1D
32d6685f55SMatthew G. Knepley 
33d6685f55SMatthew G. Knepley   Not collective
34d6685f55SMatthew G. Knepley 
35d6685f55SMatthew G. Knepley   Input Parameter:
36d6685f55SMatthew G. Knepley . x - Speed in [0, \infty]
37d6685f55SMatthew G. Knepley 
38d6685f55SMatthew G. Knepley   Output Parameter:
39d6685f55SMatthew G. Knepley . p - The probability density at x
40d6685f55SMatthew G. Knepley 
41d6685f55SMatthew G. Knepley   Level: beginner
42d6685f55SMatthew G. Knepley 
43*dce8aebaSBarry Smith .seealso: `PetscPDFMaxwellBoltzmann1D()`
44d6685f55SMatthew G. Knepley @*/
45d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
46d71ae5a4SJacob Faibussowitsch {
47d6685f55SMatthew G. Knepley   p[0] = PetscErfReal(x[0] / PETSC_SQRT2);
48d6685f55SMatthew G. Knepley   return 0;
49d6685f55SMatthew G. Knepley }
50d6685f55SMatthew G. Knepley 
51d6685f55SMatthew G. Knepley /*@
52d6685f55SMatthew G. Knepley   PetscPDFMaxwellBoltzmann2D - PDF for the Maxwell-Boltzmann distribution in 2D
53d6685f55SMatthew G. Knepley 
54d6685f55SMatthew G. Knepley   Not collective
55d6685f55SMatthew G. Knepley 
56d6685f55SMatthew G. Knepley   Input Parameter:
57d6685f55SMatthew G. Knepley . x - Speed in [0, \infty]
58d6685f55SMatthew G. Knepley 
59d6685f55SMatthew G. Knepley   Output Parameter:
60d6685f55SMatthew G. Knepley . p - The probability density at x
61d6685f55SMatthew G. Knepley 
62d6685f55SMatthew G. Knepley   Level: beginner
63d6685f55SMatthew G. Knepley 
64*dce8aebaSBarry Smith .seealso: `PetscCDFMaxwellBoltzmann2D()`
65d6685f55SMatthew G. Knepley @*/
66d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
67d71ae5a4SJacob Faibussowitsch {
68d6685f55SMatthew G. Knepley   p[0] = x[0] * PetscExpReal(-0.5 * PetscSqr(x[0]));
69d6685f55SMatthew G. Knepley   return 0;
70d6685f55SMatthew G. Knepley }
71d6685f55SMatthew G. Knepley 
72d6685f55SMatthew G. Knepley /*@
73d6685f55SMatthew G. Knepley   PetscCDFMaxwellBoltzmann2D - CDF for the Maxwell-Boltzmann distribution in 2D
74d6685f55SMatthew G. Knepley 
75d6685f55SMatthew G. Knepley   Not collective
76d6685f55SMatthew G. Knepley 
77d6685f55SMatthew G. Knepley   Input Parameter:
78d6685f55SMatthew G. Knepley . x - Speed in [0, \infty]
79d6685f55SMatthew G. Knepley 
80d6685f55SMatthew G. Knepley   Output Parameter:
81d6685f55SMatthew G. Knepley . p - The probability density at x
82d6685f55SMatthew G. Knepley 
83d6685f55SMatthew G. Knepley   Level: beginner
84d6685f55SMatthew G. Knepley 
85*dce8aebaSBarry Smith .seealso: `PetscPDFMaxwellBoltzmann2D()`
86d6685f55SMatthew G. Knepley @*/
87d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
88d71ae5a4SJacob Faibussowitsch {
89d6685f55SMatthew G. Knepley   p[0] = 1. - PetscExpReal(-0.5 * PetscSqr(x[0]));
90d6685f55SMatthew G. Knepley   return 0;
91d6685f55SMatthew G. Knepley }
92d6685f55SMatthew G. Knepley 
93d6685f55SMatthew G. Knepley /*@
94d6685f55SMatthew G. Knepley   PetscPDFMaxwellBoltzmann3D - PDF for the Maxwell-Boltzmann distribution in 3D
95d6685f55SMatthew G. Knepley 
96d6685f55SMatthew G. Knepley   Not collective
97d6685f55SMatthew G. Knepley 
98d6685f55SMatthew G. Knepley   Input Parameter:
99d6685f55SMatthew G. Knepley . x - Speed in [0, \infty]
100d6685f55SMatthew G. Knepley 
101d6685f55SMatthew G. Knepley   Output Parameter:
102d6685f55SMatthew G. Knepley . p - The probability density at x
103d6685f55SMatthew G. Knepley 
104d6685f55SMatthew G. Knepley   Level: beginner
105d6685f55SMatthew G. Knepley 
106*dce8aebaSBarry Smith .seealso: `PetscCDFMaxwellBoltzmann3D()`
107d6685f55SMatthew G. Knepley @*/
108d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
109d71ae5a4SJacob Faibussowitsch {
110d6685f55SMatthew G. Knepley   p[0] = PetscSqrtReal(2. / PETSC_PI) * PetscSqr(x[0]) * PetscExpReal(-0.5 * PetscSqr(x[0]));
111d6685f55SMatthew G. Knepley   return 0;
112d6685f55SMatthew G. Knepley }
113d6685f55SMatthew G. Knepley 
114d6685f55SMatthew G. Knepley /*@
115d6685f55SMatthew G. Knepley   PetscCDFMaxwellBoltzmann3D - CDF for the Maxwell-Boltzmann distribution in 3D
116d6685f55SMatthew G. Knepley 
117d6685f55SMatthew G. Knepley   Not collective
118d6685f55SMatthew G. Knepley 
119d6685f55SMatthew G. Knepley   Input Parameter:
120d6685f55SMatthew G. Knepley . x - Speed in [0, \infty]
121d6685f55SMatthew G. Knepley 
122d6685f55SMatthew G. Knepley   Output Parameter:
123d6685f55SMatthew G. Knepley . p - The probability density at x
124d6685f55SMatthew G. Knepley 
125d6685f55SMatthew G. Knepley   Level: beginner
126d6685f55SMatthew G. Knepley 
127*dce8aebaSBarry Smith .seealso: `PetscPDFMaxwellBoltzmann3D()`
128d6685f55SMatthew G. Knepley @*/
129d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
130d71ae5a4SJacob Faibussowitsch {
131d6685f55SMatthew G. Knepley   p[0] = PetscErfReal(x[0] / PETSC_SQRT2) - PetscSqrtReal(2. / PETSC_PI) * x[0] * PetscExpReal(-0.5 * PetscSqr(x[0]));
132d6685f55SMatthew G. Knepley   return 0;
133d6685f55SMatthew G. Knepley }
134d6685f55SMatthew G. Knepley 
135d6685f55SMatthew G. Knepley /*@
136d6685f55SMatthew G. Knepley   PetscPDFGaussian1D - PDF for the Gaussian distribution in 1D
137d6685f55SMatthew G. Knepley 
138d6685f55SMatthew G. Knepley   Not collective
139d6685f55SMatthew G. Knepley 
140d6685f55SMatthew G. Knepley   Input Parameter:
141d6685f55SMatthew G. Knepley . x - Coordinate in [-\infty, \infty]
142d6685f55SMatthew G. Knepley 
143d6685f55SMatthew G. Knepley   Output Parameter:
144d6685f55SMatthew G. Knepley . p - The probability density at x
145d6685f55SMatthew G. Knepley 
146d6685f55SMatthew G. Knepley   Level: beginner
147d6685f55SMatthew G. Knepley 
148*dce8aebaSBarry Smith .seealso: `PetscPDFMaxwellBoltzmann3D()`
149d6685f55SMatthew G. Knepley @*/
150d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFGaussian1D(const PetscReal x[], const PetscReal scale[], PetscReal p[])
151d71ae5a4SJacob Faibussowitsch {
152d6685f55SMatthew G. Knepley   const PetscReal sigma = scale ? scale[0] : 1.;
153d6685f55SMatthew G. Knepley   p[0]                  = PetscSqrtReal(1. / (2. * PETSC_PI)) * PetscExpReal(-0.5 * PetscSqr(x[0] / sigma)) / sigma;
154d6685f55SMatthew G. Knepley   return 0;
155d6685f55SMatthew G. Knepley }
156d6685f55SMatthew G. Knepley 
157d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCDFGaussian1D(const PetscReal x[], const PetscReal scale[], PetscReal p[])
158d71ae5a4SJacob Faibussowitsch {
159d6685f55SMatthew G. Knepley   const PetscReal sigma = scale ? scale[0] : 1.;
160d6685f55SMatthew G. Knepley   p[0]                  = 0.5 * (1. + PetscErfReal(x[0] / PETSC_SQRT2 / sigma));
161d6685f55SMatthew G. Knepley   return 0;
162d6685f55SMatthew G. Knepley }
163d6685f55SMatthew G. Knepley 
164d6685f55SMatthew G. Knepley /*@
165d6685f55SMatthew G. Knepley   PetscPDFSampleGaussian1D - Sample uniformly from a Gaussian distribution in 1D
166d6685f55SMatthew G. Knepley 
167d6685f55SMatthew G. Knepley   Not collective
168d6685f55SMatthew G. Knepley 
169817da375SSatish Balay   Input Parameters:
170817da375SSatish Balay + p - A uniform variable on [0, 1]
171817da375SSatish Balay - dummy - ignored
172d6685f55SMatthew G. Knepley 
173d6685f55SMatthew G. Knepley   Output Parameter:
174d6685f55SMatthew G. Knepley . x - Coordinate in [-\infty, \infty]
175d6685f55SMatthew G. Knepley 
176d6685f55SMatthew G. Knepley   Level: beginner
177d6685f55SMatthew G. Knepley 
178d6685f55SMatthew G. Knepley   Note: http://www.mimirgames.com/articles/programming/approximations-of-the-inverse-error-function/
179d6685f55SMatthew G. Knepley   https://stackoverflow.com/questions/27229371/inverse-error-function-in-c
180d6685f55SMatthew G. Knepley @*/
181d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
182d71ae5a4SJacob Faibussowitsch {
183d6685f55SMatthew G. Knepley   const PetscReal q       = 2 * p[0] - 1.;
184d6685f55SMatthew G. Knepley   const PetscInt  maxIter = 100;
185d6685f55SMatthew G. Knepley   PetscReal       ck[100], r = 0.;
186d6685f55SMatthew G. Knepley   PetscInt        k, m;
187d6685f55SMatthew G. Knepley 
188d6685f55SMatthew G. Knepley   PetscFunctionBeginHot;
189d6685f55SMatthew G. Knepley   /* Transform input to [-1, 1] since the code below computes the inverse error function */
190d6685f55SMatthew G. Knepley   for (k = 0; k < maxIter; ++k) ck[k] = 0.;
191d6685f55SMatthew G. Knepley   ck[0] = 1;
192d6685f55SMatthew G. Knepley   r     = ck[0] * (PetscSqrtReal(PETSC_PI) / 2.) * q;
193d6685f55SMatthew G. Knepley   for (k = 1; k < maxIter; ++k) {
194d6685f55SMatthew G. Knepley     const PetscReal temp = 2. * k + 1.;
195d6685f55SMatthew G. Knepley 
196d6685f55SMatthew G. Knepley     for (m = 0; m <= k - 1; ++m) {
197d6685f55SMatthew G. Knepley       PetscReal denom = (m + 1.) * (2. * m + 1.);
198d6685f55SMatthew G. Knepley 
199d6685f55SMatthew G. Knepley       ck[k] += (ck[m] * ck[k - 1 - m]) / denom;
200d6685f55SMatthew G. Knepley     }
201d6685f55SMatthew G. Knepley     r += (ck[k] / temp) * PetscPowReal((PetscSqrtReal(PETSC_PI) / 2.) * q, 2. * k + 1.);
202d6685f55SMatthew G. Knepley   }
203d6685f55SMatthew G. Knepley   /* Scale erfinv() by \sqrt{\pi/2} */
204d6685f55SMatthew G. Knepley   x[0] = PetscSqrtReal(PETSC_PI * 0.5) * r;
205d6685f55SMatthew G. Knepley   PetscFunctionReturn(0);
206d6685f55SMatthew G. Knepley }
207d6685f55SMatthew G. Knepley 
208d6685f55SMatthew G. Knepley /*@
209d6685f55SMatthew G. Knepley   PetscPDFGaussian2D - PDF for the Gaussian distribution in 2D
210d6685f55SMatthew G. Knepley 
211d6685f55SMatthew G. Knepley   Not collective
212d6685f55SMatthew G. Knepley 
213817da375SSatish Balay   Input Parameters:
214817da375SSatish Balay + x - Coordinate in [-\infty, \infty]^2
215817da375SSatish Balay - dummy - ignored
216d6685f55SMatthew G. Knepley 
217d6685f55SMatthew G. Knepley   Output Parameter:
218d6685f55SMatthew G. Knepley . p - The probability density at x
219d6685f55SMatthew G. Knepley 
220d6685f55SMatthew G. Knepley   Level: beginner
221d6685f55SMatthew G. Knepley 
222db781477SPatrick Sanan .seealso: `PetscPDFSampleGaussian2D()`, `PetscPDFMaxwellBoltzmann3D()`
223d6685f55SMatthew G. Knepley @*/
224d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFGaussian2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
225d71ae5a4SJacob Faibussowitsch {
226d6685f55SMatthew G. Knepley   p[0] = (1. / PETSC_PI) * PetscExpReal(-0.5 * (PetscSqr(x[0]) + PetscSqr(x[1])));
227d6685f55SMatthew G. Knepley   return 0;
228d6685f55SMatthew G. Knepley }
229d6685f55SMatthew G. Knepley 
230d6685f55SMatthew G. Knepley /*@
231d6685f55SMatthew G. Knepley   PetscPDFSampleGaussian2D - Sample uniformly from a Gaussian distribution in 2D
232d6685f55SMatthew G. Knepley 
233d6685f55SMatthew G. Knepley   Not collective
234d6685f55SMatthew G. Knepley 
235817da375SSatish Balay   Input Parameters:
236817da375SSatish Balay + p - A uniform variable on [0, 1]^2
237817da375SSatish Balay - dummy - ignored
238d6685f55SMatthew G. Knepley 
239d6685f55SMatthew G. Knepley   Output Parameter:
240d6685f55SMatthew G. Knepley . x - Coordinate in [-\infty, \infty]^2
241d6685f55SMatthew G. Knepley 
242d6685f55SMatthew G. Knepley   Level: beginner
243d6685f55SMatthew G. Knepley 
244d6685f55SMatthew G. Knepley   Note: https://en.wikipedia.org/wiki/Box%E2%80%93Muller_transform
245d6685f55SMatthew G. Knepley 
246db781477SPatrick Sanan .seealso: `PetscPDFGaussian2D()`, `PetscPDFMaxwellBoltzmann3D()`
247d6685f55SMatthew G. Knepley @*/
248d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
249d71ae5a4SJacob Faibussowitsch {
250d6685f55SMatthew G. Knepley   const PetscReal mag = PetscSqrtReal(-2.0 * PetscLogReal(p[0]));
251d6685f55SMatthew G. Knepley   x[0]                = mag * PetscCosReal(2.0 * PETSC_PI * p[1]);
252d6685f55SMatthew G. Knepley   x[1]                = mag * PetscSinReal(2.0 * PETSC_PI * p[1]);
253d6685f55SMatthew G. Knepley   return 0;
254d6685f55SMatthew G. Knepley }
255d6685f55SMatthew G. Knepley 
256d6685f55SMatthew G. Knepley /*@
257ea1b28ebSMatthew G. Knepley   PetscPDFGaussian3D - PDF for the Gaussian distribution in 3D
258ea1b28ebSMatthew G. Knepley 
259ea1b28ebSMatthew G. Knepley   Not collective
260ea1b28ebSMatthew G. Knepley 
261ea1b28ebSMatthew G. Knepley   Input Parameters:
262ea1b28ebSMatthew G. Knepley + x - Coordinate in [-\infty, \infty]^3
263ea1b28ebSMatthew G. Knepley - dummy - ignored
264ea1b28ebSMatthew G. Knepley 
265ea1b28ebSMatthew G. Knepley   Output Parameter:
266ea1b28ebSMatthew G. Knepley . p - The probability density at x
267ea1b28ebSMatthew G. Knepley 
268ea1b28ebSMatthew G. Knepley   Level: beginner
269ea1b28ebSMatthew G. Knepley 
270db781477SPatrick Sanan .seealso: `PetscPDFSampleGaussian3D()`, `PetscPDFMaxwellBoltzmann3D()`
271ea1b28ebSMatthew G. Knepley @*/
272d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFGaussian3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
273d71ae5a4SJacob Faibussowitsch {
274ea1b28ebSMatthew G. Knepley   p[0] = (1. / PETSC_PI * PetscSqrtReal(PETSC_PI)) * PetscExpReal(-0.5 * (PetscSqr(x[0]) + PetscSqr(x[1]) + PetscSqr(x[2])));
275ea1b28ebSMatthew G. Knepley   return 0;
276ea1b28ebSMatthew G. Knepley }
277ea1b28ebSMatthew G. Knepley 
278ea1b28ebSMatthew G. Knepley /*@
279ea1b28ebSMatthew G. Knepley   PetscPDFSampleGaussian3D - Sample uniformly from a Gaussian distribution in 3D
280ea1b28ebSMatthew G. Knepley 
281ea1b28ebSMatthew G. Knepley   Not collective
282ea1b28ebSMatthew G. Knepley 
283ea1b28ebSMatthew G. Knepley   Input Parameters:
284ea1b28ebSMatthew G. Knepley + p - A uniform variable on [0, 1]^3
285ea1b28ebSMatthew G. Knepley - dummy - ignored
286ea1b28ebSMatthew G. Knepley 
287ea1b28ebSMatthew G. Knepley   Output Parameter:
288ea1b28ebSMatthew G. Knepley . x - Coordinate in [-\infty, \infty]^3
289ea1b28ebSMatthew G. Knepley 
290ea1b28ebSMatthew G. Knepley   Level: beginner
291ea1b28ebSMatthew G. Knepley 
292*dce8aebaSBarry Smith   Reference:
293*dce8aebaSBarry Smith   https://en.wikipedia.org/wiki/Box%E2%80%93Muller_transform
294ea1b28ebSMatthew G. Knepley 
295db781477SPatrick Sanan .seealso: `PetscPDFGaussian3D()`, `PetscPDFMaxwellBoltzmann3D()`
296ea1b28ebSMatthew G. Knepley @*/
297d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
298d71ae5a4SJacob Faibussowitsch {
299ea1b28ebSMatthew G. Knepley   PetscCall(PetscPDFSampleGaussian1D(p, dummy, x));
300ea1b28ebSMatthew G. Knepley   PetscCall(PetscPDFSampleGaussian2D(&p[1], dummy, &x[1]));
301ea1b28ebSMatthew G. Knepley   return 0;
302ea1b28ebSMatthew G. Knepley }
303ea1b28ebSMatthew G. Knepley 
304ea1b28ebSMatthew G. Knepley /*@
305d6685f55SMatthew G. Knepley   PetscPDFConstant1D - PDF for the uniform distribution in 1D
306d6685f55SMatthew G. Knepley 
307d6685f55SMatthew G. Knepley   Not collective
308d6685f55SMatthew G. Knepley 
309d6685f55SMatthew G. Knepley   Input Parameter:
310d6685f55SMatthew G. Knepley . x - Coordinate in [-1, 1]
311d6685f55SMatthew G. Knepley 
312d6685f55SMatthew G. Knepley   Output Parameter:
313d6685f55SMatthew G. Knepley . p - The probability density at x
314d6685f55SMatthew G. Knepley 
315d6685f55SMatthew G. Knepley   Level: beginner
316d6685f55SMatthew G. Knepley 
317db781477SPatrick Sanan .seealso: `PetscCDFConstant1D()`, `PetscPDFSampleConstant1D()`, `PetscPDFConstant2D()`, `PetscPDFConstant3D()`
318d6685f55SMatthew G. Knepley @*/
319d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFConstant1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
320d71ae5a4SJacob Faibussowitsch {
321d6685f55SMatthew G. Knepley   p[0] = x[0] > -1. && x[0] <= 1. ? 0.5 : 0.;
322d6685f55SMatthew G. Knepley   return 0;
323d6685f55SMatthew G. Knepley }
324d6685f55SMatthew G. Knepley 
325ea1b28ebSMatthew G. Knepley /*@
326ea1b28ebSMatthew G. Knepley   PetscCDFConstant1D - CDF for the uniform distribution in 1D
327ea1b28ebSMatthew G. Knepley 
328ea1b28ebSMatthew G. Knepley   Not collective
329ea1b28ebSMatthew G. Knepley 
330ea1b28ebSMatthew G. Knepley   Input Parameter:
331ea1b28ebSMatthew G. Knepley . x - Coordinate in [-1, 1]
332ea1b28ebSMatthew G. Knepley 
333ea1b28ebSMatthew G. Knepley   Output Parameter:
334ea1b28ebSMatthew G. Knepley . p - The cumulative probability at x
335ea1b28ebSMatthew G. Knepley 
336ea1b28ebSMatthew G. Knepley   Level: beginner
337ea1b28ebSMatthew G. Knepley 
338db781477SPatrick Sanan .seealso: `PetscPDFConstant1D()`, `PetscPDFSampleConstant1D()`, `PetscCDFConstant2D()`, `PetscCDFConstant3D()`
339ea1b28ebSMatthew G. Knepley @*/
340d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCDFConstant1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
341d71ae5a4SJacob Faibussowitsch {
342d6685f55SMatthew G. Knepley   p[0] = x[0] <= -1. ? 0. : (x[0] > 1. ? 1. : 0.5 * (x[0] + 1.));
343d6685f55SMatthew G. Knepley   return 0;
344d6685f55SMatthew G. Knepley }
345d6685f55SMatthew G. Knepley 
346d6685f55SMatthew G. Knepley /*@
347d6685f55SMatthew G. Knepley   PetscPDFSampleConstant1D - Sample uniformly from a uniform distribution on [-1, 1] in 1D
348d6685f55SMatthew G. Knepley 
349d6685f55SMatthew G. Knepley   Not collective
350d6685f55SMatthew G. Knepley 
351d6685f55SMatthew G. Knepley   Input Parameter:
352d6685f55SMatthew G. Knepley . p - A uniform variable on [0, 1]
353d6685f55SMatthew G. Knepley 
354d6685f55SMatthew G. Knepley   Output Parameter:
355d6685f55SMatthew G. Knepley . x - Coordinate in [-1, 1]
356d6685f55SMatthew G. Knepley 
357d6685f55SMatthew G. Knepley   Level: beginner
358d6685f55SMatthew G. Knepley 
359db781477SPatrick Sanan .seealso: `PetscPDFConstant1D()`, `PetscCDFConstant1D()`, `PetscPDFSampleConstant2D()`, `PetscPDFSampleConstant3D()`
360d6685f55SMatthew G. Knepley @*/
361d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFSampleConstant1D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
362d71ae5a4SJacob Faibussowitsch {
363ea1b28ebSMatthew G. Knepley   x[0] = 2. * p[0] - 1.;
364ea1b28ebSMatthew G. Knepley   return 0;
365ea1b28ebSMatthew G. Knepley }
366ea1b28ebSMatthew G. Knepley 
367ea1b28ebSMatthew G. Knepley /*@
368ea1b28ebSMatthew G. Knepley   PetscPDFConstant2D - PDF for the uniform distribution in 2D
369ea1b28ebSMatthew G. Knepley 
370ea1b28ebSMatthew G. Knepley   Not collective
371ea1b28ebSMatthew G. Knepley 
372ea1b28ebSMatthew G. Knepley   Input Parameter:
373ea1b28ebSMatthew G. Knepley . x - Coordinate in [-1, 1] x [-1, 1]
374ea1b28ebSMatthew G. Knepley 
375ea1b28ebSMatthew G. Knepley   Output Parameter:
376ea1b28ebSMatthew G. Knepley . p - The probability density at x
377ea1b28ebSMatthew G. Knepley 
378ea1b28ebSMatthew G. Knepley   Level: beginner
379ea1b28ebSMatthew G. Knepley 
380db781477SPatrick Sanan .seealso: `PetscCDFConstant2D()`, `PetscPDFSampleConstant2D()`, `PetscPDFConstant1D()`, `PetscPDFConstant3D()`
381ea1b28ebSMatthew G. Knepley @*/
382d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFConstant2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
383d71ae5a4SJacob Faibussowitsch {
384ea1b28ebSMatthew G. Knepley   p[0] = x[0] > -1. && x[0] <= 1. && x[1] > -1. && x[1] <= 1. ? 0.25 : 0.;
385ea1b28ebSMatthew G. Knepley   return 0;
386ea1b28ebSMatthew G. Knepley }
387ea1b28ebSMatthew G. Knepley 
388ea1b28ebSMatthew G. Knepley /*@
389ea1b28ebSMatthew G. Knepley   PetscCDFConstant2D - CDF for the uniform distribution in 2D
390ea1b28ebSMatthew G. Knepley 
391ea1b28ebSMatthew G. Knepley   Not collective
392ea1b28ebSMatthew G. Knepley 
393ea1b28ebSMatthew G. Knepley   Input Parameter:
394ea1b28ebSMatthew G. Knepley . x - Coordinate in [-1, 1] x [-1, 1]
395ea1b28ebSMatthew G. Knepley 
396ea1b28ebSMatthew G. Knepley   Output Parameter:
397ea1b28ebSMatthew G. Knepley . p - The cumulative probability at x
398ea1b28ebSMatthew G. Knepley 
399ea1b28ebSMatthew G. Knepley   Level: beginner
400ea1b28ebSMatthew G. Knepley 
401db781477SPatrick Sanan .seealso: `PetscPDFConstant2D()`, `PetscPDFSampleConstant2D()`, `PetscCDFConstant1D()`, `PetscCDFConstant3D()`
402ea1b28ebSMatthew G. Knepley @*/
403d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCDFConstant2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
404d71ae5a4SJacob Faibussowitsch {
405ea1b28ebSMatthew G. Knepley   p[0] = x[0] <= -1. || x[1] <= -1. ? 0. : (x[0] > 1. ? 1. : 0.5 * (x[0] + 1.)) * (x[1] > 1. ? 1. : 0.5 * (x[1] + 1.));
406ea1b28ebSMatthew G. Knepley   return 0;
407ea1b28ebSMatthew G. Knepley }
408ea1b28ebSMatthew G. Knepley 
409ea1b28ebSMatthew G. Knepley /*@
410ea1b28ebSMatthew G. Knepley   PetscPDFSampleConstant2D - Sample uniformly from a uniform distribution on [-1, 1] x [-1, 1] in 2D
411ea1b28ebSMatthew G. Knepley 
412ea1b28ebSMatthew G. Knepley   Not collective
413ea1b28ebSMatthew G. Knepley 
414ea1b28ebSMatthew G. Knepley   Input Parameter:
415ea1b28ebSMatthew G. Knepley . p - Two uniform variables on [0, 1]
416ea1b28ebSMatthew G. Knepley 
417ea1b28ebSMatthew G. Knepley   Output Parameter:
418ea1b28ebSMatthew G. Knepley . x - Coordinate in [-1, 1] x [-1, 1]
419ea1b28ebSMatthew G. Knepley 
420ea1b28ebSMatthew G. Knepley   Level: beginner
421ea1b28ebSMatthew G. Knepley 
422db781477SPatrick Sanan .seealso: `PetscPDFConstant2D()`, `PetscCDFConstant2D()`, `PetscPDFSampleConstant1D()`, `PetscPDFSampleConstant3D()`
423ea1b28ebSMatthew G. Knepley @*/
424d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFSampleConstant2D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
425d71ae5a4SJacob Faibussowitsch {
426ea1b28ebSMatthew G. Knepley   x[0] = 2. * p[0] - 1.;
427ea1b28ebSMatthew G. Knepley   x[1] = 2. * p[1] - 1.;
428ea1b28ebSMatthew G. Knepley   return 0;
429ea1b28ebSMatthew G. Knepley }
430ea1b28ebSMatthew G. Knepley 
431ea1b28ebSMatthew G. Knepley /*@
432ea1b28ebSMatthew G. Knepley   PetscPDFConstant3D - PDF for the uniform distribution in 3D
433ea1b28ebSMatthew G. Knepley 
434ea1b28ebSMatthew G. Knepley   Not collective
435ea1b28ebSMatthew G. Knepley 
436ea1b28ebSMatthew G. Knepley   Input Parameter:
437ea1b28ebSMatthew G. Knepley . x - Coordinate in [-1, 1] x [-1, 1] x [-1, 1]
438ea1b28ebSMatthew G. Knepley 
439ea1b28ebSMatthew G. Knepley   Output Parameter:
440ea1b28ebSMatthew G. Knepley . p - The probability density at x
441ea1b28ebSMatthew G. Knepley 
442ea1b28ebSMatthew G. Knepley   Level: beginner
443ea1b28ebSMatthew G. Knepley 
444db781477SPatrick Sanan .seealso: `PetscCDFConstant3D()`, `PetscPDFSampleConstant3D()`, `PetscPDFSampleConstant1D()`, `PetscPDFSampleConstant2D()`
445ea1b28ebSMatthew G. Knepley @*/
446d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFConstant3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
447d71ae5a4SJacob Faibussowitsch {
448ea1b28ebSMatthew G. Knepley   p[0] = x[0] > -1. && x[0] <= 1. && x[1] > -1. && x[1] <= 1. && x[2] > -1. && x[2] <= 1. ? 0.125 : 0.;
449ea1b28ebSMatthew G. Knepley   return 0;
450ea1b28ebSMatthew G. Knepley }
451ea1b28ebSMatthew G. Knepley 
452ea1b28ebSMatthew G. Knepley /*@
453ea1b28ebSMatthew G. Knepley   PetscCDFConstant3D - CDF for the uniform distribution in 3D
454ea1b28ebSMatthew G. Knepley 
455ea1b28ebSMatthew G. Knepley   Not collective
456ea1b28ebSMatthew G. Knepley 
457ea1b28ebSMatthew G. Knepley   Input Parameter:
458ea1b28ebSMatthew G. Knepley . x - Coordinate in [-1, 1] x [-1, 1] x [-1, 1]
459ea1b28ebSMatthew G. Knepley 
460ea1b28ebSMatthew G. Knepley   Output Parameter:
461ea1b28ebSMatthew G. Knepley . p - The cumulative probability at x
462ea1b28ebSMatthew G. Knepley 
463ea1b28ebSMatthew G. Knepley   Level: beginner
464ea1b28ebSMatthew G. Knepley 
465db781477SPatrick Sanan .seealso: `PetscPDFConstant3D()`, `PetscPDFSampleConstant3D()`, `PetscCDFConstant1D()`, `PetscCDFConstant2D()`
466ea1b28ebSMatthew G. Knepley @*/
467d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCDFConstant3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
468d71ae5a4SJacob Faibussowitsch {
469ea1b28ebSMatthew G. Knepley   p[0] = x[0] <= -1. || x[1] <= -1. || x[2] <= -1. ? 0. : (x[0] > 1. ? 1. : 0.5 * (x[0] + 1.)) * (x[1] > 1. ? 1. : 0.5 * (x[1] + 1.)) * (x[2] > 1. ? 1. : 0.5 * (x[2] + 1.));
470ea1b28ebSMatthew G. Knepley   return 0;
471ea1b28ebSMatthew G. Knepley }
472ea1b28ebSMatthew G. Knepley 
473ea1b28ebSMatthew G. Knepley /*@
474ea1b28ebSMatthew G. Knepley   PetscPDFSampleConstant3D - Sample uniformly from a uniform distribution on [-1, 1] x [-1, 1] in 3D
475ea1b28ebSMatthew G. Knepley 
476ea1b28ebSMatthew G. Knepley   Not collective
477ea1b28ebSMatthew G. Knepley 
478ea1b28ebSMatthew G. Knepley   Input Parameter:
479ea1b28ebSMatthew G. Knepley . p - Three uniform variables on [0, 1]
480ea1b28ebSMatthew G. Knepley 
481ea1b28ebSMatthew G. Knepley   Output Parameter:
482ea1b28ebSMatthew G. Knepley . x - Coordinate in [-1, 1] x [-1, 1] x [-1, 1]
483ea1b28ebSMatthew G. Knepley 
484ea1b28ebSMatthew G. Knepley   Level: beginner
485ea1b28ebSMatthew G. Knepley 
486db781477SPatrick Sanan .seealso: `PetscPDFConstant3D()`, `PetscCDFConstant3D()`, `PetscPDFSampleConstant1D()`, `PetscPDFSampleConstant2D()`
487ea1b28ebSMatthew G. Knepley @*/
488d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFSampleConstant3D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
489d71ae5a4SJacob Faibussowitsch {
490ea1b28ebSMatthew G. Knepley   x[0] = 2. * p[0] - 1.;
491ea1b28ebSMatthew G. Knepley   x[1] = 2. * p[1] - 1.;
492ea1b28ebSMatthew G. Knepley   x[2] = 2. * p[2] - 1.;
493d6685f55SMatthew G. Knepley   return 0;
494d6685f55SMatthew G. Knepley }
495d6685f55SMatthew G. Knepley 
496d6685f55SMatthew G. Knepley /*@C
497d6685f55SMatthew G. Knepley   PetscProbCreateFromOptions - Return the probability distribution specified by the argumetns and options
498d6685f55SMatthew G. Knepley 
499d6685f55SMatthew G. Knepley   Not collective
500d6685f55SMatthew G. Knepley 
501d6685f55SMatthew G. Knepley   Input Parameters:
502d6685f55SMatthew G. Knepley + dim    - The dimension of sample points
503d6685f55SMatthew G. Knepley . prefix - The options prefix, or NULL
504d6685f55SMatthew G. Knepley - name   - The option name for the probility distribution type
505d6685f55SMatthew G. Knepley 
506d6685f55SMatthew G. Knepley   Output Parameters:
507d6685f55SMatthew G. Knepley + pdf - The PDF of this type
508d6685f55SMatthew G. Knepley . cdf - The CDF of this type
509d6685f55SMatthew G. Knepley - sampler - The PDF sampler of this type
510d6685f55SMatthew G. Knepley 
511d6685f55SMatthew G. Knepley   Level: intermediate
512d6685f55SMatthew G. Knepley 
513*dce8aebaSBarry Smith .seealso: `PetscProbFunc`, `PetscPDFMaxwellBoltzmann1D()`, `PetscPDFGaussian1D()`, `PetscPDFConstant1D()`
514d6685f55SMatthew G. Knepley @*/
515d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscProbCreateFromOptions(PetscInt dim, const char prefix[], const char name[], PetscProbFunc *pdf, PetscProbFunc *cdf, PetscProbFunc *sampler)
516d71ae5a4SJacob Faibussowitsch {
517d6685f55SMatthew G. Knepley   DTProbDensityType den = DTPROB_DENSITY_GAUSSIAN;
518d6685f55SMatthew G. Knepley 
519d6685f55SMatthew G. Knepley   PetscFunctionBegin;
520d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF, prefix, "PetscProb Options", "DT");
5219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum(name, "Method to compute PDF <constant, gaussian>", "", DTProbDensityTypes, (PetscEnum)den, (PetscEnum *)&den, NULL));
522d0609cedSBarry Smith   PetscOptionsEnd();
523d6685f55SMatthew G. Knepley 
5249371c9d4SSatish Balay   if (pdf) {
5259371c9d4SSatish Balay     PetscValidPointer(pdf, 4);
5269371c9d4SSatish Balay     *pdf = NULL;
5279371c9d4SSatish Balay   }
5289371c9d4SSatish Balay   if (cdf) {
5299371c9d4SSatish Balay     PetscValidPointer(cdf, 5);
5309371c9d4SSatish Balay     *cdf = NULL;
5319371c9d4SSatish Balay   }
5329371c9d4SSatish Balay   if (sampler) {
5339371c9d4SSatish Balay     PetscValidPointer(sampler, 6);
5349371c9d4SSatish Balay     *sampler = NULL;
5359371c9d4SSatish Balay   }
536d6685f55SMatthew G. Knepley   switch (den) {
537d6685f55SMatthew G. Knepley   case DTPROB_DENSITY_CONSTANT:
538d6685f55SMatthew G. Knepley     switch (dim) {
539d6685f55SMatthew G. Knepley     case 1:
540d6685f55SMatthew G. Knepley       if (pdf) *pdf = PetscPDFConstant1D;
541d6685f55SMatthew G. Knepley       if (cdf) *cdf = PetscCDFConstant1D;
542d6685f55SMatthew G. Knepley       if (sampler) *sampler = PetscPDFSampleConstant1D;
543d6685f55SMatthew G. Knepley       break;
544ea1b28ebSMatthew G. Knepley     case 2:
545ea1b28ebSMatthew G. Knepley       if (pdf) *pdf = PetscPDFConstant2D;
546ea1b28ebSMatthew G. Knepley       if (cdf) *cdf = PetscCDFConstant2D;
547ea1b28ebSMatthew G. Knepley       if (sampler) *sampler = PetscPDFSampleConstant2D;
548ea1b28ebSMatthew G. Knepley       break;
549ea1b28ebSMatthew G. Knepley     case 3:
550ea1b28ebSMatthew G. Knepley       if (pdf) *pdf = PetscPDFConstant3D;
551ea1b28ebSMatthew G. Knepley       if (cdf) *cdf = PetscCDFConstant3D;
552ea1b28ebSMatthew G. Knepley       if (sampler) *sampler = PetscPDFSampleConstant3D;
553ea1b28ebSMatthew G. Knepley       break;
554d71ae5a4SJacob Faibussowitsch     default:
555d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]);
556d6685f55SMatthew G. Knepley     }
557d6685f55SMatthew G. Knepley     break;
558d6685f55SMatthew G. Knepley   case DTPROB_DENSITY_GAUSSIAN:
559d6685f55SMatthew G. Knepley     switch (dim) {
560d6685f55SMatthew G. Knepley     case 1:
561d6685f55SMatthew G. Knepley       if (pdf) *pdf = PetscPDFGaussian1D;
562d6685f55SMatthew G. Knepley       if (cdf) *cdf = PetscCDFGaussian1D;
563d6685f55SMatthew G. Knepley       if (sampler) *sampler = PetscPDFSampleGaussian1D;
564d6685f55SMatthew G. Knepley       break;
565d6685f55SMatthew G. Knepley     case 2:
566d6685f55SMatthew G. Knepley       if (pdf) *pdf = PetscPDFGaussian2D;
567d6685f55SMatthew G. Knepley       if (sampler) *sampler = PetscPDFSampleGaussian2D;
568d6685f55SMatthew G. Knepley       break;
569ea1b28ebSMatthew G. Knepley     case 3:
570ea1b28ebSMatthew G. Knepley       if (pdf) *pdf = PetscPDFGaussian3D;
571ea1b28ebSMatthew G. Knepley       if (sampler) *sampler = PetscPDFSampleGaussian3D;
572ea1b28ebSMatthew G. Knepley       break;
573d71ae5a4SJacob Faibussowitsch     default:
574d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]);
575d6685f55SMatthew G. Knepley     }
576d6685f55SMatthew G. Knepley     break;
577d6685f55SMatthew G. Knepley   case DTPROB_DENSITY_MAXWELL_BOLTZMANN:
578d6685f55SMatthew G. Knepley     switch (dim) {
579d6685f55SMatthew G. Knepley     case 1:
580d6685f55SMatthew G. Knepley       if (pdf) *pdf = PetscPDFMaxwellBoltzmann1D;
581d6685f55SMatthew G. Knepley       if (cdf) *cdf = PetscCDFMaxwellBoltzmann1D;
582d6685f55SMatthew G. Knepley       break;
583d6685f55SMatthew G. Knepley     case 2:
584d6685f55SMatthew G. Knepley       if (pdf) *pdf = PetscPDFMaxwellBoltzmann2D;
585d6685f55SMatthew G. Knepley       if (cdf) *cdf = PetscCDFMaxwellBoltzmann2D;
586d6685f55SMatthew G. Knepley       break;
587d6685f55SMatthew G. Knepley     case 3:
588d6685f55SMatthew G. Knepley       if (pdf) *pdf = PetscPDFMaxwellBoltzmann3D;
589d6685f55SMatthew G. Knepley       if (cdf) *cdf = PetscCDFMaxwellBoltzmann3D;
590d6685f55SMatthew G. Knepley       break;
591d71ae5a4SJacob Faibussowitsch     default:
592d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]);
593d6685f55SMatthew G. Knepley     }
594d6685f55SMatthew G. Knepley     break;
595d71ae5a4SJacob Faibussowitsch   default:
596d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Density type %s is not supported", DTProbDensityTypes[PetscMax(0, PetscMin(den, DTPROB_NUM_DENSITY))]);
597d6685f55SMatthew G. Knepley   }
598d6685f55SMatthew G. Knepley   PetscFunctionReturn(0);
599d6685f55SMatthew G. Knepley }
600d6685f55SMatthew G. Knepley 
601d6685f55SMatthew G. Knepley #ifdef PETSC_HAVE_KS
602d6685f55SMatthew G. Knepley EXTERN_C_BEGIN
603d6685f55SMatthew G. Knepley   #include <KolmogorovSmirnovDist.h>
604d6685f55SMatthew G. Knepley EXTERN_C_END
605d6685f55SMatthew G. Knepley #endif
606d6685f55SMatthew G. Knepley 
607d6685f55SMatthew G. Knepley /*@C
608d6685f55SMatthew G. Knepley   PetscProbComputeKSStatistic - Compute the Kolmogorov-Smirnov statistic for the empirical distribution for an input vector, compared to an analytic CDF.
609d6685f55SMatthew G. Knepley 
610d6685f55SMatthew G. Knepley   Collective on v
611d6685f55SMatthew G. Knepley 
612d6685f55SMatthew G. Knepley   Input Parameters:
613d6685f55SMatthew G. Knepley + v   - The data vector, blocksize is the sample dimension
614d6685f55SMatthew G. Knepley - cdf - The analytic CDF
615d6685f55SMatthew G. Knepley 
616d6685f55SMatthew G. Knepley   Output Parameter:
617d6685f55SMatthew G. Knepley . alpha - The KS statisic
618d6685f55SMatthew G. Knepley 
619d6685f55SMatthew G. Knepley   Level: advanced
620d6685f55SMatthew G. Knepley 
621*dce8aebaSBarry Smith   Notes:
622*dce8aebaSBarry Smith   The Kolmogorov-Smirnov statistic for a given cumulative distribution function $F(x)$ is
623*dce8aebaSBarry Smith .vb
624*dce8aebaSBarry Smith     D_n = \sup_x \left| F_n(x) - F(x) \right|
625*dce8aebaSBarry Smith 
626*dce8aebaSBarry Smith   where $\sup_x$ is the supremum of the set of distances, and the empirical distribution function $F_n(x)$ is discrete, and given by
627*dce8aebaSBarry Smith 
628*dce8aebaSBarry Smith     F_n =  # of samples <= x / n
629*dce8aebaSBarry Smith .ve
630*dce8aebaSBarry Smith 
631d6685f55SMatthew G. Knepley   The empirical distribution function $F_n(x)$ is discrete, and thus had a ``stairstep''
632d6685f55SMatthew G. Knepley   cumulative distribution, making $n$ the number of stairs. Intuitively, the statistic takes
633d6685f55SMatthew G. Knepley   the largest absolute difference between the two distribution functions across all $x$ values.
634d6685f55SMatthew G. Knepley 
635db781477SPatrick Sanan .seealso: `PetscProbFunc`
636d6685f55SMatthew G. Knepley @*/
637d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscProbComputeKSStatistic(Vec v, PetscProbFunc cdf, PetscReal *alpha)
638d71ae5a4SJacob Faibussowitsch {
639d6685f55SMatthew G. Knepley #if !defined(PETSC_HAVE_KS)
640d6685f55SMatthew G. Knepley   SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "No support for Kolmogorov-Smirnov test.\nReconfigure using --download-ks");
641d6685f55SMatthew G. Knepley #else
642d6685f55SMatthew G. Knepley   PetscViewer        viewer = NULL;
643d6685f55SMatthew G. Knepley   PetscViewerFormat  format;
644d6685f55SMatthew G. Knepley   PetscDraw          draw;
645d6685f55SMatthew G. Knepley   PetscDrawLG        lg;
646d6685f55SMatthew G. Knepley   PetscDrawAxis      axis;
647d6685f55SMatthew G. Knepley   const PetscScalar *a;
648d6685f55SMatthew G. Knepley   PetscReal         *speed, Dn = PETSC_MIN_REAL;
649d6685f55SMatthew G. Knepley   PetscBool          isascii = PETSC_FALSE, isdraw = PETSC_FALSE, flg;
650d6685f55SMatthew G. Knepley   PetscInt           n, p, dim, d;
651d6685f55SMatthew G. Knepley   PetscMPIInt        size;
652d6685f55SMatthew G. Knepley   const char        *names[2] = {"Analytic", "Empirical"};
653d6685f55SMatthew G. Knepley   char               title[PETSC_MAX_PATH_LEN];
654d6685f55SMatthew G. Knepley   PetscOptions       options;
655d6685f55SMatthew G. Knepley   const char        *prefix;
656d6685f55SMatthew G. Knepley   MPI_Comm           comm;
657d6685f55SMatthew G. Knepley 
658d6685f55SMatthew G. Knepley   PetscFunctionBegin;
6599566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
6609566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)v, &prefix));
6619566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptions((PetscObject)v, &options));
6629566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(comm, options, prefix, "-ks_monitor", &viewer, &format, &flg));
663d6685f55SMatthew G. Knepley   if (flg) {
6649566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
6659566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
666d6685f55SMatthew G. Knepley   }
667d6685f55SMatthew G. Knepley   if (isascii) {
6689566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
669d6685f55SMatthew G. Knepley   } else if (isdraw) {
6709566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
6719566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
6729566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGCreate(draw, 2, &lg));
6739566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSetLegend(lg, names));
674d6685f55SMatthew G. Knepley   }
675d6685f55SMatthew G. Knepley 
6769566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size));
6779566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(v, &n));
6789566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(v, &dim));
679d6685f55SMatthew G. Knepley   n /= dim;
680d6685f55SMatthew G. Knepley   /* TODO Parallel algorithm is harder */
681d6685f55SMatthew G. Knepley   if (size == 1) {
6829566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &speed));
6839566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(v, &a));
684d6685f55SMatthew G. Knepley     for (p = 0; p < n; ++p) {
685d6685f55SMatthew G. Knepley       PetscReal mag = 0.;
686d6685f55SMatthew G. Knepley 
687d6685f55SMatthew G. Knepley       for (d = 0; d < dim; ++d) mag += PetscSqr(PetscRealPart(a[p * dim + d]));
688d6685f55SMatthew G. Knepley       speed[p] = PetscSqrtReal(mag);
689d6685f55SMatthew G. Knepley     }
6909566063dSJacob Faibussowitsch     PetscCall(PetscSortReal(n, speed));
6919566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(v, &a));
692d6685f55SMatthew G. Knepley     for (p = 0; p < n; ++p) {
693d6685f55SMatthew G. Knepley       const PetscReal x = speed[p], Fn = ((PetscReal)p) / n;
694d6685f55SMatthew G. Knepley       PetscReal       F, vals[2];
695d6685f55SMatthew G. Knepley 
6969566063dSJacob Faibussowitsch       PetscCall(cdf(&x, NULL, &F));
697d6685f55SMatthew G. Knepley       Dn = PetscMax(PetscAbsReal(Fn - F), Dn);
6989566063dSJacob Faibussowitsch       if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "x: %g F: %g Fn: %g Dn: %.2g\n", x, F, Fn, Dn));
6999371c9d4SSatish Balay       if (isdraw) {
7009371c9d4SSatish Balay         vals[0] = F;
7019371c9d4SSatish Balay         vals[1] = Fn;
7029371c9d4SSatish Balay         PetscCall(PetscDrawLGAddCommonPoint(lg, x, vals));
7039371c9d4SSatish Balay       }
704d6685f55SMatthew G. Knepley     }
705d6685f55SMatthew G. Knepley     if (speed[n - 1] < 6.) {
706d6685f55SMatthew G. Knepley       const PetscReal k = (PetscInt)(6. - speed[n - 1]) + 1, dx = (6. - speed[n - 1]) / k;
707d6685f55SMatthew G. Knepley       for (p = 0; p <= k; ++p) {
708d6685f55SMatthew G. Knepley         const PetscReal x = speed[n - 1] + p * dx, Fn = 1.0;
709d6685f55SMatthew G. Knepley         PetscReal       F, vals[2];
710d6685f55SMatthew G. Knepley 
7119566063dSJacob Faibussowitsch         PetscCall(cdf(&x, NULL, &F));
712d6685f55SMatthew G. Knepley         Dn = PetscMax(PetscAbsReal(Fn - F), Dn);
7139566063dSJacob Faibussowitsch         if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "x: %g F: %g Fn: %g Dn: %.2g\n", x, F, Fn, Dn));
7149371c9d4SSatish Balay         if (isdraw) {
7159371c9d4SSatish Balay           vals[0] = F;
7169371c9d4SSatish Balay           vals[1] = Fn;
7179371c9d4SSatish Balay           PetscCall(PetscDrawLGAddCommonPoint(lg, x, vals));
7189371c9d4SSatish Balay         }
719d6685f55SMatthew G. Knepley       }
720d6685f55SMatthew G. Knepley     }
7219566063dSJacob Faibussowitsch     PetscCall(PetscFree(speed));
722d6685f55SMatthew G. Knepley   }
723d6685f55SMatthew G. Knepley   if (isdraw) {
7249566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(lg, &axis));
7259566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(title, PETSC_MAX_PATH_LEN, "Kolmogorov-Smirnov Test (Dn %.2g)", Dn));
7269566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, title, "x", "CDF(x)"));
7279566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(lg));
7289566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDestroy(&lg));
729d6685f55SMatthew G. Knepley   }
730d6685f55SMatthew G. Knepley   if (viewer) {
7319566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
7329566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewer));
7339566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
734d6685f55SMatthew G. Knepley   }
735d6685f55SMatthew G. Knepley   *alpha = KSfbar((int)n, (double)Dn);
736d6685f55SMatthew G. Knepley   PetscFunctionReturn(0);
737d6685f55SMatthew G. Knepley #endif
738d6685f55SMatthew G. Knepley }
739