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 22d6685f55SMatthew G. Knepley .seealso: PetscCDFMaxwellBoltzmann1D 23d6685f55SMatthew G. Knepley @*/ 24d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 25d6685f55SMatthew G. Knepley { 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 43d6685f55SMatthew G. Knepley .seealso: PetscPDFMaxwellBoltzmann1D 44d6685f55SMatthew G. Knepley @*/ 45d6685f55SMatthew G. Knepley PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 46d6685f55SMatthew G. Knepley { 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 64d6685f55SMatthew G. Knepley .seealso: PetscCDFMaxwellBoltzmann2D 65d6685f55SMatthew G. Knepley @*/ 66d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 67d6685f55SMatthew G. Knepley { 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 85d6685f55SMatthew G. Knepley .seealso: PetscPDFMaxwellBoltzmann2D 86d6685f55SMatthew G. Knepley @*/ 87d6685f55SMatthew G. Knepley PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 88d6685f55SMatthew G. Knepley { 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 106d6685f55SMatthew G. Knepley .seealso: PetscCDFMaxwellBoltzmann3D 107d6685f55SMatthew G. Knepley @*/ 108d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 109d6685f55SMatthew G. Knepley { 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 127d6685f55SMatthew G. Knepley .seealso: PetscPDFMaxwellBoltzmann3D 128d6685f55SMatthew G. Knepley @*/ 129d6685f55SMatthew G. Knepley PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 130d6685f55SMatthew G. Knepley { 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 148d6685f55SMatthew G. Knepley .seealso: PetscPDFMaxwellBoltzmann3D 149d6685f55SMatthew G. Knepley @*/ 150d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFGaussian1D(const PetscReal x[], const PetscReal scale[], PetscReal p[]) 151d6685f55SMatthew G. Knepley { 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 157d6685f55SMatthew G. Knepley PetscErrorCode PetscCDFGaussian1D(const PetscReal x[], const PetscReal scale[], PetscReal p[]) 158d6685f55SMatthew G. Knepley { 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 @*/ 181d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal p[], const PetscReal dummy[], PetscReal x[]) 182d6685f55SMatthew G. Knepley { 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 222d6685f55SMatthew G. Knepley .seealso: PetscPDFSampleGaussian2D(), PetscPDFMaxwellBoltzmann3D() 223d6685f55SMatthew G. Knepley @*/ 224d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFGaussian2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 225d6685f55SMatthew G. Knepley { 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 246d6685f55SMatthew G. Knepley .seealso: PetscPDFGaussian2D(), PetscPDFMaxwellBoltzmann3D() 247d6685f55SMatthew G. Knepley @*/ 248d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal p[], const PetscReal dummy[], PetscReal x[]) 249d6685f55SMatthew G. Knepley { 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 /*@ 257*ea1b28ebSMatthew G. Knepley PetscPDFGaussian3D - PDF for the Gaussian distribution in 3D 258*ea1b28ebSMatthew G. Knepley 259*ea1b28ebSMatthew G. Knepley Not collective 260*ea1b28ebSMatthew G. Knepley 261*ea1b28ebSMatthew G. Knepley Input Parameters: 262*ea1b28ebSMatthew G. Knepley + x - Coordinate in [-\infty, \infty]^3 263*ea1b28ebSMatthew G. Knepley - dummy - ignored 264*ea1b28ebSMatthew G. Knepley 265*ea1b28ebSMatthew G. Knepley Output Parameter: 266*ea1b28ebSMatthew G. Knepley . p - The probability density at x 267*ea1b28ebSMatthew G. Knepley 268*ea1b28ebSMatthew G. Knepley Level: beginner 269*ea1b28ebSMatthew G. Knepley 270*ea1b28ebSMatthew G. Knepley .seealso: PetscPDFSampleGaussian3D(), PetscPDFMaxwellBoltzmann3D() 271*ea1b28ebSMatthew G. Knepley @*/ 272*ea1b28ebSMatthew G. Knepley PetscErrorCode PetscPDFGaussian3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 273*ea1b28ebSMatthew G. Knepley { 274*ea1b28ebSMatthew G. Knepley p[0] = (1. / PETSC_PI*PetscSqrtReal(PETSC_PI)) * PetscExpReal(-0.5 * (PetscSqr(x[0]) + PetscSqr(x[1]) + PetscSqr(x[2]))); 275*ea1b28ebSMatthew G. Knepley return 0; 276*ea1b28ebSMatthew G. Knepley } 277*ea1b28ebSMatthew G. Knepley 278*ea1b28ebSMatthew G. Knepley /*@ 279*ea1b28ebSMatthew G. Knepley PetscPDFSampleGaussian3D - Sample uniformly from a Gaussian distribution in 3D 280*ea1b28ebSMatthew G. Knepley 281*ea1b28ebSMatthew G. Knepley Not collective 282*ea1b28ebSMatthew G. Knepley 283*ea1b28ebSMatthew G. Knepley Input Parameters: 284*ea1b28ebSMatthew G. Knepley + p - A uniform variable on [0, 1]^3 285*ea1b28ebSMatthew G. Knepley - dummy - ignored 286*ea1b28ebSMatthew G. Knepley 287*ea1b28ebSMatthew G. Knepley Output Parameter: 288*ea1b28ebSMatthew G. Knepley . x - Coordinate in [-\infty, \infty]^3 289*ea1b28ebSMatthew G. Knepley 290*ea1b28ebSMatthew G. Knepley Level: beginner 291*ea1b28ebSMatthew G. Knepley 292*ea1b28ebSMatthew G. Knepley Note: https://en.wikipedia.org/wiki/Box%E2%80%93Muller_transform 293*ea1b28ebSMatthew G. Knepley 294*ea1b28ebSMatthew G. Knepley .seealso: PetscPDFGaussian3D(), PetscPDFMaxwellBoltzmann3D() 295*ea1b28ebSMatthew G. Knepley @*/ 296*ea1b28ebSMatthew G. Knepley PetscErrorCode PetscPDFSampleGaussian3D(const PetscReal p[], const PetscReal dummy[], PetscReal x[]) 297*ea1b28ebSMatthew G. Knepley { 298*ea1b28ebSMatthew G. Knepley PetscCall(PetscPDFSampleGaussian1D(p, dummy, x)); 299*ea1b28ebSMatthew G. Knepley PetscCall(PetscPDFSampleGaussian2D(&p[1], dummy, &x[1])); 300*ea1b28ebSMatthew G. Knepley return 0; 301*ea1b28ebSMatthew G. Knepley } 302*ea1b28ebSMatthew G. Knepley 303*ea1b28ebSMatthew 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 316*ea1b28ebSMatthew G. Knepley .seealso: PetscCDFConstant1D(), PetscPDFSampleConstant1D(), PetscPDFConstant2D(), PetscPDFConstant3D() 317d6685f55SMatthew G. Knepley @*/ 318d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFConstant1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 319d6685f55SMatthew G. Knepley { 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 324*ea1b28ebSMatthew G. Knepley /*@ 325*ea1b28ebSMatthew G. Knepley PetscCDFConstant1D - CDF for the uniform distribution in 1D 326*ea1b28ebSMatthew G. Knepley 327*ea1b28ebSMatthew G. Knepley Not collective 328*ea1b28ebSMatthew G. Knepley 329*ea1b28ebSMatthew G. Knepley Input Parameter: 330*ea1b28ebSMatthew G. Knepley . x - Coordinate in [-1, 1] 331*ea1b28ebSMatthew G. Knepley 332*ea1b28ebSMatthew G. Knepley Output Parameter: 333*ea1b28ebSMatthew G. Knepley . p - The cumulative probability at x 334*ea1b28ebSMatthew G. Knepley 335*ea1b28ebSMatthew G. Knepley Level: beginner 336*ea1b28ebSMatthew G. Knepley 337*ea1b28ebSMatthew G. Knepley .seealso: PetscPDFConstant1D(), PetscPDFSampleConstant1D(), PetscCDFConstant2D(), PetscCDFConstant3D() 338*ea1b28ebSMatthew G. Knepley @*/ 339d6685f55SMatthew G. Knepley PetscErrorCode PetscCDFConstant1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 340d6685f55SMatthew G. Knepley { 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 358*ea1b28ebSMatthew G. Knepley .seealso: PetscPDFConstant1D(), PetscCDFConstant1D(), PetscPDFSampleConstant2D(), PetscPDFSampleConstant3D() 359d6685f55SMatthew G. Knepley @*/ 360d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFSampleConstant1D(const PetscReal p[], const PetscReal dummy[], PetscReal x[]) 361d6685f55SMatthew G. Knepley { 362*ea1b28ebSMatthew G. Knepley x[0] = 2.*p[0] - 1.; 363*ea1b28ebSMatthew G. Knepley return 0; 364*ea1b28ebSMatthew G. Knepley } 365*ea1b28ebSMatthew G. Knepley 366*ea1b28ebSMatthew G. Knepley /*@ 367*ea1b28ebSMatthew G. Knepley PetscPDFConstant2D - PDF for the uniform distribution in 2D 368*ea1b28ebSMatthew G. Knepley 369*ea1b28ebSMatthew G. Knepley Not collective 370*ea1b28ebSMatthew G. Knepley 371*ea1b28ebSMatthew G. Knepley Input Parameter: 372*ea1b28ebSMatthew G. Knepley . x - Coordinate in [-1, 1] x [-1, 1] 373*ea1b28ebSMatthew G. Knepley 374*ea1b28ebSMatthew G. Knepley Output Parameter: 375*ea1b28ebSMatthew G. Knepley . p - The probability density at x 376*ea1b28ebSMatthew G. Knepley 377*ea1b28ebSMatthew G. Knepley Level: beginner 378*ea1b28ebSMatthew G. Knepley 379*ea1b28ebSMatthew G. Knepley .seealso: PetscCDFConstant2D(), PetscPDFSampleConstant2D(), PetscPDFConstant1D(), PetscPDFConstant3D() 380*ea1b28ebSMatthew G. Knepley @*/ 381*ea1b28ebSMatthew G. Knepley PetscErrorCode PetscPDFConstant2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 382*ea1b28ebSMatthew G. Knepley { 383*ea1b28ebSMatthew G. Knepley p[0] = x[0] > -1. && x[0] <= 1. && x[1] > -1. && x[1] <= 1. ? 0.25 : 0.; 384*ea1b28ebSMatthew G. Knepley return 0; 385*ea1b28ebSMatthew G. Knepley } 386*ea1b28ebSMatthew G. Knepley 387*ea1b28ebSMatthew G. Knepley /*@ 388*ea1b28ebSMatthew G. Knepley PetscCDFConstant2D - CDF for the uniform distribution in 2D 389*ea1b28ebSMatthew G. Knepley 390*ea1b28ebSMatthew G. Knepley Not collective 391*ea1b28ebSMatthew G. Knepley 392*ea1b28ebSMatthew G. Knepley Input Parameter: 393*ea1b28ebSMatthew G. Knepley . x - Coordinate in [-1, 1] x [-1, 1] 394*ea1b28ebSMatthew G. Knepley 395*ea1b28ebSMatthew G. Knepley Output Parameter: 396*ea1b28ebSMatthew G. Knepley . p - The cumulative probability at x 397*ea1b28ebSMatthew G. Knepley 398*ea1b28ebSMatthew G. Knepley Level: beginner 399*ea1b28ebSMatthew G. Knepley 400*ea1b28ebSMatthew G. Knepley .seealso: PetscPDFConstant2D(), PetscPDFSampleConstant2D(), PetscCDFConstant1D(), PetscCDFConstant3D() 401*ea1b28ebSMatthew G. Knepley @*/ 402*ea1b28ebSMatthew G. Knepley PetscErrorCode PetscCDFConstant2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 403*ea1b28ebSMatthew G. Knepley { 404*ea1b28ebSMatthew 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.)); 405*ea1b28ebSMatthew G. Knepley return 0; 406*ea1b28ebSMatthew G. Knepley } 407*ea1b28ebSMatthew G. Knepley 408*ea1b28ebSMatthew G. Knepley /*@ 409*ea1b28ebSMatthew G. Knepley PetscPDFSampleConstant2D - Sample uniformly from a uniform distribution on [-1, 1] x [-1, 1] in 2D 410*ea1b28ebSMatthew G. Knepley 411*ea1b28ebSMatthew G. Knepley Not collective 412*ea1b28ebSMatthew G. Knepley 413*ea1b28ebSMatthew G. Knepley Input Parameter: 414*ea1b28ebSMatthew G. Knepley . p - Two uniform variables on [0, 1] 415*ea1b28ebSMatthew G. Knepley 416*ea1b28ebSMatthew G. Knepley Output Parameter: 417*ea1b28ebSMatthew G. Knepley . x - Coordinate in [-1, 1] x [-1, 1] 418*ea1b28ebSMatthew G. Knepley 419*ea1b28ebSMatthew G. Knepley Level: beginner 420*ea1b28ebSMatthew G. Knepley 421*ea1b28ebSMatthew G. Knepley .seealso: PetscPDFConstant2D(), PetscCDFConstant2D(), PetscPDFSampleConstant1D(), PetscPDFSampleConstant3D() 422*ea1b28ebSMatthew G. Knepley @*/ 423*ea1b28ebSMatthew G. Knepley PetscErrorCode PetscPDFSampleConstant2D(const PetscReal p[], const PetscReal dummy[], PetscReal x[]) 424*ea1b28ebSMatthew G. Knepley { 425*ea1b28ebSMatthew G. Knepley x[0] = 2.*p[0] - 1.; 426*ea1b28ebSMatthew G. Knepley x[1] = 2.*p[1] - 1.; 427*ea1b28ebSMatthew G. Knepley return 0; 428*ea1b28ebSMatthew G. Knepley } 429*ea1b28ebSMatthew G. Knepley 430*ea1b28ebSMatthew G. Knepley /*@ 431*ea1b28ebSMatthew G. Knepley PetscPDFConstant3D - PDF for the uniform distribution in 3D 432*ea1b28ebSMatthew G. Knepley 433*ea1b28ebSMatthew G. Knepley Not collective 434*ea1b28ebSMatthew G. Knepley 435*ea1b28ebSMatthew G. Knepley Input Parameter: 436*ea1b28ebSMatthew G. Knepley . x - Coordinate in [-1, 1] x [-1, 1] x [-1, 1] 437*ea1b28ebSMatthew G. Knepley 438*ea1b28ebSMatthew G. Knepley Output Parameter: 439*ea1b28ebSMatthew G. Knepley . p - The probability density at x 440*ea1b28ebSMatthew G. Knepley 441*ea1b28ebSMatthew G. Knepley Level: beginner 442*ea1b28ebSMatthew G. Knepley 443*ea1b28ebSMatthew G. Knepley .seealso: PetscCDFConstant3D(), PetscPDFSampleConstant3D(), PetscPDFSampleConstant1D(), PetscPDFSampleConstant2D() 444*ea1b28ebSMatthew G. Knepley @*/ 445*ea1b28ebSMatthew G. Knepley PetscErrorCode PetscPDFConstant3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 446*ea1b28ebSMatthew G. Knepley { 447*ea1b28ebSMatthew 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.; 448*ea1b28ebSMatthew G. Knepley return 0; 449*ea1b28ebSMatthew G. Knepley } 450*ea1b28ebSMatthew G. Knepley 451*ea1b28ebSMatthew G. Knepley /*@ 452*ea1b28ebSMatthew G. Knepley PetscCDFConstant3D - CDF for the uniform distribution in 3D 453*ea1b28ebSMatthew G. Knepley 454*ea1b28ebSMatthew G. Knepley Not collective 455*ea1b28ebSMatthew G. Knepley 456*ea1b28ebSMatthew G. Knepley Input Parameter: 457*ea1b28ebSMatthew G. Knepley . x - Coordinate in [-1, 1] x [-1, 1] x [-1, 1] 458*ea1b28ebSMatthew G. Knepley 459*ea1b28ebSMatthew G. Knepley Output Parameter: 460*ea1b28ebSMatthew G. Knepley . p - The cumulative probability at x 461*ea1b28ebSMatthew G. Knepley 462*ea1b28ebSMatthew G. Knepley Level: beginner 463*ea1b28ebSMatthew G. Knepley 464*ea1b28ebSMatthew G. Knepley .seealso: PetscPDFConstant3D(), PetscPDFSampleConstant3D(), PetscCDFConstant1D(), PetscCDFConstant2D() 465*ea1b28ebSMatthew G. Knepley @*/ 466*ea1b28ebSMatthew G. Knepley PetscErrorCode PetscCDFConstant3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 467*ea1b28ebSMatthew G. Knepley { 468*ea1b28ebSMatthew 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.)); 469*ea1b28ebSMatthew G. Knepley return 0; 470*ea1b28ebSMatthew G. Knepley } 471*ea1b28ebSMatthew G. Knepley 472*ea1b28ebSMatthew G. Knepley /*@ 473*ea1b28ebSMatthew G. Knepley PetscPDFSampleConstant3D - Sample uniformly from a uniform distribution on [-1, 1] x [-1, 1] in 3D 474*ea1b28ebSMatthew G. Knepley 475*ea1b28ebSMatthew G. Knepley Not collective 476*ea1b28ebSMatthew G. Knepley 477*ea1b28ebSMatthew G. Knepley Input Parameter: 478*ea1b28ebSMatthew G. Knepley . p - Three uniform variables on [0, 1] 479*ea1b28ebSMatthew G. Knepley 480*ea1b28ebSMatthew G. Knepley Output Parameter: 481*ea1b28ebSMatthew G. Knepley . x - Coordinate in [-1, 1] x [-1, 1] x [-1, 1] 482*ea1b28ebSMatthew G. Knepley 483*ea1b28ebSMatthew G. Knepley Level: beginner 484*ea1b28ebSMatthew G. Knepley 485*ea1b28ebSMatthew G. Knepley .seealso: PetscPDFConstant3D(), PetscCDFConstant3D(), PetscPDFSampleConstant1D(), PetscPDFSampleConstant2D() 486*ea1b28ebSMatthew G. Knepley @*/ 487*ea1b28ebSMatthew G. Knepley PetscErrorCode PetscPDFSampleConstant3D(const PetscReal p[], const PetscReal dummy[], PetscReal x[]) 488*ea1b28ebSMatthew G. Knepley { 489*ea1b28ebSMatthew G. Knepley x[0] = 2.*p[0] - 1.; 490*ea1b28ebSMatthew G. Knepley x[1] = 2.*p[1] - 1.; 491*ea1b28ebSMatthew 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 512d6685f55SMatthew G. Knepley .seealso: PetscPDFMaxwellBoltzmann1D(), PetscPDFGaussian1D(), PetscPDFConstant1D() 513d6685f55SMatthew G. Knepley @*/ 514d6685f55SMatthew G. Knepley PetscErrorCode PetscProbCreateFromOptions(PetscInt dim, const char prefix[], const char name[], PetscProbFunc *pdf, PetscProbFunc *cdf, PetscProbFunc *sampler) 515d6685f55SMatthew G. Knepley { 516d6685f55SMatthew G. Knepley DTProbDensityType den = DTPROB_DENSITY_GAUSSIAN; 517d6685f55SMatthew G. Knepley PetscErrorCode ierr; 518d6685f55SMatthew G. Knepley 519d6685f55SMatthew G. Knepley PetscFunctionBegin; 5209566063dSJacob Faibussowitsch ierr = PetscOptionsBegin(PETSC_COMM_SELF, prefix, "PetscProb Options", "DT");PetscCall(ierr); 5219566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum(name, "Method to compute PDF <constant, gaussian>", "", DTProbDensityTypes, (PetscEnum) den, (PetscEnum *) &den, NULL)); 5229566063dSJacob Faibussowitsch ierr = PetscOptionsEnd();PetscCall(ierr); 523d6685f55SMatthew G. Knepley 524d6685f55SMatthew G. Knepley if (pdf) {PetscValidPointer(pdf, 4); *pdf = NULL;} 525d6685f55SMatthew G. Knepley if (cdf) {PetscValidPointer(cdf, 5); *cdf = NULL;} 526d6685f55SMatthew G. Knepley if (sampler) {PetscValidPointer(sampler, 6); *sampler = NULL;} 527d6685f55SMatthew G. Knepley switch (den) { 528d6685f55SMatthew G. Knepley case DTPROB_DENSITY_CONSTANT: 529d6685f55SMatthew G. Knepley switch (dim) { 530d6685f55SMatthew G. Knepley case 1: 531d6685f55SMatthew G. Knepley if (pdf) *pdf = PetscPDFConstant1D; 532d6685f55SMatthew G. Knepley if (cdf) *cdf = PetscCDFConstant1D; 533d6685f55SMatthew G. Knepley if (sampler) *sampler = PetscPDFSampleConstant1D; 534d6685f55SMatthew G. Knepley break; 535*ea1b28ebSMatthew G. Knepley case 2: 536*ea1b28ebSMatthew G. Knepley if (pdf) *pdf = PetscPDFConstant2D; 537*ea1b28ebSMatthew G. Knepley if (cdf) *cdf = PetscCDFConstant2D; 538*ea1b28ebSMatthew G. Knepley if (sampler) *sampler = PetscPDFSampleConstant2D; 539*ea1b28ebSMatthew G. Knepley break; 540*ea1b28ebSMatthew G. Knepley case 3: 541*ea1b28ebSMatthew G. Knepley if (pdf) *pdf = PetscPDFConstant3D; 542*ea1b28ebSMatthew G. Knepley if (cdf) *cdf = PetscCDFConstant3D; 543*ea1b28ebSMatthew G. Knepley if (sampler) *sampler = PetscPDFSampleConstant3D; 544*ea1b28ebSMatthew G. Knepley break; 545d6685f55SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]); 546d6685f55SMatthew G. Knepley } 547d6685f55SMatthew G. Knepley break; 548d6685f55SMatthew G. Knepley case DTPROB_DENSITY_GAUSSIAN: 549d6685f55SMatthew G. Knepley switch (dim) { 550d6685f55SMatthew G. Knepley case 1: 551d6685f55SMatthew G. Knepley if (pdf) *pdf = PetscPDFGaussian1D; 552d6685f55SMatthew G. Knepley if (cdf) *cdf = PetscCDFGaussian1D; 553d6685f55SMatthew G. Knepley if (sampler) *sampler = PetscPDFSampleGaussian1D; 554d6685f55SMatthew G. Knepley break; 555d6685f55SMatthew G. Knepley case 2: 556d6685f55SMatthew G. Knepley if (pdf) *pdf = PetscPDFGaussian2D; 557d6685f55SMatthew G. Knepley if (sampler) *sampler = PetscPDFSampleGaussian2D; 558d6685f55SMatthew G. Knepley break; 559*ea1b28ebSMatthew G. Knepley case 3: 560*ea1b28ebSMatthew G. Knepley if (pdf) *pdf = PetscPDFGaussian3D; 561*ea1b28ebSMatthew G. Knepley if (sampler) *sampler = PetscPDFSampleGaussian3D; 562*ea1b28ebSMatthew G. Knepley break; 563d6685f55SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]); 564d6685f55SMatthew G. Knepley } 565d6685f55SMatthew G. Knepley break; 566d6685f55SMatthew G. Knepley case DTPROB_DENSITY_MAXWELL_BOLTZMANN: 567d6685f55SMatthew G. Knepley switch (dim) { 568d6685f55SMatthew G. Knepley case 1: 569d6685f55SMatthew G. Knepley if (pdf) *pdf = PetscPDFMaxwellBoltzmann1D; 570d6685f55SMatthew G. Knepley if (cdf) *cdf = PetscCDFMaxwellBoltzmann1D; 571d6685f55SMatthew G. Knepley break; 572d6685f55SMatthew G. Knepley case 2: 573d6685f55SMatthew G. Knepley if (pdf) *pdf = PetscPDFMaxwellBoltzmann2D; 574d6685f55SMatthew G. Knepley if (cdf) *cdf = PetscCDFMaxwellBoltzmann2D; 575d6685f55SMatthew G. Knepley break; 576d6685f55SMatthew G. Knepley case 3: 577d6685f55SMatthew G. Knepley if (pdf) *pdf = PetscPDFMaxwellBoltzmann3D; 578d6685f55SMatthew G. Knepley if (cdf) *cdf = PetscCDFMaxwellBoltzmann3D; 579d6685f55SMatthew G. Knepley break; 580d6685f55SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]); 581d6685f55SMatthew G. Knepley } 582d6685f55SMatthew G. Knepley break; 583d6685f55SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Density type %s is not supported", DTProbDensityTypes[PetscMax(0, PetscMin(den, DTPROB_NUM_DENSITY))]); 584d6685f55SMatthew G. Knepley } 585d6685f55SMatthew G. Knepley PetscFunctionReturn(0); 586d6685f55SMatthew G. Knepley } 587d6685f55SMatthew G. Knepley 588d6685f55SMatthew G. Knepley #ifdef PETSC_HAVE_KS 589d6685f55SMatthew G. Knepley EXTERN_C_BEGIN 590d6685f55SMatthew G. Knepley #include <KolmogorovSmirnovDist.h> 591d6685f55SMatthew G. Knepley EXTERN_C_END 592d6685f55SMatthew G. Knepley #endif 593d6685f55SMatthew G. Knepley 594d6685f55SMatthew G. Knepley /*@C 595d6685f55SMatthew G. Knepley PetscProbComputeKSStatistic - Compute the Kolmogorov-Smirnov statistic for the empirical distribution for an input vector, compared to an analytic CDF. 596d6685f55SMatthew G. Knepley 597d6685f55SMatthew G. Knepley Collective on v 598d6685f55SMatthew G. Knepley 599d6685f55SMatthew G. Knepley Input Parameters: 600d6685f55SMatthew G. Knepley + v - The data vector, blocksize is the sample dimension 601d6685f55SMatthew G. Knepley - cdf - The analytic CDF 602d6685f55SMatthew G. Knepley 603d6685f55SMatthew G. Knepley Output Parameter: 604d6685f55SMatthew G. Knepley . alpha - The KS statisic 605d6685f55SMatthew G. Knepley 606d6685f55SMatthew G. Knepley Level: advanced 607d6685f55SMatthew G. Knepley 608d6685f55SMatthew G. Knepley Note: The Kolmogorov-Smirnov statistic for a given cumulative distribution function $F(x)$ is 609d6685f55SMatthew G. Knepley $ 610d6685f55SMatthew G. Knepley $ D_n = \sup_x \left| F_n(x) - F(x) \right| 611d6685f55SMatthew G. Knepley $ 612d6685f55SMatthew G. Knepley where $\sup_x$ is the supremum of the set of distances, and the empirical distribution 613d6685f55SMatthew G. Knepley function $F_n(x)$ is discrete, and given by 614d6685f55SMatthew G. Knepley $ 615d6685f55SMatthew G. Knepley $ F_n = # of samples <= x / n 616d6685f55SMatthew G. Knepley $ 617d6685f55SMatthew G. Knepley The empirical distribution function $F_n(x)$ is discrete, and thus had a ``stairstep'' 618d6685f55SMatthew G. Knepley cumulative distribution, making $n$ the number of stairs. Intuitively, the statistic takes 619d6685f55SMatthew G. Knepley the largest absolute difference between the two distribution functions across all $x$ values. 620d6685f55SMatthew G. Knepley 621d6685f55SMatthew G. Knepley .seealso: PetscProbFunc 622d6685f55SMatthew G. Knepley @*/ 623d6685f55SMatthew G. Knepley PetscErrorCode PetscProbComputeKSStatistic(Vec v, PetscProbFunc cdf, PetscReal *alpha) 624d6685f55SMatthew G. Knepley { 625d6685f55SMatthew G. Knepley #if !defined(PETSC_HAVE_KS) 626d6685f55SMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject) v), PETSC_ERR_SUP, "No support for Kolmogorov-Smirnov test.\nReconfigure using --download-ks"); 627d6685f55SMatthew G. Knepley #else 628d6685f55SMatthew G. Knepley PetscViewer viewer = NULL; 629d6685f55SMatthew G. Knepley PetscViewerFormat format; 630d6685f55SMatthew G. Knepley PetscDraw draw; 631d6685f55SMatthew G. Knepley PetscDrawLG lg; 632d6685f55SMatthew G. Knepley PetscDrawAxis axis; 633d6685f55SMatthew G. Knepley const PetscScalar *a; 634d6685f55SMatthew G. Knepley PetscReal *speed, Dn = PETSC_MIN_REAL; 635d6685f55SMatthew G. Knepley PetscBool isascii = PETSC_FALSE, isdraw = PETSC_FALSE, flg; 636d6685f55SMatthew G. Knepley PetscInt n, p, dim, d; 637d6685f55SMatthew G. Knepley PetscMPIInt size; 638d6685f55SMatthew G. Knepley const char *names[2] = {"Analytic", "Empirical"}; 639d6685f55SMatthew G. Knepley char title[PETSC_MAX_PATH_LEN]; 640d6685f55SMatthew G. Knepley PetscOptions options; 641d6685f55SMatthew G. Knepley const char *prefix; 642d6685f55SMatthew G. Knepley MPI_Comm comm; 643d6685f55SMatthew G. Knepley 644d6685f55SMatthew G. Knepley PetscFunctionBegin; 6459566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject) v, &comm)); 6469566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject) v, &prefix)); 6479566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptions((PetscObject) v, &options)); 6489566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(comm, options, prefix, "-ks_monitor", &viewer, &format, &flg)); 649d6685f55SMatthew G. Knepley if (flg) { 6509566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii)); 6519566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERDRAW, &isdraw)); 652d6685f55SMatthew G. Knepley } 653d6685f55SMatthew G. Knepley if (isascii) { 6549566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, format)); 655d6685f55SMatthew G. Knepley } else if (isdraw) { 6569566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, format)); 6579566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 6589566063dSJacob Faibussowitsch PetscCall(PetscDrawLGCreate(draw, 2, &lg)); 6599566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetLegend(lg, names)); 660d6685f55SMatthew G. Knepley } 661d6685f55SMatthew G. Knepley 6629566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject) v), &size)); 6639566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(v, &n)); 6649566063dSJacob Faibussowitsch PetscCall(VecGetBlockSize(v, &dim)); 665d6685f55SMatthew G. Knepley n /= dim; 666d6685f55SMatthew G. Knepley /* TODO Parallel algorithm is harder */ 667d6685f55SMatthew G. Knepley if (size == 1) { 6689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &speed)); 6699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(v, &a)); 670d6685f55SMatthew G. Knepley for (p = 0; p < n; ++p) { 671d6685f55SMatthew G. Knepley PetscReal mag = 0.; 672d6685f55SMatthew G. Knepley 673d6685f55SMatthew G. Knepley for (d = 0; d < dim; ++d) mag += PetscSqr(PetscRealPart(a[p*dim+d])); 674d6685f55SMatthew G. Knepley speed[p] = PetscSqrtReal(mag); 675d6685f55SMatthew G. Knepley } 6769566063dSJacob Faibussowitsch PetscCall(PetscSortReal(n, speed)); 6779566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(v, &a)); 678d6685f55SMatthew G. Knepley for (p = 0; p < n; ++p) { 679d6685f55SMatthew G. Knepley const PetscReal x = speed[p], Fn = ((PetscReal) p) / n; 680d6685f55SMatthew G. Knepley PetscReal F, vals[2]; 681d6685f55SMatthew G. Knepley 6829566063dSJacob Faibussowitsch PetscCall(cdf(&x, NULL, &F)); 683d6685f55SMatthew G. Knepley Dn = PetscMax(PetscAbsReal(Fn - F), Dn); 6849566063dSJacob Faibussowitsch if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "x: %g F: %g Fn: %g Dn: %.2g\n", x, F, Fn, Dn)); 6859566063dSJacob Faibussowitsch if (isdraw) {vals[0] = F; vals[1] = Fn; PetscCall(PetscDrawLGAddCommonPoint(lg, x, vals));} 686d6685f55SMatthew G. Knepley } 687d6685f55SMatthew G. Knepley if (speed[n-1] < 6.) { 688d6685f55SMatthew G. Knepley const PetscReal k = (PetscInt) (6. - speed[n-1]) + 1, dx = (6. - speed[n-1])/k; 689d6685f55SMatthew G. Knepley for (p = 0; p <= k; ++p) { 690d6685f55SMatthew G. Knepley const PetscReal x = speed[n-1] + p*dx, Fn = 1.0; 691d6685f55SMatthew G. Knepley PetscReal F, vals[2]; 692d6685f55SMatthew G. Knepley 6939566063dSJacob Faibussowitsch PetscCall(cdf(&x, NULL, &F)); 694d6685f55SMatthew G. Knepley Dn = PetscMax(PetscAbsReal(Fn - F), Dn); 6959566063dSJacob Faibussowitsch if (isascii) PetscCall(PetscViewerASCIIPrintf(viewer, "x: %g F: %g Fn: %g Dn: %.2g\n", x, F, Fn, Dn)); 6969566063dSJacob Faibussowitsch if (isdraw) {vals[0] = F; vals[1] = Fn; PetscCall(PetscDrawLGAddCommonPoint(lg, x, vals));} 697d6685f55SMatthew G. Knepley } 698d6685f55SMatthew G. Knepley } 6999566063dSJacob Faibussowitsch PetscCall(PetscFree(speed)); 700d6685f55SMatthew G. Knepley } 701d6685f55SMatthew G. Knepley if (isdraw) { 7029566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(lg, &axis)); 7039566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(title, PETSC_MAX_PATH_LEN, "Kolmogorov-Smirnov Test (Dn %.2g)", Dn)); 7049566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, title, "x", "CDF(x)")); 7059566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(lg)); 7069566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDestroy(&lg)); 707d6685f55SMatthew G. Knepley } 708d6685f55SMatthew G. Knepley if (viewer) { 7099566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 7109566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 7119566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 712d6685f55SMatthew G. Knepley } 713d6685f55SMatthew G. Knepley *alpha = KSfbar((int) n, (double) Dn); 714d6685f55SMatthew G. Knepley PetscFunctionReturn(0); 715d6685f55SMatthew G. Knepley #endif 716d6685f55SMatthew G. Knepley } 717