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