1*6af0ca60SMatthew G. Knepley #include <petscconvest.h> /*I "petscconvest.h" I*/ 2*6af0ca60SMatthew G. Knepley #include <petscdmplex.h> 3*6af0ca60SMatthew G. Knepley #include <petscds.h> 4*6af0ca60SMatthew G. Knepley #include <petscblaslapack.h> 5*6af0ca60SMatthew G. Knepley 6*6af0ca60SMatthew G. Knepley #include <petsc/private/petscconvestimpl.h> 7*6af0ca60SMatthew G. Knepley 8*6af0ca60SMatthew G. Knepley static PetscErrorCode zero_private(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 9*6af0ca60SMatthew G. Knepley { 10*6af0ca60SMatthew G. Knepley PetscInt c; 11*6af0ca60SMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = 0.0; 12*6af0ca60SMatthew G. Knepley return 0; 13*6af0ca60SMatthew G. Knepley } 14*6af0ca60SMatthew G. Knepley 15*6af0ca60SMatthew G. Knepley 16*6af0ca60SMatthew G. Knepley /*@ 17*6af0ca60SMatthew G. Knepley PetscConvEstCreate - Create a PetscConvEst object 18*6af0ca60SMatthew G. Knepley 19*6af0ca60SMatthew G. Knepley Collective on MPI_Comm 20*6af0ca60SMatthew G. Knepley 21*6af0ca60SMatthew G. Knepley Input Parameter: 22*6af0ca60SMatthew G. Knepley . comm - The communicator for the PetscConvEst object 23*6af0ca60SMatthew G. Knepley 24*6af0ca60SMatthew G. Knepley Output Parameter: 25*6af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 26*6af0ca60SMatthew G. Knepley 27*6af0ca60SMatthew G. Knepley Level: beginner 28*6af0ca60SMatthew G. Knepley 29*6af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence, create 30*6af0ca60SMatthew G. Knepley .seealso: PetscConvEstDestroy(), PetscConvEstGetConvRate() 31*6af0ca60SMatthew G. Knepley @*/ 32*6af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce) 33*6af0ca60SMatthew G. Knepley { 34*6af0ca60SMatthew G. Knepley PetscErrorCode ierr; 35*6af0ca60SMatthew G. Knepley 36*6af0ca60SMatthew G. Knepley PetscFunctionBegin; 37*6af0ca60SMatthew G. Knepley PetscValidPointer(ce, 2); 38*6af0ca60SMatthew G. Knepley ierr = PetscSysInitializePackage();CHKERRQ(ierr); 39*6af0ca60SMatthew G. Knepley ierr = PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView);CHKERRQ(ierr); 40*6af0ca60SMatthew G. Knepley (*ce)->monitor = PETSC_FALSE; 41*6af0ca60SMatthew G. Knepley (*ce)->Nr = 4; 42*6af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 43*6af0ca60SMatthew G. Knepley } 44*6af0ca60SMatthew G. Knepley 45*6af0ca60SMatthew G. Knepley /*@ 46*6af0ca60SMatthew G. Knepley PetscConvEstDestroy - Destroys a PetscConvEst object 47*6af0ca60SMatthew G. Knepley 48*6af0ca60SMatthew G. Knepley Collective on PetscConvEst 49*6af0ca60SMatthew G. Knepley 50*6af0ca60SMatthew G. Knepley Input Parameter: 51*6af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 52*6af0ca60SMatthew G. Knepley 53*6af0ca60SMatthew G. Knepley Level: beginner 54*6af0ca60SMatthew G. Knepley 55*6af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence, destroy 56*6af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate() 57*6af0ca60SMatthew G. Knepley @*/ 58*6af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstDestroy(PetscConvEst *ce) 59*6af0ca60SMatthew G. Knepley { 60*6af0ca60SMatthew G. Knepley PetscErrorCode ierr; 61*6af0ca60SMatthew G. Knepley 62*6af0ca60SMatthew G. Knepley PetscFunctionBegin; 63*6af0ca60SMatthew G. Knepley if (!*ce) PetscFunctionReturn(0); 64*6af0ca60SMatthew G. Knepley PetscValidHeaderSpecific((*ce),PETSC_OBJECT_CLASSID,1); 65*6af0ca60SMatthew G. Knepley if (--((PetscObject)(*ce))->refct > 0) { 66*6af0ca60SMatthew G. Knepley *ce = NULL; 67*6af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 68*6af0ca60SMatthew G. Knepley } 69*6af0ca60SMatthew G. Knepley ierr = PetscFree2((*ce)->initGuess, (*ce)->exactSol);CHKERRQ(ierr); 70*6af0ca60SMatthew G. Knepley ierr = PetscFree((*ce)->errors);CHKERRQ(ierr); 71*6af0ca60SMatthew G. Knepley ierr = PetscHeaderDestroy(ce);CHKERRQ(ierr); 72*6af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 73*6af0ca60SMatthew G. Knepley } 74*6af0ca60SMatthew G. Knepley 75*6af0ca60SMatthew G. Knepley /*@ 76*6af0ca60SMatthew G. Knepley PetscConvEstSetFromOptions - Sets a PetscConvEst object from options 77*6af0ca60SMatthew G. Knepley 78*6af0ca60SMatthew G. Knepley Collective on PetscConvEst 79*6af0ca60SMatthew G. Knepley 80*6af0ca60SMatthew G. Knepley Input Parameters: 81*6af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 82*6af0ca60SMatthew G. Knepley 83*6af0ca60SMatthew G. Knepley Level: beginner 84*6af0ca60SMatthew G. Knepley 85*6af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence, options 86*6af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate() 87*6af0ca60SMatthew G. Knepley @*/ 88*6af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstSetFromOptions(PetscConvEst ce) 89*6af0ca60SMatthew G. Knepley { 90*6af0ca60SMatthew G. Knepley PetscErrorCode ierr; 91*6af0ca60SMatthew G. Knepley 92*6af0ca60SMatthew G. Knepley PetscFunctionBegin; 93*6af0ca60SMatthew G. Knepley ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) ce), "", "Convergence Estimator Options", "PetscConvEst");CHKERRQ(ierr); 94*6af0ca60SMatthew G. Knepley ierr = PetscOptionsInt("-num_refine", "The number of refinements for the convergence check", "PetscConvEst", ce->Nr, &ce->Nr, NULL);CHKERRQ(ierr); 95*6af0ca60SMatthew G. Knepley ierr = PetscOptionsEnd(); 96*6af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 97*6af0ca60SMatthew G. Knepley } 98*6af0ca60SMatthew G. Knepley 99*6af0ca60SMatthew G. Knepley /*@ 100*6af0ca60SMatthew G. Knepley PetscConvEstView - Views a PetscConvEst object 101*6af0ca60SMatthew G. Knepley 102*6af0ca60SMatthew G. Knepley Collective on PetscConvEst 103*6af0ca60SMatthew G. Knepley 104*6af0ca60SMatthew G. Knepley Input Parameters: 105*6af0ca60SMatthew G. Knepley + ce - The PetscConvEst object 106*6af0ca60SMatthew G. Knepley - viewer - The PetscViewer object 107*6af0ca60SMatthew G. Knepley 108*6af0ca60SMatthew G. Knepley Level: beginner 109*6af0ca60SMatthew G. Knepley 110*6af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence, view 111*6af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate() 112*6af0ca60SMatthew G. Knepley @*/ 113*6af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstView(PetscConvEst ce, PetscViewer viewer) 114*6af0ca60SMatthew G. Knepley { 115*6af0ca60SMatthew G. Knepley PetscErrorCode ierr; 116*6af0ca60SMatthew G. Knepley 117*6af0ca60SMatthew G. Knepley PetscFunctionBegin; 118*6af0ca60SMatthew G. Knepley ierr = PetscObjectPrintClassNamePrefixType((PetscObject) ce, viewer);CHKERRQ(ierr); 119*6af0ca60SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "ConvEst with %D levels\n", ce->Nr+1);CHKERRQ(ierr); 120*6af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 121*6af0ca60SMatthew G. Knepley } 122*6af0ca60SMatthew G. Knepley 123*6af0ca60SMatthew G. Knepley /*@ 124*6af0ca60SMatthew G. Knepley PetscConvEstGetSolver - Gets the solver used to produce discrete solutions 125*6af0ca60SMatthew G. Knepley 126*6af0ca60SMatthew G. Knepley Not collective 127*6af0ca60SMatthew G. Knepley 128*6af0ca60SMatthew G. Knepley Input Parameter: 129*6af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 130*6af0ca60SMatthew G. Knepley 131*6af0ca60SMatthew G. Knepley Output Parameter: 132*6af0ca60SMatthew G. Knepley . snes - The solver 133*6af0ca60SMatthew G. Knepley 134*6af0ca60SMatthew G. Knepley Level: intermediate 135*6af0ca60SMatthew G. Knepley 136*6af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence 137*6af0ca60SMatthew G. Knepley .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate() 138*6af0ca60SMatthew G. Knepley @*/ 139*6af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstGetSolver(PetscConvEst ce, SNES *snes) 140*6af0ca60SMatthew G. Knepley { 141*6af0ca60SMatthew G. Knepley PetscFunctionBegin; 142*6af0ca60SMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 143*6af0ca60SMatthew G. Knepley PetscValidPointer(snes, 2); 144*6af0ca60SMatthew G. Knepley *snes = ce->snes; 145*6af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 146*6af0ca60SMatthew G. Knepley } 147*6af0ca60SMatthew G. Knepley 148*6af0ca60SMatthew G. Knepley /*@ 149*6af0ca60SMatthew G. Knepley PetscConvEstSetSolver - Sets the solver used to produce discrete solutions 150*6af0ca60SMatthew G. Knepley 151*6af0ca60SMatthew G. Knepley Not collective 152*6af0ca60SMatthew G. Knepley 153*6af0ca60SMatthew G. Knepley Input Parameters: 154*6af0ca60SMatthew G. Knepley + ce - The PetscConvEst object 155*6af0ca60SMatthew G. Knepley - snes - The solver 156*6af0ca60SMatthew G. Knepley 157*6af0ca60SMatthew G. Knepley Level: intermediate 158*6af0ca60SMatthew G. Knepley 159*6af0ca60SMatthew G. Knepley Note: The solver MUST have an attached DM/DS, so that we know the exact solution 160*6af0ca60SMatthew G. Knepley 161*6af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence 162*6af0ca60SMatthew G. Knepley .seealso: PetscConvEstGetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate() 163*6af0ca60SMatthew G. Knepley @*/ 164*6af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstSetSolver(PetscConvEst ce, SNES snes) 165*6af0ca60SMatthew G. Knepley { 166*6af0ca60SMatthew G. Knepley PetscErrorCode ierr; 167*6af0ca60SMatthew G. Knepley 168*6af0ca60SMatthew G. Knepley PetscFunctionBegin; 169*6af0ca60SMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 170*6af0ca60SMatthew G. Knepley PetscValidHeaderSpecific(snes, SNES_CLASSID, 2); 171*6af0ca60SMatthew G. Knepley ce->snes = snes; 172*6af0ca60SMatthew G. Knepley ierr = SNESGetDM(ce->snes, &ce->idm);CHKERRQ(ierr); 173*6af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 174*6af0ca60SMatthew G. Knepley } 175*6af0ca60SMatthew G. Knepley 176*6af0ca60SMatthew G. Knepley /*@ 177*6af0ca60SMatthew G. Knepley PetscConvEstSetUp - After the solver is specified, we create structures for estimating convergence 178*6af0ca60SMatthew G. Knepley 179*6af0ca60SMatthew G. Knepley Collective on PetscConvEst 180*6af0ca60SMatthew G. Knepley 181*6af0ca60SMatthew G. Knepley Input Parameters: 182*6af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 183*6af0ca60SMatthew G. Knepley 184*6af0ca60SMatthew G. Knepley Level: beginner 185*6af0ca60SMatthew G. Knepley 186*6af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence, setup 187*6af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate() 188*6af0ca60SMatthew G. Knepley @*/ 189*6af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstSetup(PetscConvEst ce) 190*6af0ca60SMatthew G. Knepley { 191*6af0ca60SMatthew G. Knepley PetscDS prob; 192*6af0ca60SMatthew G. Knepley PetscInt f; 193*6af0ca60SMatthew G. Knepley PetscErrorCode ierr; 194*6af0ca60SMatthew G. Knepley 195*6af0ca60SMatthew G. Knepley PetscFunctionBegin; 196*6af0ca60SMatthew G. Knepley ierr = DMGetDS(ce->idm, &prob);CHKERRQ(ierr); 197*6af0ca60SMatthew G. Knepley ierr = PetscDSGetNumFields(prob, &ce->Nf);CHKERRQ(ierr); 198*6af0ca60SMatthew G. Knepley ierr = PetscMalloc1((ce->Nr+1)*ce->Nf, &ce->errors);CHKERRQ(ierr); 199*6af0ca60SMatthew G. Knepley ierr = PetscMalloc2(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol);CHKERRQ(ierr); 200*6af0ca60SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) ce->initGuess[f] = zero_private; 201*6af0ca60SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 202*6af0ca60SMatthew G. Knepley ierr = PetscDSGetExactSolution(prob, f, &ce->exactSol[f]);CHKERRQ(ierr); 203*6af0ca60SMatthew G. Knepley if (!ce->exactSol[f]) SETERRQ1(PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "DS must contain exact solution functions in order to estimate convergence, missing for field %D", f); 204*6af0ca60SMatthew G. Knepley } 205*6af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 206*6af0ca60SMatthew G. Knepley } 207*6af0ca60SMatthew G. Knepley 208*6af0ca60SMatthew G. Knepley static PetscErrorCode PetscConvEstLinearRegression_Private(PetscConvEst ce, PetscInt n, const PetscReal x[], const PetscReal y[], PetscReal *slope, PetscReal *intercept) 209*6af0ca60SMatthew G. Knepley { 210*6af0ca60SMatthew G. Knepley PetscReal *X, *Y, H[4], beta[2]; 211*6af0ca60SMatthew G. Knepley PetscInt i, j, k; 212*6af0ca60SMatthew G. Knepley PetscErrorCode ierr; 213*6af0ca60SMatthew G. Knepley 214*6af0ca60SMatthew G. Knepley PetscFunctionBegin; 215*6af0ca60SMatthew G. Knepley ierr = PetscMalloc2(n*2, &X, n*2, &Y);CHKERRQ(ierr); 216*6af0ca60SMatthew G. Knepley for (k = 0; k < n; ++k) { 217*6af0ca60SMatthew G. Knepley /* X[n,2] = [1, x] */ 218*6af0ca60SMatthew G. Knepley X[k*2+0] = 1.0; 219*6af0ca60SMatthew G. Knepley X[k*2+1] = x[k]; 220*6af0ca60SMatthew G. Knepley } 221*6af0ca60SMatthew G. Knepley /* H = X^T X */ 222*6af0ca60SMatthew G. Knepley for (i = 0; i < 2; ++i) { 223*6af0ca60SMatthew G. Knepley for (j = 0; j < 2; ++j) { 224*6af0ca60SMatthew G. Knepley H[i*2+j] = 0.0; 225*6af0ca60SMatthew G. Knepley for (k = 0; k < n; ++k) { 226*6af0ca60SMatthew G. Knepley H[i*2+j] += X[k*2+i] * X[k*2+j]; 227*6af0ca60SMatthew G. Knepley } 228*6af0ca60SMatthew G. Knepley } 229*6af0ca60SMatthew G. Knepley } 230*6af0ca60SMatthew G. Knepley /* H = (X^T X)^{-1} */ 231*6af0ca60SMatthew G. Knepley { 232*6af0ca60SMatthew G. Knepley PetscBLASInt two = 2, ipiv[2], info; 233*6af0ca60SMatthew G. Knepley PetscReal work[2]; 234*6af0ca60SMatthew G. Knepley 235*6af0ca60SMatthew G. Knepley ierr = PetscFPTrapPush(PETSC_FP_TRAP_OFF);CHKERRQ(ierr); 236*6af0ca60SMatthew G. Knepley PetscStackCallBLAS("LAPACKgetrf", LAPACKgetrf_(&two, &two, H, &two, ipiv, &info)); 237*6af0ca60SMatthew G. Knepley PetscStackCallBLAS("LAPACKgetri", LAPACKgetri_(&two, H, &two, ipiv, work, &two, &info)); 238*6af0ca60SMatthew G. Knepley ierr = PetscFPTrapPop();CHKERRQ(ierr); 239*6af0ca60SMatthew G. Knepley } 240*6af0ca60SMatthew G. Knepley /* Y = H X^T */ 241*6af0ca60SMatthew G. Knepley for (i = 0; i < 2; ++i) { 242*6af0ca60SMatthew G. Knepley for (k = 0; k < n; ++k) { 243*6af0ca60SMatthew G. Knepley Y[i*n+k] = 0.0; 244*6af0ca60SMatthew G. Knepley for (j = 0; j < 2; ++j) { 245*6af0ca60SMatthew G. Knepley Y[i*n+k] += H[i*2+j] * X[k*2+j]; 246*6af0ca60SMatthew G. Knepley } 247*6af0ca60SMatthew G. Knepley } 248*6af0ca60SMatthew G. Knepley } 249*6af0ca60SMatthew G. Knepley /* beta = Y error = [y-intercept, slope] */ 250*6af0ca60SMatthew G. Knepley for (i = 0; i < 2; ++i) { 251*6af0ca60SMatthew G. Knepley beta[i] = 0.0; 252*6af0ca60SMatthew G. Knepley for (k = 0; k < n; ++k) { 253*6af0ca60SMatthew G. Knepley beta[i] += Y[i*n+k] * y[k]; 254*6af0ca60SMatthew G. Knepley } 255*6af0ca60SMatthew G. Knepley } 256*6af0ca60SMatthew G. Knepley ierr = PetscFree2(X, Y);CHKERRQ(ierr); 257*6af0ca60SMatthew G. Knepley *intercept = beta[0]; 258*6af0ca60SMatthew G. Knepley *slope = beta[1]; 259*6af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 260*6af0ca60SMatthew G. Knepley } 261*6af0ca60SMatthew G. Knepley 262*6af0ca60SMatthew G. Knepley /*@ 263*6af0ca60SMatthew G. Knepley PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization 264*6af0ca60SMatthew G. Knepley 265*6af0ca60SMatthew G. Knepley Not collective 266*6af0ca60SMatthew G. Knepley 267*6af0ca60SMatthew G. Knepley Input Parameter: 268*6af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 269*6af0ca60SMatthew G. Knepley 270*6af0ca60SMatthew G. Knepley Output Parameter: 271*6af0ca60SMatthew G. Knepley . alpha - The convergence rate 272*6af0ca60SMatthew G. Knepley 273*6af0ca60SMatthew G. Knepley Note: The convergence rate alpha is defined by 274*6af0ca60SMatthew G. Knepley $ || u_h - u_exact || < C h^alpha 275*6af0ca60SMatthew G. Knepley where u_h is the discrete solution, and h is a measure of the discretization size. 276*6af0ca60SMatthew G. Knepley 277*6af0ca60SMatthew G. Knepley We solve a series of problems on refined meshes, calculate an error based upon the exact solution in the DS, 278*6af0ca60SMatthew G. Knepley and then fit the result to our model above using linear regression. 279*6af0ca60SMatthew G. Knepley 280*6af0ca60SMatthew G. Knepley Options database keys: 281*6af0ca60SMatthew G. Knepley . -snes_convergence_estimate : Execute convergence estimation and print out the rate 282*6af0ca60SMatthew G. Knepley 283*6af0ca60SMatthew G. Knepley Level: intermediate 284*6af0ca60SMatthew G. Knepley 285*6af0ca60SMatthew G. Knepley .keywords: PetscConvEst, convergence 286*6af0ca60SMatthew G. Knepley .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate() 287*6af0ca60SMatthew G. Knepley @*/ 288*6af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal *alpha) 289*6af0ca60SMatthew G. Knepley { 290*6af0ca60SMatthew G. Knepley DM *dm; 291*6af0ca60SMatthew G. Knepley PetscDS prob; 292*6af0ca60SMatthew G. Knepley PetscObject disc; 293*6af0ca60SMatthew G. Knepley MPI_Comm comm; 294*6af0ca60SMatthew G. Knepley const char *uname, *dmname; 295*6af0ca60SMatthew G. Knepley void *ctx; 296*6af0ca60SMatthew G. Knepley Vec u; 297*6af0ca60SMatthew G. Knepley PetscReal t = 0.0, *x, *y, slope, intercept; 298*6af0ca60SMatthew G. Knepley PetscInt *dof, dim, Nr = ce->Nr, r; 299*6af0ca60SMatthew G. Knepley PetscErrorCode ierr; 300*6af0ca60SMatthew G. Knepley 301*6af0ca60SMatthew G. Knepley PetscFunctionBegin; 302*6af0ca60SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) ce, &comm);CHKERRQ(ierr); 303*6af0ca60SMatthew G. Knepley ierr = DMGetDimension(ce->idm, &dim);CHKERRQ(ierr); 304*6af0ca60SMatthew G. Knepley ierr = DMGetApplicationContext(ce->idm, &ctx);CHKERRQ(ierr); 305*6af0ca60SMatthew G. Knepley ierr = DMGetDS(ce->idm, &prob);CHKERRQ(ierr); 306*6af0ca60SMatthew G. Knepley ierr = DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE);CHKERRQ(ierr); 307*6af0ca60SMatthew G. Knepley ierr = PetscDSGetDiscretization(prob, 0, &disc);CHKERRQ(ierr); 308*6af0ca60SMatthew G. Knepley ierr = PetscObjectGetName(disc, &uname);CHKERRQ(ierr); 309*6af0ca60SMatthew G. Knepley ierr = PetscMalloc2((Nr+1), &dm, (Nr+1), &dof);CHKERRQ(ierr); 310*6af0ca60SMatthew G. Knepley dm[0] = ce->idm; 311*6af0ca60SMatthew G. Knepley *alpha = 0.0; 312*6af0ca60SMatthew G. Knepley /* Loop over meshes */ 313*6af0ca60SMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 314*6af0ca60SMatthew G. Knepley if (r > 0) { 315*6af0ca60SMatthew G. Knepley ierr = DMRefine(dm[r-1], MPI_COMM_NULL, &dm[r]);CHKERRQ(ierr); 316*6af0ca60SMatthew G. Knepley ierr = DMSetCoarseDM(dm[r], dm[r-1]);CHKERRQ(ierr); 317*6af0ca60SMatthew G. Knepley ierr = DMSetDS(dm[r], prob);CHKERRQ(ierr); 318*6af0ca60SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) dm[r-1], &dmname);CHKERRQ(ierr); 319*6af0ca60SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) dm[r], dmname);CHKERRQ(ierr); 320*6af0ca60SMatthew G. Knepley } 321*6af0ca60SMatthew G. Knepley ierr = DMViewFromOptions(dm[r], NULL, "-conv_dm_view");CHKERRQ(ierr); 322*6af0ca60SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm[r], 0, NULL, &dof[r]);CHKERRQ(ierr); 323*6af0ca60SMatthew G. Knepley /* Create solution */ 324*6af0ca60SMatthew G. Knepley ierr = DMCreateGlobalVector(dm[r], &u);CHKERRQ(ierr); 325*6af0ca60SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) u, uname);CHKERRQ(ierr); 326*6af0ca60SMatthew G. Knepley /* Setup solver */ 327*6af0ca60SMatthew G. Knepley ierr = SNESReset(ce->snes);CHKERRQ(ierr); 328*6af0ca60SMatthew G. Knepley ierr = SNESSetDM(ce->snes, dm[r]);CHKERRQ(ierr); 329*6af0ca60SMatthew G. Knepley ierr = DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx);CHKERRQ(ierr); 330*6af0ca60SMatthew G. Knepley ierr = SNESSetFromOptions(ce->snes);CHKERRQ(ierr); 331*6af0ca60SMatthew G. Knepley /* Create initial guess */ 332*6af0ca60SMatthew G. Knepley ierr = DMProjectFunction(dm[r], t, ce->initGuess, NULL, INSERT_VALUES, u);CHKERRQ(ierr); 333*6af0ca60SMatthew G. Knepley ierr = SNESSolve(ce->snes, NULL, u);CHKERRQ(ierr); 334*6af0ca60SMatthew G. Knepley ierr = DMComputeL2FieldDiff(dm[r], t, ce->exactSol, NULL, u, &ce->errors[r*ce->Nf]);CHKERRQ(ierr); 335*6af0ca60SMatthew G. Knepley /* Monitor */ 336*6af0ca60SMatthew G. Knepley if (ce->monitor) { 337*6af0ca60SMatthew G. Knepley PetscReal *errors = &ce->errors[r*ce->Nf]; 338*6af0ca60SMatthew G. Knepley PetscInt f; 339*6af0ca60SMatthew G. Knepley 340*6af0ca60SMatthew G. Knepley ierr = PetscPrintf(comm, "L_2 Error: [");CHKERRQ(ierr); 341*6af0ca60SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 342*6af0ca60SMatthew G. Knepley if (f > 0) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);} 343*6af0ca60SMatthew G. Knepley if (errors[f] < 1.0e-11) {ierr = PetscPrintf(comm, "< 1e-11");CHKERRQ(ierr);} 344*6af0ca60SMatthew G. Knepley else {ierr = PetscPrintf(comm, "%g", errors[f]);CHKERRQ(ierr);} 345*6af0ca60SMatthew G. Knepley } 346*6af0ca60SMatthew G. Knepley ierr = PetscPrintf(comm, "]\n");CHKERRQ(ierr); 347*6af0ca60SMatthew G. Knepley } 348*6af0ca60SMatthew G. Knepley /* Cleanup */ 349*6af0ca60SMatthew G. Knepley ierr = VecDestroy(&u);CHKERRQ(ierr); 350*6af0ca60SMatthew G. Knepley } 351*6af0ca60SMatthew G. Knepley for (r = 1; r <= Nr; ++r) { 352*6af0ca60SMatthew G. Knepley ierr = DMDestroy(&dm[r]);CHKERRQ(ierr); 353*6af0ca60SMatthew G. Knepley } 354*6af0ca60SMatthew G. Knepley /* Fit convergence rate */ 355*6af0ca60SMatthew G. Knepley ierr = PetscMalloc2(Nr+1, &x, Nr+1, &y);CHKERRQ(ierr); 356*6af0ca60SMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 357*6af0ca60SMatthew G. Knepley x[r] = PetscLog10Real(dof[r]); 358*6af0ca60SMatthew G. Knepley y[r] = PetscLog10Real(ce->errors[r*ce->Nf+0]); 359*6af0ca60SMatthew G. Knepley } 360*6af0ca60SMatthew G. Knepley ierr = PetscConvEstLinearRegression_Private(ce, Nr+1, x, y, &slope, &intercept);CHKERRQ(ierr); 361*6af0ca60SMatthew G. Knepley ierr = PetscFree2(x, y);CHKERRQ(ierr); 362*6af0ca60SMatthew G. Knepley /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */ 363*6af0ca60SMatthew G. Knepley *alpha = -slope * dim; 364*6af0ca60SMatthew G. Knepley ierr = PetscFree2(dm, dof);CHKERRQ(ierr); 365*6af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 366*6af0ca60SMatthew G. Knepley } 367*6af0ca60SMatthew G. Knepley 368*6af0ca60SMatthew G. Knepley /*@ 369*6af0ca60SMatthew G. Knepley PetscConvEstRateView - Displays the convergence rate to a viewer 370*6af0ca60SMatthew G. Knepley 371*6af0ca60SMatthew G. Knepley Collective on SNES 372*6af0ca60SMatthew G. Knepley 373*6af0ca60SMatthew G. Knepley Parameter: 374*6af0ca60SMatthew G. Knepley + snes - iterative context obtained from SNESCreate() 375*6af0ca60SMatthew G. Knepley . alpha - the convergence rate 376*6af0ca60SMatthew G. Knepley - viewer - the viewer to display the reason 377*6af0ca60SMatthew G. Knepley 378*6af0ca60SMatthew G. Knepley Options Database Keys: 379*6af0ca60SMatthew G. Knepley . -snes_convergence_estimate - print the convergence rate 380*6af0ca60SMatthew G. Knepley 381*6af0ca60SMatthew G. Knepley Level: developer 382*6af0ca60SMatthew G. Knepley 383*6af0ca60SMatthew G. Knepley .seealso: PetscConvEstGetRate() 384*6af0ca60SMatthew G. Knepley @*/ 385*6af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstRateView(PetscConvEst ce, PetscReal alpha, PetscViewer viewer) 386*6af0ca60SMatthew G. Knepley { 387*6af0ca60SMatthew G. Knepley PetscBool isAscii; 388*6af0ca60SMatthew G. Knepley PetscErrorCode ierr; 389*6af0ca60SMatthew G. Knepley 390*6af0ca60SMatthew G. Knepley PetscFunctionBegin; 391*6af0ca60SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr); 392*6af0ca60SMatthew G. Knepley if (isAscii) { 393*6af0ca60SMatthew G. Knepley ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr); 394*6af0ca60SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: %0.2g\n", (double) alpha);CHKERRQ(ierr); 395*6af0ca60SMatthew G. Knepley ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr); 396*6af0ca60SMatthew G. Knepley } 397*6af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 398*6af0ca60SMatthew G. Knepley } 399