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 7d71ae5a4SJacob Faibussowitsch static PetscErrorCode zero_private(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt Nc, PetscScalar *u, void *ctx) 8d71ae5a4SJacob Faibussowitsch { 96af0ca60SMatthew G. Knepley PetscInt c; 106af0ca60SMatthew G. Knepley for (c = 0; c < Nc; ++c) u[c] = 0.0; 113ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 126af0ca60SMatthew G. Knepley } 136af0ca60SMatthew G. Knepley 146af0ca60SMatthew G. Knepley /*@ 156af0ca60SMatthew G. Knepley PetscConvEstDestroy - Destroys a PetscConvEst object 166af0ca60SMatthew G. Knepley 17c3339decSBarry Smith Collective 186af0ca60SMatthew G. Knepley 196af0ca60SMatthew G. Knepley Input Parameter: 2020f4b53cSBarry Smith . ce - The `PetscConvEst` object 216af0ca60SMatthew G. Knepley 226af0ca60SMatthew G. Knepley Level: beginner 236af0ca60SMatthew G. Knepley 2420f4b53cSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()` 256af0ca60SMatthew G. Knepley @*/ 26d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstDestroy(PetscConvEst *ce) 27d71ae5a4SJacob Faibussowitsch { 286af0ca60SMatthew G. Knepley PetscFunctionBegin; 293ba16761SJacob Faibussowitsch if (!*ce) PetscFunctionReturn(PETSC_SUCCESS); 306af0ca60SMatthew G. Knepley PetscValidHeaderSpecific((*ce), PETSC_OBJECT_CLASSID, 1); 316af0ca60SMatthew G. Knepley if (--((PetscObject)(*ce))->refct > 0) { 326af0ca60SMatthew G. Knepley *ce = NULL; 333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 346af0ca60SMatthew G. Knepley } 359566063dSJacob Faibussowitsch PetscCall(PetscFree3((*ce)->initGuess, (*ce)->exactSol, (*ce)->ctxs)); 369566063dSJacob Faibussowitsch PetscCall(PetscFree2((*ce)->dofs, (*ce)->errors)); 379566063dSJacob Faibussowitsch PetscCall(PetscHeaderDestroy(ce)); 383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 396af0ca60SMatthew G. Knepley } 406af0ca60SMatthew G. Knepley 416af0ca60SMatthew G. Knepley /*@ 4220f4b53cSBarry Smith PetscConvEstSetFromOptions - Sets a `PetscConvEst` object based on values in the options database 436af0ca60SMatthew G. Knepley 44c3339decSBarry Smith Collective 456af0ca60SMatthew G. Knepley 462fe279fdSBarry Smith Input Parameter: 4720f4b53cSBarry Smith . ce - The `PetscConvEst` object 486af0ca60SMatthew G. Knepley 496af0ca60SMatthew G. Knepley Level: beginner 506af0ca60SMatthew G. Knepley 5120f4b53cSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()` 526af0ca60SMatthew G. Knepley @*/ 53d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstSetFromOptions(PetscConvEst ce) 54d71ae5a4SJacob Faibussowitsch { 556af0ca60SMatthew G. Knepley PetscFunctionBegin; 56d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)ce), "", "Convergence Estimator Options", "PetscConvEst"); 579566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-convest_num_refine", "The number of refinements for the convergence check", "PetscConvEst", ce->Nr, &ce->Nr, NULL)); 589566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-convest_refine_factor", "The increase in resolution in each dimension", "PetscConvEst", ce->r, &ce->r, NULL)); 599566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-convest_monitor", "Monitor the error for each convergence check", "PetscConvEst", ce->monitor, &ce->monitor, NULL)); 609566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-convest_no_refine", "Debugging flag to run on the same mesh each time", "PetscConvEst", ce->noRefine, &ce->noRefine, NULL)); 61d0609cedSBarry Smith PetscOptionsEnd(); 623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 636af0ca60SMatthew G. Knepley } 646af0ca60SMatthew G. Knepley 656af0ca60SMatthew G. Knepley /*@ 6620f4b53cSBarry Smith PetscConvEstView - Views a `PetscConvEst` object 676af0ca60SMatthew G. Knepley 68c3339decSBarry Smith Collective 696af0ca60SMatthew G. Knepley 706af0ca60SMatthew G. Knepley Input Parameters: 7120f4b53cSBarry Smith + ce - The `PetscConvEst` object 7220f4b53cSBarry Smith - viewer - The `PetscViewer` object 736af0ca60SMatthew G. Knepley 746af0ca60SMatthew G. Knepley Level: beginner 756af0ca60SMatthew G. Knepley 7620f4b53cSBarry Smith .seealso: `PetscConvEst`, `PetscViewer`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()` 776af0ca60SMatthew G. Knepley @*/ 78d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstView(PetscConvEst ce, PetscViewer viewer) 79d71ae5a4SJacob Faibussowitsch { 806af0ca60SMatthew G. Knepley PetscFunctionBegin; 819566063dSJacob Faibussowitsch PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)ce, viewer)); 8263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "ConvEst with %" PetscInt_FMT " levels\n", ce->Nr + 1)); 833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 846af0ca60SMatthew G. Knepley } 856af0ca60SMatthew G. Knepley 866af0ca60SMatthew G. Knepley /*@ 876af0ca60SMatthew G. Knepley PetscConvEstGetSolver - Gets the solver used to produce discrete solutions 886af0ca60SMatthew G. Knepley 8920f4b53cSBarry Smith Not Collective 906af0ca60SMatthew G. Knepley 916af0ca60SMatthew G. Knepley Input Parameter: 9220f4b53cSBarry Smith . ce - The `PetscConvEst` object 936af0ca60SMatthew G. Knepley 946af0ca60SMatthew G. Knepley Output Parameter: 95900f6b5bSMatthew G. Knepley . solver - The solver 966af0ca60SMatthew G. Knepley 976af0ca60SMatthew G. Knepley Level: intermediate 986af0ca60SMatthew G. Knepley 9920f4b53cSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstSetSolver()`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()` 1006af0ca60SMatthew G. Knepley @*/ 101d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstGetSolver(PetscConvEst ce, PetscObject *solver) 102d71ae5a4SJacob Faibussowitsch { 1036af0ca60SMatthew G. Knepley PetscFunctionBegin; 1046af0ca60SMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 105900f6b5bSMatthew G. Knepley PetscValidPointer(solver, 2); 106900f6b5bSMatthew G. Knepley *solver = ce->solver; 1073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1086af0ca60SMatthew G. Knepley } 1096af0ca60SMatthew G. Knepley 1106af0ca60SMatthew G. Knepley /*@ 1116af0ca60SMatthew G. Knepley PetscConvEstSetSolver - Sets the solver used to produce discrete solutions 1126af0ca60SMatthew G. Knepley 11320f4b53cSBarry Smith Not Collective 1146af0ca60SMatthew G. Knepley 1156af0ca60SMatthew G. Knepley Input Parameters: 11620f4b53cSBarry Smith + ce - The `PetscConvEst` object 117900f6b5bSMatthew G. Knepley - solver - The solver 1186af0ca60SMatthew G. Knepley 1196af0ca60SMatthew G. Knepley Level: intermediate 1206af0ca60SMatthew G. Knepley 12120f4b53cSBarry Smith Note: 12220f4b53cSBarry Smith The solver MUST have an attached `DM`/`DS`, so that we know the exact solution 1236af0ca60SMatthew G. Knepley 12420f4b53cSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstGetSNES()`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()` 1256af0ca60SMatthew G. Knepley @*/ 126d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstSetSolver(PetscConvEst ce, PetscObject solver) 127d71ae5a4SJacob Faibussowitsch { 1286af0ca60SMatthew G. Knepley PetscFunctionBegin; 1296af0ca60SMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 130900f6b5bSMatthew G. Knepley PetscValidHeader(solver, 2); 131900f6b5bSMatthew G. Knepley ce->solver = solver; 132dbbe0bcdSBarry Smith PetscUseTypeMethod(ce, setsolver, solver); 1333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1346af0ca60SMatthew G. Knepley } 1356af0ca60SMatthew G. Knepley 1366af0ca60SMatthew G. Knepley /*@ 13720f4b53cSBarry Smith PetscConvEstSetUp - After the solver is specified, create structures for estimating convergence 1386af0ca60SMatthew G. Knepley 139c3339decSBarry Smith Collective 1406af0ca60SMatthew G. Knepley 1412fe279fdSBarry Smith Input Parameter: 14220f4b53cSBarry Smith . ce - The `PetscConvEst` object 1436af0ca60SMatthew G. Knepley 1446af0ca60SMatthew G. Knepley Level: beginner 1456af0ca60SMatthew G. Knepley 14620f4b53cSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()` 1476af0ca60SMatthew G. Knepley @*/ 148d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstSetUp(PetscConvEst ce) 149d71ae5a4SJacob Faibussowitsch { 150083401c6SMatthew G. Knepley PetscInt Nf, f, Nds, s; 1516af0ca60SMatthew G. Knepley 1526af0ca60SMatthew G. Knepley PetscFunctionBegin; 1539566063dSJacob Faibussowitsch PetscCall(DMGetNumFields(ce->idm, &Nf)); 154900f6b5bSMatthew G. Knepley ce->Nf = PetscMax(Nf, 1); 1559566063dSJacob Faibussowitsch PetscCall(PetscMalloc2((ce->Nr + 1) * ce->Nf, &ce->dofs, (ce->Nr + 1) * ce->Nf, &ce->errors)); 1569566063dSJacob Faibussowitsch PetscCall(PetscCalloc3(ce->Nf, &ce->initGuess, ce->Nf, &ce->exactSol, ce->Nf, &ce->ctxs)); 157900f6b5bSMatthew G. Knepley for (f = 0; f < Nf; ++f) ce->initGuess[f] = zero_private; 1589566063dSJacob Faibussowitsch PetscCall(DMGetNumDS(ce->idm, &Nds)); 159083401c6SMatthew G. Knepley for (s = 0; s < Nds; ++s) { 160083401c6SMatthew G. Knepley PetscDS ds; 161083401c6SMatthew G. Knepley DMLabel label; 162083401c6SMatthew G. Knepley IS fieldIS; 163083401c6SMatthew G. Knepley const PetscInt *fields; 164083401c6SMatthew G. Knepley PetscInt dsNf; 165083401c6SMatthew G. Knepley 16607218a29SMatthew G. Knepley PetscCall(DMGetRegionNumDS(ce->idm, s, &label, &fieldIS, &ds, NULL)); 1679566063dSJacob Faibussowitsch PetscCall(PetscDSGetNumFields(ds, &dsNf)); 1689566063dSJacob Faibussowitsch if (fieldIS) PetscCall(ISGetIndices(fieldIS, &fields)); 169083401c6SMatthew G. Knepley for (f = 0; f < dsNf; ++f) { 170083401c6SMatthew G. Knepley const PetscInt field = fields[f]; 1719566063dSJacob Faibussowitsch PetscCall(PetscDSGetExactSolution(ds, field, &ce->exactSol[field], &ce->ctxs[field])); 172083401c6SMatthew G. Knepley } 1739566063dSJacob Faibussowitsch if (fieldIS) PetscCall(ISRestoreIndices(fieldIS, &fields)); 174083401c6SMatthew G. Knepley } 175ad540459SPierre 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); 1763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1776af0ca60SMatthew G. Knepley } 1786af0ca60SMatthew G. Knepley 179d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstComputeInitialGuess(PetscConvEst ce, PetscInt r, DM dm, Vec u) 180d71ae5a4SJacob Faibussowitsch { 1816af0ca60SMatthew G. Knepley PetscFunctionBegin; 182900f6b5bSMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 183900f6b5bSMatthew G. Knepley if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3); 184900f6b5bSMatthew G. Knepley PetscValidHeaderSpecific(u, VEC_CLASSID, 4); 185dbbe0bcdSBarry Smith PetscUseTypeMethod(ce, initguess, r, dm, u); 1863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 187900f6b5bSMatthew G. Knepley } 188900f6b5bSMatthew G. Knepley 189d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstComputeError(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[]) 190d71ae5a4SJacob Faibussowitsch { 191900f6b5bSMatthew G. Knepley PetscFunctionBegin; 192900f6b5bSMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 193900f6b5bSMatthew G. Knepley if (dm) PetscValidHeaderSpecific(dm, DM_CLASSID, 3); 194900f6b5bSMatthew G. Knepley PetscValidHeaderSpecific(u, VEC_CLASSID, 4); 195900f6b5bSMatthew G. Knepley PetscValidRealPointer(errors, 5); 196dbbe0bcdSBarry Smith PetscUseTypeMethod(ce, computeerror, r, dm, u, errors); 1973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 198900f6b5bSMatthew G. Knepley } 199900f6b5bSMatthew G. Knepley 200f2cacb80SMatthew G. Knepley /*@ 201f2cacb80SMatthew G. Knepley PetscConvEstMonitorDefault - Monitors the convergence estimation loop 202f2cacb80SMatthew G. Knepley 203c3339decSBarry Smith Collective 204f2cacb80SMatthew G. Knepley 205d8d19677SJose E. Roman Input Parameters: 206f6dfbefdSBarry Smith + ce - The `PetscConvEst` object 207f2cacb80SMatthew G. Knepley - r - The refinement level 208f2cacb80SMatthew G. Knepley 209*e4094ef1SJacob Faibussowitsch Options Database Keys: 210ee300463SSatish Balay . -convest_monitor - Activate the monitor 211f2cacb80SMatthew G. Knepley 212f2cacb80SMatthew G. Knepley Level: intermediate 213f2cacb80SMatthew G. Knepley 214f6dfbefdSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`, `SNESSolve()`, `TSSolve()` 215f2cacb80SMatthew G. Knepley @*/ 216d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstMonitorDefault(PetscConvEst ce, PetscInt r) 217d71ae5a4SJacob Faibussowitsch { 218900f6b5bSMatthew G. Knepley MPI_Comm comm; 219900f6b5bSMatthew G. Knepley PetscInt f; 220900f6b5bSMatthew G. Knepley 221900f6b5bSMatthew G. Knepley PetscFunctionBegin; 222900f6b5bSMatthew G. Knepley if (ce->monitor) { 2237809adefSMatthew G. Knepley PetscInt *dofs = &ce->dofs[r * ce->Nf]; 224900f6b5bSMatthew G. Knepley PetscReal *errors = &ce->errors[r * ce->Nf]; 225900f6b5bSMatthew G. Knepley 2269566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)ce, &comm)); 2279566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "N: ")); 2289566063dSJacob Faibussowitsch if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "[")); 2297809adefSMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 2309566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscPrintf(comm, ", ")); 23163a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "%7" PetscInt_FMT, dofs[f])); 2327809adefSMatthew G. Knepley } 2339566063dSJacob Faibussowitsch if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "]")); 2349566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, " ")); 2359566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "L_2 Error: ")); 2369566063dSJacob Faibussowitsch if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "[")); 237900f6b5bSMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 2389566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscPrintf(comm, ", ")); 2399566063dSJacob Faibussowitsch if (errors[f] < 1.0e-11) PetscCall(PetscPrintf(comm, "< 1e-11")); 2409566063dSJacob Faibussowitsch else PetscCall(PetscPrintf(comm, "%g", (double)errors[f])); 241900f6b5bSMatthew G. Knepley } 2429566063dSJacob Faibussowitsch if (ce->Nf > 1) PetscCall(PetscPrintf(comm, "]")); 2439566063dSJacob Faibussowitsch PetscCall(PetscPrintf(comm, "\n")); 244900f6b5bSMatthew G. Knepley } 2453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 246900f6b5bSMatthew G. Knepley } 247900f6b5bSMatthew G. Knepley 248d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstSetSNES_Private(PetscConvEst ce, PetscObject solver) 249d71ae5a4SJacob Faibussowitsch { 250900f6b5bSMatthew G. Knepley PetscClassId id; 251900f6b5bSMatthew G. Knepley 252900f6b5bSMatthew G. Knepley PetscFunctionBegin; 2539566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(ce->solver, &id)); 25408401ef6SPierre Jolivet PetscCheck(id == SNES_CLASSID, PetscObjectComm((PetscObject)ce), PETSC_ERR_ARG_WRONG, "Solver was not a SNES"); 2559566063dSJacob Faibussowitsch PetscCall(SNESGetDM((SNES)ce->solver, &ce->idm)); 2563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 257900f6b5bSMatthew G. Knepley } 258900f6b5bSMatthew G. Knepley 259d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstInitGuessSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u) 260d71ae5a4SJacob Faibussowitsch { 261900f6b5bSMatthew G. Knepley PetscFunctionBegin; 2629566063dSJacob Faibussowitsch PetscCall(DMProjectFunction(dm, 0.0, ce->initGuess, ce->ctxs, INSERT_VALUES, u)); 2633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 264900f6b5bSMatthew G. Knepley } 265900f6b5bSMatthew G. Knepley 266d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstComputeErrorSNES_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[]) 267d71ae5a4SJacob Faibussowitsch { 268900f6b5bSMatthew G. Knepley PetscFunctionBegin; 2699566063dSJacob Faibussowitsch PetscCall(DMComputeL2FieldDiff(dm, 0.0, ce->exactSol, ce->ctxs, u, errors)); 2703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 271900f6b5bSMatthew G. Knepley } 272900f6b5bSMatthew G. Knepley 2738c401167SJed Brown static PetscErrorCode PetscConvEstSetJacobianNullSpace_Private(PetscConvEst ce, SNES snes) 274d71ae5a4SJacob Faibussowitsch { 275478db826SMatthew G. Knepley DM dm; 276478db826SMatthew G. Knepley PetscInt f; 277478db826SMatthew G. Knepley 278478db826SMatthew G. Knepley PetscFunctionBegin; 2799566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 280478db826SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 281478db826SMatthew G. Knepley PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *); 282478db826SMatthew G. Knepley 2839566063dSJacob Faibussowitsch PetscCall(DMGetNullSpaceConstructor(dm, f, &nspconstr)); 284478db826SMatthew G. Knepley if (nspconstr) { 285478db826SMatthew G. Knepley MatNullSpace nullsp; 286478db826SMatthew G. Knepley Mat J; 287478db826SMatthew G. Knepley 2889566063dSJacob Faibussowitsch PetscCall((*nspconstr)(dm, f, f, &nullsp)); 2899566063dSJacob Faibussowitsch PetscCall(SNESSetUp(snes)); 2909566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, &J, NULL, NULL, NULL)); 2919566063dSJacob Faibussowitsch PetscCall(MatSetNullSpace(J, nullsp)); 2929566063dSJacob Faibussowitsch PetscCall(MatNullSpaceDestroy(&nullsp)); 293478db826SMatthew G. Knepley break; 294478db826SMatthew G. Knepley } 295478db826SMatthew G. Knepley } 2963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 297478db826SMatthew G. Knepley } 298478db826SMatthew G. Knepley 299d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstGetConvRateSNES_Private(PetscConvEst ce, PetscReal alpha[]) 300d71ae5a4SJacob Faibussowitsch { 301900f6b5bSMatthew G. Knepley SNES snes = (SNES)ce->solver; 302900f6b5bSMatthew G. Knepley DM *dm; 303900f6b5bSMatthew G. Knepley PetscObject disc; 304900f6b5bSMatthew G. Knepley PetscReal *x, *y, slope, intercept; 3057809adefSMatthew G. Knepley PetscInt Nr = ce->Nr, r, f, dim, oldlevel, oldnlev; 306900f6b5bSMatthew G. Knepley void *ctx; 307900f6b5bSMatthew G. Knepley 308900f6b5bSMatthew G. Knepley PetscFunctionBegin; 30908401ef6SPierre Jolivet PetscCheck(ce->r == 2.0, PetscObjectComm((PetscObject)ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double)ce->r); 3109566063dSJacob Faibussowitsch PetscCall(DMGetDimension(ce->idm, &dim)); 3119566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(ce->idm, &ctx)); 3129566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE)); 3139566063dSJacob Faibussowitsch PetscCall(DMGetRefineLevel(ce->idm, &oldlevel)); 3149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((Nr + 1), &dm)); 3156af0ca60SMatthew G. Knepley /* Loop over meshes */ 316900f6b5bSMatthew G. Knepley dm[0] = ce->idm; 3176af0ca60SMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 318900f6b5bSMatthew G. Knepley Vec u; 319e5ed2c37SJose E. Roman #if defined(PETSC_USE_LOG) 320adfa7136SMatthew G. Knepley PetscLogStage stage; 321e5ed2c37SJose E. Roman #endif 322adfa7136SMatthew G. Knepley char stageName[PETSC_MAX_PATH_LEN]; 323900f6b5bSMatthew G. Knepley const char *dmname, *uname; 324adfa7136SMatthew G. Knepley 32563a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN - 1, "ConvEst Refinement Level %" PetscInt_FMT, r)); 326608e5a7aSMatthew G. Knepley #if defined(PETSC_USE_LOG) 3279566063dSJacob Faibussowitsch PetscCall(PetscLogStageGetId(stageName, &stage)); 3289566063dSJacob Faibussowitsch if (stage < 0) PetscCall(PetscLogStageRegister(stageName, &stage)); 329608e5a7aSMatthew G. Knepley #endif 3309566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 3316af0ca60SMatthew G. Knepley if (r > 0) { 332b2df8587SMatthew G. Knepley if (!ce->noRefine) { 3339566063dSJacob Faibussowitsch PetscCall(DMRefine(dm[r - 1], MPI_COMM_NULL, &dm[r])); 3349566063dSJacob Faibussowitsch PetscCall(DMSetCoarseDM(dm[r], dm[r - 1])); 335b2df8587SMatthew G. Knepley } else { 336b2df8587SMatthew G. Knepley DM cdm, rcdm; 337b2df8587SMatthew G. Knepley 3389566063dSJacob Faibussowitsch PetscCall(DMClone(dm[r - 1], &dm[r])); 3399566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm[r - 1], dm[r])); 3409566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm[r - 1], &cdm)); 3419566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm[r], &rcdm)); 3429566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(cdm, rcdm)); 343b2df8587SMatthew G. Knepley } 3449566063dSJacob Faibussowitsch PetscCall(DMCopyTransform(ce->idm, dm[r])); 3459566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)dm[r - 1], &dmname)); 3469566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dm[r], dmname)); 347478db826SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 3488cda7954SMatthew G. Knepley PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *); 349900f6b5bSMatthew G. Knepley 3509566063dSJacob Faibussowitsch PetscCall(DMGetNullSpaceConstructor(dm[r - 1], f, &nspconstr)); 3519566063dSJacob Faibussowitsch PetscCall(DMSetNullSpaceConstructor(dm[r], f, nspconstr)); 3528f892730SMatthew G. Knepley } 3536af0ca60SMatthew G. Knepley } 3549566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm[r], NULL, "-conv_dm_view")); 3556af0ca60SMatthew G. Knepley /* Create solution */ 3569566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm[r], &u)); 3579566063dSJacob Faibussowitsch PetscCall(DMGetField(dm[r], 0, NULL, &disc)); 3589566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName(disc, &uname)); 3599566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, uname)); 3606af0ca60SMatthew G. Knepley /* Setup solver */ 3619566063dSJacob Faibussowitsch PetscCall(SNESReset(snes)); 3629566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, dm[r])); 3639566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(dm[r], ctx, ctx, ctx)); 3649566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 365478db826SMatthew G. Knepley /* Set nullspace for Jacobian */ 3668c401167SJed Brown PetscCall(PetscConvEstSetJacobianNullSpace_Private(ce, snes)); 3676af0ca60SMatthew G. Knepley /* Create initial guess */ 3689566063dSJacob Faibussowitsch PetscCall(PetscConvEstComputeInitialGuess(ce, r, dm[r], u)); 3699566063dSJacob Faibussowitsch PetscCall(SNESSolve(snes, NULL, u)); 3709566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ce->event, ce, 0, 0, 0)); 3719566063dSJacob Faibussowitsch PetscCall(PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r * ce->Nf])); 3729566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ce->event, ce, 0, 0, 0)); 373adfa7136SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 374a56d3bf6SMatthew G. Knepley PetscSection s, fs; 375a56d3bf6SMatthew G. Knepley PetscInt lsize; 376a56d3bf6SMatthew G. Knepley 377a56d3bf6SMatthew G. Knepley /* Could use DMGetOutputDM() to add in Dirichlet dofs */ 3789566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm[r], &s)); 3799566063dSJacob Faibussowitsch PetscCall(PetscSectionGetField(s, f, &fs)); 3809566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstrainedStorageSize(fs, &lsize)); 381712fec58SPierre Jolivet PetscCall(MPIU_Allreduce(&lsize, &ce->dofs[r * ce->Nf + f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)snes))); 3829566063dSJacob Faibussowitsch PetscCall(PetscLogEventSetDof(ce->event, f, ce->dofs[r * ce->Nf + f])); 3839566063dSJacob Faibussowitsch PetscCall(PetscLogEventSetError(ce->event, f, ce->errors[r * ce->Nf + f])); 384adfa7136SMatthew G. Knepley } 3856af0ca60SMatthew G. Knepley /* Monitor */ 3869566063dSJacob Faibussowitsch PetscCall(PetscConvEstMonitorDefault(ce, r)); 38773269098SMatthew G. Knepley if (!r) { 38873269098SMatthew G. Knepley /* PCReset() does not wipe out the level structure */ 38973269098SMatthew G. Knepley KSP ksp; 39073269098SMatthew G. Knepley PC pc; 39173269098SMatthew G. Knepley 3929566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 3939566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 3949566063dSJacob Faibussowitsch PetscCall(PCMGGetLevels(pc, &oldnlev)); 39573269098SMatthew G. Knepley } 3966af0ca60SMatthew G. Knepley /* Cleanup */ 3979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&u)); 3989566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 3996af0ca60SMatthew G. Knepley } 40048a46eb9SPierre Jolivet for (r = 1; r <= Nr; ++r) PetscCall(DMDestroy(&dm[r])); 4016af0ca60SMatthew G. Knepley /* Fit convergence rate */ 4029566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nr + 1, &x, Nr + 1, &y)); 40346079b62SMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) { 4046af0ca60SMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 4057809adefSMatthew G. Knepley x[r] = PetscLog10Real(ce->dofs[r * ce->Nf + f]); 40646079b62SMatthew G. Knepley y[r] = PetscLog10Real(ce->errors[r * ce->Nf + f]); 4076af0ca60SMatthew G. Knepley } 4089566063dSJacob Faibussowitsch PetscCall(PetscLinearRegression(Nr + 1, x, y, &slope, &intercept)); 4096af0ca60SMatthew G. Knepley /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */ 41046079b62SMatthew G. Knepley alpha[f] = -slope * dim; 41146079b62SMatthew G. Knepley } 4129566063dSJacob Faibussowitsch PetscCall(PetscFree2(x, y)); 4139566063dSJacob Faibussowitsch PetscCall(PetscFree(dm)); 4142cae373cSMatthew G. Knepley /* Restore solver */ 4159566063dSJacob Faibussowitsch PetscCall(SNESReset(snes)); 41673269098SMatthew G. Knepley { 41773269098SMatthew G. Knepley /* PCReset() does not wipe out the level structure */ 41873269098SMatthew G. Knepley KSP ksp; 41973269098SMatthew G. Knepley PC pc; 42073269098SMatthew G. Knepley 4219566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 4229566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 4239566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels(pc, oldnlev, NULL)); 4249566063dSJacob Faibussowitsch PetscCall(DMSetRefineLevel(ce->idm, oldlevel)); /* The damn DMCoarsen() calls in PCMG can reset this */ 42573269098SMatthew G. Knepley } 4269566063dSJacob Faibussowitsch PetscCall(SNESSetDM(snes, ce->idm)); 4279566063dSJacob Faibussowitsch PetscCall(DMPlexSetSNESLocalFEM(ce->idm, ctx, ctx, ctx)); 4289566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 4298c401167SJed Brown PetscCall(PetscConvEstSetJacobianNullSpace_Private(ce, snes)); 4303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 431900f6b5bSMatthew G. Knepley } 432900f6b5bSMatthew G. Knepley 433900f6b5bSMatthew G. Knepley /*@ 434900f6b5bSMatthew G. Knepley PetscConvEstGetConvRate - Returns an estimate of the convergence rate for the discretization 435900f6b5bSMatthew G. Knepley 43620f4b53cSBarry Smith Not Collective 437900f6b5bSMatthew G. Knepley 438900f6b5bSMatthew G. Knepley Input Parameter: 439f6dfbefdSBarry Smith . ce - The `PetscConvEst` object 440900f6b5bSMatthew G. Knepley 441900f6b5bSMatthew G. Knepley Output Parameter: 442900f6b5bSMatthew G. Knepley . alpha - The convergence rate for each field 443900f6b5bSMatthew G. Knepley 44420f4b53cSBarry Smith Options Database Keys: 445f6dfbefdSBarry Smith + -snes_convergence_estimate - Execute convergence estimation inside `SNESSolve()` and print out the rate 446f6dfbefdSBarry Smith - -ts_convergence_estimate - Execute convergence estimation inside `TSSolve()` and print out the rate 447f6dfbefdSBarry Smith 44820f4b53cSBarry Smith Level: intermediate 44920f4b53cSBarry Smith 450f6dfbefdSBarry Smith Notes: 451f6dfbefdSBarry Smith The convergence rate alpha is defined by 452900f6b5bSMatthew G. Knepley $ || u_\Delta - u_exact || < C \Delta^alpha 453900f6b5bSMatthew G. Knepley where u_\Delta is the discrete solution, and Delta is a measure of the discretization size. We usually use h for the 454900f6b5bSMatthew G. Knepley spatial resolution and \Delta t for the temporal resolution. 455900f6b5bSMatthew G. Knepley 456900f6b5bSMatthew G. Knepley We solve a series of problems using increasing resolution (refined meshes or decreased timesteps), calculate an error 45720f4b53cSBarry Smith based upon the exact solution in the `PetscDS`, and then fit the result to our model above using linear regression. 458900f6b5bSMatthew G. Knepley 459db781477SPatrick Sanan .seealso: `PetscConvEstSetSolver()`, `PetscConvEstCreate()`, `PetscConvEstGetConvRate()`, `SNESSolve()`, `TSSolve()` 460900f6b5bSMatthew G. Knepley @*/ 461d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstGetConvRate(PetscConvEst ce, PetscReal alpha[]) 462d71ae5a4SJacob Faibussowitsch { 463900f6b5bSMatthew G. Knepley PetscInt f; 464900f6b5bSMatthew G. Knepley 465900f6b5bSMatthew G. Knepley PetscFunctionBegin; 4669566063dSJacob Faibussowitsch if (ce->event < 0) PetscCall(PetscLogEventRegister("ConvEst Error", PETSC_OBJECT_CLASSID, &ce->event)); 467900f6b5bSMatthew G. Knepley for (f = 0; f < ce->Nf; ++f) alpha[f] = 0.0; 468dbbe0bcdSBarry Smith PetscUseTypeMethod(ce, getconvrate, alpha); 4693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4706af0ca60SMatthew G. Knepley } 4716af0ca60SMatthew G. Knepley 4726af0ca60SMatthew G. Knepley /*@ 4736af0ca60SMatthew G. Knepley PetscConvEstRateView - Displays the convergence rate to a viewer 4746af0ca60SMatthew G. Knepley 475c3339decSBarry Smith Collective 4766af0ca60SMatthew G. Knepley 477*e4094ef1SJacob Faibussowitsch Input Parameters: 478*e4094ef1SJacob Faibussowitsch + ce - iterative context obtained from `SNESCreate()` 47946079b62SMatthew G. Knepley . alpha - the convergence rate for each field 4806af0ca60SMatthew G. Knepley - viewer - the viewer to display the reason 4816af0ca60SMatthew G. Knepley 48220f4b53cSBarry Smith Options Database Key: 4836af0ca60SMatthew G. Knepley . -snes_convergence_estimate - print the convergence rate 4846af0ca60SMatthew G. Knepley 4856af0ca60SMatthew G. Knepley Level: developer 4866af0ca60SMatthew G. Knepley 487f6dfbefdSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstGetRate()` 4886af0ca60SMatthew G. Knepley @*/ 489d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstRateView(PetscConvEst ce, const PetscReal alpha[], PetscViewer viewer) 490d71ae5a4SJacob Faibussowitsch { 4916af0ca60SMatthew G. Knepley PetscBool isAscii; 4926af0ca60SMatthew G. Knepley 4936af0ca60SMatthew G. Knepley PetscFunctionBegin; 4949566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isAscii)); 4956af0ca60SMatthew G. Knepley if (isAscii) { 496900f6b5bSMatthew G. Knepley PetscInt Nf = ce->Nf, f; 497a56d3bf6SMatthew G. Knepley 4989566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ce)->tablevel)); 4999566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "L_2 convergence rate: ")); 5009566063dSJacob Faibussowitsch if (Nf > 1) PetscCall(PetscViewerASCIIPrintf(viewer, "[")); 501a56d3bf6SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 5029566063dSJacob Faibussowitsch if (f > 0) PetscCall(PetscViewerASCIIPrintf(viewer, ", ")); 5039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%#.2g", (double)alpha[f])); 50446079b62SMatthew G. Knepley } 5059566063dSJacob Faibussowitsch if (Nf > 1) PetscCall(PetscViewerASCIIPrintf(viewer, "]")); 5069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 5079566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ce)->tablevel)); 5086af0ca60SMatthew G. Knepley } 5093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 5106af0ca60SMatthew G. Knepley } 511900f6b5bSMatthew G. Knepley 512900f6b5bSMatthew G. Knepley /*@ 513f6dfbefdSBarry Smith PetscConvEstCreate - Create a `PetscConvEst` object 514900f6b5bSMatthew G. Knepley 515900f6b5bSMatthew G. Knepley Collective 516900f6b5bSMatthew G. Knepley 517900f6b5bSMatthew G. Knepley Input Parameter: 518f6dfbefdSBarry Smith . comm - The communicator for the `PetscConvEst` object 519900f6b5bSMatthew G. Knepley 520900f6b5bSMatthew G. Knepley Output Parameter: 521f6dfbefdSBarry Smith . ce - The `PetscConvEst` object 522900f6b5bSMatthew G. Knepley 523900f6b5bSMatthew G. Knepley Level: beginner 524900f6b5bSMatthew G. Knepley 525f6dfbefdSBarry Smith .seealso: `PetscConvEst`, `PetscConvEstDestroy()`, `PetscConvEstGetConvRate()` 526900f6b5bSMatthew G. Knepley @*/ 527d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstCreate(MPI_Comm comm, PetscConvEst *ce) 528d71ae5a4SJacob Faibussowitsch { 529900f6b5bSMatthew G. Knepley PetscFunctionBegin; 530900f6b5bSMatthew G. Knepley PetscValidPointer(ce, 2); 5319566063dSJacob Faibussowitsch PetscCall(PetscSysInitializePackage()); 5329566063dSJacob Faibussowitsch PetscCall(PetscHeaderCreate(*ce, PETSC_OBJECT_CLASSID, "PetscConvEst", "ConvergenceEstimator", "SNES", comm, PetscConvEstDestroy, PetscConvEstView)); 533900f6b5bSMatthew G. Knepley (*ce)->monitor = PETSC_FALSE; 5342e61be88SMatthew G. Knepley (*ce)->r = 2.0; 535900f6b5bSMatthew G. Knepley (*ce)->Nr = 4; 536900f6b5bSMatthew G. Knepley (*ce)->event = -1; 537900f6b5bSMatthew G. Knepley (*ce)->ops->setsolver = PetscConvEstSetSNES_Private; 538900f6b5bSMatthew G. Knepley (*ce)->ops->initguess = PetscConvEstInitGuessSNES_Private; 539900f6b5bSMatthew G. Knepley (*ce)->ops->computeerror = PetscConvEstComputeErrorSNES_Private; 540900f6b5bSMatthew G. Knepley (*ce)->ops->getconvrate = PetscConvEstGetConvRateSNES_Private; 5413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 542900f6b5bSMatthew G. Knepley } 543