12eace19aSMatthew G. Knepley #include "petscsys.h" 26af0ca60SMatthew G. Knepley #include <petscconvest.h> /*I "petscconvest.h" I*/ 36af0ca60SMatthew G. Knepley #include <petscdmplex.h> 46af0ca60SMatthew G. Knepley #include <petscds.h> 56af0ca60SMatthew G. Knepley 66af0ca60SMatthew G. Knepley #include <petsc/private/petscconvestimpl.h> 76af0ca60SMatthew G. Knepley 8d71ae5a4SJacob Faibussowitsch static PetscErrorCode zero_private(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 9d71ae5a4SJacob Faibussowitsch { 106af0ca60SMatthew G. Knepley PetscInt c; 116af0ca60SMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = 0.0; 123ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 136af0ca60SMatthew G. Knepley } 146af0ca60SMatthew G. Knepley 156af0ca60SMatthew G. Knepley /*@ 16420bcc1bSBarry Smith PetscConvEstDestroy - Destroys a PETSc convergence estimator `PetscConvEst` object 176af0ca60SMatthew G. Knepley 18c3339decSBarry Smith Collective 196af0ca60SMatthew G. Knepley 206af0ca60SMatthew G. Knepley Input Parameter: 2120f4b53cSBarry Smith . ce - The `PetscConvEst` object 226af0ca60SMatthew G. Knepley 236af0ca60SMatthew G. Knepley Level: beginner 246af0ca60SMatthew G. Knepley 2520f4b53cSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()` 266af0ca60SMatthew G. Knepley @*/ 27d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstDestroy(PetscConvEst *ce) 28d71ae5a4SJacob Faibussowitsch { 296af0ca60SMatthew G. Knepley PetscFunctionBegin; 303ba16761SJacob Faibussowitsch if (!*ce) PetscFunctionReturn(PETSC_SUCCESS); 31f4f49eeaSPierre Jolivet PetscValidHeaderSpecific(*ce, PETSC_OBJECT_CLASSID, 1); 32f4f49eeaSPierre Jolivet if (--((PetscObject)*ce)->refct > 0) { 336af0ca60SMatthew G. Knepley *ce = NULL; 343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 356af0ca60SMatthew G. Knepley } 369566063dSJacob Faibussowitsch PetscCall(PetscFree3((*ce)->initGuess, (*ce)->exactSol, (*ce)->ctxs)); 379566063dSJacob Faibussowitsch PetscCall(PetscFree2((*ce)->dofs, (*ce)->errors)); 389566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(ce)); 393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 406af0ca60SMatthew G. Knepley } 416af0ca60SMatthew G. Knepley 426af0ca60SMatthew G. Knepley /*@ 43420bcc1bSBarry Smith PetscConvEstSetFromOptions - Sets a convergence estimator `PetscConvEst` object based on values in the options database 446af0ca60SMatthew G. Knepley 45c3339decSBarry Smith Collective 466af0ca60SMatthew G. Knepley 472fe279fdSBarry Smith Input Parameter: 4820f4b53cSBarry Smith . ce - The `PetscConvEst` object 496af0ca60SMatthew G. Knepley 506af0ca60SMatthew G. Knepley Level: beginner 516af0ca60SMatthew G. Knepley 5220f4b53cSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()` 536af0ca60SMatthew G. Knepley @*/ 54d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstSetFromOptions(PetscConvEst ce) 55d71ae5a4SJacob Faibussowitsch { 566af0ca60SMatthew G. Knepley PetscFunctionBegin; 57d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)ce), "", "Convergence Estimator Options", "PetscConvEst"); 589566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-convest_num_refine", "The number of refinements for the convergence check", "PetscConvEst", ce->Nr, &ce->Nr, NULL)); 599566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-convest_refine_factor", "The increase in resolution in each dimension", "PetscConvEst", ce->r, &ce->r, NULL)); 609566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-convest_monitor", "Monitor the error for each convergence check", "PetscConvEst", ce->monitor, &ce->monitor, NULL)); 619566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-convest_no_refine", "Debugging flag to run on the same mesh each time", "PetscConvEst", ce->noRefine, &ce->noRefine, NULL)); 62d0609cedSBarry Smith PetscOptionsEnd(); 633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 646af0ca60SMatthew G. Knepley } 656af0ca60SMatthew G. Knepley 666af0ca60SMatthew G. Knepley /*@ 6720f4b53cSBarry Smith PetscConvEstView - Views a `PetscConvEst` object 686af0ca60SMatthew G. Knepley 69c3339decSBarry Smith Collective 706af0ca60SMatthew G. Knepley 716af0ca60SMatthew G. Knepley Input Parameters: 7220f4b53cSBarry Smith + ce - The `PetscConvEst` object 73420bcc1bSBarry Smith - viewer - The `PetscViewer` 746af0ca60SMatthew G. Knepley 756af0ca60SMatthew G. Knepley Level: beginner 766af0ca60SMatthew G. Knepley 7720f4b53cSBarry Smith .seealso: `PetscConvEst`, `PetscViewer`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()` 786af0ca60SMatthew G. Knepley @*/ 79d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstView(PetscConvEst ce, PetscViewer viewer) 80d71ae5a4SJacob Faibussowitsch { 816af0ca60SMatthew G. Knepley PetscFunctionBegin; 829566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)ce, viewer)); 8363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "ConvEst with %" PetscInt_FMT " levels\n", ce->Nr + 1)); 843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 856af0ca60SMatthew G. Knepley } 866af0ca60SMatthew G. Knepley 876af0ca60SMatthew G. Knepley /*@ 886af0ca60SMatthew G. Knepley PetscConvEstGetSolver - Gets the solver used to produce discrete solutions 896af0ca60SMatthew G. Knepley 9020f4b53cSBarry Smith Not Collective 916af0ca60SMatthew G. Knepley 926af0ca60SMatthew G. Knepley Input Parameter: 9320f4b53cSBarry Smith . ce - The `PetscConvEst` object 946af0ca60SMatthew G. Knepley 956af0ca60SMatthew G. Knepley Output Parameter: 96900f6b5bSMatthew G. Knepley . solver - The solver 976af0ca60SMatthew G. Knepley 986af0ca60SMatthew G. Knepley Level: intermediate 996af0ca60SMatthew G. Knepley 10020f4b53cSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstSetSolver()`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()` 1016af0ca60SMatthew G. Knepley @*/ 102d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstGetSolver(PetscConvEst ce, PetscObject *solver) 103d71ae5a4SJacob Faibussowitsch { 1046af0ca60SMatthew G. Knepley PetscFunctionBegin; 1056af0ca60SMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 1064f572ea9SToby Isaac PetscAssertPointer(solver, 2); 107900f6b5bSMatthew G. Knepley *solver = ce->solver; 1083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1096af0ca60SMatthew G. Knepley } 1106af0ca60SMatthew G. Knepley 1116af0ca60SMatthew G. Knepley /*@ 1126af0ca60SMatthew G. Knepley PetscConvEstSetSolver - Sets the solver used to produce discrete solutions 1136af0ca60SMatthew G. Knepley 11420f4b53cSBarry Smith Not Collective 1156af0ca60SMatthew G. Knepley 1166af0ca60SMatthew G. Knepley Input Parameters: 11720f4b53cSBarry Smith + ce - The `PetscConvEst` object 118420bcc1bSBarry Smith - solver - The solver, must be a `KSP`, `SNES`, or `TS` object with an attached `DM`/`DS`, that can compute an exact solution 1196af0ca60SMatthew G. Knepley 1206af0ca60SMatthew G. Knepley Level: intermediate 1216af0ca60SMatthew G. Knepley 12220f4b53cSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstGetSNES()`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()` 1236af0ca60SMatthew G. Knepley @*/ 124d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstSetSolver(PetscConvEst ce, PetscObject solver) 125d71ae5a4SJacob Faibussowitsch { 1266af0ca60SMatthew G. Knepley PetscFunctionBegin; 1276af0ca60SMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 128900f6b5bSMatthew G. Knepley PetscValidHeader(solver, 2); 129900f6b5bSMatthew G. Knepley ce->solver = solver; 130dbbe0bcdSBarry Smith PetscUseTypeMethod(ce, setsolver, solver); 1313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1326af0ca60SMatthew G. Knepley } 1336af0ca60SMatthew G. Knepley 1346af0ca60SMatthew G. Knepley /*@ 135420bcc1bSBarry Smith PetscConvEstSetUp - After the solver is specified, create data structures needed for estimating convergence 1366af0ca60SMatthew G. Knepley 137c3339decSBarry Smith Collective 1386af0ca60SMatthew G. Knepley 1392fe279fdSBarry Smith Input Parameter: 14020f4b53cSBarry Smith . ce - The `PetscConvEst` object 1416af0ca60SMatthew G. Knepley 1426af0ca60SMatthew G. Knepley Level: beginner 1436af0ca60SMatthew G. Knepley 14420f4b53cSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()` 1456af0ca60SMatthew G. Knepley @*/ 146d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstSetUp(PetscConvEst ce) 147d71ae5a4SJacob Faibussowitsch { 148083401c6SMatthew G. Knepley PetscInt Nf, f, Nds, s; 1496af0ca60SMatthew G. Knepley 1506af0ca60SMatthew G. Knepley PetscFunctionBegin; 1519566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(ce->idm, &Nf)); 152900f6b5bSMatthew G. Knepley ce->Nf = PetscMax(Nf, 1); 1539566063dSJacob Faibussowitsch PetscCall(PetscMalloc2((ce->Nr + 1) * ce->Nf, &ce->dofs, (ce->Nr + 1) * ce->Nf, &ce->errors)); 1549566063dSJacob Faibussowitsch PetscCall(PetscCalloc3(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol, ce->Nf, &ce->ctxs)); 155900f6b5bSMatthew G. Knepley for (f = 0; f < Nf; ++f) ce->initGuess[f] = zero_private; 1569566063dSJacob Faibussowitsch PetscCall(DMGetNumDS(ce->idm, &Nds)); 157083401c6SMatthew G. Knepley for (s = 0; s < Nds; ++s) { 158083401c6SMatthew G. Knepley PetscDS ds; 159083401c6SMatthew G. Knepley DMLabel label; 160083401c6SMatthew G. Knepley IS fieldIS; 161083401c6SMatthew G. Knepley const PetscInt *fields; 162083401c6SMatthew G. Knepley PetscInt dsNf; 163083401c6SMatthew G. Knepley 16407218a29SMatthew G. Knepley PetscCall(DMGetRegionNumDS(ce->idm, s, &label, &fieldIS, &ds, NULL)); 1659566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &dsNf)); 1669566063dSJacob Faibussowitsch if (fieldIS) PetscCall(ISGetIndices(fieldIS, &fields)); 167083401c6SMatthew G. Knepley for (f = 0; f < dsNf; ++f) { 168083401c6SMatthew G. Knepley const PetscInt field = fields[f]; 1699566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, field, &ce->exactSol[field], &ce->ctxs[field])); 170083401c6SMatthew G. Knepley } 1719566063dSJacob Faibussowitsch if (fieldIS) PetscCall(ISRestoreIndices(fieldIS, &fields)); 172083401c6SMatthew G. Knepley } 173ad540459SPierre Jolivet for (f = 0; f < Nf; ++f) PetscCheck(ce->exactSol[f], PetscObjectComm((PetscObject)ce), PETSC_ERR_ARG_WRONG, "DS must contain exact solution functions in order to estimate convergence, missing for field %" PetscInt_FMT, f); 1743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1756af0ca60SMatthew G. Knepley } 1766af0ca60SMatthew G. Knepley 177d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstComputeInitialGuess(PetscConvEst ce, PetscInt r, DM dm, Vec u) 178d71ae5a4SJacob Faibussowitsch { 1796af0ca60SMatthew G. Knepley PetscFunctionBegin; 180900f6b5bSMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 181900f6b5bSMatthew G. Knepley if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3); 182900f6b5bSMatthew G. Knepley PetscValidHeaderSpecific(u, VEC_CLASSID, 4); 183dbbe0bcdSBarry Smith PetscUseTypeMethod(ce, initguess, r, dm, u); 1843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 185900f6b5bSMatthew G. Knepley } 186900f6b5bSMatthew G. Knepley 187d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstComputeError(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[]) 188d71ae5a4SJacob Faibussowitsch { 189900f6b5bSMatthew G. Knepley PetscFunctionBegin; 190900f6b5bSMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 191900f6b5bSMatthew G. Knepley if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3); 192900f6b5bSMatthew G. Knepley PetscValidHeaderSpecific(u, VEC_CLASSID, 4); 1934f572ea9SToby Isaac PetscAssertPointer(errors, 5); 194dbbe0bcdSBarry Smith PetscUseTypeMethod(ce, computeerror, r, dm, u, errors); 1953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 196900f6b5bSMatthew G. Knepley } 197900f6b5bSMatthew G. Knepley 198f2cacb80SMatthew G. Knepley /*@ 199f2cacb80SMatthew G. Knepley PetscConvEstMonitorDefault - Monitors the convergence estimation loop 200f2cacb80SMatthew G. Knepley 201c3339decSBarry Smith Collective 202f2cacb80SMatthew G. Knepley 203d8d19677SJose E. Roman Input Parameters: 204f6dfbefdSBarry Smith + ce - The `PetscConvEst` object 205f2cacb80SMatthew G. Knepley - r - The refinement level 206f2cacb80SMatthew G. Knepley 207420bcc1bSBarry Smith Options Database Key: 208ee300463SSatish Balay . -convest_monitor - Activate the monitor 209f2cacb80SMatthew G. Knepley 210f2cacb80SMatthew G. Knepley Level: intermediate 211f2cacb80SMatthew G. Knepley 212f6dfbefdSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`, `SNESSolve()`, `TSSolve()` 213f2cacb80SMatthew G. Knepley @*/ 214d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstMonitorDefault(PetscConvEst ce, PetscInt r) 215d71ae5a4SJacob Faibussowitsch { 216900f6b5bSMatthew G. Knepley MPI_Comm comm; 217900f6b5bSMatthew G. Knepley PetscInt f; 218900f6b5bSMatthew G. Knepley 219900f6b5bSMatthew G. Knepley PetscFunctionBegin; 220900f6b5bSMatthew G. Knepley if (ce->monitor) { 2217809adefSMatthew G. Knepley PetscInt *dofs = &ce->dofs[r * ce->Nf]; 222900f6b5bSMatthew G. Knepley PetscReal *errors = &ce->errors[r * ce->Nf]; 223900f6b5bSMatthew G. Knepley 2249566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ce, &comm)); 2259566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "N: ")); 2269566063dSJacob Faibussowitsch if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "[")); 2277809adefSMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 2289566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscPrintf(comm, ", ")); 22963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "%7" PetscInt_FMT, dofs[f])); 2307809adefSMatthew G. Knepley } 2319566063dSJacob Faibussowitsch if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "]")); 2329566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " ")); 2339566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "L_2 Error: ")); 2349566063dSJacob Faibussowitsch if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "[")); 235900f6b5bSMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 2369566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscPrintf(comm, ", ")); 2379566063dSJacob Faibussowitsch if (errors[f] < 1.0e-11) PetscCall(PetscPrintf(comm, "< 1e-11")); 2389566063dSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "%g", (double)errors[f])); 239900f6b5bSMatthew G. Knepley } 2409566063dSJacob Faibussowitsch if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "]")); 2419566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "\n")); 242900f6b5bSMatthew G. Knepley } 2433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 244900f6b5bSMatthew G. Knepley } 245900f6b5bSMatthew G. Knepley 246d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstSetSNES_Private(PetscConvEst ce, PetscObject solver) 247d71ae5a4SJacob Faibussowitsch { 248900f6b5bSMatthew G. Knepley PetscClassId id; 249900f6b5bSMatthew G. Knepley 250900f6b5bSMatthew G. Knepley PetscFunctionBegin; 2519566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(ce->solver, &id)); 25208401ef6SPierre Jolivet PetscCheck(id == SNES_CLASSID, PetscObjectComm((PetscObject)ce), PETSC_ERR_ARG_WRONG, "Solver was not a SNES"); 2539566063dSJacob Faibussowitsch PetscCall(SNESGetDM((SNES)ce->solver, &ce->idm)); 2543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 255900f6b5bSMatthew G. Knepley } 256900f6b5bSMatthew G. Knepley 257d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstInitGuessSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u) 258d71ae5a4SJacob Faibussowitsch { 259900f6b5bSMatthew G. Knepley PetscFunctionBegin; 2609566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, ce->initGuess, ce->ctxs, INSERT_VALUES, u)); 2613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 262900f6b5bSMatthew G. Knepley } 263900f6b5bSMatthew G. Knepley 264d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstComputeErrorSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[]) 265d71ae5a4SJacob Faibussowitsch { 2662eace19aSMatthew G. Knepley const char *prefix; 2672eace19aSMatthew G. Knepley PetscBool errorView = PETSC_FALSE; 2682eace19aSMatthew G. Knepley 269900f6b5bSMatthew G. Knepley PetscFunctionBegin; 2709566063dSJacob Faibussowitsch PetscCall(DMComputeL2FieldDiff(dm, 0.0, ce->exactSol, ce->ctxs, u, errors)); 2712eace19aSMatthew G. Knepley PetscCall(PetscObjectGetOptionsPrefix((PetscObject)ce, &prefix)); 2722eace19aSMatthew G. Knepley PetscCall(PetscOptionsHasName(NULL, prefix, "-convest_error_view", &errorView)); 2732eace19aSMatthew G. Knepley if (errorView) { 2742eace19aSMatthew G. Knepley DM dmError; 2752eace19aSMatthew G. Knepley PetscFE feError, fe; 2762eace19aSMatthew G. Knepley PetscQuadrature quad; 2772eace19aSMatthew G. Knepley Vec e; 2782eace19aSMatthew G. Knepley PetscDS ds; 279*2192575eSBarry Smith PetscSimplePointFn **funcs; 2802eace19aSMatthew G. Knepley void **ctxs; 2812eace19aSMatthew G. Knepley PetscInt dim, Nf; 2822eace19aSMatthew G. Knepley 2832eace19aSMatthew G. Knepley PetscCall(DMGetDimension(dm, &dim)); 2842eace19aSMatthew G. Knepley PetscCall(DMGetDS(dm, &ds)); 2852eace19aSMatthew G. Knepley PetscCall(PetscDSGetNumFields(ds, &Nf)); 2862eace19aSMatthew G. Knepley PetscCall(PetscMalloc2(Nf, &funcs, Nf, &ctxs)); 2872eace19aSMatthew G. Knepley for (PetscInt f = 0; f < Nf; ++f) PetscCall(PetscDSGetExactSolution(ds, f, &funcs[f], &ctxs[f])); 2882eace19aSMatthew G. Knepley PetscCall(DMClone(dm, &dmError)); 2892eace19aSMatthew G. Knepley PetscCall(PetscFECreateLagrange(PETSC_COMM_SELF, dim, 1, PETSC_FALSE, 0, PETSC_DETERMINE, &feError)); 2902eace19aSMatthew G. Knepley PetscCall(DMGetField(dm, 0, NULL, (PetscObject *)&fe)); 2912eace19aSMatthew G. Knepley PetscCall(PetscFEGetQuadrature(fe, &quad)); 2922eace19aSMatthew G. Knepley PetscCall(PetscFESetQuadrature(feError, quad)); 2932eace19aSMatthew G. Knepley PetscCall(DMSetField(dmError, 0, NULL, (PetscObject)feError)); 2942eace19aSMatthew G. Knepley PetscCall(PetscFEDestroy(&feError)); 2952eace19aSMatthew G. Knepley PetscCall(DMCreateDS(dmError)); 2962eace19aSMatthew G. Knepley PetscCall(DMGetGlobalVector(dmError, &e)); 2972eace19aSMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)e, "Error")); 2982eace19aSMatthew G. Knepley PetscCall(DMPlexComputeL2DiffVec(dm, 0., funcs, ctxs, u, e)); 2992eace19aSMatthew G. Knepley PetscCall(VecViewFromOptions(e, NULL, "-convest_error_view")); 3002eace19aSMatthew G. Knepley PetscCall(DMRestoreGlobalVector(dmError, &e)); 3012eace19aSMatthew G. Knepley PetscCall(DMDestroy(&dmError)); 3022eace19aSMatthew G. Knepley PetscCall(PetscFree2(funcs, ctxs)); 3032eace19aSMatthew G. Knepley } 3043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 305900f6b5bSMatthew G. Knepley } 306900f6b5bSMatthew G. Knepley 3078c401167SJed Brown static PetscErrorCode PetscConvEstSetJacobianNullSpace_Private(PetscConvEst ce, SNES snes) 308d71ae5a4SJacob Faibussowitsch { 309478db826SMatthew G. Knepley DM dm; 310478db826SMatthew G. Knepley PetscInt f; 311478db826SMatthew G. Knepley 312478db826SMatthew G. Knepley PetscFunctionBegin; 3139566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 314478db826SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 315478db826SMatthew G. Knepley PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *); 316478db826SMatthew G. Knepley 3179566063dSJacob Faibussowitsch PetscCall(DMGetNullSpaceConstructor(dm, f, &nspconstr)); 318478db826SMatthew G. Knepley if (nspconstr) { 319478db826SMatthew G. Knepley MatNullSpace nullsp; 320478db826SMatthew G. Knepley Mat J; 321478db826SMatthew G. Knepley 3229566063dSJacob Faibussowitsch PetscCall((*nspconstr)(dm, f, f, &nullsp)); 3239566063dSJacob Faibussowitsch PetscCall(SNESSetUp(snes)); 3249566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL)); 3259566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(J, nullsp)); 3269566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nullsp)); 327478db826SMatthew G. Knepley break; 328478db826SMatthew G. Knepley } 329478db826SMatthew G. Knepley } 3303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 331478db826SMatthew G. Knepley } 332478db826SMatthew G. Knepley 333d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstGetConvRateSNES_Private(PetscConvEst ce, PetscReal alpha[]) 334d71ae5a4SJacob Faibussowitsch { 335900f6b5bSMatthew G. Knepley SNES snes = (SNES)ce->solver; 336900f6b5bSMatthew G. Knepley DM *dm; 337900f6b5bSMatthew G. Knepley PetscObject disc; 338900f6b5bSMatthew G. Knepley PetscReal *x, *y, slope, intercept; 3397809adefSMatthew G. Knepley PetscInt Nr = ce->Nr, r, f, dim, oldlevel, oldnlev; 340900f6b5bSMatthew G. Knepley void *ctx; 341900f6b5bSMatthew G. Knepley 342900f6b5bSMatthew G. Knepley PetscFunctionBegin; 34308401ef6SPierre Jolivet PetscCheck(ce->r == 2.0, PetscObjectComm((PetscObject)ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double)ce->r); 3449566063dSJacob Faibussowitsch PetscCall(DMGetDimension(ce->idm, &dim)); 3459566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(ce->idm, &ctx)); 3469566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE)); 3479566063dSJacob Faibussowitsch PetscCall(DMGetRefineLevel(ce->idm, &oldlevel)); 34857508eceSPierre Jolivet PetscCall(PetscMalloc1(Nr + 1, &dm)); 3496af0ca60SMatthew G. Knepley /* Loop over meshes */ 350900f6b5bSMatthew G. Knepley dm[0] = ce->idm; 3516af0ca60SMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 352900f6b5bSMatthew G. Knepley Vec u; 353adfa7136SMatthew G. Knepley PetscLogStage stage; 354adfa7136SMatthew G. Knepley char stageName[PETSC_MAX_PATH_LEN]; 355900f6b5bSMatthew G. Knepley const char *dmname, *uname; 356adfa7136SMatthew G. Knepley 35763a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN - 1, "ConvEst Refinement Level %" PetscInt_FMT, r)); 3589566063dSJacob Faibussowitsch PetscCall(PetscLogStageGetId(stageName, &stage)); 3599566063dSJacob Faibussowitsch if (stage < 0) PetscCall(PetscLogStageRegister(stageName, &stage)); 3609566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 3616af0ca60SMatthew G. Knepley if (r > 0) { 362b2df8587SMatthew G. Knepley if (!ce->noRefine) { 3639566063dSJacob Faibussowitsch PetscCall(DMRefine(dm[r - 1], MPI_COMM_NULL, &dm[r])); 3649566063dSJacob Faibussowitsch PetscCall(DMSetCoarseDM(dm[r], dm[r - 1])); 365b2df8587SMatthew G. Knepley } else { 366b2df8587SMatthew G. Knepley DM cdm, rcdm; 367b2df8587SMatthew G. Knepley 3689566063dSJacob Faibussowitsch PetscCall(DMClone(dm[r - 1], &dm[r])); 3699566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm[r - 1], dm[r])); 3709566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm[r - 1], &cdm)); 3719566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm[r], &rcdm)); 3729566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(cdm, rcdm)); 373b2df8587SMatthew G. Knepley } 3749566063dSJacob Faibussowitsch PetscCall(DMCopyTransform(ce->idm, dm[r])); 3759566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)dm[r - 1], &dmname)); 3769566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dm[r], dmname)); 377478db826SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 3788cda7954SMatthew G. Knepley PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *); 379900f6b5bSMatthew G. Knepley 3809566063dSJacob Faibussowitsch PetscCall(DMGetNullSpaceConstructor(dm[r - 1], f, &nspconstr)); 3819566063dSJacob Faibussowitsch PetscCall(DMSetNullSpaceConstructor(dm[r], f, nspconstr)); 3828f892730SMatthew G. Knepley } 3836af0ca60SMatthew G. Knepley } 3849566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm[r], NULL, "-conv_dm_view")); 3856af0ca60SMatthew G. Knepley /* Create solution */ 3869566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm[r], &u)); 3879566063dSJacob Faibussowitsch PetscCall(DMGetField(dm[r], 0, NULL, &disc)); 3889566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName(disc, &uname)); 3899566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, uname)); 3906af0ca60SMatthew G. Knepley /* Setup solver */ 3919566063dSJacob Faibussowitsch PetscCall(SNESReset(snes)); 3929566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm[r])); 3936493148fSStefano Zampini PetscCall(DMPlexSetSNESLocalFEM(dm[r], PETSC_FALSE, ctx)); 3942eace19aSMatthew G. Knepley PetscCall(DMPlexSetSNESVariableBounds(dm[r], snes)); 3959566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 396478db826SMatthew G. Knepley /* Set nullspace for Jacobian */ 3978c401167SJed Brown PetscCall(PetscConvEstSetJacobianNullSpace_Private(ce, snes)); 3986af0ca60SMatthew G. Knepley /* Create initial guess */ 3999566063dSJacob Faibussowitsch PetscCall(PetscConvEstComputeInitialGuess(ce, r, dm[r], u)); 4009566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 4019566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ce->event, ce, 0, 0, 0)); 4029566063dSJacob Faibussowitsch PetscCall(PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r * ce->Nf])); 4039566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ce->event, ce, 0, 0, 0)); 404adfa7136SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 405a56d3bf6SMatthew G. Knepley PetscSection s, fs; 406a56d3bf6SMatthew G. Knepley PetscInt lsize; 407a56d3bf6SMatthew G. Knepley 408a56d3bf6SMatthew G. Knepley /* Could use DMGetOutputDM() to add in Dirichlet dofs */ 4099566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm[r], &s)); 4109566063dSJacob Faibussowitsch PetscCall(PetscSectionGetField(s, f, &fs)); 4119566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstrainedStorageSize(fs, &lsize)); 412462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&lsize, &ce->dofs[r * ce->Nf + f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes))); 4139566063dSJacob Faibussowitsch PetscCall(PetscLogEventSetDof(ce->event, f, ce->dofs[r * ce->Nf + f])); 4149566063dSJacob Faibussowitsch PetscCall(PetscLogEventSetError(ce->event, f, ce->errors[r * ce->Nf + f])); 415adfa7136SMatthew G. Knepley } 4166af0ca60SMatthew G. Knepley /* Monitor */ 4179566063dSJacob Faibussowitsch PetscCall(PetscConvEstMonitorDefault(ce, r)); 41873269098SMatthew G. Knepley if (!r) { 41973269098SMatthew G. Knepley /* PCReset() does not wipe out the level structure */ 42073269098SMatthew G. Knepley KSP ksp; 42173269098SMatthew G. Knepley PC pc; 42273269098SMatthew G. Knepley 4239566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 4249566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 4259566063dSJacob Faibussowitsch PetscCall(PCMGGetLevels(pc, &oldnlev)); 42673269098SMatthew G. Knepley } 4276af0ca60SMatthew G. Knepley /* Cleanup */ 4289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 4299566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 4306af0ca60SMatthew G. Knepley } 43148a46eb9SPierre Jolivet for (r = 1; r <= Nr; ++r) PetscCall(DMDestroy(&dm[r])); 4326af0ca60SMatthew G. Knepley /* Fit convergence rate */ 4339566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nr + 1, &x, Nr + 1, &y)); 43446079b62SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 4356af0ca60SMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 4367809adefSMatthew G. Knepley x[r] = PetscLog10Real(ce->dofs[r * ce->Nf + f]); 43746079b62SMatthew G. Knepley y[r] = PetscLog10Real(ce->errors[r * ce->Nf + f]); 4386af0ca60SMatthew G. Knepley } 4399566063dSJacob Faibussowitsch PetscCall(PetscLinearRegression(Nr + 1, x, y, &slope, &intercept)); 4406af0ca60SMatthew G. Knepley /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */ 44146079b62SMatthew G. Knepley alpha[f] = -slope * dim; 44246079b62SMatthew G. Knepley } 4439566063dSJacob Faibussowitsch PetscCall(PetscFree2(x, y)); 4449566063dSJacob Faibussowitsch PetscCall(PetscFree(dm)); 4452cae373cSMatthew G. Knepley /* Restore solver */ 4469566063dSJacob Faibussowitsch PetscCall(SNESReset(snes)); 44773269098SMatthew G. Knepley { 44873269098SMatthew G. Knepley /* PCReset() does not wipe out the level structure */ 44973269098SMatthew G. Knepley KSP ksp; 45073269098SMatthew G. Knepley PC pc; 45173269098SMatthew G. Knepley 4529566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 4539566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 4549566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels(pc, oldnlev, NULL)); 4559566063dSJacob Faibussowitsch PetscCall(DMSetRefineLevel(ce->idm, oldlevel)); /* The damn DMCoarsen() calls in PCMG can reset this */ 45673269098SMatthew G. Knepley } 4579566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, ce->idm)); 4586493148fSStefano Zampini PetscCall(DMPlexSetSNESLocalFEM(ce->idm, PETSC_FALSE, ctx)); 4592eace19aSMatthew G. Knepley PetscCall(DMPlexSetSNESVariableBounds(ce->idm, snes)); 4609566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 4618c401167SJed Brown PetscCall(PetscConvEstSetJacobianNullSpace_Private(ce, snes)); 4623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 463900f6b5bSMatthew G. Knepley } 464900f6b5bSMatthew G. Knepley 465900f6b5bSMatthew G. Knepley /*@ 466900f6b5bSMatthew G. Knepley PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization 467900f6b5bSMatthew G. Knepley 46820f4b53cSBarry Smith Not Collective 469900f6b5bSMatthew G. Knepley 470900f6b5bSMatthew G. Knepley Input Parameter: 471f6dfbefdSBarry Smith . ce - The `PetscConvEst` object 472900f6b5bSMatthew G. Knepley 473900f6b5bSMatthew G. Knepley Output Parameter: 474900f6b5bSMatthew G. Knepley . alpha - The convergence rate for each field 475900f6b5bSMatthew G. Knepley 47620f4b53cSBarry Smith Options Database Keys: 477f6dfbefdSBarry Smith + -snes_convergence_estimate - Execute convergence estimation inside `SNESSolve()` and print out the rate 478f6dfbefdSBarry Smith - -ts_convergence_estimate - Execute convergence estimation inside `TSSolve()` and print out the rate 479f6dfbefdSBarry Smith 48020f4b53cSBarry Smith Level: intermediate 48120f4b53cSBarry Smith 482f6dfbefdSBarry Smith Notes: 483f6dfbefdSBarry Smith The convergence rate alpha is defined by 484b44f4de4SBarry Smith 485b44f4de4SBarry Smith $$ 486b44f4de4SBarry Smith || u_\Delta - u_{exact} || < C \Delta^\alpha 487b44f4de4SBarry Smith $$ 488b44f4de4SBarry Smith 489b44f4de4SBarry Smith where $u_{\Delta} $ is the discrete solution, and $\Delta$ is a measure of the discretization size. We usually use $h$ for the 490b44f4de4SBarry Smith spatial resolution and $\Delta t $ for the temporal resolution. 491900f6b5bSMatthew G. Knepley 492900f6b5bSMatthew G. Knepley We solve a series of problems using increasing resolution (refined meshes or decreased timesteps), calculate an error 49320f4b53cSBarry Smith based upon the exact solution in the `PetscDS`, and then fit the result to our model above using linear regression. 494900f6b5bSMatthew G. Knepley 49542747ad1SJacob Faibussowitsch .seealso: `PetscConvEstSetSolver()`, `PetscConvEstCreate()`, `SNESSolve()`, `TSSolve()` 496900f6b5bSMatthew G. Knepley @*/ 497d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[]) 498d71ae5a4SJacob Faibussowitsch { 499900f6b5bSMatthew G. Knepley PetscInt f; 500900f6b5bSMatthew G. Knepley 501900f6b5bSMatthew G. Knepley PetscFunctionBegin; 5029566063dSJacob Faibussowitsch if (ce->event < 0) PetscCall(PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &ce->event)); 503900f6b5bSMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0; 504dbbe0bcdSBarry Smith PetscUseTypeMethod(ce, getconvrate, alpha); 5053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5066af0ca60SMatthew G. Knepley } 5076af0ca60SMatthew G. Knepley 5086af0ca60SMatthew G. Knepley /*@ 509420bcc1bSBarry Smith PetscConvEstRateView - Displays the convergence rate obtained from `PetscConvEstGetConvRate()` using a `PetscViewer` 5106af0ca60SMatthew G. Knepley 511c3339decSBarry Smith Collective 5126af0ca60SMatthew G. Knepley 513e4094ef1SJacob Faibussowitsch Input Parameters: 514e4094ef1SJacob Faibussowitsch + ce - iterative context obtained from `SNESCreate()` 51546079b62SMatthew G. Knepley . alpha - the convergence rate for each field 5166af0ca60SMatthew G. Knepley - viewer - the viewer to display the reason 5176af0ca60SMatthew G. Knepley 51820f4b53cSBarry Smith Options Database Key: 5196af0ca60SMatthew G. Knepley . -snes_convergence_estimate - print the convergence rate 5206af0ca60SMatthew G. Knepley 5216af0ca60SMatthew G. Knepley Level: developer 5226af0ca60SMatthew G. Knepley 523420bcc1bSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstGetConvRate()` 5246af0ca60SMatthew G. Knepley @*/ 525d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstRateView(PetscConvEst ce, const PetscReal alpha[], PetscViewer viewer) 526d71ae5a4SJacob Faibussowitsch { 5276af0ca60SMatthew G. Knepley PetscBool isAscii; 5286af0ca60SMatthew G. Knepley 5296af0ca60SMatthew G. Knepley PetscFunctionBegin; 5309566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isAscii)); 5316af0ca60SMatthew G. Knepley if (isAscii) { 532900f6b5bSMatthew G. Knepley PetscInt Nf = ce->Nf, f; 533a56d3bf6SMatthew G. Knepley 5349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ce)->tablevel)); 5359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: ")); 5369566063dSJacob Faibussowitsch if (Nf > 1) PetscCall(PetscViewerASCIIPrintf(viewer, "[")); 537a56d3bf6SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 5389566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", ")); 5399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%#.2g", (double)alpha[f])); 54046079b62SMatthew G. Knepley } 5419566063dSJacob Faibussowitsch if (Nf > 1) PetscCall(PetscViewerASCIIPrintf(viewer, "]")); 5429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 5439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ce)->tablevel)); 5446af0ca60SMatthew G. Knepley } 5453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5466af0ca60SMatthew G. Knepley } 547900f6b5bSMatthew G. Knepley 548900f6b5bSMatthew G. Knepley /*@ 54926a11704SBarry Smith PetscConvEstCreate - Create a `PetscConvEst` object. This is used to study the convergence rate of approximations on grids to a continuum solution 550900f6b5bSMatthew G. Knepley 551900f6b5bSMatthew G. Knepley Collective 552900f6b5bSMatthew G. Knepley 553900f6b5bSMatthew G. Knepley Input Parameter: 554f6dfbefdSBarry Smith . comm - The communicator for the `PetscConvEst` object 555900f6b5bSMatthew G. Knepley 556900f6b5bSMatthew G. Knepley Output Parameter: 557f6dfbefdSBarry Smith . ce - The `PetscConvEst` object 558900f6b5bSMatthew G. Knepley 559900f6b5bSMatthew G. Knepley Level: beginner 560900f6b5bSMatthew G. Knepley 561420bcc1bSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstDestroy()`, `PetscConvEstGetConvRate()`, `DMAdaptorCreate()`, `DMAdaptor` 562900f6b5bSMatthew G. Knepley @*/ 563d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce) 564d71ae5a4SJacob Faibussowitsch { 565900f6b5bSMatthew G. Knepley PetscFunctionBegin; 5664f572ea9SToby Isaac PetscAssertPointer(ce, 2); 5679566063dSJacob Faibussowitsch PetscCall(PetscSysInitializePackage()); 5689566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView)); 569900f6b5bSMatthew G. Knepley (*ce)->monitor = PETSC_FALSE; 5702e61be88SMatthew G. Knepley (*ce)->r = 2.0; 571900f6b5bSMatthew G. Knepley (*ce)->Nr = 4; 572900f6b5bSMatthew G. Knepley (*ce)->event = -1; 573900f6b5bSMatthew G. Knepley (*ce)->ops->setsolver = PetscConvEstSetSNES_Private; 574900f6b5bSMatthew G. Knepley (*ce)->ops->initguess = PetscConvEstInitGuessSNES_Private; 575900f6b5bSMatthew G. Knepley (*ce)->ops->computeerror = PetscConvEstComputeErrorSNES_Private; 576900f6b5bSMatthew G. Knepley (*ce)->ops->getconvrate = PetscConvEstGetConvRateSNES_Private; 5773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 578900f6b5bSMatthew G. Knepley } 579