xref: /petsc/src/dm/dt/interface/dtprob.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
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]$
16*2a8381b2SBarry Smith - unused - 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 @*/
PetscPDFMaxwellBoltzmann1D(const PetscReal x[],const PetscReal unused[],PetscReal p[])25*2a8381b2SBarry Smith PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal x[], const PetscReal unused[], 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]$
38*2a8381b2SBarry Smith - unused - 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 @*/
PetscCDFMaxwellBoltzmann1D(const PetscReal x[],const PetscReal unused[],PetscReal p[])47*2a8381b2SBarry Smith PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal x[], const PetscReal unused[], 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]$
60*2a8381b2SBarry Smith - unused - 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 @*/
PetscPDFMaxwellBoltzmann2D(const PetscReal x[],const PetscReal unused[],PetscReal p[])69*2a8381b2SBarry Smith PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal x[], const PetscReal unused[], 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]$
82*2a8381b2SBarry Smith - unused - 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 @*/
PetscCDFMaxwellBoltzmann2D(const PetscReal x[],const PetscReal unused[],PetscReal p[])91*2a8381b2SBarry Smith PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal x[], const PetscReal unused[], 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]$
104*2a8381b2SBarry Smith - unused - 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 @*/
PetscPDFMaxwellBoltzmann3D(const PetscReal x[],const PetscReal unused[],PetscReal p[])113*2a8381b2SBarry Smith PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal x[], const PetscReal unused[], 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]$
126*2a8381b2SBarry Smith - unused - 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 @*/
PetscCDFMaxwellBoltzmann3D(const PetscReal x[],const PetscReal unused[],PetscReal p[])135*2a8381b2SBarry Smith PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal x[], const PetscReal unused[], 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 @*/
PetscPDFGaussian1D(const PetscReal x[],const PetscReal scale[],PetscReal p[])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 
PetscCDFGaussian1D(const PetscReal x[],const PetscReal scale[],PetscReal p[])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]$
178*2a8381b2SBarry Smith - unused - 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 @*/
PetscPDFSampleGaussian1D(const PetscReal p[],const PetscReal unused[],PetscReal x[])191*2a8381b2SBarry Smith PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal p[], const PetscReal unused[], 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$
225*2a8381b2SBarry Smith - unused - 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 @*/
PetscPDFGaussian2D(const PetscReal x[],const PetscReal unused[],PetscReal p[])234*2a8381b2SBarry Smith PetscErrorCode PetscPDFGaussian2D(const PetscReal x[], const PetscReal unused[], 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$
248*2a8381b2SBarry Smith - unused - 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 @*/
PetscPDFSampleGaussian2D(const PetscReal p[],const PetscReal unused[],PetscReal x[])257*2a8381b2SBarry Smith PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal p[], const PetscReal unused[], 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$
272*2a8381b2SBarry Smith - unused - 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 @*/
PetscPDFGaussian3D(const PetscReal x[],const PetscReal unused[],PetscReal p[])281*2a8381b2SBarry Smith PetscErrorCode PetscPDFGaussian3D(const PetscReal x[], const PetscReal unused[], 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$
295*2a8381b2SBarry Smith - unused - 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 @*/
PetscPDFSampleGaussian3D(const PetscReal p[],const PetscReal unused[],PetscReal x[])304*2a8381b2SBarry Smith PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal p[], const PetscReal unused[], PetscReal x[])
305d71ae5a4SJacob Faibussowitsch {
306*2a8381b2SBarry Smith   PetscCall(PetscPDFSampleGaussian1D(p, unused, x));
307*2a8381b2SBarry Smith   PetscCall(PetscPDFSampleGaussian2D(&p[1], unused, &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]$
318*2a8381b2SBarry Smith - unused - 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 @*/
PetscPDFConstant1D(const PetscReal x[],const PetscReal unused[],PetscReal p[])327*2a8381b2SBarry Smith PetscErrorCode PetscPDFConstant1D(const PetscReal x[], const PetscReal unused[], 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]$
340*2a8381b2SBarry Smith - unused - 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 @*/
PetscCDFConstant1D(const PetscReal x[],const PetscReal unused[],PetscReal p[])349*2a8381b2SBarry Smith PetscErrorCode PetscCDFConstant1D(const PetscReal x[], const PetscReal unused[], 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]$
362*2a8381b2SBarry Smith - unused - 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 @*/
PetscPDFSampleConstant1D(const PetscReal p[],const PetscReal unused[],PetscReal x[])371*2a8381b2SBarry Smith PetscErrorCode PetscPDFSampleConstant1D(const PetscReal p[], const PetscReal unused[], 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$
384*2a8381b2SBarry Smith - unused - 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 @*/
PetscPDFConstant2D(const PetscReal x[],const PetscReal unused[],PetscReal p[])393*2a8381b2SBarry Smith PetscErrorCode PetscPDFConstant2D(const PetscReal x[], const PetscReal unused[], 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$
406*2a8381b2SBarry Smith - unused - 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 @*/
PetscCDFConstant2D(const PetscReal x[],const PetscReal unused[],PetscReal p[])415*2a8381b2SBarry Smith PetscErrorCode PetscCDFConstant2D(const PetscReal x[], const PetscReal unused[], 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]$
428*2a8381b2SBarry Smith - unused - 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 @*/
PetscPDFSampleConstant2D(const PetscReal p[],const PetscReal unused[],PetscReal x[])437*2a8381b2SBarry Smith PetscErrorCode PetscPDFSampleConstant2D(const PetscReal p[], const PetscReal unused[], 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$
451*2a8381b2SBarry Smith - unused - 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 @*/
PetscPDFConstant3D(const PetscReal x[],const PetscReal unused[],PetscReal p[])460*2a8381b2SBarry Smith PetscErrorCode PetscPDFConstant3D(const PetscReal x[], const PetscReal unused[], 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$
473*2a8381b2SBarry Smith - unused - 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 @*/
PetscCDFConstant3D(const PetscReal x[],const PetscReal unused[],PetscReal p[])482*2a8381b2SBarry Smith PetscErrorCode PetscCDFConstant3D(const PetscReal x[], const PetscReal unused[], 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]$
495*2a8381b2SBarry Smith - unused - 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 @*/
PetscPDFSampleConstant3D(const PetscReal p[],const PetscReal unused[],PetscReal x[])504*2a8381b2SBarry Smith PetscErrorCode PetscPDFSampleConstant3D(const PetscReal p[], const PetscReal unused[], 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`
520f13dfd9eSBarry Smith - name   - The options database name for the probability distribution type
521d6685f55SMatthew G. Knepley 
522d6685f55SMatthew G. Knepley   Output Parameters:
523f13dfd9eSBarry Smith + pdf     - The PDF of this type, or `NULL`
524f13dfd9eSBarry Smith . cdf     - The CDF of this type, or `NULL`
525f13dfd9eSBarry Smith - sampler - The PDF sampler of this type, or `NULL`
526d6685f55SMatthew G. Knepley 
527d6685f55SMatthew G. Knepley   Level: intermediate
528d6685f55SMatthew G. Knepley 
529f8662bd6SBarry Smith .seealso: `PetscProbFn`, `PetscPDFMaxwellBoltzmann1D()`, `PetscPDFGaussian1D()`, `PetscPDFConstant1D()`
530d6685f55SMatthew G. Knepley @*/
PetscProbCreateFromOptions(PetscInt dim,const char prefix[],const char name[],PetscProbFn ** pdf,PetscProbFn ** cdf,PetscProbFn ** sampler)531f8662bd6SBarry Smith PetscErrorCode PetscProbCreateFromOptions(PetscInt dim, const char prefix[], const char name[], PetscProbFn **pdf, PetscProbFn **cdf, PetscProbFn **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 
623e3ce8211SMatthew G. Knepley typedef enum {
624e3ce8211SMatthew G. Knepley   NONE,
625e3ce8211SMatthew G. Knepley   ASCII,
626e3ce8211SMatthew G. Knepley   DRAW
627e3ce8211SMatthew G. Knepley } OutputType;
628e3ce8211SMatthew G. Knepley 
KSViewerCreate(PetscObject obj,OutputType * outputType,PetscViewer * viewer)629e3ce8211SMatthew G. Knepley static PetscErrorCode KSViewerCreate(PetscObject obj, OutputType *outputType, PetscViewer *viewer)
630e3ce8211SMatthew G. Knepley {
631e3ce8211SMatthew G. Knepley   PetscViewerFormat format;
632e3ce8211SMatthew G. Knepley   PetscOptions      options;
633e3ce8211SMatthew G. Knepley   const char       *prefix;
634e3ce8211SMatthew G. Knepley   PetscBool         flg;
635e3ce8211SMatthew G. Knepley   MPI_Comm          comm;
636e3ce8211SMatthew G. Knepley 
637e3ce8211SMatthew G. Knepley   PetscFunctionBegin;
638e3ce8211SMatthew G. Knepley   *outputType = NONE;
639e3ce8211SMatthew G. Knepley   PetscCall(PetscObjectGetComm(obj, &comm));
640e3ce8211SMatthew G. Knepley   PetscCall(PetscObjectGetOptionsPrefix(obj, &prefix));
641e3ce8211SMatthew G. Knepley   PetscCall(PetscObjectGetOptions(obj, &options));
642e3ce8211SMatthew G. Knepley   PetscCall(PetscOptionsCreateViewer(comm, options, prefix, "-ks_monitor", viewer, &format, &flg));
643e3ce8211SMatthew G. Knepley   if (flg) {
644e3ce8211SMatthew G. Knepley     PetscCall(PetscObjectTypeCompare((PetscObject)*viewer, PETSCVIEWERASCII, &flg));
645e3ce8211SMatthew G. Knepley     if (flg) *outputType = ASCII;
646e3ce8211SMatthew G. Knepley     PetscCall(PetscObjectTypeCompare((PetscObject)*viewer, PETSCVIEWERDRAW, &flg));
647e3ce8211SMatthew G. Knepley     if (flg) *outputType = DRAW;
648e3ce8211SMatthew G. Knepley     PetscCall(PetscViewerPushFormat(*viewer, format));
649e3ce8211SMatthew G. Knepley   }
650e3ce8211SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
651e3ce8211SMatthew G. Knepley }
652e3ce8211SMatthew G. Knepley 
KSViewerDestroy(PetscViewer * viewer)653e3ce8211SMatthew G. Knepley static PetscErrorCode KSViewerDestroy(PetscViewer *viewer)
654e3ce8211SMatthew G. Knepley {
655e3ce8211SMatthew G. Knepley   PetscFunctionBegin;
656e3ce8211SMatthew G. Knepley   if (*viewer) {
657e3ce8211SMatthew G. Knepley     PetscCall(PetscViewerFlush(*viewer));
658e3ce8211SMatthew G. Knepley     PetscCall(PetscViewerPopFormat(*viewer));
659e3ce8211SMatthew G. Knepley     PetscCall(PetscViewerDestroy(viewer));
660e3ce8211SMatthew G. Knepley   }
661e3ce8211SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
662e3ce8211SMatthew G. Knepley }
663e3ce8211SMatthew G. Knepley 
PetscProbComputeKSStatistic_Internal(MPI_Comm comm,PetscInt n,PetscReal val[],PetscReal wgt[],PetscProbFn * cdf,PetscReal * alpha,OutputType outputType,PetscViewer viewer)664f8662bd6SBarry Smith static PetscErrorCode PetscProbComputeKSStatistic_Internal(MPI_Comm comm, PetscInt n, PetscReal val[], PetscReal wgt[], PetscProbFn *cdf, PetscReal *alpha, OutputType outputType, PetscViewer viewer)
665e3ce8211SMatthew G. Knepley {
666e3ce8211SMatthew G. Knepley #if !defined(PETSC_HAVE_KS)
667e3ce8211SMatthew G. Knepley   SETERRQ(comm, PETSC_ERR_SUP, "No support for Kolmogorov-Smirnov test.\nReconfigure using --download-ks");
668e3ce8211SMatthew G. Knepley #else
669e3ce8211SMatthew G. Knepley   PetscDraw     draw;
670e3ce8211SMatthew G. Knepley   PetscDrawLG   lg;
671e3ce8211SMatthew G. Knepley   PetscDrawAxis axis;
672e3ce8211SMatthew G. Knepley   const char   *names[2] = {"Analytic", "Empirical"};
673e3ce8211SMatthew G. Knepley   char          title[PETSC_MAX_PATH_LEN];
674e3ce8211SMatthew G. Knepley   PetscReal     Fn = 0., Dn = PETSC_MIN_REAL;
675e3ce8211SMatthew G. Knepley   PetscMPIInt   size;
676e3ce8211SMatthew G. Knepley 
677e3ce8211SMatthew G. Knepley   PetscFunctionBegin;
678e3ce8211SMatthew G. Knepley   PetscCallMPI(MPI_Comm_size(comm, &size));
679e3ce8211SMatthew G. Knepley   PetscCheck(size == 1, comm, PETSC_ERR_SUP, "Parallel K-S test not yet supported");
680e3ce8211SMatthew G. Knepley   if (outputType == DRAW) {
681e3ce8211SMatthew G. Knepley     PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw));
682e3ce8211SMatthew G. Knepley     PetscCall(PetscDrawLGCreate(draw, 2, &lg));
683e3ce8211SMatthew G. Knepley     PetscCall(PetscDrawLGSetLegend(lg, names));
684e3ce8211SMatthew G. Knepley   }
685e3ce8211SMatthew G. Knepley   if (wgt) {
686e3ce8211SMatthew G. Knepley     PetscReal *tmpv, *tmpw;
687e3ce8211SMatthew G. Knepley     PetscInt  *perm;
688e3ce8211SMatthew G. Knepley 
689e3ce8211SMatthew G. Knepley     PetscCall(PetscMalloc3(n, &perm, n, &tmpv, n, &tmpw));
690e3ce8211SMatthew G. Knepley     for (PetscInt i = 0; i < n; ++i) perm[i] = i;
691e3ce8211SMatthew G. Knepley     PetscCall(PetscSortRealWithPermutation(n, val, perm));
692e3ce8211SMatthew G. Knepley     for (PetscInt i = 0; i < n; ++i) {
693e3ce8211SMatthew G. Knepley       tmpv[i] = val[perm[i]];
694e3ce8211SMatthew G. Knepley       tmpw[i] = wgt[perm[i]];
695e3ce8211SMatthew G. Knepley     }
696e3ce8211SMatthew G. Knepley     for (PetscInt i = 0; i < n; ++i) {
697e3ce8211SMatthew G. Knepley       val[i] = tmpv[i];
698e3ce8211SMatthew G. Knepley       wgt[i] = tmpw[i];
699e3ce8211SMatthew G. Knepley     }
700e3ce8211SMatthew G. Knepley     PetscCall(PetscFree3(perm, tmpv, tmpw));
701e3ce8211SMatthew G. Knepley   } else PetscCall(PetscSortReal(n, val));
702e3ce8211SMatthew G. Knepley   // Compute empirical cumulative distribution F_n and discrepancy D_n
703e3ce8211SMatthew G. Knepley   for (PetscInt p = 0; p < n; ++p) {
704e3ce8211SMatthew G. Knepley     const PetscReal x = val[p];
705e3ce8211SMatthew G. Knepley     const PetscReal w = wgt ? wgt[p] : 1. / n;
706e3ce8211SMatthew G. Knepley     PetscReal       F, vals[2];
707e3ce8211SMatthew G. Knepley 
708e3ce8211SMatthew G. Knepley     Fn += w;
709e3ce8211SMatthew G. Knepley     PetscCall(cdf(&x, NULL, &F));
710e3ce8211SMatthew G. Knepley     Dn = PetscMax(PetscAbsReal(Fn - F), Dn);
711e3ce8211SMatthew G. Knepley     switch (outputType) {
712e3ce8211SMatthew G. Knepley     case ASCII:
713e3ce8211SMatthew G. Knepley       PetscCall(PetscViewerASCIIPrintf(viewer, "x: %g F: %g Fn: %g Dn: %.2g\n", x, F, Fn, Dn));
714e3ce8211SMatthew G. Knepley       break;
715e3ce8211SMatthew G. Knepley     case DRAW:
716e3ce8211SMatthew G. Knepley       vals[0] = F;
717e3ce8211SMatthew G. Knepley       vals[1] = Fn;
718e3ce8211SMatthew G. Knepley       PetscCall(PetscDrawLGAddCommonPoint(lg, x, vals));
719e3ce8211SMatthew G. Knepley       break;
720e3ce8211SMatthew G. Knepley     case NONE:
721e3ce8211SMatthew G. Knepley       break;
722e3ce8211SMatthew G. Knepley     }
723e3ce8211SMatthew G. Knepley   }
724e3ce8211SMatthew G. Knepley   if (outputType == DRAW) {
725e3ce8211SMatthew G. Knepley     PetscCall(PetscDrawLGGetAxis(lg, &axis));
726e3ce8211SMatthew G. Knepley     PetscCall(PetscSNPrintf(title, PETSC_MAX_PATH_LEN, "Kolmogorov-Smirnov Test (Dn %.2g)", Dn));
727e3ce8211SMatthew G. Knepley     PetscCall(PetscDrawAxisSetLabels(axis, title, "x", "CDF(x)"));
728e3ce8211SMatthew G. Knepley     PetscCall(PetscDrawLGDraw(lg));
729e3ce8211SMatthew G. Knepley     PetscCall(PetscDrawLGDestroy(&lg));
730e3ce8211SMatthew G. Knepley   }
731e3ce8211SMatthew G. Knepley   *alpha = KSfbar((int)n, (double)Dn);
732e3ce8211SMatthew G. Knepley   if (outputType == ASCII) PetscCall(PetscViewerASCIIPrintf(viewer, "KSfbar(%" PetscInt_FMT ", %.2g): %g\n", n, Dn, *alpha));
733e3ce8211SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
734e3ce8211SMatthew G. Knepley #endif
735e3ce8211SMatthew G. Knepley }
736e3ce8211SMatthew G. Knepley 
737d6685f55SMatthew G. Knepley /*@C
738d6685f55SMatthew G. Knepley   PetscProbComputeKSStatistic - Compute the Kolmogorov-Smirnov statistic for the empirical distribution for an input vector, compared to an analytic CDF.
739d6685f55SMatthew G. Knepley 
74020f4b53cSBarry Smith   Collective
741d6685f55SMatthew G. Knepley 
742d6685f55SMatthew G. Knepley   Input Parameters:
743d6685f55SMatthew G. Knepley + v   - The data vector, blocksize is the sample dimension
744d6685f55SMatthew G. Knepley - cdf - The analytic CDF
745d6685f55SMatthew G. Knepley 
746d6685f55SMatthew G. Knepley   Output Parameter:
747d7c1f440SPierre Jolivet . alpha - The KS statistic
748d6685f55SMatthew G. Knepley 
749d6685f55SMatthew G. Knepley   Level: advanced
750d6685f55SMatthew G. Knepley 
751dce8aebaSBarry Smith   Notes:
752dce8aebaSBarry Smith   The Kolmogorov-Smirnov statistic for a given cumulative distribution function $F(x)$ is
7531d27aa22SBarry Smith 
7541d27aa22SBarry Smith   $$
755dce8aebaSBarry Smith   D_n = \sup_x \left| F_n(x) - F(x) \right|
7561d27aa22SBarry Smith   $$
757dce8aebaSBarry Smith 
758dce8aebaSBarry 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
759dce8aebaSBarry Smith 
7601d27aa22SBarry Smith   $$
761dce8aebaSBarry Smith   F_n =  # of samples <= x / n
7621d27aa22SBarry Smith   $$
763dce8aebaSBarry Smith 
764d6685f55SMatthew G. Knepley   The empirical distribution function $F_n(x)$ is discrete, and thus had a ``stairstep''
765d6685f55SMatthew G. Knepley   cumulative distribution, making $n$ the number of stairs. Intuitively, the statistic takes
766d6685f55SMatthew G. Knepley   the largest absolute difference between the two distribution functions across all $x$ values.
767d6685f55SMatthew G. Knepley 
768e3ce8211SMatthew G. Knepley   The goodness-of-fit test, or Kolmogorov-Smirnov test, is constructed using the Kolmogorov
769e3ce8211SMatthew G. Knepley   distribution. It rejects the null hypothesis at level $\alpha$ if
770e3ce8211SMatthew G. Knepley 
771e3ce8211SMatthew G. Knepley   $$
772e3ce8211SMatthew G. Knepley   \sqrt{n} D_{n} > K_{\alpha},
773e3ce8211SMatthew G. Knepley   $$
774e3ce8211SMatthew G. Knepley 
775e3ce8211SMatthew G. Knepley   where $K_\alpha$ is found from
776e3ce8211SMatthew G. Knepley 
777e3ce8211SMatthew G. Knepley   $$
778e3ce8211SMatthew G. Knepley   \operatorname{Pr}(K \leq K_{\alpha}) = 1 - \alpha.
779e3ce8211SMatthew G. Knepley   $$
780e3ce8211SMatthew G. Knepley 
781e3ce8211SMatthew G. Knepley   This means that getting a small alpha says that we have high confidence that the data did not come
782e3ce8211SMatthew G. Knepley   from the input distribution, so we say that it rejects the null hypothesis.
783e3ce8211SMatthew G. Knepley 
784f8662bd6SBarry Smith .seealso: `PetscProbComputeKSStatisticWeighted()`, `PetscProbComputeKSStatisticMagnitude()`, `PetscProbFn`
785d6685f55SMatthew G. Knepley @*/
PetscProbComputeKSStatistic(Vec v,PetscProbFn * cdf,PetscReal * alpha)786f8662bd6SBarry Smith PetscErrorCode PetscProbComputeKSStatistic(Vec v, PetscProbFn *cdf, PetscReal *alpha)
787d71ae5a4SJacob Faibussowitsch {
788d6685f55SMatthew G. Knepley   PetscViewer        viewer     = NULL;
789e3ce8211SMatthew G. Knepley   OutputType         outputType = NONE;
790e3ce8211SMatthew G. Knepley   const PetscScalar *val;
791e3ce8211SMatthew G. Knepley   PetscInt           n;
792d6685f55SMatthew G. Knepley 
793d6685f55SMatthew G. Knepley   PetscFunctionBegin;
794e3ce8211SMatthew G. Knepley   PetscCall(KSViewerCreate((PetscObject)v, &outputType, &viewer));
795e3ce8211SMatthew G. Knepley   PetscCall(VecGetLocalSize(v, &n));
796e3ce8211SMatthew G. Knepley   PetscCall(VecGetArrayRead(v, &val));
797e3ce8211SMatthew G. Knepley   PetscCheck(!PetscDefined(USE_COMPLEX), PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "K-S test does not support complex");
798e3ce8211SMatthew G. Knepley   PetscCall(PetscProbComputeKSStatistic_Internal(PetscObjectComm((PetscObject)v), n, (PetscReal *)val, NULL, cdf, alpha, outputType, viewer));
799e3ce8211SMatthew G. Knepley   PetscCall(VecRestoreArrayRead(v, &val));
800e3ce8211SMatthew G. Knepley   PetscCall(KSViewerDestroy(&viewer));
801e3ce8211SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
802d6685f55SMatthew G. Knepley }
803d6685f55SMatthew G. Knepley 
804e3ce8211SMatthew G. Knepley /*@C
805e3ce8211SMatthew G. Knepley   PetscProbComputeKSStatisticWeighted - Compute the Kolmogorov-Smirnov statistic for the weighted empirical distribution for an input vector, compared to an analytic CDF.
806e3ce8211SMatthew G. Knepley 
807e3ce8211SMatthew G. Knepley   Collective
808e3ce8211SMatthew G. Knepley 
809e3ce8211SMatthew G. Knepley   Input Parameters:
810e3ce8211SMatthew G. Knepley + v   - The data vector, blocksize is the sample dimension
811e3ce8211SMatthew G. Knepley . w   - The vector of weights for each sample, instead of the default 1/n
812e3ce8211SMatthew G. Knepley - cdf - The analytic CDF
813e3ce8211SMatthew G. Knepley 
814e3ce8211SMatthew G. Knepley   Output Parameter:
815e3ce8211SMatthew G. Knepley . alpha - The KS statistic
816e3ce8211SMatthew G. Knepley 
817e3ce8211SMatthew G. Knepley   Level: advanced
818e3ce8211SMatthew G. Knepley 
819e3ce8211SMatthew G. Knepley   Notes:
820e3ce8211SMatthew G. Knepley   The Kolmogorov-Smirnov statistic for a given cumulative distribution function $F(x)$ is
821e3ce8211SMatthew G. Knepley 
822e3ce8211SMatthew G. Knepley   $$
823e3ce8211SMatthew G. Knepley   D_n = \sup_x \left| F_n(x) - F(x) \right|
824e3ce8211SMatthew G. Knepley   $$
825e3ce8211SMatthew G. Knepley 
826e3ce8211SMatthew G. Knepley   where $\sup_x$ is the supremum of the set of distances, and the empirical distribution function $F_n(x)$ is discrete, and given by
827e3ce8211SMatthew G. Knepley 
828e3ce8211SMatthew G. Knepley   $$
829e3ce8211SMatthew G. Knepley   F_n =  # of samples <= x / n
830e3ce8211SMatthew G. Knepley   $$
831e3ce8211SMatthew G. Knepley 
832e3ce8211SMatthew G. Knepley   The empirical distribution function $F_n(x)$ is discrete, and thus had a ``stairstep''
833e3ce8211SMatthew G. Knepley   cumulative distribution, making $n$ the number of stairs. Intuitively, the statistic takes
834e3ce8211SMatthew G. Knepley   the largest absolute difference between the two distribution functions across all $x$ values.
835e3ce8211SMatthew G. Knepley 
836e3ce8211SMatthew G. Knepley   The goodness-of-fit test, or Kolmogorov-Smirnov test, is constructed using the Kolmogorov
837e3ce8211SMatthew G. Knepley   distribution. It rejects the null hypothesis at level $\alpha$ if
838e3ce8211SMatthew G. Knepley 
839e3ce8211SMatthew G. Knepley   $$
840e3ce8211SMatthew G. Knepley   \sqrt{n} D_{n} > K_{\alpha},
841e3ce8211SMatthew G. Knepley   $$
842e3ce8211SMatthew G. Knepley 
843e3ce8211SMatthew G. Knepley   where $K_\alpha$ is found from
844e3ce8211SMatthew G. Knepley 
845e3ce8211SMatthew G. Knepley   $$
846e3ce8211SMatthew G. Knepley   \operatorname{Pr}(K \leq K_{\alpha}) = 1 - \alpha.
847e3ce8211SMatthew G. Knepley   $$
848e3ce8211SMatthew G. Knepley 
849e3ce8211SMatthew G. Knepley   This means that getting a small alpha says that we have high confidence that the data did not come
850e3ce8211SMatthew G. Knepley   from the input distribution, so we say that it rejects the null hypothesis.
851e3ce8211SMatthew G. Knepley 
852f8662bd6SBarry Smith .seealso: `PetscProbComputeKSStatistic()`, `PetscProbComputeKSStatisticMagnitude()`, `PetscProbFn`
853e3ce8211SMatthew G. Knepley @*/
PetscProbComputeKSStatisticWeighted(Vec v,Vec w,PetscProbFn * cdf,PetscReal * alpha)854f8662bd6SBarry Smith PetscErrorCode PetscProbComputeKSStatisticWeighted(Vec v, Vec w, PetscProbFn *cdf, PetscReal *alpha)
855e3ce8211SMatthew G. Knepley {
856e3ce8211SMatthew G. Knepley   PetscViewer        viewer     = NULL;
857e3ce8211SMatthew G. Knepley   OutputType         outputType = NONE;
858e3ce8211SMatthew G. Knepley   const PetscScalar *val, *wgt;
859e3ce8211SMatthew G. Knepley   PetscInt           n;
860e3ce8211SMatthew G. Knepley 
861e3ce8211SMatthew G. Knepley   PetscFunctionBegin;
862e3ce8211SMatthew G. Knepley   PetscCall(KSViewerCreate((PetscObject)v, &outputType, &viewer));
863e3ce8211SMatthew G. Knepley   PetscCall(VecGetLocalSize(v, &n));
864e3ce8211SMatthew G. Knepley   PetscCall(VecGetArrayRead(v, &val));
865e3ce8211SMatthew G. Knepley   PetscCall(VecGetArrayRead(w, &wgt));
866e3ce8211SMatthew G. Knepley   PetscCheck(!PetscDefined(USE_COMPLEX), PetscObjectComm((PetscObject)v), PETSC_ERR_SUP, "K-S test does not support complex");
867e3ce8211SMatthew G. Knepley   PetscCall(PetscProbComputeKSStatistic_Internal(PetscObjectComm((PetscObject)v), n, (PetscReal *)val, (PetscReal *)wgt, cdf, alpha, outputType, viewer));
868e3ce8211SMatthew G. Knepley   PetscCall(VecRestoreArrayRead(v, &val));
869e3ce8211SMatthew G. Knepley   PetscCall(VecRestoreArrayRead(w, &wgt));
870e3ce8211SMatthew G. Knepley   PetscCall(KSViewerDestroy(&viewer));
871e3ce8211SMatthew G. Knepley   PetscFunctionReturn(PETSC_SUCCESS);
872e3ce8211SMatthew G. Knepley }
873e3ce8211SMatthew G. Knepley 
874e3ce8211SMatthew G. Knepley /*@C
875e3ce8211SMatthew G. Knepley   PetscProbComputeKSStatisticMagnitude - Compute the Kolmogorov-Smirnov statistic for the empirical distribution for the magnitude over each block of an input vector, compared to an analytic CDF.
876e3ce8211SMatthew G. Knepley 
877e3ce8211SMatthew G. Knepley   Collective
878e3ce8211SMatthew G. Knepley 
879e3ce8211SMatthew G. Knepley   Input Parameters:
880e3ce8211SMatthew G. Knepley + v   - The data vector, blocksize is the sample dimension
881e3ce8211SMatthew G. Knepley - cdf - The analytic CDF
882e3ce8211SMatthew G. Knepley 
883e3ce8211SMatthew G. Knepley   Output Parameter:
884e3ce8211SMatthew G. Knepley . alpha - The KS statistic
885e3ce8211SMatthew G. Knepley 
886e3ce8211SMatthew G. Knepley   Level: advanced
887e3ce8211SMatthew G. Knepley 
888e3ce8211SMatthew G. Knepley   Notes:
889e3ce8211SMatthew G. Knepley   The Kolmogorov-Smirnov statistic for a given cumulative distribution function $F(x)$ is
890e3ce8211SMatthew G. Knepley 
891e3ce8211SMatthew G. Knepley   $$
892e3ce8211SMatthew G. Knepley   D_n = \sup_x \left| F_n(x) - F(x) \right|
893e3ce8211SMatthew G. Knepley   $$
894e3ce8211SMatthew G. Knepley 
895e3ce8211SMatthew G. Knepley   where $\sup_x$ is the supremum of the set of distances, and the empirical distribution function $F_n(x)$ is discrete, and given by
896e3ce8211SMatthew G. Knepley 
897e3ce8211SMatthew G. Knepley   $$
898e3ce8211SMatthew G. Knepley   F_n =  # of samples <= x / n
899e3ce8211SMatthew G. Knepley   $$
900e3ce8211SMatthew G. Knepley 
901e3ce8211SMatthew G. Knepley   The empirical distribution function $F_n(x)$ is discrete, and thus had a ``stairstep''
902e3ce8211SMatthew G. Knepley   cumulative distribution, making $n$ the number of stairs. Intuitively, the statistic takes
903e3ce8211SMatthew G. Knepley   the largest absolute difference between the two distribution functions across all $x$ values.
904e3ce8211SMatthew G. Knepley 
905e3ce8211SMatthew G. Knepley   The goodness-of-fit test, or Kolmogorov-Smirnov test, is constructed using the Kolmogorov
906e3ce8211SMatthew G. Knepley   distribution. It rejects the null hypothesis at level $\alpha$ if
907e3ce8211SMatthew G. Knepley 
908e3ce8211SMatthew G. Knepley   $$
909e3ce8211SMatthew G. Knepley   \sqrt{n} D_{n} > K_{\alpha},
910e3ce8211SMatthew G. Knepley   $$
911e3ce8211SMatthew G. Knepley 
912e3ce8211SMatthew G. Knepley   where $K_\alpha$ is found from
913e3ce8211SMatthew G. Knepley 
914e3ce8211SMatthew G. Knepley   $$
915e3ce8211SMatthew G. Knepley   \operatorname{Pr}(K \leq K_{\alpha}) = 1 - \alpha.
916e3ce8211SMatthew G. Knepley   $$
917e3ce8211SMatthew G. Knepley 
918e3ce8211SMatthew G. Knepley   This means that getting a small alpha says that we have high confidence that the data did not come
919e3ce8211SMatthew G. Knepley   from the input distribution, so we say that it rejects the null hypothesis.
920e3ce8211SMatthew G. Knepley 
921f8662bd6SBarry Smith .seealso: `PetscProbComputeKSStatistic()`, `PetscProbComputeKSStatisticWeighted()`, `PetscProbFn`
922e3ce8211SMatthew G. Knepley @*/
PetscProbComputeKSStatisticMagnitude(Vec v,PetscProbFn * cdf,PetscReal * alpha)923f8662bd6SBarry Smith PetscErrorCode PetscProbComputeKSStatisticMagnitude(Vec v, PetscProbFn *cdf, PetscReal *alpha)
924e3ce8211SMatthew G. Knepley {
925e3ce8211SMatthew G. Knepley   PetscViewer        viewer     = NULL;
926e3ce8211SMatthew G. Knepley   OutputType         outputType = NONE;
927e3ce8211SMatthew G. Knepley   const PetscScalar *a;
928e3ce8211SMatthew G. Knepley   PetscReal         *speed;
929e3ce8211SMatthew G. Knepley   PetscInt           n, dim;
930e3ce8211SMatthew G. Knepley 
931e3ce8211SMatthew G. Knepley   PetscFunctionBegin;
932e3ce8211SMatthew G. Knepley   PetscCall(KSViewerCreate((PetscObject)v, &outputType, &viewer));
933e3ce8211SMatthew G. Knepley   // Convert to a scalar value
9349566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(v, &n));
9359566063dSJacob Faibussowitsch   PetscCall(VecGetBlockSize(v, &dim));
936d6685f55SMatthew G. Knepley   n /= dim;
9379566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &speed));
9389566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(v, &a));
939e3ce8211SMatthew G. Knepley   for (PetscInt p = 0; p < n; ++p) {
940d6685f55SMatthew G. Knepley     PetscReal mag = 0.;
941d6685f55SMatthew G. Knepley 
942e3ce8211SMatthew G. Knepley     for (PetscInt d = 0; d < dim; ++d) mag += PetscSqr(PetscRealPart(a[p * dim + d]));
943d6685f55SMatthew G. Knepley     speed[p] = PetscSqrtReal(mag);
944d6685f55SMatthew G. Knepley   }
9459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(v, &a));
946e3ce8211SMatthew G. Knepley   // Compute statistic
947e3ce8211SMatthew G. Knepley   PetscCall(PetscProbComputeKSStatistic_Internal(PetscObjectComm((PetscObject)v), n, speed, NULL, cdf, alpha, outputType, viewer));
9489566063dSJacob Faibussowitsch   PetscCall(PetscFree(speed));
949e3ce8211SMatthew G. Knepley   PetscCall(KSViewerDestroy(&viewer));
9503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
951d6685f55SMatthew G. Knepley }
952