xref: /petsc/src/dm/dt/interface/dtprob.c (revision a4e35b1925eceef64945ea472b84f2bf06a67b5e)
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 
1220f4b53cSBarry Smith   Not Collective
13d6685f55SMatthew G. Knepley 
14*a4e35b19SJacob Faibussowitsch   Input Parameters:
15*a4e35b19SJacob Faibussowitsch + x     - Speed in [0, \infty]
16*a4e35b19SJacob Faibussowitsch - dummy - Unused
17d6685f55SMatthew G. Knepley 
18d6685f55SMatthew G. Knepley   Output Parameter:
19d6685f55SMatthew G. Knepley . p - The probability density at x
20d6685f55SMatthew G. Knepley 
21d6685f55SMatthew G. Knepley   Level: beginner
22d6685f55SMatthew G. Knepley 
23dce8aebaSBarry Smith .seealso: `PetscCDFMaxwellBoltzmann1D()`
24d6685f55SMatthew G. Knepley @*/
25d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
26d71ae5a4SJacob Faibussowitsch {
27d6685f55SMatthew G. Knepley   p[0] = PetscSqrtReal(2. / PETSC_PI) * PetscExpReal(-0.5 * PetscSqr(x[0]));
283ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
29d6685f55SMatthew G. Knepley }
30d6685f55SMatthew G. Knepley 
31d6685f55SMatthew G. Knepley /*@
32d6685f55SMatthew G. Knepley   PetscCDFMaxwellBoltzmann1D - CDF for the Maxwell-Boltzmann distribution in 1D
33d6685f55SMatthew G. Knepley 
3420f4b53cSBarry Smith   Not Collective
35d6685f55SMatthew G. Knepley 
36*a4e35b19SJacob Faibussowitsch   Input Parameters:
37*a4e35b19SJacob Faibussowitsch + x     - Speed in [0, \infty]
38*a4e35b19SJacob Faibussowitsch - dummy - Unused
39d6685f55SMatthew G. Knepley 
40d6685f55SMatthew G. Knepley   Output Parameter:
41d6685f55SMatthew G. Knepley . p - The probability density at x
42d6685f55SMatthew G. Knepley 
43d6685f55SMatthew G. Knepley   Level: beginner
44d6685f55SMatthew G. Knepley 
45dce8aebaSBarry Smith .seealso: `PetscPDFMaxwellBoltzmann1D()`
46d6685f55SMatthew G. Knepley @*/
47d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
48d71ae5a4SJacob Faibussowitsch {
49d6685f55SMatthew G. Knepley   p[0] = PetscErfReal(x[0] / PETSC_SQRT2);
503ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
51d6685f55SMatthew G. Knepley }
52d6685f55SMatthew G. Knepley 
53d6685f55SMatthew G. Knepley /*@
54d6685f55SMatthew G. Knepley   PetscPDFMaxwellBoltzmann2D - PDF for the Maxwell-Boltzmann distribution in 2D
55d6685f55SMatthew G. Knepley 
5620f4b53cSBarry Smith   Not Collective
57d6685f55SMatthew G. Knepley 
58*a4e35b19SJacob Faibussowitsch   Input Parameters:
59*a4e35b19SJacob Faibussowitsch + x     - Speed in [0, \infty]
60*a4e35b19SJacob Faibussowitsch - dummy - Unused
61d6685f55SMatthew G. Knepley 
62d6685f55SMatthew G. Knepley   Output Parameter:
63d6685f55SMatthew G. Knepley . p - The probability density at x
64d6685f55SMatthew G. Knepley 
65d6685f55SMatthew G. Knepley   Level: beginner
66d6685f55SMatthew G. Knepley 
67dce8aebaSBarry Smith .seealso: `PetscCDFMaxwellBoltzmann2D()`
68d6685f55SMatthew G. Knepley @*/
69d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
70d71ae5a4SJacob Faibussowitsch {
71d6685f55SMatthew G. Knepley   p[0] = x[0] * PetscExpReal(-0.5 * PetscSqr(x[0]));
723ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
73d6685f55SMatthew G. Knepley }
74d6685f55SMatthew G. Knepley 
75d6685f55SMatthew G. Knepley /*@
76d6685f55SMatthew G. Knepley   PetscCDFMaxwellBoltzmann2D - CDF for the Maxwell-Boltzmann distribution in 2D
77d6685f55SMatthew G. Knepley 
7820f4b53cSBarry Smith   Not Collective
79d6685f55SMatthew G. Knepley 
80*a4e35b19SJacob Faibussowitsch   Input Parameters:
81*a4e35b19SJacob Faibussowitsch + x     - Speed in [0, \infty]
82*a4e35b19SJacob Faibussowitsch - dummy - Unused
83d6685f55SMatthew G. Knepley 
84d6685f55SMatthew G. Knepley   Output Parameter:
85d6685f55SMatthew G. Knepley . p - The probability density at x
86d6685f55SMatthew G. Knepley 
87d6685f55SMatthew G. Knepley   Level: beginner
88d6685f55SMatthew G. Knepley 
89dce8aebaSBarry Smith .seealso: `PetscPDFMaxwellBoltzmann2D()`
90d6685f55SMatthew G. Knepley @*/
91d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
92d71ae5a4SJacob Faibussowitsch {
93d6685f55SMatthew G. Knepley   p[0] = 1. - PetscExpReal(-0.5 * PetscSqr(x[0]));
943ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
95d6685f55SMatthew G. Knepley }
96d6685f55SMatthew G. Knepley 
97d6685f55SMatthew G. Knepley /*@
98d6685f55SMatthew G. Knepley   PetscPDFMaxwellBoltzmann3D - PDF for the Maxwell-Boltzmann distribution in 3D
99d6685f55SMatthew G. Knepley 
10020f4b53cSBarry Smith   Not Collective
101d6685f55SMatthew G. Knepley 
102*a4e35b19SJacob Faibussowitsch   Input Parameters:
103*a4e35b19SJacob Faibussowitsch + x     - Speed in [0, \infty]
104*a4e35b19SJacob Faibussowitsch - dummy - Unused
105d6685f55SMatthew G. Knepley 
106d6685f55SMatthew G. Knepley   Output Parameter:
107d6685f55SMatthew G. Knepley . p - The probability density at x
108d6685f55SMatthew G. Knepley 
109d6685f55SMatthew G. Knepley   Level: beginner
110d6685f55SMatthew G. Knepley 
111dce8aebaSBarry Smith .seealso: `PetscCDFMaxwellBoltzmann3D()`
112d6685f55SMatthew G. Knepley @*/
113d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
114d71ae5a4SJacob Faibussowitsch {
115d6685f55SMatthew G. Knepley   p[0] = PetscSqrtReal(2. / PETSC_PI) * PetscSqr(x[0]) * PetscExpReal(-0.5 * PetscSqr(x[0]));
1163ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
117d6685f55SMatthew G. Knepley }
118d6685f55SMatthew G. Knepley 
119d6685f55SMatthew G. Knepley /*@
120d6685f55SMatthew G. Knepley   PetscCDFMaxwellBoltzmann3D - CDF for the Maxwell-Boltzmann distribution in 3D
121d6685f55SMatthew G. Knepley 
12220f4b53cSBarry Smith   Not Collective
123d6685f55SMatthew G. Knepley 
124*a4e35b19SJacob Faibussowitsch   Input Parameters:
125*a4e35b19SJacob Faibussowitsch + x     - Speed in [0, \infty]
126*a4e35b19SJacob Faibussowitsch - dummy - Unused
127d6685f55SMatthew G. Knepley 
128d6685f55SMatthew G. Knepley   Output Parameter:
129d6685f55SMatthew G. Knepley . p - The probability density at x
130d6685f55SMatthew G. Knepley 
131d6685f55SMatthew G. Knepley   Level: beginner
132d6685f55SMatthew G. Knepley 
133dce8aebaSBarry Smith .seealso: `PetscPDFMaxwellBoltzmann3D()`
134d6685f55SMatthew G. Knepley @*/
135d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
136d71ae5a4SJacob Faibussowitsch {
137d6685f55SMatthew G. Knepley   p[0] = PetscErfReal(x[0] / PETSC_SQRT2) - PetscSqrtReal(2. / PETSC_PI) * x[0] * PetscExpReal(-0.5 * PetscSqr(x[0]));
1383ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
139d6685f55SMatthew G. Knepley }
140d6685f55SMatthew G. Knepley 
141d6685f55SMatthew G. Knepley /*@
142d6685f55SMatthew G. Knepley   PetscPDFGaussian1D - PDF for the Gaussian distribution in 1D
143d6685f55SMatthew G. Knepley 
14420f4b53cSBarry Smith   Not Collective
145d6685f55SMatthew G. Knepley 
146*a4e35b19SJacob Faibussowitsch   Input Parameters:
147*a4e35b19SJacob Faibussowitsch + x     - Coordinate in [-\infty, \infty]
148*a4e35b19SJacob Faibussowitsch - scale - Scaling value
149d6685f55SMatthew G. Knepley 
150d6685f55SMatthew G. Knepley   Output Parameter:
151d6685f55SMatthew G. Knepley . p - The probability density at x
152d6685f55SMatthew G. Knepley 
153d6685f55SMatthew G. Knepley   Level: beginner
154d6685f55SMatthew G. Knepley 
155dce8aebaSBarry Smith .seealso: `PetscPDFMaxwellBoltzmann3D()`
156d6685f55SMatthew G. Knepley @*/
157d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFGaussian1D(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]                  = PetscSqrtReal(1. / (2. * PETSC_PI)) * PetscExpReal(-0.5 * PetscSqr(x[0] / sigma)) / sigma;
1613ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
162d6685f55SMatthew G. Knepley }
163d6685f55SMatthew G. Knepley 
164d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCDFGaussian1D(const PetscReal x[], const PetscReal scale[], PetscReal p[])
165d71ae5a4SJacob Faibussowitsch {
166d6685f55SMatthew G. Knepley   const PetscReal sigma = scale ? scale[0] : 1.;
167d6685f55SMatthew G. Knepley   p[0]                  = 0.5 * (1. + PetscErfReal(x[0] / PETSC_SQRT2 / sigma));
1683ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
169d6685f55SMatthew G. Knepley }
170d6685f55SMatthew G. Knepley 
171d6685f55SMatthew G. Knepley /*@
172d6685f55SMatthew G. Knepley   PetscPDFSampleGaussian1D - Sample uniformly from a Gaussian distribution in 1D
173d6685f55SMatthew G. Knepley 
17420f4b53cSBarry Smith   Not Collective
175d6685f55SMatthew G. Knepley 
176817da375SSatish Balay   Input Parameters:
177817da375SSatish Balay + p     - A uniform variable on [0, 1]
178817da375SSatish Balay - dummy - ignored
179d6685f55SMatthew G. Knepley 
180d6685f55SMatthew G. Knepley   Output Parameter:
181d6685f55SMatthew G. Knepley . x - Coordinate in [-\infty, \infty]
182d6685f55SMatthew G. Knepley 
183d6685f55SMatthew G. Knepley   Level: beginner
184d6685f55SMatthew G. Knepley 
185*a4e35b19SJacob Faibussowitsch   References:
186*a4e35b19SJacob Faibussowitsch + * - http://www.mimirgames.com/articles/programming/approximations-of-the-inverse-error-function/
187*a4e35b19SJacob Faibussowitsch - * - https://stackoverflow.com/questions/27229371/inverse-error-function-in-c
188*a4e35b19SJacob Faibussowitsch 
189*a4e35b19SJacob Faibussowitsch .seealso: `PetscPDFGaussian2D()`
190d6685f55SMatthew G. Knepley @*/
191d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
192d71ae5a4SJacob Faibussowitsch {
193d6685f55SMatthew G. Knepley   const PetscReal q       = 2 * p[0] - 1.;
194d6685f55SMatthew G. Knepley   const PetscInt  maxIter = 100;
195d6685f55SMatthew G. Knepley   PetscReal       ck[100], r = 0.;
196d6685f55SMatthew G. Knepley   PetscInt        k, m;
197d6685f55SMatthew G. Knepley 
198d6685f55SMatthew G. Knepley   PetscFunctionBeginHot;
199d6685f55SMatthew G. Knepley   /* Transform input to [-1, 1] since the code below computes the inverse error function */
200d6685f55SMatthew G. Knepley   for (k = 0; k < maxIter; ++k) ck[k] = 0.;
201d6685f55SMatthew G. Knepley   ck[0] = 1;
202d6685f55SMatthew G. Knepley   r     = ck[0] * (PetscSqrtReal(PETSC_PI) / 2.) * q;
203d6685f55SMatthew G. Knepley   for (k = 1; k < maxIter; ++k) {
204d6685f55SMatthew G. Knepley     const PetscReal temp = 2. * k + 1.;
205d6685f55SMatthew G. Knepley 
206d6685f55SMatthew G. Knepley     for (m = 0; m <= k - 1; ++m) {
207d6685f55SMatthew G. Knepley       PetscReal denom = (m + 1.) * (2. * m + 1.);
208d6685f55SMatthew G. Knepley 
209d6685f55SMatthew G. Knepley       ck[k] += (ck[m] * ck[k - 1 - m]) / denom;
210d6685f55SMatthew G. Knepley     }
211d6685f55SMatthew G. Knepley     r += (ck[k] / temp) * PetscPowReal((PetscSqrtReal(PETSC_PI) / 2.) * q, 2. * k + 1.);
212d6685f55SMatthew G. Knepley   }
213d6685f55SMatthew G. Knepley   /* Scale erfinv() by \sqrt{\pi/2} */
214d6685f55SMatthew G. Knepley   x[0] = PetscSqrtReal(PETSC_PI * 0.5) * r;
2153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
216d6685f55SMatthew G. Knepley }
217d6685f55SMatthew G. Knepley 
218d6685f55SMatthew G. Knepley /*@
219d6685f55SMatthew G. Knepley   PetscPDFGaussian2D - PDF for the Gaussian distribution in 2D
220d6685f55SMatthew G. Knepley 
22120f4b53cSBarry Smith   Not Collective
222d6685f55SMatthew G. Knepley 
223817da375SSatish Balay   Input Parameters:
224817da375SSatish Balay + x     - Coordinate in [-\infty, \infty]^2
225817da375SSatish Balay - dummy - ignored
226d6685f55SMatthew G. Knepley 
227d6685f55SMatthew G. Knepley   Output Parameter:
228d6685f55SMatthew G. Knepley . p - The probability density at x
229d6685f55SMatthew G. Knepley 
230d6685f55SMatthew G. Knepley   Level: beginner
231d6685f55SMatthew G. Knepley 
232db781477SPatrick Sanan .seealso: `PetscPDFSampleGaussian2D()`, `PetscPDFMaxwellBoltzmann3D()`
233d6685f55SMatthew G. Knepley @*/
234d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFGaussian2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
235d71ae5a4SJacob Faibussowitsch {
236d6685f55SMatthew G. Knepley   p[0] = (1. / PETSC_PI) * PetscExpReal(-0.5 * (PetscSqr(x[0]) + PetscSqr(x[1])));
2373ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
238d6685f55SMatthew G. Knepley }
239d6685f55SMatthew G. Knepley 
240d6685f55SMatthew G. Knepley /*@
241d6685f55SMatthew G. Knepley   PetscPDFSampleGaussian2D - Sample uniformly from a Gaussian distribution in 2D
242d6685f55SMatthew G. Knepley 
24320f4b53cSBarry Smith   Not Collective
244d6685f55SMatthew G. Knepley 
245817da375SSatish Balay   Input Parameters:
246817da375SSatish Balay + p     - A uniform variable on [0, 1]^2
247817da375SSatish Balay - dummy - ignored
248d6685f55SMatthew G. Knepley 
249d6685f55SMatthew G. Knepley   Output Parameter:
250d6685f55SMatthew G. Knepley . x - Coordinate in [-\infty, \infty]^2
251d6685f55SMatthew G. Knepley 
252d6685f55SMatthew G. Knepley   Level: beginner
253d6685f55SMatthew G. Knepley 
254*a4e35b19SJacob Faibussowitsch   References:
255*a4e35b19SJacob Faibussowitsch . * - https://en.wikipedia.org/wiki/Box%E2%80%93Muller_transform
256d6685f55SMatthew G. Knepley 
257db781477SPatrick Sanan .seealso: `PetscPDFGaussian2D()`, `PetscPDFMaxwellBoltzmann3D()`
258d6685f55SMatthew G. Knepley @*/
259d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
260d71ae5a4SJacob Faibussowitsch {
261d6685f55SMatthew G. Knepley   const PetscReal mag = PetscSqrtReal(-2.0 * PetscLogReal(p[0]));
262d6685f55SMatthew G. Knepley   x[0]                = mag * PetscCosReal(2.0 * PETSC_PI * p[1]);
263d6685f55SMatthew G. Knepley   x[1]                = mag * PetscSinReal(2.0 * PETSC_PI * p[1]);
2643ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
265d6685f55SMatthew G. Knepley }
266d6685f55SMatthew G. Knepley 
267d6685f55SMatthew G. Knepley /*@
268ea1b28ebSMatthew G. Knepley   PetscPDFGaussian3D - PDF for the Gaussian distribution in 3D
269ea1b28ebSMatthew G. Knepley 
27020f4b53cSBarry Smith   Not Collective
271ea1b28ebSMatthew G. Knepley 
272ea1b28ebSMatthew G. Knepley   Input Parameters:
273ea1b28ebSMatthew G. Knepley + x     - Coordinate in [-\infty, \infty]^3
274ea1b28ebSMatthew G. Knepley - dummy - ignored
275ea1b28ebSMatthew G. Knepley 
276ea1b28ebSMatthew G. Knepley   Output Parameter:
277ea1b28ebSMatthew G. Knepley . p - The probability density at x
278ea1b28ebSMatthew G. Knepley 
279ea1b28ebSMatthew G. Knepley   Level: beginner
280ea1b28ebSMatthew G. Knepley 
281db781477SPatrick Sanan .seealso: `PetscPDFSampleGaussian3D()`, `PetscPDFMaxwellBoltzmann3D()`
282ea1b28ebSMatthew G. Knepley @*/
283d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFGaussian3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
284d71ae5a4SJacob Faibussowitsch {
285ea1b28ebSMatthew G. Knepley   p[0] = (1. / PETSC_PI * PetscSqrtReal(PETSC_PI)) * PetscExpReal(-0.5 * (PetscSqr(x[0]) + PetscSqr(x[1]) + PetscSqr(x[2])));
2863ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
287ea1b28ebSMatthew G. Knepley }
288ea1b28ebSMatthew G. Knepley 
289ea1b28ebSMatthew G. Knepley /*@
290ea1b28ebSMatthew G. Knepley   PetscPDFSampleGaussian3D - Sample uniformly from a Gaussian distribution in 3D
291ea1b28ebSMatthew G. Knepley 
29220f4b53cSBarry Smith   Not Collective
293ea1b28ebSMatthew G. Knepley 
294ea1b28ebSMatthew G. Knepley   Input Parameters:
295ea1b28ebSMatthew G. Knepley + p     - A uniform variable on [0, 1]^3
296ea1b28ebSMatthew G. Knepley - dummy - ignored
297ea1b28ebSMatthew G. Knepley 
298ea1b28ebSMatthew G. Knepley   Output Parameter:
299ea1b28ebSMatthew G. Knepley . x - Coordinate in [-\infty, \infty]^3
300ea1b28ebSMatthew G. Knepley 
301ea1b28ebSMatthew G. Knepley   Level: beginner
302ea1b28ebSMatthew G. Knepley 
30360225df5SJacob Faibussowitsch   References:
304dce8aebaSBarry Smith   https://en.wikipedia.org/wiki/Box%E2%80%93Muller_transform
305ea1b28ebSMatthew G. Knepley 
306db781477SPatrick Sanan .seealso: `PetscPDFGaussian3D()`, `PetscPDFMaxwellBoltzmann3D()`
307ea1b28ebSMatthew G. Knepley @*/
308d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
309d71ae5a4SJacob Faibussowitsch {
310ea1b28ebSMatthew G. Knepley   PetscCall(PetscPDFSampleGaussian1D(p, dummy, x));
311ea1b28ebSMatthew G. Knepley   PetscCall(PetscPDFSampleGaussian2D(&p[1], dummy, &x[1]));
3123ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
313ea1b28ebSMatthew G. Knepley }
314ea1b28ebSMatthew G. Knepley 
315ea1b28ebSMatthew G. Knepley /*@
316d6685f55SMatthew G. Knepley   PetscPDFConstant1D - PDF for the uniform distribution in 1D
317d6685f55SMatthew G. Knepley 
31820f4b53cSBarry Smith   Not Collective
319d6685f55SMatthew G. Knepley 
320*a4e35b19SJacob Faibussowitsch   Input Parameters:
321*a4e35b19SJacob Faibussowitsch + x     - Coordinate in [-1, 1]
322*a4e35b19SJacob Faibussowitsch - dummy - Unused
323d6685f55SMatthew G. Knepley 
324d6685f55SMatthew G. Knepley   Output Parameter:
325d6685f55SMatthew G. Knepley . p - The probability density at x
326d6685f55SMatthew G. Knepley 
327d6685f55SMatthew G. Knepley   Level: beginner
328d6685f55SMatthew G. Knepley 
329db781477SPatrick Sanan .seealso: `PetscCDFConstant1D()`, `PetscPDFSampleConstant1D()`, `PetscPDFConstant2D()`, `PetscPDFConstant3D()`
330d6685f55SMatthew G. Knepley @*/
331d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFConstant1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
332d71ae5a4SJacob Faibussowitsch {
333d6685f55SMatthew G. Knepley   p[0] = x[0] > -1. && x[0] <= 1. ? 0.5 : 0.;
3343ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
335d6685f55SMatthew G. Knepley }
336d6685f55SMatthew G. Knepley 
337ea1b28ebSMatthew G. Knepley /*@
338ea1b28ebSMatthew G. Knepley   PetscCDFConstant1D - CDF for the uniform distribution in 1D
339ea1b28ebSMatthew G. Knepley 
34020f4b53cSBarry Smith   Not Collective
341ea1b28ebSMatthew G. Knepley 
342*a4e35b19SJacob Faibussowitsch   Input Parameters:
343*a4e35b19SJacob Faibussowitsch + x     - Coordinate in [-1, 1]
344*a4e35b19SJacob Faibussowitsch - dummy - Unused
345ea1b28ebSMatthew G. Knepley 
346ea1b28ebSMatthew G. Knepley   Output Parameter:
347ea1b28ebSMatthew G. Knepley . p - The cumulative probability at x
348ea1b28ebSMatthew G. Knepley 
349ea1b28ebSMatthew G. Knepley   Level: beginner
350ea1b28ebSMatthew G. Knepley 
351db781477SPatrick Sanan .seealso: `PetscPDFConstant1D()`, `PetscPDFSampleConstant1D()`, `PetscCDFConstant2D()`, `PetscCDFConstant3D()`
352ea1b28ebSMatthew G. Knepley @*/
353d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCDFConstant1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
354d71ae5a4SJacob Faibussowitsch {
355d6685f55SMatthew G. Knepley   p[0] = x[0] <= -1. ? 0. : (x[0] > 1. ? 1. : 0.5 * (x[0] + 1.));
3563ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
357d6685f55SMatthew G. Knepley }
358d6685f55SMatthew G. Knepley 
359d6685f55SMatthew G. Knepley /*@
360d6685f55SMatthew G. Knepley   PetscPDFSampleConstant1D - Sample uniformly from a uniform distribution on [-1, 1] in 1D
361d6685f55SMatthew G. Knepley 
36220f4b53cSBarry Smith   Not Collective
363d6685f55SMatthew G. Knepley 
364*a4e35b19SJacob Faibussowitsch   Input Parameters:
365*a4e35b19SJacob Faibussowitsch + p     - A uniform variable on [0, 1]
366*a4e35b19SJacob Faibussowitsch - dummy - Unused
367d6685f55SMatthew G. Knepley 
368d6685f55SMatthew G. Knepley   Output Parameter:
369d6685f55SMatthew G. Knepley . x - Coordinate in [-1, 1]
370d6685f55SMatthew G. Knepley 
371d6685f55SMatthew G. Knepley   Level: beginner
372d6685f55SMatthew G. Knepley 
373db781477SPatrick Sanan .seealso: `PetscPDFConstant1D()`, `PetscCDFConstant1D()`, `PetscPDFSampleConstant2D()`, `PetscPDFSampleConstant3D()`
374d6685f55SMatthew G. Knepley @*/
375d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFSampleConstant1D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
376d71ae5a4SJacob Faibussowitsch {
377ea1b28ebSMatthew G. Knepley   x[0] = 2. * p[0] - 1.;
3783ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
379ea1b28ebSMatthew G. Knepley }
380ea1b28ebSMatthew G. Knepley 
381ea1b28ebSMatthew G. Knepley /*@
382ea1b28ebSMatthew G. Knepley   PetscPDFConstant2D - PDF for the uniform distribution in 2D
383ea1b28ebSMatthew G. Knepley 
38420f4b53cSBarry Smith   Not Collective
385ea1b28ebSMatthew G. Knepley 
386*a4e35b19SJacob Faibussowitsch   Input Parameters:
387*a4e35b19SJacob Faibussowitsch + x     - Coordinate in [-1, 1] x [-1, 1]
388*a4e35b19SJacob Faibussowitsch - dummy - Unused
389ea1b28ebSMatthew G. Knepley 
390ea1b28ebSMatthew G. Knepley   Output Parameter:
391ea1b28ebSMatthew G. Knepley . p - The probability density at x
392ea1b28ebSMatthew G. Knepley 
393ea1b28ebSMatthew G. Knepley   Level: beginner
394ea1b28ebSMatthew G. Knepley 
395db781477SPatrick Sanan .seealso: `PetscCDFConstant2D()`, `PetscPDFSampleConstant2D()`, `PetscPDFConstant1D()`, `PetscPDFConstant3D()`
396ea1b28ebSMatthew G. Knepley @*/
397d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFConstant2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
398d71ae5a4SJacob Faibussowitsch {
399ea1b28ebSMatthew G. Knepley   p[0] = x[0] > -1. && x[0] <= 1. && x[1] > -1. && x[1] <= 1. ? 0.25 : 0.;
4003ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
401ea1b28ebSMatthew G. Knepley }
402ea1b28ebSMatthew G. Knepley 
403ea1b28ebSMatthew G. Knepley /*@
404ea1b28ebSMatthew G. Knepley   PetscCDFConstant2D - CDF for the uniform distribution in 2D
405ea1b28ebSMatthew G. Knepley 
40620f4b53cSBarry Smith   Not Collective
407ea1b28ebSMatthew G. Knepley 
408*a4e35b19SJacob Faibussowitsch   Input Parameters:
409*a4e35b19SJacob Faibussowitsch + x     - Coordinate in [-1, 1] x [-1, 1]
410*a4e35b19SJacob Faibussowitsch - dummy - Unused
411ea1b28ebSMatthew G. Knepley 
412ea1b28ebSMatthew G. Knepley   Output Parameter:
413ea1b28ebSMatthew G. Knepley . p - The cumulative probability at x
414ea1b28ebSMatthew G. Knepley 
415ea1b28ebSMatthew G. Knepley   Level: beginner
416ea1b28ebSMatthew G. Knepley 
417db781477SPatrick Sanan .seealso: `PetscPDFConstant2D()`, `PetscPDFSampleConstant2D()`, `PetscCDFConstant1D()`, `PetscCDFConstant3D()`
418ea1b28ebSMatthew G. Knepley @*/
419d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCDFConstant2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
420d71ae5a4SJacob Faibussowitsch {
421ea1b28ebSMatthew 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.));
4223ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
423ea1b28ebSMatthew G. Knepley }
424ea1b28ebSMatthew G. Knepley 
425ea1b28ebSMatthew G. Knepley /*@
426ea1b28ebSMatthew G. Knepley   PetscPDFSampleConstant2D - Sample uniformly from a uniform distribution on [-1, 1] x [-1, 1] in 2D
427ea1b28ebSMatthew G. Knepley 
42820f4b53cSBarry Smith   Not Collective
429ea1b28ebSMatthew G. Knepley 
430*a4e35b19SJacob Faibussowitsch   Input Parameters:
431*a4e35b19SJacob Faibussowitsch + p     - Two uniform variables on [0, 1]
432*a4e35b19SJacob Faibussowitsch - dummy - Unused
433ea1b28ebSMatthew G. Knepley 
434ea1b28ebSMatthew G. Knepley   Output Parameter:
435ea1b28ebSMatthew G. Knepley . x - Coordinate in [-1, 1] x [-1, 1]
436ea1b28ebSMatthew G. Knepley 
437ea1b28ebSMatthew G. Knepley   Level: beginner
438ea1b28ebSMatthew G. Knepley 
439db781477SPatrick Sanan .seealso: `PetscPDFConstant2D()`, `PetscCDFConstant2D()`, `PetscPDFSampleConstant1D()`, `PetscPDFSampleConstant3D()`
440ea1b28ebSMatthew G. Knepley @*/
441d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFSampleConstant2D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
442d71ae5a4SJacob Faibussowitsch {
443ea1b28ebSMatthew G. Knepley   x[0] = 2. * p[0] - 1.;
444ea1b28ebSMatthew G. Knepley   x[1] = 2. * p[1] - 1.;
4453ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
446ea1b28ebSMatthew G. Knepley }
447ea1b28ebSMatthew G. Knepley 
448ea1b28ebSMatthew G. Knepley /*@
449ea1b28ebSMatthew G. Knepley   PetscPDFConstant3D - PDF for the uniform distribution in 3D
450ea1b28ebSMatthew G. Knepley 
45120f4b53cSBarry Smith   Not Collective
452ea1b28ebSMatthew G. Knepley 
453*a4e35b19SJacob Faibussowitsch   Input Parameters:
454*a4e35b19SJacob Faibussowitsch + x     - Coordinate in [-1, 1] x [-1, 1] x [-1, 1]
455*a4e35b19SJacob Faibussowitsch - dummy - Unused
456ea1b28ebSMatthew G. Knepley 
457ea1b28ebSMatthew G. Knepley   Output Parameter:
458ea1b28ebSMatthew G. Knepley . p - The probability density at x
459ea1b28ebSMatthew G. Knepley 
460ea1b28ebSMatthew G. Knepley   Level: beginner
461ea1b28ebSMatthew G. Knepley 
462db781477SPatrick Sanan .seealso: `PetscCDFConstant3D()`, `PetscPDFSampleConstant3D()`, `PetscPDFSampleConstant1D()`, `PetscPDFSampleConstant2D()`
463ea1b28ebSMatthew G. Knepley @*/
464d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFConstant3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
465d71ae5a4SJacob Faibussowitsch {
466ea1b28ebSMatthew 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.;
4673ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
468ea1b28ebSMatthew G. Knepley }
469ea1b28ebSMatthew G. Knepley 
470ea1b28ebSMatthew G. Knepley /*@
471ea1b28ebSMatthew G. Knepley   PetscCDFConstant3D - CDF for the uniform distribution in 3D
472ea1b28ebSMatthew G. Knepley 
47320f4b53cSBarry Smith   Not Collective
474ea1b28ebSMatthew G. Knepley 
475*a4e35b19SJacob Faibussowitsch   Input Parameters:
476*a4e35b19SJacob Faibussowitsch + x     - Coordinate in [-1, 1] x [-1, 1] x [-1, 1]
477*a4e35b19SJacob Faibussowitsch - dummy - Unused
478ea1b28ebSMatthew G. Knepley 
479ea1b28ebSMatthew G. Knepley   Output Parameter:
480ea1b28ebSMatthew G. Knepley . p - The cumulative probability at x
481ea1b28ebSMatthew G. Knepley 
482ea1b28ebSMatthew G. Knepley   Level: beginner
483ea1b28ebSMatthew G. Knepley 
484db781477SPatrick Sanan .seealso: `PetscPDFConstant3D()`, `PetscPDFSampleConstant3D()`, `PetscCDFConstant1D()`, `PetscCDFConstant2D()`
485ea1b28ebSMatthew G. Knepley @*/
486d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCDFConstant3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
487d71ae5a4SJacob Faibussowitsch {
488ea1b28ebSMatthew 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.));
4893ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
490ea1b28ebSMatthew G. Knepley }
491ea1b28ebSMatthew G. Knepley 
492ea1b28ebSMatthew G. Knepley /*@
493ea1b28ebSMatthew G. Knepley   PetscPDFSampleConstant3D - Sample uniformly from a uniform distribution on [-1, 1] x [-1, 1] in 3D
494ea1b28ebSMatthew G. Knepley 
49520f4b53cSBarry Smith   Not Collective
496ea1b28ebSMatthew G. Knepley 
497*a4e35b19SJacob Faibussowitsch   Input Parameters:
498*a4e35b19SJacob Faibussowitsch + p     - Three uniform variables on [0, 1]
499*a4e35b19SJacob Faibussowitsch - dummy - Unused
500ea1b28ebSMatthew G. Knepley 
501ea1b28ebSMatthew G. Knepley   Output Parameter:
502ea1b28ebSMatthew G. Knepley . x - Coordinate in [-1, 1] x [-1, 1] x [-1, 1]
503ea1b28ebSMatthew G. Knepley 
504ea1b28ebSMatthew G. Knepley   Level: beginner
505ea1b28ebSMatthew G. Knepley 
506db781477SPatrick Sanan .seealso: `PetscPDFConstant3D()`, `PetscCDFConstant3D()`, `PetscPDFSampleConstant1D()`, `PetscPDFSampleConstant2D()`
507ea1b28ebSMatthew G. Knepley @*/
508d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFSampleConstant3D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
509d71ae5a4SJacob Faibussowitsch {
510ea1b28ebSMatthew G. Knepley   x[0] = 2. * p[0] - 1.;
511ea1b28ebSMatthew G. Knepley   x[1] = 2. * p[1] - 1.;
512ea1b28ebSMatthew G. Knepley   x[2] = 2. * p[2] - 1.;
5133ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
514d6685f55SMatthew G. Knepley }
515d6685f55SMatthew G. Knepley 
516d6685f55SMatthew G. Knepley /*@C
517da81f932SPierre Jolivet   PetscProbCreateFromOptions - Return the probability distribution specified by the arguments and options
518d6685f55SMatthew G. Knepley 
51920f4b53cSBarry Smith   Not Collective
520d6685f55SMatthew G. Knepley 
521d6685f55SMatthew G. Knepley   Input Parameters:
522d6685f55SMatthew G. Knepley + dim    - The dimension of sample points
523d6685f55SMatthew G. Knepley . prefix - The options prefix, or NULL
524d6685f55SMatthew G. Knepley - name   - The option name for the probility distribution type
525d6685f55SMatthew G. Knepley 
526d6685f55SMatthew G. Knepley   Output Parameters:
527d6685f55SMatthew G. Knepley + pdf     - The PDF of this type
528d6685f55SMatthew G. Knepley . cdf     - The CDF of this type
529d6685f55SMatthew G. Knepley - sampler - The PDF sampler of this type
530d6685f55SMatthew G. Knepley 
531d6685f55SMatthew G. Knepley   Level: intermediate
532d6685f55SMatthew G. Knepley 
533dce8aebaSBarry Smith .seealso: `PetscProbFunc`, `PetscPDFMaxwellBoltzmann1D()`, `PetscPDFGaussian1D()`, `PetscPDFConstant1D()`
534d6685f55SMatthew G. Knepley @*/
535d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscProbCreateFromOptions(PetscInt dim, const char prefix[], const char name[], PetscProbFunc *pdf, PetscProbFunc *cdf, PetscProbFunc *sampler)
536d71ae5a4SJacob Faibussowitsch {
537d6685f55SMatthew G. Knepley   DTProbDensityType den = DTPROB_DENSITY_GAUSSIAN;
538d6685f55SMatthew G. Knepley 
539d6685f55SMatthew G. Knepley   PetscFunctionBegin;
540d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF, prefix, "PetscProb Options", "DT");
5419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum(name, "Method to compute PDF <constant, gaussian>", "", DTProbDensityTypes, (PetscEnum)den, (PetscEnum *)&den, NULL));
542d0609cedSBarry Smith   PetscOptionsEnd();
543d6685f55SMatthew G. Knepley 
5449371c9d4SSatish Balay   if (pdf) {
5454f572ea9SToby Isaac     PetscAssertPointer(pdf, 4);
5469371c9d4SSatish Balay     *pdf = NULL;
5479371c9d4SSatish Balay   }
5489371c9d4SSatish Balay   if (cdf) {
5494f572ea9SToby Isaac     PetscAssertPointer(cdf, 5);
5509371c9d4SSatish Balay     *cdf = NULL;
5519371c9d4SSatish Balay   }
5529371c9d4SSatish Balay   if (sampler) {
5534f572ea9SToby Isaac     PetscAssertPointer(sampler, 6);
5549371c9d4SSatish Balay     *sampler = NULL;
5559371c9d4SSatish Balay   }
556d6685f55SMatthew G. Knepley   switch (den) {
557d6685f55SMatthew G. Knepley   case DTPROB_DENSITY_CONSTANT:
558d6685f55SMatthew G. Knepley     switch (dim) {
559d6685f55SMatthew G. Knepley     case 1:
560d6685f55SMatthew G. Knepley       if (pdf) *pdf = PetscPDFConstant1D;
561d6685f55SMatthew G. Knepley       if (cdf) *cdf = PetscCDFConstant1D;
562d6685f55SMatthew G. Knepley       if (sampler) *sampler = PetscPDFSampleConstant1D;
563d6685f55SMatthew G. Knepley       break;
564ea1b28ebSMatthew G. Knepley     case 2:
565ea1b28ebSMatthew G. Knepley       if (pdf) *pdf = PetscPDFConstant2D;
566ea1b28ebSMatthew G. Knepley       if (cdf) *cdf = PetscCDFConstant2D;
567ea1b28ebSMatthew G. Knepley       if (sampler) *sampler = PetscPDFSampleConstant2D;
568ea1b28ebSMatthew G. Knepley       break;
569ea1b28ebSMatthew G. Knepley     case 3:
570ea1b28ebSMatthew G. Knepley       if (pdf) *pdf = PetscPDFConstant3D;
571ea1b28ebSMatthew G. Knepley       if (cdf) *cdf = PetscCDFConstant3D;
572ea1b28ebSMatthew G. Knepley       if (sampler) *sampler = PetscPDFSampleConstant3D;
573ea1b28ebSMatthew G. Knepley       break;
574d71ae5a4SJacob Faibussowitsch     default:
575d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]);
576d6685f55SMatthew G. Knepley     }
577d6685f55SMatthew G. Knepley     break;
578d6685f55SMatthew G. Knepley   case DTPROB_DENSITY_GAUSSIAN:
579d6685f55SMatthew G. Knepley     switch (dim) {
580d6685f55SMatthew G. Knepley     case 1:
581d6685f55SMatthew G. Knepley       if (pdf) *pdf = PetscPDFGaussian1D;
582d6685f55SMatthew G. Knepley       if (cdf) *cdf = PetscCDFGaussian1D;
583d6685f55SMatthew G. Knepley       if (sampler) *sampler = PetscPDFSampleGaussian1D;
584d6685f55SMatthew G. Knepley       break;
585d6685f55SMatthew G. Knepley     case 2:
586d6685f55SMatthew G. Knepley       if (pdf) *pdf = PetscPDFGaussian2D;
587d6685f55SMatthew G. Knepley       if (sampler) *sampler = PetscPDFSampleGaussian2D;
588d6685f55SMatthew G. Knepley       break;
589ea1b28ebSMatthew G. Knepley     case 3:
590ea1b28ebSMatthew G. Knepley       if (pdf) *pdf = PetscPDFGaussian3D;
591ea1b28ebSMatthew G. Knepley       if (sampler) *sampler = PetscPDFSampleGaussian3D;
592ea1b28ebSMatthew G. Knepley       break;
593d71ae5a4SJacob Faibussowitsch     default:
594d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]);
595d6685f55SMatthew G. Knepley     }
596d6685f55SMatthew G. Knepley     break;
597d6685f55SMatthew G. Knepley   case DTPROB_DENSITY_MAXWELL_BOLTZMANN:
598d6685f55SMatthew G. Knepley     switch (dim) {
599d6685f55SMatthew G. Knepley     case 1:
600d6685f55SMatthew G. Knepley       if (pdf) *pdf = PetscPDFMaxwellBoltzmann1D;
601d6685f55SMatthew G. Knepley       if (cdf) *cdf = PetscCDFMaxwellBoltzmann1D;
602d6685f55SMatthew G. Knepley       break;
603d6685f55SMatthew G. Knepley     case 2:
604d6685f55SMatthew G. Knepley       if (pdf) *pdf = PetscPDFMaxwellBoltzmann2D;
605d6685f55SMatthew G. Knepley       if (cdf) *cdf = PetscCDFMaxwellBoltzmann2D;
606d6685f55SMatthew G. Knepley       break;
607d6685f55SMatthew G. Knepley     case 3:
608d6685f55SMatthew G. Knepley       if (pdf) *pdf = PetscPDFMaxwellBoltzmann3D;
609d6685f55SMatthew G. Knepley       if (cdf) *cdf = PetscCDFMaxwellBoltzmann3D;
610d6685f55SMatthew G. Knepley       break;
611d71ae5a4SJacob Faibussowitsch     default:
612d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]);
613d6685f55SMatthew G. Knepley     }
614d6685f55SMatthew G. Knepley     break;
615d71ae5a4SJacob Faibussowitsch   default:
616d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Density type %s is not supported", DTProbDensityTypes[PetscMax(0, PetscMin(den, DTPROB_NUM_DENSITY))]);
617d6685f55SMatthew G. Knepley   }
6183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
619d6685f55SMatthew G. Knepley }
620d6685f55SMatthew G. Knepley 
621d6685f55SMatthew G. Knepley #ifdef PETSC_HAVE_KS
622d6685f55SMatthew G. Knepley EXTERN_C_BEGIN
623d6685f55SMatthew G. Knepley   #include <KolmogorovSmirnovDist.h>
624d6685f55SMatthew G. Knepley EXTERN_C_END
625d6685f55SMatthew G. Knepley #endif
626d6685f55SMatthew G. Knepley 
627d6685f55SMatthew G. Knepley /*@C
628d6685f55SMatthew G. Knepley   PetscProbComputeKSStatistic - Compute the Kolmogorov-Smirnov statistic for the empirical distribution for an input vector, compared to an analytic CDF.
629d6685f55SMatthew G. Knepley 
63020f4b53cSBarry Smith   Collective
631d6685f55SMatthew G. Knepley 
632d6685f55SMatthew G. Knepley   Input Parameters:
633d6685f55SMatthew G. Knepley + v   - The data vector, blocksize is the sample dimension
634d6685f55SMatthew G. Knepley - cdf - The analytic CDF
635d6685f55SMatthew G. Knepley 
636d6685f55SMatthew G. Knepley   Output Parameter:
637d6685f55SMatthew G. Knepley . alpha - The KS statisic
638d6685f55SMatthew G. Knepley 
639d6685f55SMatthew G. Knepley   Level: advanced
640d6685f55SMatthew G. Knepley 
641dce8aebaSBarry Smith   Notes:
642dce8aebaSBarry Smith   The Kolmogorov-Smirnov statistic for a given cumulative distribution function $F(x)$ is
643dce8aebaSBarry Smith .vb
644dce8aebaSBarry Smith     D_n = \sup_x \left| F_n(x) - F(x) \right|
645dce8aebaSBarry Smith 
646dce8aebaSBarry 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
647dce8aebaSBarry Smith 
648dce8aebaSBarry Smith     F_n =  # of samples <= x / n
649dce8aebaSBarry Smith .ve
650dce8aebaSBarry Smith 
651d6685f55SMatthew G. Knepley   The empirical distribution function $F_n(x)$ is discrete, and thus had a ``stairstep''
652d6685f55SMatthew G. Knepley   cumulative distribution, making $n$ the number of stairs. Intuitively, the statistic takes
653d6685f55SMatthew G. Knepley   the largest absolute difference between the two distribution functions across all $x$ values.
654d6685f55SMatthew G. Knepley 
655db781477SPatrick Sanan .seealso: `PetscProbFunc`
656d6685f55SMatthew G. Knepley @*/
657d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscProbComputeKSStatistic(Vec v, PetscProbFunc cdf, PetscReal *alpha)
658d71ae5a4SJacob Faibussowitsch {
659d6685f55SMatthew G. Knepley #if !defined(PETSC_HAVE_KS)
660d6685f55SMatthew G. Knepley   SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "No support for Kolmogorov-Smirnov test.\nReconfigure using --download-ks");
661d6685f55SMatthew G. Knepley #else
662d6685f55SMatthew G. Knepley   PetscViewer        viewer = NULL;
663d6685f55SMatthew G. Knepley   PetscViewerFormat  format;
664d6685f55SMatthew G. Knepley   PetscDraw          draw;
665d6685f55SMatthew G. Knepley   PetscDrawLG        lg;
666d6685f55SMatthew G. Knepley   PetscDrawAxis      axis;
667d6685f55SMatthew G. Knepley   const PetscScalar *a;
668d6685f55SMatthew G. Knepley   PetscReal         *speed, Dn = PETSC_MIN_REAL;
669d6685f55SMatthew G. Knepley   PetscBool          isascii = PETSC_FALSE, isdraw = PETSC_FALSE, flg;
670d6685f55SMatthew G. Knepley   PetscInt           n, p, dim, d;
671d6685f55SMatthew G. Knepley   PetscMPIInt        size;
672d6685f55SMatthew G. Knepley   const char        *names[2] = {"Analytic", "Empirical"};
673d6685f55SMatthew G. Knepley   char               title[PETSC_MAX_PATH_LEN];
674d6685f55SMatthew G. Knepley   PetscOptions       options;
675d6685f55SMatthew G. Knepley   const char        *prefix;
676d6685f55SMatthew G. Knepley   MPI_Comm           comm;
677d6685f55SMatthew G. Knepley 
678d6685f55SMatthew G. Knepley   PetscFunctionBegin;
6799566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
6809566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)v, &prefix));
6819566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptions((PetscObject)v, &options));
6829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(comm, options, prefix, "-ks_monitor", &viewer, &format, &flg));
683d6685f55SMatthew G. Knepley   if (flg) {
6849566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
6859566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
686d6685f55SMatthew G. Knepley   }
687d6685f55SMatthew G. Knepley   if (isascii) {
6889566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
689d6685f55SMatthew G. Knepley   } else if (isdraw) {
6909566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
6919566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
6929566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGCreate(draw, 2, &lg));
6939566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSetLegend(lg, names));
694d6685f55SMatthew G. Knepley   }
695d6685f55SMatthew G. Knepley 
6969566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size));
6979566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(v, &n));
6989566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(v, &dim));
699d6685f55SMatthew G. Knepley   n /= dim;
700d6685f55SMatthew G. Knepley   /* TODO Parallel algorithm is harder */
701d6685f55SMatthew G. Knepley   if (size == 1) {
7029566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &speed));
7039566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(v, &a));
704d6685f55SMatthew G. Knepley     for (p = 0; p < n; ++p) {
705d6685f55SMatthew G. Knepley       PetscReal mag = 0.;
706d6685f55SMatthew G. Knepley 
707d6685f55SMatthew G. Knepley       for (d = 0; d < dim; ++d) mag += PetscSqr(PetscRealPart(a[p * dim + d]));
708d6685f55SMatthew G. Knepley       speed[p] = PetscSqrtReal(mag);
709d6685f55SMatthew G. Knepley     }
7109566063dSJacob Faibussowitsch     PetscCall(PetscSortReal(n, speed));
7119566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(v, &a));
712d6685f55SMatthew G. Knepley     for (p = 0; p < n; ++p) {
713d6685f55SMatthew G. Knepley       const PetscReal x = speed[p], Fn = ((PetscReal)p) / n;
714d6685f55SMatthew G. Knepley       PetscReal       F, vals[2];
715d6685f55SMatthew G. Knepley 
7169566063dSJacob Faibussowitsch       PetscCall(cdf(&x, NULL, &F));
717d6685f55SMatthew G. Knepley       Dn = PetscMax(PetscAbsReal(Fn - F), Dn);
7189566063dSJacob Faibussowitsch       if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "x: %g F: %g Fn: %g Dn: %.2g\n", x, F, Fn, Dn));
7199371c9d4SSatish Balay       if (isdraw) {
7209371c9d4SSatish Balay         vals[0] = F;
7219371c9d4SSatish Balay         vals[1] = Fn;
7229371c9d4SSatish Balay         PetscCall(PetscDrawLGAddCommonPoint(lg, x, vals));
7239371c9d4SSatish Balay       }
724d6685f55SMatthew G. Knepley     }
725d6685f55SMatthew G. Knepley     if (speed[n - 1] < 6.) {
726d6685f55SMatthew G. Knepley       const PetscReal k = (PetscInt)(6. - speed[n - 1]) + 1, dx = (6. - speed[n - 1]) / k;
727d6685f55SMatthew G. Knepley       for (p = 0; p <= k; ++p) {
728d6685f55SMatthew G. Knepley         const PetscReal x = speed[n - 1] + p * dx, Fn = 1.0;
729d6685f55SMatthew G. Knepley         PetscReal       F, vals[2];
730d6685f55SMatthew G. Knepley 
7319566063dSJacob Faibussowitsch         PetscCall(cdf(&x, NULL, &F));
732d6685f55SMatthew G. Knepley         Dn = PetscMax(PetscAbsReal(Fn - F), Dn);
7339566063dSJacob Faibussowitsch         if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "x: %g F: %g Fn: %g Dn: %.2g\n", x, F, Fn, Dn));
7349371c9d4SSatish Balay         if (isdraw) {
7359371c9d4SSatish Balay           vals[0] = F;
7369371c9d4SSatish Balay           vals[1] = Fn;
7379371c9d4SSatish Balay           PetscCall(PetscDrawLGAddCommonPoint(lg, x, vals));
7389371c9d4SSatish Balay         }
739d6685f55SMatthew G. Knepley       }
740d6685f55SMatthew G. Knepley     }
7419566063dSJacob Faibussowitsch     PetscCall(PetscFree(speed));
742d6685f55SMatthew G. Knepley   }
743d6685f55SMatthew G. Knepley   if (isdraw) {
7449566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(lg, &axis));
7459566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(title, PETSC_MAX_PATH_LEN, "Kolmogorov-Smirnov Test (Dn %.2g)", Dn));
7469566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, title, "x", "CDF(x)"));
7479566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(lg));
7489566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDestroy(&lg));
749d6685f55SMatthew G. Knepley   }
750d6685f55SMatthew G. Knepley   if (viewer) {
7519566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
7529566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewer));
7539566063dSJacob Faibussowitsch     PetscCall(PetscViewerDestroy(&viewer));
754d6685f55SMatthew G. Knepley   }
755d6685f55SMatthew G. Knepley   *alpha = KSfbar((int)n, (double)Dn);
7563ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
757d6685f55SMatthew G. Knepley #endif
758d6685f55SMatthew G. Knepley }
759