xref: /petsc/src/dm/dt/interface/dtprob.c (revision f13dfd9ea68e0ddeee984e65c377a1819eab8a8a)
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 
14a4e35b19SJacob Faibussowitsch   Input Parameters:
151d27aa22SBarry Smith + x     - Speed in $[0, \infty]$
16a4e35b19SJacob Faibussowitsch - dummy - Unused
17d6685f55SMatthew G. Knepley 
18d6685f55SMatthew G. Knepley   Output Parameter:
191d27aa22SBarry Smith . 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 
36a4e35b19SJacob Faibussowitsch   Input Parameters:
371d27aa22SBarry Smith + x     - Speed in $[0, \infty]$
38a4e35b19SJacob Faibussowitsch - dummy - Unused
39d6685f55SMatthew G. Knepley 
40d6685f55SMatthew G. Knepley   Output Parameter:
411d27aa22SBarry Smith . 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 
58a4e35b19SJacob Faibussowitsch   Input Parameters:
591d27aa22SBarry Smith + x     - Speed in $[0, \infty]$
60a4e35b19SJacob Faibussowitsch - dummy - Unused
61d6685f55SMatthew G. Knepley 
62d6685f55SMatthew G. Knepley   Output Parameter:
631d27aa22SBarry Smith . 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 
80a4e35b19SJacob Faibussowitsch   Input Parameters:
811d27aa22SBarry Smith + x     - Speed in $[0, \infty]$
82a4e35b19SJacob Faibussowitsch - dummy - Unused
83d6685f55SMatthew G. Knepley 
84d6685f55SMatthew G. Knepley   Output Parameter:
851d27aa22SBarry Smith . 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 
102a4e35b19SJacob Faibussowitsch   Input Parameters:
1031d27aa22SBarry Smith + x     - Speed in $[0, \infty]$
104a4e35b19SJacob Faibussowitsch - dummy - Unused
105d6685f55SMatthew G. Knepley 
106d6685f55SMatthew G. Knepley   Output Parameter:
1071d27aa22SBarry Smith . 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 
124a4e35b19SJacob Faibussowitsch   Input Parameters:
1251d27aa22SBarry Smith + x     - Speed in $[0, \infty]$
126a4e35b19SJacob Faibussowitsch - dummy - Unused
127d6685f55SMatthew G. Knepley 
128d6685f55SMatthew G. Knepley   Output Parameter:
1291d27aa22SBarry Smith . 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 
146a4e35b19SJacob Faibussowitsch   Input Parameters:
1471d27aa22SBarry Smith + x     - Coordinate in $[-\infty, \infty]$
148a4e35b19SJacob Faibussowitsch - scale - Scaling value
149d6685f55SMatthew G. Knepley 
150d6685f55SMatthew G. Knepley   Output Parameter:
1511d27aa22SBarry Smith . 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:
1771d27aa22SBarry Smith + p     - A uniform variable on $[0, 1]$
178817da375SSatish Balay - dummy - ignored
179d6685f55SMatthew G. Knepley 
180d6685f55SMatthew G. Knepley   Output Parameter:
1811d27aa22SBarry Smith . x - Coordinate in $[-\infty, \infty]$
182d6685f55SMatthew G. Knepley 
183d6685f55SMatthew G. Knepley   Level: beginner
184d6685f55SMatthew G. Knepley 
1851d27aa22SBarry Smith   Note:
1861d27aa22SBarry Smith   See <http://www.mimirgames.com/articles/programming/approximations-of-the-inverse-error-function> and
1871d27aa22SBarry Smith   <https://stackoverflow.com/questions/27229371/inverse-error-function-in-c>
188a4e35b19SJacob Faibussowitsch 
189a4e35b19SJacob 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:
2241d27aa22SBarry Smith + x     - Coordinate in $[-\infty, \infty]^2$
225817da375SSatish Balay - dummy - ignored
226d6685f55SMatthew G. Knepley 
227d6685f55SMatthew G. Knepley   Output Parameter:
2281d27aa22SBarry Smith . 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
2421d27aa22SBarry Smith   <https://en.wikipedia.org/wiki/Box%E2%80%93Muller_transform>
243d6685f55SMatthew G. Knepley 
24420f4b53cSBarry Smith   Not Collective
245d6685f55SMatthew G. Knepley 
246817da375SSatish Balay   Input Parameters:
2471d27aa22SBarry Smith + p     - A uniform variable on $[0, 1]^2$
248817da375SSatish Balay - dummy - ignored
249d6685f55SMatthew G. Knepley 
250d6685f55SMatthew G. Knepley   Output Parameter:
2511d27aa22SBarry Smith . x - Coordinate in $[-\infty, \infty]^2 $
252d6685f55SMatthew G. Knepley 
253d6685f55SMatthew G. Knepley   Level: beginner
254d6685f55SMatthew G. Knepley 
255db781477SPatrick Sanan .seealso: `PetscPDFGaussian2D()`, `PetscPDFMaxwellBoltzmann3D()`
256d6685f55SMatthew G. Knepley @*/
257d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
258d71ae5a4SJacob Faibussowitsch {
259d6685f55SMatthew G. Knepley   const PetscReal mag = PetscSqrtReal(-2.0 * PetscLogReal(p[0]));
260d6685f55SMatthew G. Knepley   x[0]                = mag * PetscCosReal(2.0 * PETSC_PI * p[1]);
261d6685f55SMatthew G. Knepley   x[1]                = mag * PetscSinReal(2.0 * PETSC_PI * p[1]);
2623ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
263d6685f55SMatthew G. Knepley }
264d6685f55SMatthew G. Knepley 
265d6685f55SMatthew G. Knepley /*@
266ea1b28ebSMatthew G. Knepley   PetscPDFGaussian3D - PDF for the Gaussian distribution in 3D
267ea1b28ebSMatthew G. Knepley 
26820f4b53cSBarry Smith   Not Collective
269ea1b28ebSMatthew G. Knepley 
270ea1b28ebSMatthew G. Knepley   Input Parameters:
2711d27aa22SBarry Smith + x     - Coordinate in $[-\infty, \infty]^3$
272ea1b28ebSMatthew G. Knepley - dummy - ignored
273ea1b28ebSMatthew G. Knepley 
274ea1b28ebSMatthew G. Knepley   Output Parameter:
2751d27aa22SBarry Smith . p - The probability density at `x`
276ea1b28ebSMatthew G. Knepley 
277ea1b28ebSMatthew G. Knepley   Level: beginner
278ea1b28ebSMatthew G. Knepley 
279db781477SPatrick Sanan .seealso: `PetscPDFSampleGaussian3D()`, `PetscPDFMaxwellBoltzmann3D()`
280ea1b28ebSMatthew G. Knepley @*/
281d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFGaussian3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
282d71ae5a4SJacob Faibussowitsch {
283ea1b28ebSMatthew G. Knepley   p[0] = (1. / PETSC_PI * PetscSqrtReal(PETSC_PI)) * PetscExpReal(-0.5 * (PetscSqr(x[0]) + PetscSqr(x[1]) + PetscSqr(x[2])));
2843ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
285ea1b28ebSMatthew G. Knepley }
286ea1b28ebSMatthew G. Knepley 
287ea1b28ebSMatthew G. Knepley /*@
288ea1b28ebSMatthew G. Knepley   PetscPDFSampleGaussian3D - Sample uniformly from a Gaussian distribution in 3D
2891d27aa22SBarry Smith   <https://en.wikipedia.org/wiki/Box%E2%80%93Muller_transform>
290ea1b28ebSMatthew G. Knepley 
29120f4b53cSBarry Smith   Not Collective
292ea1b28ebSMatthew G. Knepley 
293ea1b28ebSMatthew G. Knepley   Input Parameters:
2941d27aa22SBarry Smith + p     - A uniform variable on $[0, 1]^3$
295ea1b28ebSMatthew G. Knepley - dummy - ignored
296ea1b28ebSMatthew G. Knepley 
297ea1b28ebSMatthew G. Knepley   Output Parameter:
2981d27aa22SBarry Smith . x - Coordinate in $[-\infty, \infty]^3$
299ea1b28ebSMatthew G. Knepley 
300ea1b28ebSMatthew G. Knepley   Level: beginner
301ea1b28ebSMatthew G. Knepley 
302db781477SPatrick Sanan .seealso: `PetscPDFGaussian3D()`, `PetscPDFMaxwellBoltzmann3D()`
303ea1b28ebSMatthew G. Knepley @*/
304d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
305d71ae5a4SJacob Faibussowitsch {
306ea1b28ebSMatthew G. Knepley   PetscCall(PetscPDFSampleGaussian1D(p, dummy, x));
307ea1b28ebSMatthew G. Knepley   PetscCall(PetscPDFSampleGaussian2D(&p[1], dummy, &x[1]));
3083ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
309ea1b28ebSMatthew G. Knepley }
310ea1b28ebSMatthew G. Knepley 
311ea1b28ebSMatthew G. Knepley /*@
312d6685f55SMatthew G. Knepley   PetscPDFConstant1D - PDF for the uniform distribution in 1D
313d6685f55SMatthew G. Knepley 
31420f4b53cSBarry Smith   Not Collective
315d6685f55SMatthew G. Knepley 
316a4e35b19SJacob Faibussowitsch   Input Parameters:
3171d27aa22SBarry Smith + x     - Coordinate in $[-1, 1]$
318a4e35b19SJacob Faibussowitsch - dummy - Unused
319d6685f55SMatthew G. Knepley 
320d6685f55SMatthew G. Knepley   Output Parameter:
3211d27aa22SBarry Smith . p - The probability density at `x`
322d6685f55SMatthew G. Knepley 
323d6685f55SMatthew G. Knepley   Level: beginner
324d6685f55SMatthew G. Knepley 
325db781477SPatrick Sanan .seealso: `PetscCDFConstant1D()`, `PetscPDFSampleConstant1D()`, `PetscPDFConstant2D()`, `PetscPDFConstant3D()`
326d6685f55SMatthew G. Knepley @*/
327d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFConstant1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
328d71ae5a4SJacob Faibussowitsch {
329d6685f55SMatthew G. Knepley   p[0] = x[0] > -1. && x[0] <= 1. ? 0.5 : 0.;
3303ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
331d6685f55SMatthew G. Knepley }
332d6685f55SMatthew G. Knepley 
333ea1b28ebSMatthew G. Knepley /*@
334ea1b28ebSMatthew G. Knepley   PetscCDFConstant1D - CDF for the uniform distribution in 1D
335ea1b28ebSMatthew G. Knepley 
33620f4b53cSBarry Smith   Not Collective
337ea1b28ebSMatthew G. Knepley 
338a4e35b19SJacob Faibussowitsch   Input Parameters:
3391d27aa22SBarry Smith + x     - Coordinate in $[-1, 1]$
340a4e35b19SJacob Faibussowitsch - dummy - Unused
341ea1b28ebSMatthew G. Knepley 
342ea1b28ebSMatthew G. Knepley   Output Parameter:
3431d27aa22SBarry Smith . p - The cumulative probability at `x`
344ea1b28ebSMatthew G. Knepley 
345ea1b28ebSMatthew G. Knepley   Level: beginner
346ea1b28ebSMatthew G. Knepley 
347db781477SPatrick Sanan .seealso: `PetscPDFConstant1D()`, `PetscPDFSampleConstant1D()`, `PetscCDFConstant2D()`, `PetscCDFConstant3D()`
348ea1b28ebSMatthew G. Knepley @*/
349d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCDFConstant1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
350d71ae5a4SJacob Faibussowitsch {
351d6685f55SMatthew G. Knepley   p[0] = x[0] <= -1. ? 0. : (x[0] > 1. ? 1. : 0.5 * (x[0] + 1.));
3523ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
353d6685f55SMatthew G. Knepley }
354d6685f55SMatthew G. Knepley 
355d6685f55SMatthew G. Knepley /*@
356d6685f55SMatthew G. Knepley   PetscPDFSampleConstant1D - Sample uniformly from a uniform distribution on [-1, 1] in 1D
357d6685f55SMatthew G. Knepley 
35820f4b53cSBarry Smith   Not Collective
359d6685f55SMatthew G. Knepley 
360a4e35b19SJacob Faibussowitsch   Input Parameters:
3611d27aa22SBarry Smith + p     - A uniform variable on $[0, 1]$
362a4e35b19SJacob Faibussowitsch - dummy - Unused
363d6685f55SMatthew G. Knepley 
364d6685f55SMatthew G. Knepley   Output Parameter:
3651d27aa22SBarry Smith . x - Coordinate in $[-1, 1]$
366d6685f55SMatthew G. Knepley 
367d6685f55SMatthew G. Knepley   Level: beginner
368d6685f55SMatthew G. Knepley 
369db781477SPatrick Sanan .seealso: `PetscPDFConstant1D()`, `PetscCDFConstant1D()`, `PetscPDFSampleConstant2D()`, `PetscPDFSampleConstant3D()`
370d6685f55SMatthew G. Knepley @*/
371d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFSampleConstant1D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
372d71ae5a4SJacob Faibussowitsch {
373ea1b28ebSMatthew G. Knepley   x[0] = 2. * p[0] - 1.;
3743ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
375ea1b28ebSMatthew G. Knepley }
376ea1b28ebSMatthew G. Knepley 
377ea1b28ebSMatthew G. Knepley /*@
378ea1b28ebSMatthew G. Knepley   PetscPDFConstant2D - PDF for the uniform distribution in 2D
379ea1b28ebSMatthew G. Knepley 
38020f4b53cSBarry Smith   Not Collective
381ea1b28ebSMatthew G. Knepley 
382a4e35b19SJacob Faibussowitsch   Input Parameters:
3831d27aa22SBarry Smith + x     - Coordinate in $[-1, 1]^2$
384a4e35b19SJacob Faibussowitsch - dummy - Unused
385ea1b28ebSMatthew G. Knepley 
386ea1b28ebSMatthew G. Knepley   Output Parameter:
3871d27aa22SBarry Smith . p - The probability density at `x`
388ea1b28ebSMatthew G. Knepley 
389ea1b28ebSMatthew G. Knepley   Level: beginner
390ea1b28ebSMatthew G. Knepley 
391db781477SPatrick Sanan .seealso: `PetscCDFConstant2D()`, `PetscPDFSampleConstant2D()`, `PetscPDFConstant1D()`, `PetscPDFConstant3D()`
392ea1b28ebSMatthew G. Knepley @*/
393d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFConstant2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
394d71ae5a4SJacob Faibussowitsch {
395ea1b28ebSMatthew G. Knepley   p[0] = x[0] > -1. && x[0] <= 1. && x[1] > -1. && x[1] <= 1. ? 0.25 : 0.;
3963ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
397ea1b28ebSMatthew G. Knepley }
398ea1b28ebSMatthew G. Knepley 
399ea1b28ebSMatthew G. Knepley /*@
400ea1b28ebSMatthew G. Knepley   PetscCDFConstant2D - CDF for the uniform distribution in 2D
401ea1b28ebSMatthew G. Knepley 
40220f4b53cSBarry Smith   Not Collective
403ea1b28ebSMatthew G. Knepley 
404a4e35b19SJacob Faibussowitsch   Input Parameters:
4051d27aa22SBarry Smith + x     - Coordinate in $[-1, 1]^2$
406a4e35b19SJacob Faibussowitsch - dummy - Unused
407ea1b28ebSMatthew G. Knepley 
408ea1b28ebSMatthew G. Knepley   Output Parameter:
4091d27aa22SBarry Smith . p - The cumulative probability at `x`
410ea1b28ebSMatthew G. Knepley 
411ea1b28ebSMatthew G. Knepley   Level: beginner
412ea1b28ebSMatthew G. Knepley 
413db781477SPatrick Sanan .seealso: `PetscPDFConstant2D()`, `PetscPDFSampleConstant2D()`, `PetscCDFConstant1D()`, `PetscCDFConstant3D()`
414ea1b28ebSMatthew G. Knepley @*/
415d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCDFConstant2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
416d71ae5a4SJacob Faibussowitsch {
417ea1b28ebSMatthew 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.));
4183ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
419ea1b28ebSMatthew G. Knepley }
420ea1b28ebSMatthew G. Knepley 
421ea1b28ebSMatthew G. Knepley /*@
4221d27aa22SBarry Smith   PetscPDFSampleConstant2D - Sample uniformly from a uniform distribution on $[-1, 1]^2$ in 2D
423ea1b28ebSMatthew G. Knepley 
42420f4b53cSBarry Smith   Not Collective
425ea1b28ebSMatthew G. Knepley 
426a4e35b19SJacob Faibussowitsch   Input Parameters:
4271d27aa22SBarry Smith + p     - Two uniform variables on $[0, 1]$
428a4e35b19SJacob Faibussowitsch - dummy - Unused
429ea1b28ebSMatthew G. Knepley 
430ea1b28ebSMatthew G. Knepley   Output Parameter:
4311d27aa22SBarry Smith . x - Coordinate in $[-1, 1]^2$
432ea1b28ebSMatthew G. Knepley 
433ea1b28ebSMatthew G. Knepley   Level: beginner
434ea1b28ebSMatthew G. Knepley 
435db781477SPatrick Sanan .seealso: `PetscPDFConstant2D()`, `PetscCDFConstant2D()`, `PetscPDFSampleConstant1D()`, `PetscPDFSampleConstant3D()`
436ea1b28ebSMatthew G. Knepley @*/
437d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFSampleConstant2D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
438d71ae5a4SJacob Faibussowitsch {
439ea1b28ebSMatthew G. Knepley   x[0] = 2. * p[0] - 1.;
440ea1b28ebSMatthew G. Knepley   x[1] = 2. * p[1] - 1.;
4413ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
442ea1b28ebSMatthew G. Knepley }
443ea1b28ebSMatthew G. Knepley 
444ea1b28ebSMatthew G. Knepley /*@
445ea1b28ebSMatthew G. Knepley   PetscPDFConstant3D - PDF for the uniform distribution in 3D
446ea1b28ebSMatthew G. Knepley 
44720f4b53cSBarry Smith   Not Collective
448ea1b28ebSMatthew G. Knepley 
449a4e35b19SJacob Faibussowitsch   Input Parameters:
4501d27aa22SBarry Smith + x     - Coordinate in $[-1, 1]^3$
451a4e35b19SJacob Faibussowitsch - dummy - Unused
452ea1b28ebSMatthew G. Knepley 
453ea1b28ebSMatthew G. Knepley   Output Parameter:
4541d27aa22SBarry Smith . p - The probability density at `x`
455ea1b28ebSMatthew G. Knepley 
456ea1b28ebSMatthew G. Knepley   Level: beginner
457ea1b28ebSMatthew G. Knepley 
458db781477SPatrick Sanan .seealso: `PetscCDFConstant3D()`, `PetscPDFSampleConstant3D()`, `PetscPDFSampleConstant1D()`, `PetscPDFSampleConstant2D()`
459ea1b28ebSMatthew G. Knepley @*/
460d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFConstant3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
461d71ae5a4SJacob Faibussowitsch {
462ea1b28ebSMatthew 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.;
4633ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
464ea1b28ebSMatthew G. Knepley }
465ea1b28ebSMatthew G. Knepley 
466ea1b28ebSMatthew G. Knepley /*@
467ea1b28ebSMatthew G. Knepley   PetscCDFConstant3D - CDF for the uniform distribution in 3D
468ea1b28ebSMatthew G. Knepley 
46920f4b53cSBarry Smith   Not Collective
470ea1b28ebSMatthew G. Knepley 
471a4e35b19SJacob Faibussowitsch   Input Parameters:
4721d27aa22SBarry Smith + x     - Coordinate in $[-1, 1]^3$
473a4e35b19SJacob Faibussowitsch - dummy - Unused
474ea1b28ebSMatthew G. Knepley 
475ea1b28ebSMatthew G. Knepley   Output Parameter:
4761d27aa22SBarry Smith . p - The cumulative probability at `x`
477ea1b28ebSMatthew G. Knepley 
478ea1b28ebSMatthew G. Knepley   Level: beginner
479ea1b28ebSMatthew G. Knepley 
480db781477SPatrick Sanan .seealso: `PetscPDFConstant3D()`, `PetscPDFSampleConstant3D()`, `PetscCDFConstant1D()`, `PetscCDFConstant2D()`
481ea1b28ebSMatthew G. Knepley @*/
482d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscCDFConstant3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[])
483d71ae5a4SJacob Faibussowitsch {
484ea1b28ebSMatthew 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.));
4853ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
486ea1b28ebSMatthew G. Knepley }
487ea1b28ebSMatthew G. Knepley 
488ea1b28ebSMatthew G. Knepley /*@
4891d27aa22SBarry Smith   PetscPDFSampleConstant3D - Sample uniformly from a uniform distribution on $[-1, 1]^3$ in 3D
490ea1b28ebSMatthew G. Knepley 
49120f4b53cSBarry Smith   Not Collective
492ea1b28ebSMatthew G. Knepley 
493a4e35b19SJacob Faibussowitsch   Input Parameters:
4941d27aa22SBarry Smith + p     - Three uniform variables on $[0, 1]$
495a4e35b19SJacob Faibussowitsch - dummy - Unused
496ea1b28ebSMatthew G. Knepley 
497ea1b28ebSMatthew G. Knepley   Output Parameter:
4981d27aa22SBarry Smith . x - Coordinate in $[-1, 1]^3$
499ea1b28ebSMatthew G. Knepley 
500ea1b28ebSMatthew G. Knepley   Level: beginner
501ea1b28ebSMatthew G. Knepley 
502db781477SPatrick Sanan .seealso: `PetscPDFConstant3D()`, `PetscCDFConstant3D()`, `PetscPDFSampleConstant1D()`, `PetscPDFSampleConstant2D()`
503ea1b28ebSMatthew G. Knepley @*/
504d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscPDFSampleConstant3D(const PetscReal p[], const PetscReal dummy[], PetscReal x[])
505d71ae5a4SJacob Faibussowitsch {
506ea1b28ebSMatthew G. Knepley   x[0] = 2. * p[0] - 1.;
507ea1b28ebSMatthew G. Knepley   x[1] = 2. * p[1] - 1.;
508ea1b28ebSMatthew G. Knepley   x[2] = 2. * p[2] - 1.;
5093ba16761SJacob Faibussowitsch   return PETSC_SUCCESS;
510d6685f55SMatthew G. Knepley }
511d6685f55SMatthew G. Knepley 
512d6685f55SMatthew G. Knepley /*@C
513da81f932SPierre Jolivet   PetscProbCreateFromOptions - Return the probability distribution specified by the arguments and options
514d6685f55SMatthew G. Knepley 
51520f4b53cSBarry Smith   Not Collective
516d6685f55SMatthew G. Knepley 
517d6685f55SMatthew G. Knepley   Input Parameters:
518d6685f55SMatthew G. Knepley + dim    - The dimension of sample points
5191d27aa22SBarry Smith . prefix - The options prefix, or `NULL`
520*f13dfd9eSBarry Smith - name   - The options database name for the probability distribution type
521d6685f55SMatthew G. Knepley 
522d6685f55SMatthew G. Knepley   Output Parameters:
523*f13dfd9eSBarry Smith + pdf     - The PDF of this type, or `NULL`
524*f13dfd9eSBarry Smith . cdf     - The CDF of this type, or `NULL`
525*f13dfd9eSBarry Smith - sampler - The PDF sampler of this type, or `NULL`
526d6685f55SMatthew G. Knepley 
527d6685f55SMatthew G. Knepley   Level: intermediate
528d6685f55SMatthew G. Knepley 
529dce8aebaSBarry Smith .seealso: `PetscProbFunc`, `PetscPDFMaxwellBoltzmann1D()`, `PetscPDFGaussian1D()`, `PetscPDFConstant1D()`
530d6685f55SMatthew G. Knepley @*/
531d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscProbCreateFromOptions(PetscInt dim, const char prefix[], const char name[], PetscProbFunc *pdf, PetscProbFunc *cdf, PetscProbFunc *sampler)
532d71ae5a4SJacob Faibussowitsch {
533d6685f55SMatthew G. Knepley   DTProbDensityType den = DTPROB_DENSITY_GAUSSIAN;
534d6685f55SMatthew G. Knepley 
535d6685f55SMatthew G. Knepley   PetscFunctionBegin;
536d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_SELF, prefix, "PetscProb Options", "DT");
5379566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum(name, "Method to compute PDF <constant, gaussian>", "", DTProbDensityTypes, (PetscEnum)den, (PetscEnum *)&den, NULL));
538d0609cedSBarry Smith   PetscOptionsEnd();
539d6685f55SMatthew G. Knepley 
5409371c9d4SSatish Balay   if (pdf) {
5414f572ea9SToby Isaac     PetscAssertPointer(pdf, 4);
5429371c9d4SSatish Balay     *pdf = NULL;
5439371c9d4SSatish Balay   }
5449371c9d4SSatish Balay   if (cdf) {
5454f572ea9SToby Isaac     PetscAssertPointer(cdf, 5);
5469371c9d4SSatish Balay     *cdf = NULL;
5479371c9d4SSatish Balay   }
5489371c9d4SSatish Balay   if (sampler) {
5494f572ea9SToby Isaac     PetscAssertPointer(sampler, 6);
5509371c9d4SSatish Balay     *sampler = NULL;
5519371c9d4SSatish Balay   }
552d6685f55SMatthew G. Knepley   switch (den) {
553d6685f55SMatthew G. Knepley   case DTPROB_DENSITY_CONSTANT:
554d6685f55SMatthew G. Knepley     switch (dim) {
555d6685f55SMatthew G. Knepley     case 1:
556d6685f55SMatthew G. Knepley       if (pdf) *pdf = PetscPDFConstant1D;
557d6685f55SMatthew G. Knepley       if (cdf) *cdf = PetscCDFConstant1D;
558d6685f55SMatthew G. Knepley       if (sampler) *sampler = PetscPDFSampleConstant1D;
559d6685f55SMatthew G. Knepley       break;
560ea1b28ebSMatthew G. Knepley     case 2:
561ea1b28ebSMatthew G. Knepley       if (pdf) *pdf = PetscPDFConstant2D;
562ea1b28ebSMatthew G. Knepley       if (cdf) *cdf = PetscCDFConstant2D;
563ea1b28ebSMatthew G. Knepley       if (sampler) *sampler = PetscPDFSampleConstant2D;
564ea1b28ebSMatthew G. Knepley       break;
565ea1b28ebSMatthew G. Knepley     case 3:
566ea1b28ebSMatthew G. Knepley       if (pdf) *pdf = PetscPDFConstant3D;
567ea1b28ebSMatthew G. Knepley       if (cdf) *cdf = PetscCDFConstant3D;
568ea1b28ebSMatthew G. Knepley       if (sampler) *sampler = PetscPDFSampleConstant3D;
569ea1b28ebSMatthew G. Knepley       break;
570d71ae5a4SJacob Faibussowitsch     default:
571d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]);
572d6685f55SMatthew G. Knepley     }
573d6685f55SMatthew G. Knepley     break;
574d6685f55SMatthew G. Knepley   case DTPROB_DENSITY_GAUSSIAN:
575d6685f55SMatthew G. Knepley     switch (dim) {
576d6685f55SMatthew G. Knepley     case 1:
577d6685f55SMatthew G. Knepley       if (pdf) *pdf = PetscPDFGaussian1D;
578d6685f55SMatthew G. Knepley       if (cdf) *cdf = PetscCDFGaussian1D;
579d6685f55SMatthew G. Knepley       if (sampler) *sampler = PetscPDFSampleGaussian1D;
580d6685f55SMatthew G. Knepley       break;
581d6685f55SMatthew G. Knepley     case 2:
582d6685f55SMatthew G. Knepley       if (pdf) *pdf = PetscPDFGaussian2D;
583d6685f55SMatthew G. Knepley       if (sampler) *sampler = PetscPDFSampleGaussian2D;
584d6685f55SMatthew G. Knepley       break;
585ea1b28ebSMatthew G. Knepley     case 3:
586ea1b28ebSMatthew G. Knepley       if (pdf) *pdf = PetscPDFGaussian3D;
587ea1b28ebSMatthew G. Knepley       if (sampler) *sampler = PetscPDFSampleGaussian3D;
588ea1b28ebSMatthew G. Knepley       break;
589d71ae5a4SJacob Faibussowitsch     default:
590d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]);
591d6685f55SMatthew G. Knepley     }
592d6685f55SMatthew G. Knepley     break;
593d6685f55SMatthew G. Knepley   case DTPROB_DENSITY_MAXWELL_BOLTZMANN:
594d6685f55SMatthew G. Knepley     switch (dim) {
595d6685f55SMatthew G. Knepley     case 1:
596d6685f55SMatthew G. Knepley       if (pdf) *pdf = PetscPDFMaxwellBoltzmann1D;
597d6685f55SMatthew G. Knepley       if (cdf) *cdf = PetscCDFMaxwellBoltzmann1D;
598d6685f55SMatthew G. Knepley       break;
599d6685f55SMatthew G. Knepley     case 2:
600d6685f55SMatthew G. Knepley       if (pdf) *pdf = PetscPDFMaxwellBoltzmann2D;
601d6685f55SMatthew G. Knepley       if (cdf) *cdf = PetscCDFMaxwellBoltzmann2D;
602d6685f55SMatthew G. Knepley       break;
603d6685f55SMatthew G. Knepley     case 3:
604d6685f55SMatthew G. Knepley       if (pdf) *pdf = PetscPDFMaxwellBoltzmann3D;
605d6685f55SMatthew G. Knepley       if (cdf) *cdf = PetscCDFMaxwellBoltzmann3D;
606d6685f55SMatthew G. Knepley       break;
607d71ae5a4SJacob Faibussowitsch     default:
608d71ae5a4SJacob Faibussowitsch       SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]);
609d6685f55SMatthew G. Knepley     }
610d6685f55SMatthew G. Knepley     break;
611d71ae5a4SJacob Faibussowitsch   default:
612d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Density type %s is not supported", DTProbDensityTypes[PetscMax(0, PetscMin(den, DTPROB_NUM_DENSITY))]);
613d6685f55SMatthew G. Knepley   }
6143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
615d6685f55SMatthew G. Knepley }
616d6685f55SMatthew G. Knepley 
617d6685f55SMatthew G. Knepley #ifdef PETSC_HAVE_KS
618d6685f55SMatthew G. Knepley EXTERN_C_BEGIN
619d6685f55SMatthew G. Knepley   #include <KolmogorovSmirnovDist.h>
620d6685f55SMatthew G. Knepley EXTERN_C_END
621d6685f55SMatthew G. Knepley #endif
622d6685f55SMatthew G. Knepley 
623d6685f55SMatthew G. Knepley /*@C
624d6685f55SMatthew G. Knepley   PetscProbComputeKSStatistic - Compute the Kolmogorov-Smirnov statistic for the empirical distribution for an input vector, compared to an analytic CDF.
625d6685f55SMatthew G. Knepley 
62620f4b53cSBarry Smith   Collective
627d6685f55SMatthew G. Knepley 
628d6685f55SMatthew G. Knepley   Input Parameters:
629d6685f55SMatthew G. Knepley + v   - The data vector, blocksize is the sample dimension
630d6685f55SMatthew G. Knepley - cdf - The analytic CDF
631d6685f55SMatthew G. Knepley 
632d6685f55SMatthew G. Knepley   Output Parameter:
633d6685f55SMatthew G. Knepley . alpha - The KS statisic
634d6685f55SMatthew G. Knepley 
635d6685f55SMatthew G. Knepley   Level: advanced
636d6685f55SMatthew G. Knepley 
637dce8aebaSBarry Smith   Notes:
638dce8aebaSBarry Smith   The Kolmogorov-Smirnov statistic for a given cumulative distribution function $F(x)$ is
6391d27aa22SBarry Smith 
6401d27aa22SBarry Smith   $$
641dce8aebaSBarry Smith   D_n = \sup_x \left| F_n(x) - F(x) \right|
6421d27aa22SBarry Smith   $$
643dce8aebaSBarry Smith 
644dce8aebaSBarry 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
645dce8aebaSBarry Smith 
6461d27aa22SBarry Smith   $$
647dce8aebaSBarry Smith   F_n =  # of samples <= x / n
6481d27aa22SBarry Smith   $$
649dce8aebaSBarry Smith 
650d6685f55SMatthew G. Knepley   The empirical distribution function $F_n(x)$ is discrete, and thus had a ``stairstep''
651d6685f55SMatthew G. Knepley   cumulative distribution, making $n$ the number of stairs. Intuitively, the statistic takes
652d6685f55SMatthew G. Knepley   the largest absolute difference between the two distribution functions across all $x$ values.
653d6685f55SMatthew G. Knepley 
654db781477SPatrick Sanan .seealso: `PetscProbFunc`
655d6685f55SMatthew G. Knepley @*/
656d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscProbComputeKSStatistic(Vec v, PetscProbFunc cdf, PetscReal *alpha)
657d71ae5a4SJacob Faibussowitsch {
658d6685f55SMatthew G. Knepley #if !defined(PETSC_HAVE_KS)
659d6685f55SMatthew G. Knepley   SETERRQ(PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "No support for Kolmogorov-Smirnov test.\nReconfigure using --download-ks");
660d6685f55SMatthew G. Knepley #else
661d6685f55SMatthew G. Knepley   PetscViewer        viewer = NULL;
662d6685f55SMatthew G. Knepley   PetscViewerFormat  format;
663d6685f55SMatthew G. Knepley   PetscDraw          draw;
664d6685f55SMatthew G. Knepley   PetscDrawLG        lg;
665d6685f55SMatthew G. Knepley   PetscDrawAxis      axis;
666d6685f55SMatthew G. Knepley   const PetscScalar *a;
667d6685f55SMatthew G. Knepley   PetscReal         *speed, Dn = PETSC_MIN_REAL;
668d6685f55SMatthew G. Knepley   PetscBool          isascii = PETSC_FALSE, isdraw = PETSC_FALSE, flg;
669d6685f55SMatthew G. Knepley   PetscInt           n, p, dim, d;
670d6685f55SMatthew G. Knepley   PetscMPIInt        size;
671d6685f55SMatthew G. Knepley   const char        *names[2] = {"Analytic", "Empirical"};
672d6685f55SMatthew G. Knepley   char               title[PETSC_MAX_PATH_LEN];
673d6685f55SMatthew G. Knepley   PetscOptions       options;
674d6685f55SMatthew G. Knepley   const char        *prefix;
675d6685f55SMatthew G. Knepley   MPI_Comm           comm;
676d6685f55SMatthew G. Knepley 
677d6685f55SMatthew G. Knepley   PetscFunctionBegin;
6789566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)v, &comm));
6799566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)v, &prefix));
6809566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetOptions((PetscObject)v, &options));
6819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetViewer(comm, options, prefix, "-ks_monitor", &viewer, &format, &flg));
682d6685f55SMatthew G. Knepley   if (flg) {
6839566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
6849566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw));
685d6685f55SMatthew G. Knepley   }
686d6685f55SMatthew G. Knepley   if (isascii) {
6879566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
688d6685f55SMatthew G. Knepley   } else if (isdraw) {
6899566063dSJacob Faibussowitsch     PetscCall(PetscViewerPushFormat(viewer, format));
6909566063dSJacob Faibussowitsch     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
6919566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGCreate(draw, 2, &lg));
6929566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGSetLegend(lg, names));
693d6685f55SMatthew G. Knepley   }
694d6685f55SMatthew G. Knepley 
6959566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)v), &size));
6969566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(v, &n));
6979566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(v, &dim));
698d6685f55SMatthew G. Knepley   n /= dim;
699d6685f55SMatthew G. Knepley   /* TODO Parallel algorithm is harder */
700d6685f55SMatthew G. Knepley   if (size == 1) {
7019566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &speed));
7029566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(v, &a));
703d6685f55SMatthew G. Knepley     for (p = 0; p < n; ++p) {
704d6685f55SMatthew G. Knepley       PetscReal mag = 0.;
705d6685f55SMatthew G. Knepley 
706d6685f55SMatthew G. Knepley       for (d = 0; d < dim; ++d) mag += PetscSqr(PetscRealPart(a[p * dim + d]));
707d6685f55SMatthew G. Knepley       speed[p] = PetscSqrtReal(mag);
708d6685f55SMatthew G. Knepley     }
7099566063dSJacob Faibussowitsch     PetscCall(PetscSortReal(n, speed));
7109566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(v, &a));
711d6685f55SMatthew G. Knepley     for (p = 0; p < n; ++p) {
712d6685f55SMatthew G. Knepley       const PetscReal x = speed[p], Fn = ((PetscReal)p) / n;
713d6685f55SMatthew G. Knepley       PetscReal       F, vals[2];
714d6685f55SMatthew G. Knepley 
7159566063dSJacob Faibussowitsch       PetscCall(cdf(&x, NULL, &F));
716d6685f55SMatthew G. Knepley       Dn = PetscMax(PetscAbsReal(Fn - F), Dn);
7179566063dSJacob Faibussowitsch       if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "x: %g F: %g Fn: %g Dn: %.2g\n", x, F, Fn, Dn));
7189371c9d4SSatish Balay       if (isdraw) {
7199371c9d4SSatish Balay         vals[0] = F;
7209371c9d4SSatish Balay         vals[1] = Fn;
7219371c9d4SSatish Balay         PetscCall(PetscDrawLGAddCommonPoint(lg, x, vals));
7229371c9d4SSatish Balay       }
723d6685f55SMatthew G. Knepley     }
724d6685f55SMatthew G. Knepley     if (speed[n - 1] < 6.) {
725d6685f55SMatthew G. Knepley       const PetscReal k = (PetscInt)(6. - speed[n - 1]) + 1, dx = (6. - speed[n - 1]) / k;
726d6685f55SMatthew G. Knepley       for (p = 0; p <= k; ++p) {
727d6685f55SMatthew G. Knepley         const PetscReal x = speed[n - 1] + p * dx, Fn = 1.0;
728d6685f55SMatthew G. Knepley         PetscReal       F, vals[2];
729d6685f55SMatthew G. Knepley 
7309566063dSJacob Faibussowitsch         PetscCall(cdf(&x, NULL, &F));
731d6685f55SMatthew G. Knepley         Dn = PetscMax(PetscAbsReal(Fn - F), Dn);
7329566063dSJacob Faibussowitsch         if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "x: %g F: %g Fn: %g Dn: %.2g\n", x, F, Fn, Dn));
7339371c9d4SSatish Balay         if (isdraw) {
7349371c9d4SSatish Balay           vals[0] = F;
7359371c9d4SSatish Balay           vals[1] = Fn;
7369371c9d4SSatish Balay           PetscCall(PetscDrawLGAddCommonPoint(lg, x, vals));
7379371c9d4SSatish Balay         }
738d6685f55SMatthew G. Knepley       }
739d6685f55SMatthew G. Knepley     }
7409566063dSJacob Faibussowitsch     PetscCall(PetscFree(speed));
741d6685f55SMatthew G. Knepley   }
742d6685f55SMatthew G. Knepley   if (isdraw) {
7439566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGGetAxis(lg, &axis));
7449566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(title, PETSC_MAX_PATH_LEN, "Kolmogorov-Smirnov Test (Dn %.2g)", Dn));
7459566063dSJacob Faibussowitsch     PetscCall(PetscDrawAxisSetLabels(axis, title, "x", "CDF(x)"));
7469566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDraw(lg));
7479566063dSJacob Faibussowitsch     PetscCall(PetscDrawLGDestroy(&lg));
748d6685f55SMatthew G. Knepley   }
749d6685f55SMatthew G. Knepley   if (viewer) {
7509566063dSJacob Faibussowitsch     PetscCall(PetscViewerFlush(viewer));
7519566063dSJacob Faibussowitsch     PetscCall(PetscViewerPopFormat(viewer));
752cd791dc2SBarry Smith     PetscCall(PetscOptionsRestoreViewer(&viewer));
753d6685f55SMatthew G. Knepley   }
754d6685f55SMatthew G. Knepley   *alpha = KSfbar((int)n, (double)Dn);
7553ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
756d6685f55SMatthew G. Knepley #endif
757d6685f55SMatthew G. Knepley }
758