1*d6685f55SMatthew G. Knepley #include <petscdt.h> /*I "petscdt.h" I*/ 2*d6685f55SMatthew G. Knepley 3*d6685f55SMatthew G. Knepley #include <petscvec.h> 4*d6685f55SMatthew G. Knepley #include <petscdraw.h> 5*d6685f55SMatthew G. Knepley #include <petsc/private/petscimpl.h> 6*d6685f55SMatthew G. Knepley 7*d6685f55SMatthew G. Knepley const char * const DTProbDensityTypes[] = {"constant", "gaussian", "maxwell_boltzmann", "DTProb Density", "DTPROB_DENSITY_", NULL}; 8*d6685f55SMatthew G. Knepley 9*d6685f55SMatthew G. Knepley /*@ 10*d6685f55SMatthew G. Knepley PetscPDFMaxwellBoltzmann1D - PDF for the Maxwell-Boltzmann distribution in 1D 11*d6685f55SMatthew G. Knepley 12*d6685f55SMatthew G. Knepley Not collective 13*d6685f55SMatthew G. Knepley 14*d6685f55SMatthew G. Knepley Input Parameter: 15*d6685f55SMatthew G. Knepley . x - Speed in [0, \infty] 16*d6685f55SMatthew G. Knepley 17*d6685f55SMatthew G. Knepley Output Parameter: 18*d6685f55SMatthew G. Knepley . p - The probability density at x 19*d6685f55SMatthew G. Knepley 20*d6685f55SMatthew G. Knepley Level: beginner 21*d6685f55SMatthew G. Knepley 22*d6685f55SMatthew G. Knepley .seealso: PetscCDFMaxwellBoltzmann1D 23*d6685f55SMatthew G. Knepley @*/ 24*d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFMaxwellBoltzmann1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 25*d6685f55SMatthew G. Knepley { 26*d6685f55SMatthew G. Knepley p[0] = PetscSqrtReal(2. / PETSC_PI) * PetscExpReal(-0.5 * PetscSqr(x[0])); 27*d6685f55SMatthew G. Knepley return 0; 28*d6685f55SMatthew G. Knepley } 29*d6685f55SMatthew G. Knepley 30*d6685f55SMatthew G. Knepley /*@ 31*d6685f55SMatthew G. Knepley PetscCDFMaxwellBoltzmann1D - CDF for the Maxwell-Boltzmann distribution in 1D 32*d6685f55SMatthew G. Knepley 33*d6685f55SMatthew G. Knepley Not collective 34*d6685f55SMatthew G. Knepley 35*d6685f55SMatthew G. Knepley Input Parameter: 36*d6685f55SMatthew G. Knepley . x - Speed in [0, \infty] 37*d6685f55SMatthew G. Knepley 38*d6685f55SMatthew G. Knepley Output Parameter: 39*d6685f55SMatthew G. Knepley . p - The probability density at x 40*d6685f55SMatthew G. Knepley 41*d6685f55SMatthew G. Knepley Level: beginner 42*d6685f55SMatthew G. Knepley 43*d6685f55SMatthew G. Knepley .seealso: PetscPDFMaxwellBoltzmann1D 44*d6685f55SMatthew G. Knepley @*/ 45*d6685f55SMatthew G. Knepley PetscErrorCode PetscCDFMaxwellBoltzmann1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 46*d6685f55SMatthew G. Knepley { 47*d6685f55SMatthew G. Knepley p[0] = PetscErfReal(x[0] / PETSC_SQRT2); 48*d6685f55SMatthew G. Knepley return 0; 49*d6685f55SMatthew G. Knepley } 50*d6685f55SMatthew G. Knepley 51*d6685f55SMatthew G. Knepley /*@ 52*d6685f55SMatthew G. Knepley PetscPDFMaxwellBoltzmann2D - PDF for the Maxwell-Boltzmann distribution in 2D 53*d6685f55SMatthew G. Knepley 54*d6685f55SMatthew G. Knepley Not collective 55*d6685f55SMatthew G. Knepley 56*d6685f55SMatthew G. Knepley Input Parameter: 57*d6685f55SMatthew G. Knepley . x - Speed in [0, \infty] 58*d6685f55SMatthew G. Knepley 59*d6685f55SMatthew G. Knepley Output Parameter: 60*d6685f55SMatthew G. Knepley . p - The probability density at x 61*d6685f55SMatthew G. Knepley 62*d6685f55SMatthew G. Knepley Level: beginner 63*d6685f55SMatthew G. Knepley 64*d6685f55SMatthew G. Knepley .seealso: PetscCDFMaxwellBoltzmann2D 65*d6685f55SMatthew G. Knepley @*/ 66*d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFMaxwellBoltzmann2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 67*d6685f55SMatthew G. Knepley { 68*d6685f55SMatthew G. Knepley p[0] = x[0] * PetscExpReal(-0.5 * PetscSqr(x[0])); 69*d6685f55SMatthew G. Knepley return 0; 70*d6685f55SMatthew G. Knepley } 71*d6685f55SMatthew G. Knepley 72*d6685f55SMatthew G. Knepley /*@ 73*d6685f55SMatthew G. Knepley PetscCDFMaxwellBoltzmann2D - CDF for the Maxwell-Boltzmann distribution in 2D 74*d6685f55SMatthew G. Knepley 75*d6685f55SMatthew G. Knepley Not collective 76*d6685f55SMatthew G. Knepley 77*d6685f55SMatthew G. Knepley Input Parameter: 78*d6685f55SMatthew G. Knepley . x - Speed in [0, \infty] 79*d6685f55SMatthew G. Knepley 80*d6685f55SMatthew G. Knepley Output Parameter: 81*d6685f55SMatthew G. Knepley . p - The probability density at x 82*d6685f55SMatthew G. Knepley 83*d6685f55SMatthew G. Knepley Level: beginner 84*d6685f55SMatthew G. Knepley 85*d6685f55SMatthew G. Knepley .seealso: PetscPDFMaxwellBoltzmann2D 86*d6685f55SMatthew G. Knepley @*/ 87*d6685f55SMatthew G. Knepley PetscErrorCode PetscCDFMaxwellBoltzmann2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 88*d6685f55SMatthew G. Knepley { 89*d6685f55SMatthew G. Knepley p[0] = 1. - PetscExpReal(-0.5 * PetscSqr(x[0])); 90*d6685f55SMatthew G. Knepley return 0; 91*d6685f55SMatthew G. Knepley } 92*d6685f55SMatthew G. Knepley 93*d6685f55SMatthew G. Knepley /*@ 94*d6685f55SMatthew G. Knepley PetscPDFMaxwellBoltzmann3D - PDF for the Maxwell-Boltzmann distribution in 3D 95*d6685f55SMatthew G. Knepley 96*d6685f55SMatthew G. Knepley Not collective 97*d6685f55SMatthew G. Knepley 98*d6685f55SMatthew G. Knepley Input Parameter: 99*d6685f55SMatthew G. Knepley . x - Speed in [0, \infty] 100*d6685f55SMatthew G. Knepley 101*d6685f55SMatthew G. Knepley Output Parameter: 102*d6685f55SMatthew G. Knepley . p - The probability density at x 103*d6685f55SMatthew G. Knepley 104*d6685f55SMatthew G. Knepley Level: beginner 105*d6685f55SMatthew G. Knepley 106*d6685f55SMatthew G. Knepley .seealso: PetscCDFMaxwellBoltzmann3D 107*d6685f55SMatthew G. Knepley @*/ 108*d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFMaxwellBoltzmann3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 109*d6685f55SMatthew G. Knepley { 110*d6685f55SMatthew G. Knepley p[0] = PetscSqrtReal(2. / PETSC_PI) * PetscSqr(x[0]) * PetscExpReal(-0.5 * PetscSqr(x[0])); 111*d6685f55SMatthew G. Knepley return 0; 112*d6685f55SMatthew G. Knepley } 113*d6685f55SMatthew G. Knepley 114*d6685f55SMatthew G. Knepley /*@ 115*d6685f55SMatthew G. Knepley PetscCDFMaxwellBoltzmann3D - CDF for the Maxwell-Boltzmann distribution in 3D 116*d6685f55SMatthew G. Knepley 117*d6685f55SMatthew G. Knepley Not collective 118*d6685f55SMatthew G. Knepley 119*d6685f55SMatthew G. Knepley Input Parameter: 120*d6685f55SMatthew G. Knepley . x - Speed in [0, \infty] 121*d6685f55SMatthew G. Knepley 122*d6685f55SMatthew G. Knepley Output Parameter: 123*d6685f55SMatthew G. Knepley . p - The probability density at x 124*d6685f55SMatthew G. Knepley 125*d6685f55SMatthew G. Knepley Level: beginner 126*d6685f55SMatthew G. Knepley 127*d6685f55SMatthew G. Knepley .seealso: PetscPDFMaxwellBoltzmann3D 128*d6685f55SMatthew G. Knepley @*/ 129*d6685f55SMatthew G. Knepley PetscErrorCode PetscCDFMaxwellBoltzmann3D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 130*d6685f55SMatthew G. Knepley { 131*d6685f55SMatthew G. Knepley p[0] = PetscErfReal(x[0] / PETSC_SQRT2) - PetscSqrtReal(2. / PETSC_PI) * x[0] * PetscExpReal(-0.5 * PetscSqr(x[0])); 132*d6685f55SMatthew G. Knepley return 0; 133*d6685f55SMatthew G. Knepley } 134*d6685f55SMatthew G. Knepley 135*d6685f55SMatthew G. Knepley /*@ 136*d6685f55SMatthew G. Knepley PetscPDFGaussian1D - PDF for the Gaussian distribution in 1D 137*d6685f55SMatthew G. Knepley 138*d6685f55SMatthew G. Knepley Not collective 139*d6685f55SMatthew G. Knepley 140*d6685f55SMatthew G. Knepley Input Parameter: 141*d6685f55SMatthew G. Knepley . x - Coordinate in [-\infty, \infty] 142*d6685f55SMatthew G. Knepley 143*d6685f55SMatthew G. Knepley Output Parameter: 144*d6685f55SMatthew G. Knepley . p - The probability density at x 145*d6685f55SMatthew G. Knepley 146*d6685f55SMatthew G. Knepley Level: beginner 147*d6685f55SMatthew G. Knepley 148*d6685f55SMatthew G. Knepley .seealso: PetscPDFMaxwellBoltzmann3D 149*d6685f55SMatthew G. Knepley @*/ 150*d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFGaussian1D(const PetscReal x[], const PetscReal scale[], PetscReal p[]) 151*d6685f55SMatthew G. Knepley { 152*d6685f55SMatthew G. Knepley const PetscReal sigma = scale ? scale[0] : 1.; 153*d6685f55SMatthew G. Knepley p[0] = PetscSqrtReal(1. / (2.*PETSC_PI)) * PetscExpReal(-0.5 * PetscSqr(x[0] / sigma)) / sigma; 154*d6685f55SMatthew G. Knepley return 0; 155*d6685f55SMatthew G. Knepley } 156*d6685f55SMatthew G. Knepley 157*d6685f55SMatthew G. Knepley PetscErrorCode PetscCDFGaussian1D(const PetscReal x[], const PetscReal scale[], PetscReal p[]) 158*d6685f55SMatthew G. Knepley { 159*d6685f55SMatthew G. Knepley const PetscReal sigma = scale ? scale[0] : 1.; 160*d6685f55SMatthew G. Knepley p[0] = 0.5 * (1. + PetscErfReal(x[0] / PETSC_SQRT2 / sigma)); 161*d6685f55SMatthew G. Knepley return 0; 162*d6685f55SMatthew G. Knepley } 163*d6685f55SMatthew G. Knepley 164*d6685f55SMatthew G. Knepley /*@ 165*d6685f55SMatthew G. Knepley PetscPDFSampleGaussian1D - Sample uniformly from a Gaussian distribution in 1D 166*d6685f55SMatthew G. Knepley 167*d6685f55SMatthew G. Knepley Not collective 168*d6685f55SMatthew G. Knepley 169*d6685f55SMatthew G. Knepley Input Parameter: 170*d6685f55SMatthew G. Knepley . p - A uniform variable on [0, 1] 171*d6685f55SMatthew G. Knepley 172*d6685f55SMatthew G. Knepley Output Parameter: 173*d6685f55SMatthew G. Knepley . x - Coordinate in [-\infty, \infty] 174*d6685f55SMatthew G. Knepley 175*d6685f55SMatthew G. Knepley Level: beginner 176*d6685f55SMatthew G. Knepley 177*d6685f55SMatthew G. Knepley Note: http://www.mimirgames.com/articles/programming/approximations-of-the-inverse-error-function/ 178*d6685f55SMatthew G. Knepley https://stackoverflow.com/questions/27229371/inverse-error-function-in-c 179*d6685f55SMatthew G. Knepley @*/ 180*d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFSampleGaussian1D(const PetscReal p[], const PetscReal dummy[], PetscReal x[]) 181*d6685f55SMatthew G. Knepley { 182*d6685f55SMatthew G. Knepley const PetscReal q = 2*p[0] - 1.; 183*d6685f55SMatthew G. Knepley const PetscInt maxIter = 100; 184*d6685f55SMatthew G. Knepley PetscReal ck[100], r = 0.; 185*d6685f55SMatthew G. Knepley PetscInt k, m; 186*d6685f55SMatthew G. Knepley 187*d6685f55SMatthew G. Knepley PetscFunctionBeginHot; 188*d6685f55SMatthew G. Knepley /* Transform input to [-1, 1] since the code below computes the inverse error function */ 189*d6685f55SMatthew G. Knepley for (k = 0; k < maxIter; ++k) ck[k] = 0.; 190*d6685f55SMatthew G. Knepley ck[0] = 1; 191*d6685f55SMatthew G. Knepley r = ck[0]* (PetscSqrtReal(PETSC_PI) / 2.) * q; 192*d6685f55SMatthew G. Knepley for (k = 1; k < maxIter; ++k) { 193*d6685f55SMatthew G. Knepley const PetscReal temp = 2.*k + 1.; 194*d6685f55SMatthew G. Knepley 195*d6685f55SMatthew G. Knepley for (m = 0; m <= k-1; ++m) { 196*d6685f55SMatthew G. Knepley PetscReal denom = (m + 1.)*(2.*m + 1.); 197*d6685f55SMatthew G. Knepley 198*d6685f55SMatthew G. Knepley ck[k] += (ck[m]*ck[k-1-m])/denom; 199*d6685f55SMatthew G. Knepley } 200*d6685f55SMatthew G. Knepley r += (ck[k]/temp) * PetscPowReal((PetscSqrtReal(PETSC_PI)/2.) * q, 2.*k+1.); 201*d6685f55SMatthew G. Knepley } 202*d6685f55SMatthew G. Knepley /* Scale erfinv() by \sqrt{\pi/2} */ 203*d6685f55SMatthew G. Knepley x[0] = PetscSqrtReal(PETSC_PI * 0.5) * r; 204*d6685f55SMatthew G. Knepley PetscFunctionReturn(0); 205*d6685f55SMatthew G. Knepley } 206*d6685f55SMatthew G. Knepley 207*d6685f55SMatthew G. Knepley /*@ 208*d6685f55SMatthew G. Knepley PetscPDFGaussian2D - PDF for the Gaussian distribution in 2D 209*d6685f55SMatthew G. Knepley 210*d6685f55SMatthew G. Knepley Not collective 211*d6685f55SMatthew G. Knepley 212*d6685f55SMatthew G. Knepley Input Parameter: 213*d6685f55SMatthew G. Knepley . x - Coordinate in [-\infty, \infty]^2 214*d6685f55SMatthew G. Knepley 215*d6685f55SMatthew G. Knepley Output Parameter: 216*d6685f55SMatthew G. Knepley . p - The probability density at x 217*d6685f55SMatthew G. Knepley 218*d6685f55SMatthew G. Knepley Level: beginner 219*d6685f55SMatthew G. Knepley 220*d6685f55SMatthew G. Knepley .seealso: PetscPDFSampleGaussian2D(), PetscPDFMaxwellBoltzmann3D() 221*d6685f55SMatthew G. Knepley @*/ 222*d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFGaussian2D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 223*d6685f55SMatthew G. Knepley { 224*d6685f55SMatthew G. Knepley p[0] = (1. / PETSC_PI) * PetscExpReal(-0.5 * (PetscSqr(x[0]) + PetscSqr(x[1]))); 225*d6685f55SMatthew G. Knepley return 0; 226*d6685f55SMatthew G. Knepley } 227*d6685f55SMatthew G. Knepley 228*d6685f55SMatthew G. Knepley /*@ 229*d6685f55SMatthew G. Knepley PetscPDFSampleGaussian2D - Sample uniformly from a Gaussian distribution in 2D 230*d6685f55SMatthew G. Knepley 231*d6685f55SMatthew G. Knepley Not collective 232*d6685f55SMatthew G. Knepley 233*d6685f55SMatthew G. Knepley Input Parameter: 234*d6685f55SMatthew G. Knepley . p - A uniform variable on [0, 1]^2 235*d6685f55SMatthew G. Knepley 236*d6685f55SMatthew G. Knepley Output Parameter: 237*d6685f55SMatthew G. Knepley . x - Coordinate in [-\infty, \infty]^2 238*d6685f55SMatthew G. Knepley 239*d6685f55SMatthew G. Knepley Level: beginner 240*d6685f55SMatthew G. Knepley 241*d6685f55SMatthew G. Knepley Note: https://en.wikipedia.org/wiki/Box%E2%80%93Muller_transform 242*d6685f55SMatthew G. Knepley 243*d6685f55SMatthew G. Knepley .seealso: PetscPDFGaussian2D(), PetscPDFMaxwellBoltzmann3D() 244*d6685f55SMatthew G. Knepley @*/ 245*d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFSampleGaussian2D(const PetscReal p[], const PetscReal dummy[], PetscReal x[]) 246*d6685f55SMatthew G. Knepley { 247*d6685f55SMatthew G. Knepley const PetscReal mag = PetscSqrtReal(-2.0 * PetscLogReal(p[0])); 248*d6685f55SMatthew G. Knepley x[0] = mag * PetscCosReal(2.0 * PETSC_PI * p[1]); 249*d6685f55SMatthew G. Knepley x[1] = mag * PetscSinReal(2.0 * PETSC_PI * p[1]); 250*d6685f55SMatthew G. Knepley return 0; 251*d6685f55SMatthew G. Knepley } 252*d6685f55SMatthew G. Knepley 253*d6685f55SMatthew G. Knepley /*@ 254*d6685f55SMatthew G. Knepley PetscPDFConstant1D - PDF for the uniform distribution in 1D 255*d6685f55SMatthew G. Knepley 256*d6685f55SMatthew G. Knepley Not collective 257*d6685f55SMatthew G. Knepley 258*d6685f55SMatthew G. Knepley Input Parameter: 259*d6685f55SMatthew G. Knepley . x - Coordinate in [-1, 1] 260*d6685f55SMatthew G. Knepley 261*d6685f55SMatthew G. Knepley Output Parameter: 262*d6685f55SMatthew G. Knepley . p - The probability density at x 263*d6685f55SMatthew G. Knepley 264*d6685f55SMatthew G. Knepley Level: beginner 265*d6685f55SMatthew G. Knepley 266*d6685f55SMatthew G. Knepley .seealso: PetscCDFConstant1D() 267*d6685f55SMatthew G. Knepley @*/ 268*d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFConstant1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 269*d6685f55SMatthew G. Knepley { 270*d6685f55SMatthew G. Knepley p[0] = x[0] > -1. && x[0] <= 1. ? 0.5 : 0.; 271*d6685f55SMatthew G. Knepley return 0; 272*d6685f55SMatthew G. Knepley } 273*d6685f55SMatthew G. Knepley 274*d6685f55SMatthew G. Knepley PetscErrorCode PetscCDFConstant1D(const PetscReal x[], const PetscReal dummy[], PetscReal p[]) 275*d6685f55SMatthew G. Knepley { 276*d6685f55SMatthew G. Knepley p[0] = x[0] <= -1. ? 0. : (x[0] > 1. ? 1. : 0.5*(x[0] + 1.)); 277*d6685f55SMatthew G. Knepley return 0; 278*d6685f55SMatthew G. Knepley } 279*d6685f55SMatthew G. Knepley 280*d6685f55SMatthew G. Knepley /*@ 281*d6685f55SMatthew G. Knepley PetscPDFSampleConstant1D - Sample uniformly from a uniform distribution on [-1, 1] in 1D 282*d6685f55SMatthew G. Knepley 283*d6685f55SMatthew G. Knepley Not collective 284*d6685f55SMatthew G. Knepley 285*d6685f55SMatthew G. Knepley Input Parameter: 286*d6685f55SMatthew G. Knepley . p - A uniform variable on [0, 1] 287*d6685f55SMatthew G. Knepley 288*d6685f55SMatthew G. Knepley Output Parameter: 289*d6685f55SMatthew G. Knepley . x - Coordinate in [-1, 1] 290*d6685f55SMatthew G. Knepley 291*d6685f55SMatthew G. Knepley Level: beginner 292*d6685f55SMatthew G. Knepley 293*d6685f55SMatthew G. Knepley .seealso: PetscPDFConstant1D(), PetscCDFConstant1D() 294*d6685f55SMatthew G. Knepley @*/ 295*d6685f55SMatthew G. Knepley PetscErrorCode PetscPDFSampleConstant1D(const PetscReal p[], const PetscReal dummy[], PetscReal x[]) 296*d6685f55SMatthew G. Knepley { 297*d6685f55SMatthew G. Knepley x[0] = 2.*p[1] - 1.; 298*d6685f55SMatthew G. Knepley return 0; 299*d6685f55SMatthew G. Knepley } 300*d6685f55SMatthew G. Knepley 301*d6685f55SMatthew G. Knepley /*@C 302*d6685f55SMatthew G. Knepley PetscProbCreateFromOptions - Return the probability distribution specified by the argumetns and options 303*d6685f55SMatthew G. Knepley 304*d6685f55SMatthew G. Knepley Not collective 305*d6685f55SMatthew G. Knepley 306*d6685f55SMatthew G. Knepley Input Parameters: 307*d6685f55SMatthew G. Knepley + dim - The dimension of sample points 308*d6685f55SMatthew G. Knepley . prefix - The options prefix, or NULL 309*d6685f55SMatthew G. Knepley - name - The option name for the probility distribution type 310*d6685f55SMatthew G. Knepley 311*d6685f55SMatthew G. Knepley Output Parameters: 312*d6685f55SMatthew G. Knepley + pdf - The PDF of this type 313*d6685f55SMatthew G. Knepley . cdf - The CDF of this type 314*d6685f55SMatthew G. Knepley - sampler - The PDF sampler of this type 315*d6685f55SMatthew G. Knepley 316*d6685f55SMatthew G. Knepley Level: intermediate 317*d6685f55SMatthew G. Knepley 318*d6685f55SMatthew G. Knepley .seealso: PetscPDFMaxwellBoltzmann1D(), PetscPDFGaussian1D(), PetscPDFConstant1D() 319*d6685f55SMatthew G. Knepley @*/ 320*d6685f55SMatthew G. Knepley PetscErrorCode PetscProbCreateFromOptions(PetscInt dim, const char prefix[], const char name[], PetscProbFunc *pdf, PetscProbFunc *cdf, PetscProbFunc *sampler) 321*d6685f55SMatthew G. Knepley { 322*d6685f55SMatthew G. Knepley DTProbDensityType den = DTPROB_DENSITY_GAUSSIAN; 323*d6685f55SMatthew G. Knepley PetscErrorCode ierr; 324*d6685f55SMatthew G. Knepley 325*d6685f55SMatthew G. Knepley PetscFunctionBegin; 326*d6685f55SMatthew G. Knepley ierr = PetscOptionsBegin(PETSC_COMM_SELF, prefix, "PetscProb Options", "DT");CHKERRQ(ierr); 327*d6685f55SMatthew G. Knepley ierr = PetscOptionsEnum(name, "Method to compute PDF <constant, gaussian>", "", DTProbDensityTypes, (PetscEnum) den, (PetscEnum *) &den, NULL);CHKERRQ(ierr); 328*d6685f55SMatthew G. Knepley ierr = PetscOptionsEnd();CHKERRQ(ierr); 329*d6685f55SMatthew G. Knepley 330*d6685f55SMatthew G. Knepley if (pdf) {PetscValidPointer(pdf, 4); *pdf = NULL;} 331*d6685f55SMatthew G. Knepley if (cdf) {PetscValidPointer(cdf, 5); *cdf = NULL;} 332*d6685f55SMatthew G. Knepley if (sampler) {PetscValidPointer(sampler, 6); *sampler = NULL;} 333*d6685f55SMatthew G. Knepley switch (den) { 334*d6685f55SMatthew G. Knepley case DTPROB_DENSITY_CONSTANT: 335*d6685f55SMatthew G. Knepley switch (dim) { 336*d6685f55SMatthew G. Knepley case 1: 337*d6685f55SMatthew G. Knepley if (pdf) *pdf = PetscPDFConstant1D; 338*d6685f55SMatthew G. Knepley if (cdf) *cdf = PetscCDFConstant1D; 339*d6685f55SMatthew G. Knepley if (sampler) *sampler = PetscPDFSampleConstant1D; 340*d6685f55SMatthew G. Knepley break; 341*d6685f55SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]); 342*d6685f55SMatthew G. Knepley } 343*d6685f55SMatthew G. Knepley break; 344*d6685f55SMatthew G. Knepley case DTPROB_DENSITY_GAUSSIAN: 345*d6685f55SMatthew G. Knepley switch (dim) { 346*d6685f55SMatthew G. Knepley case 1: 347*d6685f55SMatthew G. Knepley if (pdf) *pdf = PetscPDFGaussian1D; 348*d6685f55SMatthew G. Knepley if (cdf) *cdf = PetscCDFGaussian1D; 349*d6685f55SMatthew G. Knepley if (sampler) *sampler = PetscPDFSampleGaussian1D; 350*d6685f55SMatthew G. Knepley break; 351*d6685f55SMatthew G. Knepley case 2: 352*d6685f55SMatthew G. Knepley if (pdf) *pdf = PetscPDFGaussian2D; 353*d6685f55SMatthew G. Knepley if (sampler) *sampler = PetscPDFSampleGaussian2D; 354*d6685f55SMatthew G. Knepley break; 355*d6685f55SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]); 356*d6685f55SMatthew G. Knepley } 357*d6685f55SMatthew G. Knepley break; 358*d6685f55SMatthew G. Knepley case DTPROB_DENSITY_MAXWELL_BOLTZMANN: 359*d6685f55SMatthew G. Knepley switch (dim) { 360*d6685f55SMatthew G. Knepley case 1: 361*d6685f55SMatthew G. Knepley if (pdf) *pdf = PetscPDFMaxwellBoltzmann1D; 362*d6685f55SMatthew G. Knepley if (cdf) *cdf = PetscCDFMaxwellBoltzmann1D; 363*d6685f55SMatthew G. Knepley break; 364*d6685f55SMatthew G. Knepley case 2: 365*d6685f55SMatthew G. Knepley if (pdf) *pdf = PetscPDFMaxwellBoltzmann2D; 366*d6685f55SMatthew G. Knepley if (cdf) *cdf = PetscCDFMaxwellBoltzmann2D; 367*d6685f55SMatthew G. Knepley break; 368*d6685f55SMatthew G. Knepley case 3: 369*d6685f55SMatthew G. Knepley if (pdf) *pdf = PetscPDFMaxwellBoltzmann3D; 370*d6685f55SMatthew G. Knepley if (cdf) *cdf = PetscCDFMaxwellBoltzmann3D; 371*d6685f55SMatthew G. Knepley break; 372*d6685f55SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %" PetscInt_FMT " not supported for density type %s", dim, DTProbDensityTypes[den]); 373*d6685f55SMatthew G. Knepley } 374*d6685f55SMatthew G. Knepley break; 375*d6685f55SMatthew G. Knepley default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Density type %s is not supported", DTProbDensityTypes[PetscMax(0, PetscMin(den, DTPROB_NUM_DENSITY))]); 376*d6685f55SMatthew G. Knepley } 377*d6685f55SMatthew G. Knepley PetscFunctionReturn(0); 378*d6685f55SMatthew G. Knepley } 379*d6685f55SMatthew G. Knepley 380*d6685f55SMatthew G. Knepley #ifdef PETSC_HAVE_KS 381*d6685f55SMatthew G. Knepley EXTERN_C_BEGIN 382*d6685f55SMatthew G. Knepley #include <KolmogorovSmirnovDist.h> 383*d6685f55SMatthew G. Knepley EXTERN_C_END 384*d6685f55SMatthew G. Knepley #endif 385*d6685f55SMatthew G. Knepley 386*d6685f55SMatthew G. Knepley /*@C 387*d6685f55SMatthew G. Knepley PetscProbComputeKSStatistic - Compute the Kolmogorov-Smirnov statistic for the empirical distribution for an input vector, compared to an analytic CDF. 388*d6685f55SMatthew G. Knepley 389*d6685f55SMatthew G. Knepley Collective on v 390*d6685f55SMatthew G. Knepley 391*d6685f55SMatthew G. Knepley Input Parameters: 392*d6685f55SMatthew G. Knepley + v - The data vector, blocksize is the sample dimension 393*d6685f55SMatthew G. Knepley - cdf - The analytic CDF 394*d6685f55SMatthew G. Knepley 395*d6685f55SMatthew G. Knepley Output Parameter: 396*d6685f55SMatthew G. Knepley . alpha - The KS statisic 397*d6685f55SMatthew G. Knepley 398*d6685f55SMatthew G. Knepley Level: advanced 399*d6685f55SMatthew G. Knepley 400*d6685f55SMatthew G. Knepley Note: The Kolmogorov-Smirnov statistic for a given cumulative distribution function $F(x)$ is 401*d6685f55SMatthew G. Knepley $ 402*d6685f55SMatthew G. Knepley $ D_n = \sup_x \left| F_n(x) - F(x) \right| 403*d6685f55SMatthew G. Knepley $ 404*d6685f55SMatthew G. Knepley where $\sup_x$ is the supremum of the set of distances, and the empirical distribution 405*d6685f55SMatthew G. Knepley function $F_n(x)$ is discrete, and given by 406*d6685f55SMatthew G. Knepley $ 407*d6685f55SMatthew G. Knepley $ F_n = # of samples <= x / n 408*d6685f55SMatthew G. Knepley $ 409*d6685f55SMatthew G. Knepley The empirical distribution function $F_n(x)$ is discrete, and thus had a ``stairstep'' 410*d6685f55SMatthew G. Knepley cumulative distribution, making $n$ the number of stairs. Intuitively, the statistic takes 411*d6685f55SMatthew G. Knepley the largest absolute difference between the two distribution functions across all $x$ values. 412*d6685f55SMatthew G. Knepley 413*d6685f55SMatthew G. Knepley .seealso: PetscProbFunc 414*d6685f55SMatthew G. Knepley @*/ 415*d6685f55SMatthew G. Knepley PetscErrorCode PetscProbComputeKSStatistic(Vec v, PetscProbFunc cdf, PetscReal *alpha) 416*d6685f55SMatthew G. Knepley { 417*d6685f55SMatthew G. Knepley #if !defined(PETSC_HAVE_KS) 418*d6685f55SMatthew G. Knepley SETERRQ(PetscObjectComm((PetscObject) v), PETSC_ERR_SUP, "No support for Kolmogorov-Smirnov test.\nReconfigure using --download-ks"); 419*d6685f55SMatthew G. Knepley #else 420*d6685f55SMatthew G. Knepley PetscViewer viewer = NULL; 421*d6685f55SMatthew G. Knepley PetscViewerFormat format; 422*d6685f55SMatthew G. Knepley PetscDraw draw; 423*d6685f55SMatthew G. Knepley PetscDrawLG lg; 424*d6685f55SMatthew G. Knepley PetscDrawAxis axis; 425*d6685f55SMatthew G. Knepley const PetscScalar *a; 426*d6685f55SMatthew G. Knepley PetscReal *speed, Dn = PETSC_MIN_REAL; 427*d6685f55SMatthew G. Knepley PetscBool isascii = PETSC_FALSE, isdraw = PETSC_FALSE, flg; 428*d6685f55SMatthew G. Knepley PetscInt n, p, dim, d; 429*d6685f55SMatthew G. Knepley PetscMPIInt size; 430*d6685f55SMatthew G. Knepley const char *names[2] = {"Analytic", "Empirical"}; 431*d6685f55SMatthew G. Knepley char title[PETSC_MAX_PATH_LEN]; 432*d6685f55SMatthew G. Knepley PetscOptions options; 433*d6685f55SMatthew G. Knepley const char *prefix; 434*d6685f55SMatthew G. Knepley MPI_Comm comm; 435*d6685f55SMatthew G. Knepley PetscErrorCode ierr; 436*d6685f55SMatthew G. Knepley 437*d6685f55SMatthew G. Knepley PetscFunctionBegin; 438*d6685f55SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) v, &comm);CHKERRQ(ierr); 439*d6685f55SMatthew G. Knepley ierr = PetscObjectGetOptionsPrefix((PetscObject) v, &prefix);CHKERRQ(ierr); 440*d6685f55SMatthew G. Knepley ierr = PetscObjectGetOptions((PetscObject) v, &options);CHKERRQ(ierr); 441*d6685f55SMatthew G. Knepley ierr = PetscOptionsGetViewer(comm, options, prefix, "-ks_monitor", &viewer, &format, &flg);CHKERRQ(ierr); 442*d6685f55SMatthew G. Knepley if (flg) { 443*d6685f55SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);CHKERRQ(ierr); 444*d6685f55SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERDRAW, &isdraw);CHKERRQ(ierr); 445*d6685f55SMatthew G. Knepley } 446*d6685f55SMatthew G. Knepley if (isascii) { 447*d6685f55SMatthew G. Knepley ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr); 448*d6685f55SMatthew G. Knepley } else if (isdraw) { 449*d6685f55SMatthew G. Knepley ierr = PetscViewerPushFormat(viewer, format);CHKERRQ(ierr); 450*d6685f55SMatthew G. Knepley ierr = PetscViewerDrawGetDraw(viewer, 0, &draw);CHKERRQ(ierr); 451*d6685f55SMatthew G. Knepley ierr = PetscDrawLGCreate(draw, 2, &lg);CHKERRQ(ierr); 452*d6685f55SMatthew G. Knepley ierr = PetscDrawLGSetLegend(lg, names);CHKERRQ(ierr); 453*d6685f55SMatthew G. Knepley } 454*d6685f55SMatthew G. Knepley 455*d6685f55SMatthew G. Knepley ierr = MPI_Comm_size(PetscObjectComm((PetscObject) v), &size);CHKERRMPI(ierr); 456*d6685f55SMatthew G. Knepley ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr); 457*d6685f55SMatthew G. Knepley ierr = VecGetBlockSize(v, &dim);CHKERRQ(ierr); 458*d6685f55SMatthew G. Knepley n /= dim; 459*d6685f55SMatthew G. Knepley /* TODO Parallel algorithm is harder */ 460*d6685f55SMatthew G. Knepley if (size == 1) { 461*d6685f55SMatthew G. Knepley ierr = PetscMalloc1(n, &speed);CHKERRQ(ierr); 462*d6685f55SMatthew G. Knepley ierr = VecGetArrayRead(v, &a);CHKERRQ(ierr); 463*d6685f55SMatthew G. Knepley for (p = 0; p < n; ++p) { 464*d6685f55SMatthew G. Knepley PetscReal mag = 0.; 465*d6685f55SMatthew G. Knepley 466*d6685f55SMatthew G. Knepley for (d = 0; d < dim; ++d) mag += PetscSqr(PetscRealPart(a[p*dim+d])); 467*d6685f55SMatthew G. Knepley speed[p] = PetscSqrtReal(mag); 468*d6685f55SMatthew G. Knepley } 469*d6685f55SMatthew G. Knepley ierr = PetscSortReal(n, speed);CHKERRQ(ierr); 470*d6685f55SMatthew G. Knepley ierr = VecRestoreArrayRead(v, &a);CHKERRQ(ierr); 471*d6685f55SMatthew G. Knepley for (p = 0; p < n; ++p) { 472*d6685f55SMatthew G. Knepley const PetscReal x = speed[p], Fn = ((PetscReal) p) / n; 473*d6685f55SMatthew G. Knepley PetscReal F, vals[2]; 474*d6685f55SMatthew G. Knepley 475*d6685f55SMatthew G. Knepley ierr = cdf(&x, NULL, &F);CHKERRQ(ierr); 476*d6685f55SMatthew G. Knepley Dn = PetscMax(PetscAbsReal(Fn - F), Dn); 477*d6685f55SMatthew G. Knepley if (isascii) {ierr = PetscViewerASCIIPrintf(viewer, "x: %g F: %g Fn: %g Dn: %.2g\n", x, F, Fn, Dn);CHKERRQ(ierr);} 478*d6685f55SMatthew G. Knepley if (isdraw) {vals[0] = F; vals[1] = Fn; ierr = PetscDrawLGAddCommonPoint(lg, x, vals);CHKERRQ(ierr);} 479*d6685f55SMatthew G. Knepley } 480*d6685f55SMatthew G. Knepley if (speed[n-1] < 6.) { 481*d6685f55SMatthew G. Knepley const PetscReal k = (PetscInt) (6. - speed[n-1]) + 1, dx = (6. - speed[n-1])/k; 482*d6685f55SMatthew G. Knepley for (p = 0; p <= k; ++p) { 483*d6685f55SMatthew G. Knepley const PetscReal x = speed[n-1] + p*dx, Fn = 1.0; 484*d6685f55SMatthew G. Knepley PetscReal F, vals[2]; 485*d6685f55SMatthew G. Knepley 486*d6685f55SMatthew G. Knepley ierr = cdf(&x, NULL, &F);CHKERRQ(ierr); 487*d6685f55SMatthew G. Knepley Dn = PetscMax(PetscAbsReal(Fn - F), Dn); 488*d6685f55SMatthew G. Knepley if (isascii) {ierr = PetscViewerASCIIPrintf(viewer, "x: %g F: %g Fn: %g Dn: %.2g\n", x, F, Fn, Dn);CHKERRQ(ierr);} 489*d6685f55SMatthew G. Knepley if (isdraw) {vals[0] = F; vals[1] = Fn; ierr = PetscDrawLGAddCommonPoint(lg, x, vals);CHKERRQ(ierr);} 490*d6685f55SMatthew G. Knepley } 491*d6685f55SMatthew G. Knepley } 492*d6685f55SMatthew G. Knepley ierr = PetscFree(speed);CHKERRQ(ierr); 493*d6685f55SMatthew G. Knepley } 494*d6685f55SMatthew G. Knepley if (isdraw) { 495*d6685f55SMatthew G. Knepley ierr = PetscDrawLGGetAxis(lg, &axis);CHKERRQ(ierr); 496*d6685f55SMatthew G. Knepley ierr = PetscSNPrintf(title, PETSC_MAX_PATH_LEN, "Kolmogorov-Smirnov Test (Dn %.2g)", Dn);CHKERRQ(ierr); 497*d6685f55SMatthew G. Knepley ierr = PetscDrawAxisSetLabels(axis, title, "x", "CDF(x)");CHKERRQ(ierr); 498*d6685f55SMatthew G. Knepley ierr = PetscDrawLGDraw(lg);CHKERRQ(ierr); 499*d6685f55SMatthew G. Knepley ierr = PetscDrawLGDestroy(&lg);CHKERRQ(ierr); 500*d6685f55SMatthew G. Knepley } 501*d6685f55SMatthew G. Knepley if (viewer) { 502*d6685f55SMatthew G. Knepley ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 503*d6685f55SMatthew G. Knepley ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 504*d6685f55SMatthew G. Knepley ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 505*d6685f55SMatthew G. Knepley } 506*d6685f55SMatthew G. Knepley *alpha = KSfbar((int) n, (double) Dn); 507*d6685f55SMatthew G. Knepley PetscFunctionReturn(0); 508*d6685f55SMatthew G. Knepley #endif 509*d6685f55SMatthew G. Knepley } 510