16af0ca60SMatthew G. Knepley #include <petscconvest.h> /*I "petscconvest.h" I*/ 26af0ca60SMatthew G. Knepley #include <petscdmplex.h> 36af0ca60SMatthew G. Knepley #include <petscds.h> 46af0ca60SMatthew G. Knepley 56af0ca60SMatthew G. Knepley #include <petsc/private/petscconvestimpl.h> 66af0ca60SMatthew G. Knepley 76af0ca60SMatthew G. Knepley static PetscErrorCode zero_private(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 86af0ca60SMatthew G. Knepley { 96af0ca60SMatthew G. Knepley PetscInt c; 106af0ca60SMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = 0.0; 116af0ca60SMatthew G. Knepley return 0; 126af0ca60SMatthew G. Knepley } 136af0ca60SMatthew G. Knepley 146af0ca60SMatthew G. Knepley /*@ 156af0ca60SMatthew G. Knepley PetscConvEstDestroy - Destroys a PetscConvEst object 166af0ca60SMatthew G. Knepley 176af0ca60SMatthew G. Knepley Collective on PetscConvEst 186af0ca60SMatthew G. Knepley 196af0ca60SMatthew G. Knepley Input Parameter: 206af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 216af0ca60SMatthew G. Knepley 226af0ca60SMatthew G. Knepley Level: beginner 236af0ca60SMatthew G. Knepley 246af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate() 256af0ca60SMatthew G. Knepley @*/ 266af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstDestroy(PetscConvEst *ce) 276af0ca60SMatthew G. Knepley { 286af0ca60SMatthew G. Knepley PetscErrorCode ierr; 296af0ca60SMatthew G. Knepley 306af0ca60SMatthew G. Knepley PetscFunctionBegin; 316af0ca60SMatthew G. Knepley if (!*ce) PetscFunctionReturn(0); 326af0ca60SMatthew G. Knepley PetscValidHeaderSpecific((*ce),PETSC_OBJECT_CLASSID,1); 336af0ca60SMatthew G. Knepley if (--((PetscObject)(*ce))->refct > 0) { 346af0ca60SMatthew G. Knepley *ce = NULL; 356af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 366af0ca60SMatthew G. Knepley } 3795cbbfd3SMatthew G. Knepley ierr = PetscFree3((*ce)->initGuess, (*ce)->exactSol, (*ce)->ctxs);CHKERRQ(ierr); 387809adefSMatthew G. Knepley ierr = PetscFree2((*ce)->dofs, (*ce)->errors);CHKERRQ(ierr); 396af0ca60SMatthew G. Knepley ierr = PetscHeaderDestroy(ce);CHKERRQ(ierr); 406af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 416af0ca60SMatthew G. Knepley } 426af0ca60SMatthew G. Knepley 436af0ca60SMatthew G. Knepley /*@ 446af0ca60SMatthew G. Knepley PetscConvEstSetFromOptions - Sets a PetscConvEst object from options 456af0ca60SMatthew G. Knepley 466af0ca60SMatthew G. Knepley Collective on PetscConvEst 476af0ca60SMatthew G. Knepley 486af0ca60SMatthew G. Knepley Input Parameters: 496af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 506af0ca60SMatthew G. Knepley 516af0ca60SMatthew G. Knepley Level: beginner 526af0ca60SMatthew G. Knepley 536af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate() 546af0ca60SMatthew G. Knepley @*/ 556af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstSetFromOptions(PetscConvEst ce) 566af0ca60SMatthew G. Knepley { 576af0ca60SMatthew G. Knepley PetscErrorCode ierr; 586af0ca60SMatthew G. Knepley 596af0ca60SMatthew G. Knepley PetscFunctionBegin; 606af0ca60SMatthew G. Knepley ierr = PetscOptionsBegin(PetscObjectComm((PetscObject) ce), "", "Convergence Estimator Options", "PetscConvEst");CHKERRQ(ierr); 6146079b62SMatthew G. Knepley ierr = PetscOptionsInt("-convest_num_refine", "The number of refinements for the convergence check", "PetscConvEst", ce->Nr, &ce->Nr, NULL);CHKERRQ(ierr); 622e61be88SMatthew G. Knepley ierr = PetscOptionsReal("-convest_refine_factor", "The increase in resolution in each dimension", "PetscConvEst", ce->r, &ce->r, NULL);CHKERRQ(ierr); 6346079b62SMatthew G. Knepley ierr = PetscOptionsBool("-convest_monitor", "Monitor the error for each convergence check", "PetscConvEst", ce->monitor, &ce->monitor, NULL);CHKERRQ(ierr); 64b2df8587SMatthew G. Knepley ierr = PetscOptionsBool("-convest_no_refine", "Debugging flag to run on the same mesh each time", "PetscConvEst", ce->noRefine, &ce->noRefine, NULL);CHKERRQ(ierr); 651e1ea65dSPierre Jolivet ierr = PetscOptionsEnd();CHKERRQ(ierr); 666af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 676af0ca60SMatthew G. Knepley } 686af0ca60SMatthew G. Knepley 696af0ca60SMatthew G. Knepley /*@ 706af0ca60SMatthew G. Knepley PetscConvEstView - Views a PetscConvEst object 716af0ca60SMatthew G. Knepley 726af0ca60SMatthew G. Knepley Collective on PetscConvEst 736af0ca60SMatthew G. Knepley 746af0ca60SMatthew G. Knepley Input Parameters: 756af0ca60SMatthew G. Knepley + ce - The PetscConvEst object 766af0ca60SMatthew G. Knepley - viewer - The PetscViewer object 776af0ca60SMatthew G. Knepley 786af0ca60SMatthew G. Knepley Level: beginner 796af0ca60SMatthew G. Knepley 806af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate() 816af0ca60SMatthew G. Knepley @*/ 826af0ca60SMatthew G. Knepley PetscErrorCode PetscConvEstView(PetscConvEst ce, PetscViewer viewer) 836af0ca60SMatthew G. Knepley { 846af0ca60SMatthew G. Knepley PetscErrorCode ierr; 856af0ca60SMatthew G. Knepley 866af0ca60SMatthew G. Knepley PetscFunctionBegin; 876af0ca60SMatthew G. Knepley ierr = PetscObjectPrintClassNamePrefixType((PetscObject) ce, viewer);CHKERRQ(ierr); 886af0ca60SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "ConvEst with %D levels\n", ce->Nr+1);CHKERRQ(ierr); 896af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 906af0ca60SMatthew G. Knepley } 916af0ca60SMatthew G. Knepley 926af0ca60SMatthew G. Knepley /*@ 936af0ca60SMatthew G. Knepley PetscConvEstGetSolver - Gets the solver used to produce discrete solutions 946af0ca60SMatthew G. Knepley 956af0ca60SMatthew G. Knepley Not collective 966af0ca60SMatthew G. Knepley 976af0ca60SMatthew G. Knepley Input Parameter: 986af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 996af0ca60SMatthew G. Knepley 1006af0ca60SMatthew G. Knepley Output Parameter: 101900f6b5bSMatthew G. Knepley . solver - The solver 1026af0ca60SMatthew G. Knepley 1036af0ca60SMatthew G. Knepley Level: intermediate 1046af0ca60SMatthew G. Knepley 1056af0ca60SMatthew G. Knepley .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate() 1066af0ca60SMatthew G. Knepley @*/ 107900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstGetSolver(PetscConvEst ce, PetscObject *solver) 1086af0ca60SMatthew G. Knepley { 1096af0ca60SMatthew G. Knepley PetscFunctionBegin; 1106af0ca60SMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 111900f6b5bSMatthew G. Knepley PetscValidPointer(solver, 2); 112900f6b5bSMatthew G. Knepley *solver = ce->solver; 1136af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 1146af0ca60SMatthew G. Knepley } 1156af0ca60SMatthew G. Knepley 1166af0ca60SMatthew G. Knepley /*@ 1176af0ca60SMatthew G. Knepley PetscConvEstSetSolver - Sets the solver used to produce discrete solutions 1186af0ca60SMatthew G. Knepley 1196af0ca60SMatthew G. Knepley Not collective 1206af0ca60SMatthew G. Knepley 1216af0ca60SMatthew G. Knepley Input Parameters: 1226af0ca60SMatthew G. Knepley + ce - The PetscConvEst object 123900f6b5bSMatthew G. Knepley - solver - The solver 1246af0ca60SMatthew G. Knepley 1256af0ca60SMatthew G. Knepley Level: intermediate 1266af0ca60SMatthew G. Knepley 1276af0ca60SMatthew G. Knepley Note: The solver MUST have an attached DM/DS, so that we know the exact solution 1286af0ca60SMatthew G. Knepley 129900f6b5bSMatthew G. Knepley .seealso: PetscConvEstGetSNES(), PetscConvEstCreate(), PetscConvEstGetConvRate() 1306af0ca60SMatthew G. Knepley @*/ 131900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstSetSolver(PetscConvEst ce, PetscObject solver) 1326af0ca60SMatthew G. Knepley { 1336af0ca60SMatthew G. Knepley PetscErrorCode ierr; 1346af0ca60SMatthew G. Knepley 1356af0ca60SMatthew G. Knepley PetscFunctionBegin; 1366af0ca60SMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 137900f6b5bSMatthew G. Knepley PetscValidHeader(solver, 2); 138900f6b5bSMatthew G. Knepley ce->solver = solver; 139900f6b5bSMatthew G. Knepley ierr = (*ce->ops->setsolver)(ce, solver);CHKERRQ(ierr); 1406af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 1416af0ca60SMatthew G. Knepley } 1426af0ca60SMatthew G. Knepley 1436af0ca60SMatthew G. Knepley /*@ 1446af0ca60SMatthew G. Knepley PetscConvEstSetUp - After the solver is specified, we create structures for estimating convergence 1456af0ca60SMatthew G. Knepley 1466af0ca60SMatthew G. Knepley Collective on PetscConvEst 1476af0ca60SMatthew G. Knepley 1486af0ca60SMatthew G. Knepley Input Parameters: 1496af0ca60SMatthew G. Knepley . ce - The PetscConvEst object 1506af0ca60SMatthew G. Knepley 1516af0ca60SMatthew G. Knepley Level: beginner 1526af0ca60SMatthew G. Knepley 1536af0ca60SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate() 1546af0ca60SMatthew G. Knepley @*/ 1550955ed61SMatthew G. Knepley PetscErrorCode PetscConvEstSetUp(PetscConvEst ce) 1566af0ca60SMatthew G. Knepley { 157083401c6SMatthew G. Knepley PetscInt Nf, f, Nds, s; 1586af0ca60SMatthew G. Knepley PetscErrorCode ierr; 1596af0ca60SMatthew G. Knepley 1606af0ca60SMatthew G. Knepley PetscFunctionBegin; 161083401c6SMatthew G. Knepley ierr = DMGetNumFields(ce->idm, &Nf);CHKERRQ(ierr); 162900f6b5bSMatthew G. Knepley ce->Nf = PetscMax(Nf, 1); 1637809adefSMatthew G. Knepley ierr = PetscMalloc2((ce->Nr+1)*ce->Nf, &ce->dofs, (ce->Nr+1)*ce->Nf, &ce->errors);CHKERRQ(ierr); 164083401c6SMatthew G. Knepley ierr = PetscCalloc3(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol, ce->Nf, &ce->ctxs);CHKERRQ(ierr); 165900f6b5bSMatthew G. Knepley for (f = 0; f < Nf; ++f) ce->initGuess[f] = zero_private; 166083401c6SMatthew G. Knepley ierr = DMGetNumDS(ce->idm, &Nds);CHKERRQ(ierr); 167083401c6SMatthew G. Knepley for (s = 0; s < Nds; ++s) { 168083401c6SMatthew G. Knepley PetscDS ds; 169083401c6SMatthew G. Knepley DMLabel label; 170083401c6SMatthew G. Knepley IS fieldIS; 171083401c6SMatthew G. Knepley const PetscInt *fields; 172083401c6SMatthew G. Knepley PetscInt dsNf; 173083401c6SMatthew G. Knepley 174083401c6SMatthew G. Knepley ierr = DMGetRegionNumDS(ce->idm, s, &label, &fieldIS, &ds);CHKERRQ(ierr); 175083401c6SMatthew G. Knepley ierr = PetscDSGetNumFields(ds, &dsNf);CHKERRQ(ierr); 176083401c6SMatthew G. Knepley if (fieldIS) {ierr = ISGetIndices(fieldIS, &fields);CHKERRQ(ierr);} 177083401c6SMatthew G. Knepley for (f = 0; f < dsNf; ++f) { 178083401c6SMatthew G. Knepley const PetscInt field = fields[f]; 179083401c6SMatthew G. Knepley ierr = PetscDSGetExactSolution(ds, field, &ce->exactSol[field], &ce->ctxs[field]);CHKERRQ(ierr); 180083401c6SMatthew G. Knepley } 181083401c6SMatthew G. Knepley if (fieldIS) {ierr = ISRestoreIndices(fieldIS, &fields);CHKERRQ(ierr);} 182083401c6SMatthew G. Knepley } 183900f6b5bSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 1846af0ca60SMatthew 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); 1856af0ca60SMatthew G. Knepley } 1866af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 1876af0ca60SMatthew G. Knepley } 1886af0ca60SMatthew G. Knepley 189900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstComputeInitialGuess(PetscConvEst ce, PetscInt r, DM dm, Vec u) 1906af0ca60SMatthew G. Knepley { 1916af0ca60SMatthew G. Knepley PetscErrorCode ierr; 1926af0ca60SMatthew G. Knepley 1936af0ca60SMatthew G. Knepley PetscFunctionBegin; 194900f6b5bSMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 195900f6b5bSMatthew G. Knepley if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3); 196900f6b5bSMatthew G. Knepley PetscValidHeaderSpecific(u, VEC_CLASSID, 4); 197900f6b5bSMatthew G. Knepley ierr = (*ce->ops->initguess)(ce, r, dm, u);CHKERRQ(ierr); 198900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 199900f6b5bSMatthew G. Knepley } 200900f6b5bSMatthew G. Knepley 201900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstComputeError(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[]) 202900f6b5bSMatthew G. Knepley { 203900f6b5bSMatthew G. Knepley PetscErrorCode ierr; 204900f6b5bSMatthew G. Knepley 205900f6b5bSMatthew G. Knepley PetscFunctionBegin; 206900f6b5bSMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 207900f6b5bSMatthew G. Knepley if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3); 208900f6b5bSMatthew G. Knepley PetscValidHeaderSpecific(u, VEC_CLASSID, 4); 209900f6b5bSMatthew G. Knepley PetscValidRealPointer(errors, 5); 210900f6b5bSMatthew G. Knepley ierr = (*ce->ops->computeerror)(ce, r, dm, u, errors);CHKERRQ(ierr); 211900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 212900f6b5bSMatthew G. Knepley } 213900f6b5bSMatthew G. Knepley 214f2cacb80SMatthew G. Knepley /*@ 215f2cacb80SMatthew G. Knepley PetscConvEstMonitorDefault - Monitors the convergence estimation loop 216f2cacb80SMatthew G. Knepley 217f2cacb80SMatthew G. Knepley Collective on PetscConvEst 218f2cacb80SMatthew G. Knepley 219*d8d19677SJose E. Roman Input Parameters: 220f2cacb80SMatthew G. Knepley + ce - The PetscConvEst object 221f2cacb80SMatthew G. Knepley - r - The refinement level 222f2cacb80SMatthew G. Knepley 223f2cacb80SMatthew G. Knepley Options database keys: 224ee300463SSatish Balay . -convest_monitor - Activate the monitor 225f2cacb80SMatthew G. Knepley 226f2cacb80SMatthew G. Knepley Level: intermediate 227f2cacb80SMatthew G. Knepley 228f2cacb80SMatthew G. Knepley .seealso: PetscConvEstCreate(), PetscConvEstGetConvRate(), SNESSolve(), TSSolve() 229f2cacb80SMatthew G. Knepley @*/ 230f2cacb80SMatthew G. Knepley PetscErrorCode PetscConvEstMonitorDefault(PetscConvEst ce, PetscInt r) 231900f6b5bSMatthew G. Knepley { 232900f6b5bSMatthew G. Knepley MPI_Comm comm; 233900f6b5bSMatthew G. Knepley PetscInt f; 234900f6b5bSMatthew G. Knepley PetscErrorCode ierr; 235900f6b5bSMatthew G. Knepley 236900f6b5bSMatthew G. Knepley PetscFunctionBegin; 237900f6b5bSMatthew G. Knepley if (ce->monitor) { 2387809adefSMatthew G. Knepley PetscInt *dofs = &ce->dofs[r*ce->Nf]; 239900f6b5bSMatthew G. Knepley PetscReal *errors = &ce->errors[r*ce->Nf]; 240900f6b5bSMatthew G. Knepley 2416af0ca60SMatthew G. Knepley ierr = PetscObjectGetComm((PetscObject) ce, &comm);CHKERRQ(ierr); 2427809adefSMatthew G. Knepley ierr = PetscPrintf(comm, "N: ");CHKERRQ(ierr); 2437809adefSMatthew G. Knepley if (ce->Nf > 1) {ierr = PetscPrintf(comm, "[");CHKERRQ(ierr);} 2447809adefSMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 2457809adefSMatthew G. Knepley if (f > 0) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);} 2467809adefSMatthew G. Knepley ierr = PetscPrintf(comm, "%7D", dofs[f]);CHKERRQ(ierr); 2477809adefSMatthew G. Knepley } 2487809adefSMatthew G. Knepley if (ce->Nf > 1) {ierr = PetscPrintf(comm, "]");CHKERRQ(ierr);} 2497809adefSMatthew G. Knepley ierr = PetscPrintf(comm, " ");CHKERRQ(ierr); 250900f6b5bSMatthew G. Knepley ierr = PetscPrintf(comm, "L_2 Error: ");CHKERRQ(ierr); 251900f6b5bSMatthew G. Knepley if (ce->Nf > 1) {ierr = PetscPrintf(comm, "[");CHKERRQ(ierr);} 252900f6b5bSMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 253900f6b5bSMatthew G. Knepley if (f > 0) {ierr = PetscPrintf(comm, ", ");CHKERRQ(ierr);} 254900f6b5bSMatthew G. Knepley if (errors[f] < 1.0e-11) {ierr = PetscPrintf(comm, "< 1e-11");CHKERRQ(ierr);} 255900f6b5bSMatthew G. Knepley else {ierr = PetscPrintf(comm, "%g", (double) errors[f]);CHKERRQ(ierr);} 256900f6b5bSMatthew G. Knepley } 257900f6b5bSMatthew G. Knepley if (ce->Nf > 1) {ierr = PetscPrintf(comm, "]");CHKERRQ(ierr);} 258900f6b5bSMatthew G. Knepley ierr = PetscPrintf(comm, "\n");CHKERRQ(ierr); 259900f6b5bSMatthew G. Knepley } 260900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 261900f6b5bSMatthew G. Knepley } 262900f6b5bSMatthew G. Knepley 263900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstSetSNES_Private(PetscConvEst ce, PetscObject solver) 264900f6b5bSMatthew G. Knepley { 265900f6b5bSMatthew G. Knepley PetscClassId id; 266900f6b5bSMatthew G. Knepley PetscErrorCode ierr; 267900f6b5bSMatthew G. Knepley 268900f6b5bSMatthew G. Knepley PetscFunctionBegin; 269900f6b5bSMatthew G. Knepley ierr = PetscObjectGetClassId(ce->solver, &id);CHKERRQ(ierr); 270900f6b5bSMatthew G. Knepley if (id != SNES_CLASSID) SETERRQ(PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "Solver was not a SNES"); 271900f6b5bSMatthew G. Knepley ierr = SNESGetDM((SNES) ce->solver, &ce->idm);CHKERRQ(ierr); 272900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 273900f6b5bSMatthew G. Knepley } 274900f6b5bSMatthew G. Knepley 275900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstInitGuessSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u) 276900f6b5bSMatthew G. Knepley { 277900f6b5bSMatthew G. Knepley PetscErrorCode ierr; 278900f6b5bSMatthew G. Knepley 279900f6b5bSMatthew G. Knepley PetscFunctionBegin; 280900f6b5bSMatthew G. Knepley ierr = DMProjectFunction(dm, 0.0, ce->initGuess, ce->ctxs, INSERT_VALUES, u);CHKERRQ(ierr); 281900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 282900f6b5bSMatthew G. Knepley } 283900f6b5bSMatthew G. Knepley 284900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstComputeErrorSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[]) 285900f6b5bSMatthew G. Knepley { 286900f6b5bSMatthew G. Knepley PetscErrorCode ierr; 287900f6b5bSMatthew G. Knepley 288900f6b5bSMatthew G. Knepley PetscFunctionBegin; 289900f6b5bSMatthew G. Knepley ierr = DMComputeL2FieldDiff(dm, 0.0, ce->exactSol, ce->ctxs, u, errors);CHKERRQ(ierr); 290900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 291900f6b5bSMatthew G. Knepley } 292900f6b5bSMatthew G. Knepley 293478db826SMatthew G. Knepley static PetscErrorCode PetscConvEstSetJacobianNullspace_Private(PetscConvEst ce, SNES snes) 294478db826SMatthew G. Knepley { 295478db826SMatthew G. Knepley DM dm; 296478db826SMatthew G. Knepley PetscInt f; 297478db826SMatthew G. Knepley PetscErrorCode ierr; 298478db826SMatthew G. Knepley 299478db826SMatthew G. Knepley PetscFunctionBegin; 300478db826SMatthew G. Knepley ierr = SNESGetDM(snes, &dm);CHKERRQ(ierr); 301478db826SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 302478db826SMatthew G. Knepley PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *); 303478db826SMatthew G. Knepley 304478db826SMatthew G. Knepley ierr = DMGetNullSpaceConstructor(dm, f, &nspconstr);CHKERRQ(ierr); 305478db826SMatthew G. Knepley if (nspconstr) { 306478db826SMatthew G. Knepley MatNullSpace nullsp; 307478db826SMatthew G. Knepley Mat J; 308478db826SMatthew G. Knepley 309478db826SMatthew G. Knepley ierr = (*nspconstr)(dm, f, f,&nullsp);CHKERRQ(ierr); 310478db826SMatthew G. Knepley ierr = SNESSetUp(snes);CHKERRQ(ierr); 311478db826SMatthew G. Knepley ierr = SNESGetJacobian(snes, &J, NULL, NULL, NULL);CHKERRQ(ierr); 312478db826SMatthew G. Knepley ierr = MatSetNullSpace(J, nullsp);CHKERRQ(ierr); 313478db826SMatthew G. Knepley ierr = MatNullSpaceDestroy(&nullsp);CHKERRQ(ierr); 314478db826SMatthew G. Knepley break; 315478db826SMatthew G. Knepley } 316478db826SMatthew G. Knepley } 317478db826SMatthew G. Knepley PetscFunctionReturn(0); 318478db826SMatthew G. Knepley } 319478db826SMatthew G. Knepley 320900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstGetConvRateSNES_Private(PetscConvEst ce, PetscReal alpha[]) 321900f6b5bSMatthew G. Knepley { 322900f6b5bSMatthew G. Knepley SNES snes = (SNES) ce->solver; 323900f6b5bSMatthew G. Knepley DM *dm; 324900f6b5bSMatthew G. Knepley PetscObject disc; 325900f6b5bSMatthew G. Knepley PetscReal *x, *y, slope, intercept; 3267809adefSMatthew G. Knepley PetscInt Nr = ce->Nr, r, f, dim, oldlevel, oldnlev; 327900f6b5bSMatthew G. Knepley void *ctx; 328900f6b5bSMatthew G. Knepley PetscErrorCode ierr; 329900f6b5bSMatthew G. Knepley 330900f6b5bSMatthew G. Knepley PetscFunctionBegin; 3312e61be88SMatthew G. Knepley if (ce->r != 2.0) SETERRQ1(PetscObjectComm((PetscObject) ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double) ce->r); 3326af0ca60SMatthew G. Knepley ierr = DMGetDimension(ce->idm, &dim);CHKERRQ(ierr); 3336af0ca60SMatthew G. Knepley ierr = DMGetApplicationContext(ce->idm, &ctx);CHKERRQ(ierr); 3346af0ca60SMatthew G. Knepley ierr = DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE);CHKERRQ(ierr); 33573269098SMatthew G. Knepley ierr = DMGetRefineLevel(ce->idm, &oldlevel);CHKERRQ(ierr); 3367809adefSMatthew G. Knepley ierr = PetscMalloc1((Nr+1), &dm);CHKERRQ(ierr); 3376af0ca60SMatthew G. Knepley /* Loop over meshes */ 338900f6b5bSMatthew G. Knepley dm[0] = ce->idm; 3396af0ca60SMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 340900f6b5bSMatthew G. Knepley Vec u; 341e5ed2c37SJose E. Roman #if defined(PETSC_USE_LOG) 342adfa7136SMatthew G. Knepley PetscLogStage stage; 343e5ed2c37SJose E. Roman #endif 344adfa7136SMatthew G. Knepley char stageName[PETSC_MAX_PATH_LEN]; 345900f6b5bSMatthew G. Knepley const char *dmname, *uname; 346adfa7136SMatthew G. Knepley 347adfa7136SMatthew G. Knepley ierr = PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN-1, "ConvEst Refinement Level %D", r);CHKERRQ(ierr); 348608e5a7aSMatthew G. Knepley #if defined(PETSC_USE_LOG) 349608e5a7aSMatthew G. Knepley ierr = PetscLogStageGetId(stageName, &stage);CHKERRQ(ierr); 350608e5a7aSMatthew G. Knepley if (stage < 0) {ierr = PetscLogStageRegister(stageName, &stage);CHKERRQ(ierr);} 351608e5a7aSMatthew G. Knepley #endif 352adfa7136SMatthew G. Knepley ierr = PetscLogStagePush(stage);CHKERRQ(ierr); 3536af0ca60SMatthew G. Knepley if (r > 0) { 354b2df8587SMatthew G. Knepley if (!ce->noRefine) { 3556af0ca60SMatthew G. Knepley ierr = DMRefine(dm[r-1], MPI_COMM_NULL, &dm[r]);CHKERRQ(ierr); 3566af0ca60SMatthew G. Knepley ierr = DMSetCoarseDM(dm[r], dm[r-1]);CHKERRQ(ierr); 357b2df8587SMatthew G. Knepley } else { 358b2df8587SMatthew G. Knepley DM cdm, rcdm; 359b2df8587SMatthew G. Knepley 360b2df8587SMatthew G. Knepley ierr = DMClone(dm[r-1], &dm[r]);CHKERRQ(ierr); 361b2df8587SMatthew G. Knepley ierr = DMCopyDisc(dm[r-1], dm[r]);CHKERRQ(ierr); 362b2df8587SMatthew G. Knepley ierr = DMGetCoordinateDM(dm[r-1], &cdm);CHKERRQ(ierr); 363b2df8587SMatthew G. Knepley ierr = DMGetCoordinateDM(dm[r], &rcdm);CHKERRQ(ierr); 364b2df8587SMatthew G. Knepley ierr = DMCopyDisc(cdm, rcdm);CHKERRQ(ierr); 365b2df8587SMatthew G. Knepley } 366ca3d3a14SMatthew G. Knepley ierr = DMCopyTransform(ce->idm, dm[r]);CHKERRQ(ierr); 3676af0ca60SMatthew G. Knepley ierr = PetscObjectGetName((PetscObject) dm[r-1], &dmname);CHKERRQ(ierr); 3686af0ca60SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) dm[r], dmname);CHKERRQ(ierr); 369478db826SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 3708cda7954SMatthew G. Knepley PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *); 371900f6b5bSMatthew G. Knepley 3728f892730SMatthew G. Knepley ierr = DMGetNullSpaceConstructor(dm[r-1], f, &nspconstr);CHKERRQ(ierr); 3738f892730SMatthew G. Knepley ierr = DMSetNullSpaceConstructor(dm[r], f, nspconstr);CHKERRQ(ierr); 3748f892730SMatthew G. Knepley } 3756af0ca60SMatthew G. Knepley } 3766af0ca60SMatthew G. Knepley ierr = DMViewFromOptions(dm[r], NULL, "-conv_dm_view");CHKERRQ(ierr); 3776af0ca60SMatthew G. Knepley /* Create solution */ 3786af0ca60SMatthew G. Knepley ierr = DMCreateGlobalVector(dm[r], &u);CHKERRQ(ierr); 379e5e52638SMatthew G. Knepley ierr = DMGetField(dm[r], 0, NULL, &disc);CHKERRQ(ierr); 380239a0881SMatthew G. Knepley ierr = PetscObjectGetName(disc, &uname);CHKERRQ(ierr); 3816af0ca60SMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) u, uname);CHKERRQ(ierr); 3826af0ca60SMatthew G. Knepley /* Setup solver */ 383900f6b5bSMatthew G. Knepley ierr = SNESReset(snes);CHKERRQ(ierr); 384900f6b5bSMatthew G. Knepley ierr = SNESSetDM(snes, dm[r]);CHKERRQ(ierr); 3856af0ca60SMatthew G. Knepley ierr = DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx);CHKERRQ(ierr); 386900f6b5bSMatthew G. Knepley ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 387478db826SMatthew G. Knepley /* Set nullspace for Jacobian */ 388478db826SMatthew G. Knepley ierr = PetscConvEstSetJacobianNullspace_Private(ce, snes);CHKERRQ(ierr); 3896af0ca60SMatthew G. Knepley /* Create initial guess */ 390900f6b5bSMatthew G. Knepley ierr = PetscConvEstComputeInitialGuess(ce, r, dm[r], u);CHKERRQ(ierr); 391900f6b5bSMatthew G. Knepley ierr = SNESSolve(snes, NULL, u);CHKERRQ(ierr); 392900f6b5bSMatthew G. Knepley ierr = PetscLogEventBegin(ce->event, ce, 0, 0, 0);CHKERRQ(ierr); 393900f6b5bSMatthew G. Knepley ierr = PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r*ce->Nf]);CHKERRQ(ierr); 394900f6b5bSMatthew G. Knepley ierr = PetscLogEventEnd(ce->event, ce, 0, 0, 0);CHKERRQ(ierr); 395adfa7136SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 396a56d3bf6SMatthew G. Knepley PetscSection s, fs; 397a56d3bf6SMatthew G. Knepley PetscInt lsize; 398a56d3bf6SMatthew G. Knepley 399a56d3bf6SMatthew G. Knepley /* Could use DMGetOutputDM() to add in Dirichlet dofs */ 40092fd8e1eSJed Brown ierr = DMGetLocalSection(dm[r], &s);CHKERRQ(ierr); 401a56d3bf6SMatthew G. Knepley ierr = PetscSectionGetField(s, f, &fs);CHKERRQ(ierr); 402a56d3bf6SMatthew G. Knepley ierr = PetscSectionGetConstrainedStorageSize(fs, &lsize);CHKERRQ(ierr); 403ffc4695bSBarry Smith ierr = MPI_Allreduce(&lsize, &ce->dofs[r*ce->Nf+f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) snes));CHKERRMPI(ierr); 4047809adefSMatthew G. Knepley ierr = PetscLogEventSetDof(ce->event, f, ce->dofs[r*ce->Nf+f]);CHKERRQ(ierr); 405900f6b5bSMatthew G. Knepley ierr = PetscLogEventSetError(ce->event, f, ce->errors[r*ce->Nf+f]);CHKERRQ(ierr); 406adfa7136SMatthew G. Knepley } 4076af0ca60SMatthew G. Knepley /* Monitor */ 408f2cacb80SMatthew G. Knepley ierr = PetscConvEstMonitorDefault(ce, r);CHKERRQ(ierr); 40973269098SMatthew G. Knepley if (!r) { 41073269098SMatthew G. Knepley /* PCReset() does not wipe out the level structure */ 41173269098SMatthew G. Knepley KSP ksp; 41273269098SMatthew G. Knepley PC pc; 41373269098SMatthew G. Knepley 414900f6b5bSMatthew G. Knepley ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 41573269098SMatthew G. Knepley ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 41673269098SMatthew G. Knepley ierr = PCMGGetLevels(pc, &oldnlev);CHKERRQ(ierr); 41773269098SMatthew G. Knepley } 4186af0ca60SMatthew G. Knepley /* Cleanup */ 4196af0ca60SMatthew G. Knepley ierr = VecDestroy(&u);CHKERRQ(ierr); 420adfa7136SMatthew G. Knepley ierr = PetscLogStagePop();CHKERRQ(ierr); 4216af0ca60SMatthew G. Knepley } 4226af0ca60SMatthew G. Knepley for (r = 1; r <= Nr; ++r) { 4236af0ca60SMatthew G. Knepley ierr = DMDestroy(&dm[r]);CHKERRQ(ierr); 4246af0ca60SMatthew G. Knepley } 4256af0ca60SMatthew G. Knepley /* Fit convergence rate */ 4266af0ca60SMatthew G. Knepley ierr = PetscMalloc2(Nr+1, &x, Nr+1, &y);CHKERRQ(ierr); 42746079b62SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 4286af0ca60SMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 4297809adefSMatthew G. Knepley x[r] = PetscLog10Real(ce->dofs[r*ce->Nf+f]); 43046079b62SMatthew G. Knepley y[r] = PetscLog10Real(ce->errors[r*ce->Nf+f]); 4316af0ca60SMatthew G. Knepley } 432bebf13c0SMatthew G. Knepley ierr = PetscLinearRegression(Nr+1, x, y, &slope, &intercept);CHKERRQ(ierr); 4336af0ca60SMatthew G. Knepley /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */ 43446079b62SMatthew G. Knepley alpha[f] = -slope * dim; 43546079b62SMatthew G. Knepley } 43646079b62SMatthew G. Knepley ierr = PetscFree2(x, y);CHKERRQ(ierr); 4377809adefSMatthew G. Knepley ierr = PetscFree(dm);CHKERRQ(ierr); 4382cae373cSMatthew G. Knepley /* Restore solver */ 439900f6b5bSMatthew G. Knepley ierr = SNESReset(snes);CHKERRQ(ierr); 44073269098SMatthew G. Knepley { 44173269098SMatthew G. Knepley /* PCReset() does not wipe out the level structure */ 44273269098SMatthew G. Knepley KSP ksp; 44373269098SMatthew G. Knepley PC pc; 44473269098SMatthew G. Knepley 445900f6b5bSMatthew G. Knepley ierr = SNESGetKSP(snes, &ksp);CHKERRQ(ierr); 44673269098SMatthew G. Knepley ierr = KSPGetPC(ksp, &pc);CHKERRQ(ierr); 44773269098SMatthew G. Knepley ierr = PCMGSetLevels(pc, oldnlev, NULL);CHKERRQ(ierr); 44873269098SMatthew G. Knepley ierr = DMSetRefineLevel(ce->idm, oldlevel);CHKERRQ(ierr); /* The damn DMCoarsen() calls in PCMG can reset this */ 44973269098SMatthew G. Knepley } 450900f6b5bSMatthew G. Knepley ierr = SNESSetDM(snes, ce->idm);CHKERRQ(ierr); 4512cae373cSMatthew G. Knepley ierr = DMPlexSetSNESLocalFEM(ce->idm, ctx, ctx, ctx);CHKERRQ(ierr); 452900f6b5bSMatthew G. Knepley ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 453478db826SMatthew G. Knepley ierr = PetscConvEstSetJacobianNullspace_Private(ce, snes);CHKERRQ(ierr); 454900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 455900f6b5bSMatthew G. Knepley } 456900f6b5bSMatthew G. Knepley 457900f6b5bSMatthew G. Knepley /*@ 458900f6b5bSMatthew G. Knepley PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization 459900f6b5bSMatthew G. Knepley 460900f6b5bSMatthew G. Knepley Not collective 461900f6b5bSMatthew G. Knepley 462900f6b5bSMatthew G. Knepley Input Parameter: 463900f6b5bSMatthew G. Knepley . ce - The PetscConvEst object 464900f6b5bSMatthew G. Knepley 465900f6b5bSMatthew G. Knepley Output Parameter: 466900f6b5bSMatthew G. Knepley . alpha - The convergence rate for each field 467900f6b5bSMatthew G. Knepley 468900f6b5bSMatthew G. Knepley Note: The convergence rate alpha is defined by 469900f6b5bSMatthew G. Knepley $ || u_\Delta - u_exact || < C \Delta^alpha 470900f6b5bSMatthew G. Knepley where u_\Delta is the discrete solution, and Delta is a measure of the discretization size. We usually use h for the 471900f6b5bSMatthew G. Knepley spatial resolution and \Delta t for the temporal resolution. 472900f6b5bSMatthew G. Knepley 473900f6b5bSMatthew G. Knepley We solve a series of problems using increasing resolution (refined meshes or decreased timesteps), calculate an error 474900f6b5bSMatthew G. Knepley based upon the exact solution in the DS, and then fit the result to our model above using linear regression. 475900f6b5bSMatthew G. Knepley 476900f6b5bSMatthew G. Knepley Options database keys: 477900f6b5bSMatthew G. Knepley + -snes_convergence_estimate : Execute convergence estimation inside SNESSolve() and print out the rate 478900f6b5bSMatthew G. Knepley - -ts_convergence_estimate : Execute convergence estimation inside TSSolve() and print out the rate 479900f6b5bSMatthew G. Knepley 480900f6b5bSMatthew G. Knepley Level: intermediate 481900f6b5bSMatthew G. Knepley 482900f6b5bSMatthew G. Knepley .seealso: PetscConvEstSetSolver(), PetscConvEstCreate(), PetscConvEstGetConvRate(), SNESSolve(), TSSolve() 483900f6b5bSMatthew G. Knepley @*/ 484900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[]) 485900f6b5bSMatthew G. Knepley { 486900f6b5bSMatthew G. Knepley PetscInt f; 487900f6b5bSMatthew G. Knepley PetscErrorCode ierr; 488900f6b5bSMatthew G. Knepley 489900f6b5bSMatthew G. Knepley PetscFunctionBegin; 490900f6b5bSMatthew G. Knepley if (ce->event < 0) {ierr = PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &ce->event);CHKERRQ(ierr);} 491900f6b5bSMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0; 492900f6b5bSMatthew G. Knepley ierr = (*ce->ops->getconvrate)(ce, alpha);CHKERRQ(ierr); 4936af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 4946af0ca60SMatthew G. Knepley } 4956af0ca60SMatthew G. Knepley 4966af0ca60SMatthew G. Knepley /*@ 4976af0ca60SMatthew G. Knepley PetscConvEstRateView - Displays the convergence rate to a viewer 4986af0ca60SMatthew G. Knepley 4996af0ca60SMatthew G. Knepley Collective on SNES 5006af0ca60SMatthew G. Knepley 5016af0ca60SMatthew G. Knepley Parameter: 5026af0ca60SMatthew G. Knepley + snes - iterative context obtained from SNESCreate() 50346079b62SMatthew G. Knepley . alpha - the convergence rate for each field 5046af0ca60SMatthew G. Knepley - viewer - the viewer to display the reason 5056af0ca60SMatthew G. Knepley 5066af0ca60SMatthew G. Knepley Options Database Keys: 5076af0ca60SMatthew G. Knepley . -snes_convergence_estimate - print the convergence rate 5086af0ca60SMatthew G. Knepley 5096af0ca60SMatthew G. Knepley Level: developer 5106af0ca60SMatthew G. Knepley 5116af0ca60SMatthew G. Knepley .seealso: PetscConvEstGetRate() 5126af0ca60SMatthew G. Knepley @*/ 513a56d3bf6SMatthew G. Knepley PetscErrorCode PetscConvEstRateView(PetscConvEst ce, const PetscReal alpha[], PetscViewer viewer) 5146af0ca60SMatthew G. Knepley { 5156af0ca60SMatthew G. Knepley PetscBool isAscii; 5166af0ca60SMatthew G. Knepley PetscErrorCode ierr; 5176af0ca60SMatthew G. Knepley 5186af0ca60SMatthew G. Knepley PetscFunctionBegin; 5196af0ca60SMatthew G. Knepley ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isAscii);CHKERRQ(ierr); 5206af0ca60SMatthew G. Knepley if (isAscii) { 521900f6b5bSMatthew G. Knepley PetscInt Nf = ce->Nf, f; 522a56d3bf6SMatthew G. Knepley 5236af0ca60SMatthew G. Knepley ierr = PetscViewerASCIIAddTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr); 524adfa7136SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: ");CHKERRQ(ierr); 525a56d3bf6SMatthew G. Knepley if (Nf > 1) {ierr = PetscViewerASCIIPrintf(viewer, "[");CHKERRQ(ierr);} 526a56d3bf6SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 52746079b62SMatthew G. Knepley if (f > 0) {ierr = PetscViewerASCIIPrintf(viewer, ", ");CHKERRQ(ierr);} 5285f9b3039SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "%#.2g", (double) alpha[f]);CHKERRQ(ierr); 52946079b62SMatthew G. Knepley } 530a56d3bf6SMatthew G. Knepley if (Nf > 1) {ierr = PetscViewerASCIIPrintf(viewer, "]");CHKERRQ(ierr);} 531adfa7136SMatthew G. Knepley ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 5326af0ca60SMatthew G. Knepley ierr = PetscViewerASCIISubtractTab(viewer, ((PetscObject) ce)->tablevel);CHKERRQ(ierr); 5336af0ca60SMatthew G. Knepley } 5346af0ca60SMatthew G. Knepley PetscFunctionReturn(0); 5356af0ca60SMatthew G. Knepley } 536900f6b5bSMatthew G. Knepley 537900f6b5bSMatthew G. Knepley /*@ 538900f6b5bSMatthew G. Knepley PetscConvEstCreate - Create a PetscConvEst object 539900f6b5bSMatthew G. Knepley 540900f6b5bSMatthew G. Knepley Collective 541900f6b5bSMatthew G. Knepley 542900f6b5bSMatthew G. Knepley Input Parameter: 543900f6b5bSMatthew G. Knepley . comm - The communicator for the PetscConvEst object 544900f6b5bSMatthew G. Knepley 545900f6b5bSMatthew G. Knepley Output Parameter: 546900f6b5bSMatthew G. Knepley . ce - The PetscConvEst object 547900f6b5bSMatthew G. Knepley 548900f6b5bSMatthew G. Knepley Level: beginner 549900f6b5bSMatthew G. Knepley 550900f6b5bSMatthew G. Knepley .seealso: PetscConvEstDestroy(), PetscConvEstGetConvRate() 551900f6b5bSMatthew G. Knepley @*/ 552900f6b5bSMatthew G. Knepley PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce) 553900f6b5bSMatthew G. Knepley { 554900f6b5bSMatthew G. Knepley PetscErrorCode ierr; 555900f6b5bSMatthew G. Knepley 556900f6b5bSMatthew G. Knepley PetscFunctionBegin; 557900f6b5bSMatthew G. Knepley PetscValidPointer(ce, 2); 558900f6b5bSMatthew G. Knepley ierr = PetscSysInitializePackage();CHKERRQ(ierr); 559900f6b5bSMatthew G. Knepley ierr = PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView);CHKERRQ(ierr); 560900f6b5bSMatthew G. Knepley (*ce)->monitor = PETSC_FALSE; 5612e61be88SMatthew G. Knepley (*ce)->r = 2.0; 562900f6b5bSMatthew G. Knepley (*ce)->Nr = 4; 563900f6b5bSMatthew G. Knepley (*ce)->event = -1; 564900f6b5bSMatthew G. Knepley (*ce)->ops->setsolver = PetscConvEstSetSNES_Private; 565900f6b5bSMatthew G. Knepley (*ce)->ops->initguess = PetscConvEstInitGuessSNES_Private; 566900f6b5bSMatthew G. Knepley (*ce)->ops->computeerror = PetscConvEstComputeErrorSNES_Private; 567900f6b5bSMatthew G. Knepley (*ce)->ops->getconvrate = PetscConvEstGetConvRateSNES_Private; 568900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 569900f6b5bSMatthew G. Knepley } 570