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