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 169*817da375SSatish Balay Input Parameters: 170*817da375SSatish Balay + p - A uniform variable on [0, 1] 171*817da375SSatish 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 213*817da375SSatish Balay Input Parameters: 214*817da375SSatish Balay + x - Coordinate in [-\infty, \infty]^2 215*817da375SSatish 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 235*817da375SSatish Balay Input Parameters: 236*817da375SSatish Balay + p - A uniform variable on [0, 1]^2 237*817da375SSatish 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 /*@ 257d6685f55SMatthew G. Knepley PetscPDFConstant1D - PDF for the uniform distribution in 1D 258d6685f55SMatthew G. Knepley 259d6685f55SMatthew G. Knepley Not collective 260d6685f55SMatthew G. Knepley 261d6685f55SMatthew G. Knepley Input Parameter: 262d6685f55SMatthew G. Knepley . x - Coordinate in [-1, 1] 263d6685f55SMatthew G. Knepley 264d6685f55SMatthew G. Knepley Output Parameter: 265d6685f55SMatthew G. Knepley . p - The probability density at x 266d6685f55SMatthew G. Knepley 267d6685f55SMatthew G. Knepley Level: beginner 268d6685f55SMatthew G. Knepley 269d6685f55SMatthew G. Knepley .seealso: PetscCDFConstant1D() 270d6685f55SMatthew G. Knepley @*/ 271d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFConstant1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 272d6685f55SMatthew G. Knepley { 273d6685f55SMatthew G. Knepley p[0] = x[0] > -1. && x[0] <= 1. ? 0.5 : 0.; 274d6685f55SMatthew G. Knepley return 0; 275d6685f55SMatthew G. Knepley } 276d6685f55SMatthew G. Knepley 277d6685f55SMatthew G. Knepley PetscErrorCode PetscCDFConstant1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 278d6685f55SMatthew G. Knepley { 279d6685f55SMatthew G. Knepley p[0] = x[0] <= -1. ? 0. : (x[0] > 1. ? 1. : 0.5*(x[0] + 1.)); 280d6685f55SMatthew G. Knepley return 0; 281d6685f55SMatthew G. Knepley } 282d6685f55SMatthew G. Knepley 283d6685f55SMatthew G. Knepley /*@ 284d6685f55SMatthew G. Knepley PetscPDFSampleConstant1D - Sample uniformly from a uniform distribution on [-1, 1] in 1D 285d6685f55SMatthew G. Knepley 286d6685f55SMatthew G. Knepley Not collective 287d6685f55SMatthew G. Knepley 288d6685f55SMatthew G. Knepley Input Parameter: 289d6685f55SMatthew G. Knepley . p - A uniform variable on [0, 1] 290d6685f55SMatthew G. Knepley 291d6685f55SMatthew G. Knepley Output Parameter: 292d6685f55SMatthew G. Knepley . x - Coordinate in [-1, 1] 293d6685f55SMatthew G. Knepley 294d6685f55SMatthew G. Knepley Level: beginner 295d6685f55SMatthew G. Knepley 296d6685f55SMatthew G. Knepley .seealso: PetscPDFConstant1D(), PetscCDFConstant1D() 297d6685f55SMatthew G. Knepley @*/ 298d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFSampleConstant1D(const PetscReal p[], const PetscReal dummy[], PetscReal x[]) 299d6685f55SMatthew G. Knepley { 300d6685f55SMatthew G. Knepley x[0] = 2.*p[1] - 1.; 301d6685f55SMatthew G. Knepley return 0; 302d6685f55SMatthew G. Knepley } 303d6685f55SMatthew G. Knepley 304d6685f55SMatthew G. Knepley /*@C 305d6685f55SMatthew G. Knepley PetscProbCreateFromOptions - Return the probability distribution specified by the argumetns and options 306d6685f55SMatthew G. Knepley 307d6685f55SMatthew G. Knepley Not collective 308d6685f55SMatthew G. Knepley 309d6685f55SMatthew G. Knepley Input Parameters: 310d6685f55SMatthew G. Knepley + dim - The dimension of sample points 311d6685f55SMatthew G. Knepley . prefix - The options prefix, or NULL 312d6685f55SMatthew G. Knepley - name - The option name for the probility distribution type 313d6685f55SMatthew G. Knepley 314d6685f55SMatthew G. Knepley Output Parameters: 315d6685f55SMatthew G. Knepley + pdf - The PDF of this type 316d6685f55SMatthew G. Knepley . cdf - The CDF of this type 317d6685f55SMatthew G. Knepley - sampler - The PDF sampler of this type 318d6685f55SMatthew G. Knepley 319d6685f55SMatthew G. Knepley Level: intermediate 320d6685f55SMatthew G. Knepley 321d6685f55SMatthew G. Knepley .seealso: PetscPDFMaxwellBoltzmann1D(), PetscPDFGaussian1D(), PetscPDFConstant1D() 322d6685f55SMatthew G. Knepley @*/ 323d6685f55SMatthew G. Knepley PetscErrorCode PetscProbCreateFromOptions(PetscInt dim, const char prefix[], const char name[], PetscProbFunc *pdf, PetscProbFunc *cdf, PetscProbFunc *sampler) 324d6685f55SMatthew G. Knepley { 325d6685f55SMatthew G. Knepley DTProbDensityType den = DTPROB_DENSITY_GAUSSIAN; 326d6685f55SMatthew G. Knepley PetscErrorCode ierr; 327d6685f55SMatthew G. Knepley 328d6685f55SMatthew G. Knepley PetscFunctionBegin; 329d6685f55SMatthew G. Knepley ierr = PetscOptionsBegin(PETSC_COMM_SELF, prefix, "PetscProb Options", "DT");CHKERRQ(ierr); 330d6685f55SMatthew G. Knepley ierr = PetscOptionsEnum(name, "Method to compute PDF <constant, gaussian>", "", DTProbDensityTypes, (PetscEnum) den, (PetscEnum *) &den, NULL);CHKERRQ(ierr); 331d6685f55SMatthew G. Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 332d6685f55SMatthew G. Knepley 333d6685f55SMatthew G. Knepley if (pdf) {PetscValidPointer(pdf, 4); *pdf = NULL;} 334d6685f55SMatthew G. Knepley if (cdf) {PetscValidPointer(cdf, 5); *cdf = NULL;} 335d6685f55SMatthew G. Knepley if (sampler) {PetscValidPointer(sampler, 6); *sampler = NULL;} 336d6685f55SMatthew G. Knepley switch (den) { 337d6685f55SMatthew G. Knepley case DTPROB_DENSITY_CONSTANT: 338d6685f55SMatthew G. Knepley switch (dim) { 339d6685f55SMatthew G. Knepley case 1: 340d6685f55SMatthew G. Knepley if (pdf) *pdf = PetscPDFConstant1D; 341d6685f55SMatthew G. Knepley if (cdf) *cdf = PetscCDFConstant1D; 342d6685f55SMatthew G. Knepley if (sampler) *sampler = PetscPDFSampleConstant1D; 343d6685f55SMatthew G. Knepley break; 344d6685f55SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]); 345d6685f55SMatthew G. Knepley } 346d6685f55SMatthew G. Knepley break; 347d6685f55SMatthew G. Knepley case DTPROB_DENSITY_GAUSSIAN: 348d6685f55SMatthew G. Knepley switch (dim) { 349d6685f55SMatthew G. Knepley case 1: 350d6685f55SMatthew G. Knepley if (pdf) *pdf = PetscPDFGaussian1D; 351d6685f55SMatthew G. Knepley if (cdf) *cdf = PetscCDFGaussian1D; 352d6685f55SMatthew G. Knepley if (sampler) *sampler = PetscPDFSampleGaussian1D; 353d6685f55SMatthew G. Knepley break; 354d6685f55SMatthew G. Knepley case 2: 355d6685f55SMatthew G. Knepley if (pdf) *pdf = PetscPDFGaussian2D; 356d6685f55SMatthew G. Knepley if (sampler) *sampler = PetscPDFSampleGaussian2D; 357d6685f55SMatthew G. Knepley break; 358d6685f55SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]); 359d6685f55SMatthew G. Knepley } 360d6685f55SMatthew G. Knepley break; 361d6685f55SMatthew G. Knepley case DTPROB_DENSITY_MAXWELL_BOLTZMANN: 362d6685f55SMatthew G. Knepley switch (dim) { 363d6685f55SMatthew G. Knepley case 1: 364d6685f55SMatthew G. Knepley if (pdf) *pdf = PetscPDFMaxwellBoltzmann1D; 365d6685f55SMatthew G. Knepley if (cdf) *cdf = PetscCDFMaxwellBoltzmann1D; 366d6685f55SMatthew G. Knepley break; 367d6685f55SMatthew G. Knepley case 2: 368d6685f55SMatthew G. Knepley if (pdf) *pdf = PetscPDFMaxwellBoltzmann2D; 369d6685f55SMatthew G. Knepley if (cdf) *cdf = PetscCDFMaxwellBoltzmann2D; 370d6685f55SMatthew G. Knepley break; 371d6685f55SMatthew G. Knepley case 3: 372d6685f55SMatthew G. Knepley if (pdf) *pdf = PetscPDFMaxwellBoltzmann3D; 373d6685f55SMatthew G. Knepley if (cdf) *cdf = PetscCDFMaxwellBoltzmann3D; 374d6685f55SMatthew G. Knepley break; 375d6685f55SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]); 376d6685f55SMatthew G. Knepley } 377d6685f55SMatthew G. Knepley break; 378d6685f55SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Density type %s is not supported", DTProbDensityTypes[PetscMax(0, PetscMin(den, DTPROB_NUM_DENSITY))]); 379d6685f55SMatthew G. Knepley } 380d6685f55SMatthew G. Knepley PetscFunctionReturn(0); 381d6685f55SMatthew G. Knepley } 382d6685f55SMatthew G. Knepley 383d6685f55SMatthew G. Knepley #ifdef PETSC_HAVE_KS 384d6685f55SMatthew G. Knepley EXTERN_C_BEGIN 385d6685f55SMatthew G. Knepley #include <KolmogorovSmirnovDist.h> 386d6685f55SMatthew G. Knepley EXTERN_C_END 387d6685f55SMatthew G. Knepley #endif 388d6685f55SMatthew G. Knepley 389d6685f55SMatthew G. Knepley /*@C 390d6685f55SMatthew G. Knepley PetscProbComputeKSStatistic - Compute the Kolmogorov-Smirnov statistic for the empirical distribution for an input vector, compared to an analytic CDF. 391d6685f55SMatthew G. Knepley 392d6685f55SMatthew G. Knepley Collective on v 393d6685f55SMatthew G. Knepley 394d6685f55SMatthew G. Knepley Input Parameters: 395d6685f55SMatthew G. Knepley + v - The data vector, blocksize is the sample dimension 396d6685f55SMatthew G. Knepley - cdf - The analytic CDF 397d6685f55SMatthew G. Knepley 398d6685f55SMatthew G. Knepley Output Parameter: 399d6685f55SMatthew G. Knepley . alpha - The KS statisic 400d6685f55SMatthew G. Knepley 401d6685f55SMatthew G. Knepley Level: advanced 402d6685f55SMatthew G. Knepley 403d6685f55SMatthew G. Knepley Note: The Kolmogorov-Smirnov statistic for a given cumulative distribution function $F(x)$ is 404d6685f55SMatthew G. Knepley $ 405d6685f55SMatthew G. Knepley $ D_n = \sup_x \left| F_n(x) - F(x) \right| 406d6685f55SMatthew G. Knepley $ 407d6685f55SMatthew G. Knepley where $\sup_x$ is the supremum of the set of distances, and the empirical distribution 408d6685f55SMatthew G. Knepley function $F_n(x)$ is discrete, and given by 409d6685f55SMatthew G. Knepley $ 410d6685f55SMatthew G. Knepley $ F_n = # of samples <= x / n 411d6685f55SMatthew G. Knepley $ 412d6685f55SMatthew G. Knepley The empirical distribution function $F_n(x)$ is discrete, and thus had a ``stairstep'' 413d6685f55SMatthew G. Knepley cumulative distribution, making $n$ the number of stairs. Intuitively, the statistic takes 414d6685f55SMatthew G. Knepley the largest absolute difference between the two distribution functions across all $x$ values. 415d6685f55SMatthew G. Knepley 416d6685f55SMatthew G. Knepley .seealso: PetscProbFunc 417d6685f55SMatthew G. Knepley @*/ 418d6685f55SMatthew G. Knepley PetscErrorCode PetscProbComputeKSStatistic(Vec v, PetscProbFunc cdf, PetscReal *alpha) 419d6685f55SMatthew G. Knepley { 420d6685f55SMatthew G. Knepley #if !defined(PETSC_HAVE_KS) 421d6685f55SMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject) v), PETSC_ERR_SUP, "No support for Kolmogorov-Smirnov test.\nReconfigure using --download-ks"); 422d6685f55SMatthew G. Knepley #else 423d6685f55SMatthew G. Knepley PetscViewer viewer = NULL; 424d6685f55SMatthew G. Knepley PetscViewerFormat format; 425d6685f55SMatthew G. Knepley PetscDraw draw; 426d6685f55SMatthew G. Knepley PetscDrawLG lg; 427d6685f55SMatthew G. Knepley PetscDrawAxis axis; 428d6685f55SMatthew G. Knepley const PetscScalar *a; 429d6685f55SMatthew G. Knepley PetscReal *speed, Dn = PETSC_MIN_REAL; 430d6685f55SMatthew G. Knepley PetscBool isascii = PETSC_FALSE, isdraw = PETSC_FALSE, flg; 431d6685f55SMatthew G. Knepley PetscInt n, p, dim, d; 432d6685f55SMatthew G. Knepley PetscMPIInt size; 433d6685f55SMatthew G. Knepley const char *names[2] = {"Analytic", "Empirical"}; 434d6685f55SMatthew G. Knepley char title[PETSC_MAX_PATH_LEN]; 435d6685f55SMatthew G. Knepley PetscOptions options; 436d6685f55SMatthew G. Knepley const char *prefix; 437d6685f55SMatthew G. Knepley MPI_Comm comm; 438d6685f55SMatthew G. Knepley PetscErrorCode ierr; 439d6685f55SMatthew G. Knepley 440d6685f55SMatthew G. Knepley PetscFunctionBegin; 441d6685f55SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr); 442d6685f55SMatthew G. Knepley ierr = PetscObjectGetOptionsPrefix((PetscObject) v, &prefix);CHKERRQ(ierr); 443d6685f55SMatthew G. Knepley ierr = PetscObjectGetOptions((PetscObject) v, &options);CHKERRQ(ierr); 444d6685f55SMatthew G. Knepley ierr = PetscOptionsGetViewer(comm, options, prefix, "-ks_monitor", &viewer, &format, &flg);CHKERRQ(ierr); 445d6685f55SMatthew G. Knepley if (flg) { 446d6685f55SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 447d6685f55SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERDRAW, &isdraw);CHKERRQ(ierr); 448d6685f55SMatthew G. Knepley } 449d6685f55SMatthew G. Knepley if (isascii) { 450d6685f55SMatthew G. Knepley ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr); 451d6685f55SMatthew G. Knepley } else if (isdraw) { 452d6685f55SMatthew G. Knepley ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr); 453d6685f55SMatthew G. Knepley ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr); 454d6685f55SMatthew G. Knepley ierr = PetscDrawLGCreate(draw, 2, &lg);CHKERRQ(ierr); 455d6685f55SMatthew G. Knepley ierr = PetscDrawLGSetLegend(lg, names);CHKERRQ(ierr); 456d6685f55SMatthew G. Knepley } 457d6685f55SMatthew G. Knepley 458d6685f55SMatthew G. Knepley ierr = MPI_Comm_size(PetscObjectComm((PetscObject) v), &size);CHKERRMPI(ierr); 459d6685f55SMatthew G. Knepley ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr); 460d6685f55SMatthew G. Knepley ierr = VecGetBlockSize(v, &dim);CHKERRQ(ierr); 461d6685f55SMatthew G. Knepley n /= dim; 462d6685f55SMatthew G. Knepley /* TODO Parallel algorithm is harder */ 463d6685f55SMatthew G. Knepley if (size == 1) { 464d6685f55SMatthew G. Knepley ierr = PetscMalloc1(n, &speed);CHKERRQ(ierr); 465d6685f55SMatthew G. Knepley ierr = VecGetArrayRead(v, &a);CHKERRQ(ierr); 466d6685f55SMatthew G. Knepley for (p = 0; p < n; ++p) { 467d6685f55SMatthew G. Knepley PetscReal mag = 0.; 468d6685f55SMatthew G. Knepley 469d6685f55SMatthew G. Knepley for (d = 0; d < dim; ++d) mag += PetscSqr(PetscRealPart(a[p*dim+d])); 470d6685f55SMatthew G. Knepley speed[p] = PetscSqrtReal(mag); 471d6685f55SMatthew G. Knepley } 472d6685f55SMatthew G. Knepley ierr = PetscSortReal(n, speed);CHKERRQ(ierr); 473d6685f55SMatthew G. Knepley ierr = VecRestoreArrayRead(v, &a);CHKERRQ(ierr); 474d6685f55SMatthew G. Knepley for (p = 0; p < n; ++p) { 475d6685f55SMatthew G. Knepley const PetscReal x = speed[p], Fn = ((PetscReal) p) / n; 476d6685f55SMatthew G. Knepley PetscReal F, vals[2]; 477d6685f55SMatthew G. Knepley 478d6685f55SMatthew G. Knepley ierr = cdf(&x, NULL, &F);CHKERRQ(ierr); 479d6685f55SMatthew G. Knepley Dn = PetscMax(PetscAbsReal(Fn - F), Dn); 480d6685f55SMatthew G. Knepley if (isascii) {ierr = PetscViewerASCIIPrintf(viewer, "x: %g F: %g Fn: %g Dn: %.2g\n", x, F, Fn, Dn);CHKERRQ(ierr);} 481d6685f55SMatthew G. Knepley if (isdraw) {vals[0] = F; vals[1] = Fn; ierr = PetscDrawLGAddCommonPoint(lg, x, vals);CHKERRQ(ierr);} 482d6685f55SMatthew G. Knepley } 483d6685f55SMatthew G. Knepley if (speed[n-1] < 6.) { 484d6685f55SMatthew G. Knepley const PetscReal k = (PetscInt) (6. - speed[n-1]) + 1, dx = (6. - speed[n-1])/k; 485d6685f55SMatthew G. Knepley for (p = 0; p <= k; ++p) { 486d6685f55SMatthew G. Knepley const PetscReal x = speed[n-1] + p*dx, Fn = 1.0; 487d6685f55SMatthew G. Knepley PetscReal F, vals[2]; 488d6685f55SMatthew G. Knepley 489d6685f55SMatthew G. Knepley ierr = cdf(&x, NULL, &F);CHKERRQ(ierr); 490d6685f55SMatthew G. Knepley Dn = PetscMax(PetscAbsReal(Fn - F), Dn); 491d6685f55SMatthew G. Knepley if (isascii) {ierr = PetscViewerASCIIPrintf(viewer, "x: %g F: %g Fn: %g Dn: %.2g\n", x, F, Fn, Dn);CHKERRQ(ierr);} 492d6685f55SMatthew G. Knepley if (isdraw) {vals[0] = F; vals[1] = Fn; ierr = PetscDrawLGAddCommonPoint(lg, x, vals);CHKERRQ(ierr);} 493d6685f55SMatthew G. Knepley } 494d6685f55SMatthew G. Knepley } 495d6685f55SMatthew G. Knepley ierr = PetscFree(speed);CHKERRQ(ierr); 496d6685f55SMatthew G. Knepley } 497d6685f55SMatthew G. Knepley if (isdraw) { 498d6685f55SMatthew G. Knepley ierr = PetscDrawLGGetAxis(lg, &axis);CHKERRQ(ierr); 499d6685f55SMatthew G. Knepley ierr = PetscSNPrintf(title, PETSC_MAX_PATH_LEN, "Kolmogorov-Smirnov Test (Dn %.2g)", Dn);CHKERRQ(ierr); 500d6685f55SMatthew G. Knepley ierr = PetscDrawAxisSetLabels(axis, title, "x", "CDF(x)");CHKERRQ(ierr); 501d6685f55SMatthew G. Knepley ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 502d6685f55SMatthew G. Knepley ierr = PetscDrawLGDestroy(&lg);CHKERRQ(ierr); 503d6685f55SMatthew G. Knepley } 504d6685f55SMatthew G. Knepley if (viewer) { 505d6685f55SMatthew G. Knepley ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 506d6685f55SMatthew G. Knepley ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 507d6685f55SMatthew G. Knepley ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 508d6685f55SMatthew G. Knepley } 509d6685f55SMatthew G. Knepley *alpha = KSfbar((int) n, (double) Dn); 510d6685f55SMatthew G. Knepley PetscFunctionReturn(0); 511d6685f55SMatthew G. Knepley #endif 512d6685f55SMatthew G. Knepley } 513