1900f6b5bSMatthew G. Knepley #include <petscconvest.h> /*I "petscconvest.h" I*/ 2900f6b5bSMatthew G. Knepley #include <petscts.h> 3f2ed2dc7SMatthew G. Knepley #include <petscdmplex.h> 4900f6b5bSMatthew G. Knepley 5900f6b5bSMatthew G. Knepley #include <petsc/private/petscconvestimpl.h> 6900f6b5bSMatthew G. Knepley 7d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstSetTS_Private(PetscConvEst ce, PetscObject solver) 8d71ae5a4SJacob Faibussowitsch { 9900f6b5bSMatthew G. Knepley PetscClassId id; 10900f6b5bSMatthew G. Knepley 11900f6b5bSMatthew G. Knepley PetscFunctionBegin; 129566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(ce->solver, &id)); 133c633725SBarry Smith PetscCheck(id == TS_CLASSID, PetscObjectComm((PetscObject)ce), PETSC_ERR_ARG_WRONG, "Solver was not a TS"); 149566063dSJacob Faibussowitsch PetscCall(TSGetDM((TS)ce->solver, &ce->idm)); 153ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16900f6b5bSMatthew G. Knepley } 17900f6b5bSMatthew G. Knepley 18d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstInitGuessTS_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u) 19d71ae5a4SJacob Faibussowitsch { 20900f6b5bSMatthew G. Knepley PetscFunctionBegin; 219566063dSJacob Faibussowitsch PetscCall(TSComputeInitialCondition((TS)ce->solver, u)); 223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 23900f6b5bSMatthew G. Knepley } 24900f6b5bSMatthew G. Knepley 25d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstComputeErrorTS_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[]) 26d71ae5a4SJacob Faibussowitsch { 27f2cacb80SMatthew G. Knepley TS ts = (TS)ce->solver; 28f2cacb80SMatthew G. Knepley PetscErrorCode (*exactError)(TS, Vec, Vec); 29900f6b5bSMatthew G. Knepley 30900f6b5bSMatthew G. Knepley PetscFunctionBegin; 319566063dSJacob Faibussowitsch PetscCall(TSGetComputeExactError(ts, &exactError)); 32f2cacb80SMatthew G. Knepley if (exactError) { 33f2cacb80SMatthew G. Knepley Vec e; 34f2cacb80SMatthew G. Knepley PetscInt f; 35f2cacb80SMatthew G. Knepley 369566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &e)); 379566063dSJacob Faibussowitsch PetscCall(TSComputeExactError(ts, u, e)); 389566063dSJacob Faibussowitsch PetscCall(VecNorm(e, NORM_2, errors)); 39f2cacb80SMatthew G. Knepley for (f = 1; f < ce->Nf; ++f) errors[f] = errors[0]; 409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&e)); 41f2cacb80SMatthew G. Knepley } else { 42f2cacb80SMatthew G. Knepley PetscReal t; 43f2cacb80SMatthew G. Knepley 449566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &t)); 459566063dSJacob Faibussowitsch PetscCall(DMComputeL2FieldDiff(dm, t, ce->exactSol, ce->ctxs, u, errors)); 46f2cacb80SMatthew G. Knepley } 473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 48900f6b5bSMatthew G. Knepley } 49900f6b5bSMatthew G. Knepley 50d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstGetConvRateTS_Temporal_Private(PetscConvEst ce, PetscReal alpha[]) 51d71ae5a4SJacob Faibussowitsch { 52900f6b5bSMatthew G. Knepley TS ts = (TS)ce->solver; 5349979724SMatthew G. Knepley Vec u, u0; 54900f6b5bSMatthew G. Knepley PetscReal *dt, *x, *y, slope, intercept; 55900f6b5bSMatthew G. Knepley PetscInt Ns, oNs, Nf = ce->Nf, f, Nr = ce->Nr, r; 56900f6b5bSMatthew G. Knepley 57900f6b5bSMatthew G. Knepley PetscFunctionBegin; 589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nr + 1, &dt)); 599566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt[0])); 609566063dSJacob Faibussowitsch PetscCall(TSGetMaxSteps(ts, &oNs)); 6149979724SMatthew G. Knepley PetscCall(TSGetSolution(ts, &u0)); 6249979724SMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)u0)); 63900f6b5bSMatthew G. Knepley Ns = oNs; 64900f6b5bSMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 65900f6b5bSMatthew G. Knepley if (r > 0) { 662e61be88SMatthew G. Knepley dt[r] = dt[r - 1] / ce->r; 672e61be88SMatthew G. Knepley Ns = PetscCeilReal(Ns * ce->r); 68900f6b5bSMatthew G. Knepley } 699566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 709566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 719566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt[r])); 729566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, Ns)); 7349979724SMatthew G. Knepley PetscCall(TSGetSolution(ts, &u)); 749566063dSJacob Faibussowitsch PetscCall(PetscConvEstComputeInitialGuess(ce, r, NULL, u)); 7549979724SMatthew G. Knepley PetscCall(TSSolve(ts, NULL)); 7649979724SMatthew G. Knepley PetscCall(TSGetSolution(ts, &u)); 779566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ce->event, ce, 0, 0, 0)); 789566063dSJacob Faibussowitsch PetscCall(PetscConvEstComputeError(ce, r, ce->idm, u, &ce->errors[r * Nf])); 799566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ce->event, ce, 0, 0, 0)); 80f2ed2dc7SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 819327cc11SMatthew G. Knepley ce->dofs[r * Nf + f] = 1.0 / dt[r]; 829566063dSJacob Faibussowitsch PetscCall(PetscLogEventSetDof(ce->event, f, ce->dofs[r * Nf + f])); 839566063dSJacob Faibussowitsch PetscCall(PetscLogEventSetError(ce->event, f, ce->errors[r * Nf + f])); 84900f6b5bSMatthew G. Knepley } 85f2cacb80SMatthew G. Knepley /* Monitor */ 869566063dSJacob Faibussowitsch PetscCall(PetscConvEstMonitorDefault(ce, r)); 87900f6b5bSMatthew G. Knepley } 88900f6b5bSMatthew G. Knepley /* Fit convergence rate */ 89900f6b5bSMatthew G. Knepley if (Nr) { 909566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nr + 1, &x, Nr + 1, &y)); 91900f6b5bSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 92900f6b5bSMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 93f2cacb80SMatthew G. Knepley x[r] = PetscLog10Real(dt[r]); 94900f6b5bSMatthew G. Knepley y[r] = PetscLog10Real(ce->errors[r * Nf + f]); 95900f6b5bSMatthew G. Knepley } 969566063dSJacob Faibussowitsch PetscCall(PetscLinearRegression(Nr + 1, x, y, &slope, &intercept)); 97900f6b5bSMatthew G. Knepley /* Since lg err = s lg dt + b */ 98900f6b5bSMatthew G. Knepley alpha[f] = slope; 99900f6b5bSMatthew G. Knepley } 1009566063dSJacob Faibussowitsch PetscCall(PetscFree2(x, y)); 101900f6b5bSMatthew G. Knepley } 102900f6b5bSMatthew G. Knepley /* Reset solver */ 10349979724SMatthew G. Knepley PetscCall(TSReset(ts)); 1049566063dSJacob Faibussowitsch PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_ITERATING)); 1059566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 1069566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 1079566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt[0])); 1089566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, oNs)); 10949979724SMatthew G. Knepley PetscCall(TSSetSolution(ts, u0)); 11049979724SMatthew G. Knepley PetscCall(PetscConvEstComputeInitialGuess(ce, 0, NULL, u0)); 11149979724SMatthew G. Knepley PetscCall(VecDestroy(&u0)); 1129566063dSJacob Faibussowitsch PetscCall(PetscFree(dt)); 1133ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 114900f6b5bSMatthew G. Knepley } 115900f6b5bSMatthew G. Knepley 116d71ae5a4SJacob Faibussowitsch static PetscErrorCode PetscConvEstGetConvRateTS_Spatial_Private(PetscConvEst ce, PetscReal alpha[]) 117d71ae5a4SJacob Faibussowitsch { 118f2ed2dc7SMatthew G. Knepley TS ts = (TS)ce->solver; 119f2ed2dc7SMatthew G. Knepley Vec uInitial; 120f2ed2dc7SMatthew G. Knepley DM *dm; 121f2ed2dc7SMatthew G. Knepley PetscObject disc; 122f2ed2dc7SMatthew G. Knepley PetscReal *x, *y, slope, intercept; 1237809adefSMatthew G. Knepley PetscInt Nr = ce->Nr, r, Nf = ce->Nf, f, dim, oldlevel, oldnlev; 124ae11cfe1SMatthew G. Knepley PetscErrorCode (*ifunc)(DM, PetscReal, Vec, Vec, Vec, void *); 125ae11cfe1SMatthew G. Knepley PetscErrorCode (*ijac)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 126ae11cfe1SMatthew G. Knepley PetscErrorCode (*rhsfunc)(DM, PetscReal, Vec, Vec, void *); 127*2a8381b2SBarry Smith PetscCtx fctx, jctx, rctx, ctx; 128f2ed2dc7SMatthew G. Knepley 129f2ed2dc7SMatthew G. Knepley PetscFunctionBegin; 1303c633725SBarry Smith PetscCheck(ce->r == 2.0, PetscObjectComm((PetscObject)ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double)ce->r); 1319566063dSJacob Faibussowitsch PetscCall(DMGetDimension(ce->idm, &dim)); 1329566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(ce->idm, &ctx)); 1339566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE)); 1349566063dSJacob Faibussowitsch PetscCall(DMGetRefineLevel(ce->idm, &oldlevel)); 13557508eceSPierre Jolivet PetscCall(PetscMalloc1(Nr + 1, &dm)); 1369566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &uInitial)); 137a7e8c07fSStefano Zampini PetscCall(PetscObjectReference((PetscObject)uInitial)); 138a7e8c07fSStefano Zampini 139f2ed2dc7SMatthew G. Knepley /* Loop over meshes */ 140f2ed2dc7SMatthew G. Knepley dm[0] = ce->idm; 141f2ed2dc7SMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 142f2ed2dc7SMatthew G. Knepley Vec u; 143f2ed2dc7SMatthew G. Knepley PetscLogStage stage; 144f2ed2dc7SMatthew G. Knepley char stageName[PETSC_MAX_PATH_LEN]; 145f2ed2dc7SMatthew G. Knepley const char *dmname, *uname; 146f2ed2dc7SMatthew G. Knepley 14763a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN - 1, "ConvEst Refinement Level %" PetscInt_FMT, r)); 1489566063dSJacob Faibussowitsch PetscCall(PetscLogStageGetId(stageName, &stage)); 1499566063dSJacob Faibussowitsch if (stage < 0) PetscCall(PetscLogStageRegister(stageName, &stage)); 1509566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 151f2ed2dc7SMatthew G. Knepley if (r > 0) { 152b2df8587SMatthew G. Knepley if (!ce->noRefine) { 1539566063dSJacob Faibussowitsch PetscCall(DMRefine(dm[r - 1], MPI_COMM_NULL, &dm[r])); 1549566063dSJacob Faibussowitsch PetscCall(DMSetCoarseDM(dm[r], dm[r - 1])); 155b2df8587SMatthew G. Knepley } else { 156b2df8587SMatthew G. Knepley DM cdm, rcdm; 157b2df8587SMatthew G. Knepley 1589566063dSJacob Faibussowitsch PetscCall(DMClone(dm[r - 1], &dm[r])); 1599566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm[r - 1], dm[r])); 1609566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm[r - 1], &cdm)); 1619566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm[r], &rcdm)); 1629566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(cdm, rcdm)); 163b2df8587SMatthew G. Knepley } 1649566063dSJacob Faibussowitsch PetscCall(DMCopyTransform(ce->idm, dm[r])); 1659566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)dm[r - 1], &dmname)); 1669566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dm[r], dmname)); 1674758e3ceSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 1688cda7954SMatthew G. Knepley PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *); 169f2ed2dc7SMatthew G. Knepley 1709566063dSJacob Faibussowitsch PetscCall(DMGetNullSpaceConstructor(dm[r - 1], f, &nspconstr)); 1719566063dSJacob Faibussowitsch PetscCall(DMSetNullSpaceConstructor(dm[r], f, nspconstr)); 172f2ed2dc7SMatthew G. Knepley } 173f2ed2dc7SMatthew G. Knepley } 1749566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm[r], NULL, "-conv_dm_view")); 175f2ed2dc7SMatthew G. Knepley /* Create solution */ 1769566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm[r], &u)); 1779566063dSJacob Faibussowitsch PetscCall(DMGetField(dm[r], 0, NULL, &disc)); 1789566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName(disc, &uname)); 1799566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, uname)); 180f2ed2dc7SMatthew G. Knepley /* Setup solver */ 1819566063dSJacob Faibussowitsch PetscCall(TSReset(ts)); 1829566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm[r])); 1839566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm[r], DMPlexTSComputeBoundary, ctx)); 184ae11cfe1SMatthew G. Knepley if (r > 0) { 185ae11cfe1SMatthew G. Knepley PetscCall(DMTSGetIFunctionLocal(dm[r - 1], &ifunc, &fctx)); 186ae11cfe1SMatthew G. Knepley PetscCall(DMTSGetIJacobianLocal(dm[r - 1], &ijac, &jctx)); 187ae11cfe1SMatthew G. Knepley PetscCall(DMTSGetRHSFunctionLocal(dm[r - 1], &rhsfunc, &rctx)); 18848a46eb9SPierre Jolivet if (ifunc) PetscCall(DMTSSetIFunctionLocal(dm[r], ifunc, fctx)); 18948a46eb9SPierre Jolivet if (ijac) PetscCall(DMTSSetIJacobianLocal(dm[r], ijac, jctx)); 19048a46eb9SPierre Jolivet if (rhsfunc) PetscCall(DMTSSetRHSFunctionLocal(dm[r], rhsfunc, rctx)); 191ae11cfe1SMatthew G. Knepley } 1929566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 1939566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 1949566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 19549979724SMatthew G. Knepley PetscCall(TSSetSolution(ts, u)); 19649979724SMatthew G. Knepley PetscCall(VecDestroy(&u)); 197f2ed2dc7SMatthew G. Knepley /* Create initial guess */ 19849979724SMatthew G. Knepley PetscCall(TSGetSolution(ts, &u)); 1999566063dSJacob Faibussowitsch PetscCall(PetscConvEstComputeInitialGuess(ce, r, dm[r], u)); 20049979724SMatthew G. Knepley PetscCall(TSSolve(ts, NULL)); 20149979724SMatthew G. Knepley PetscCall(TSGetSolution(ts, &u)); 2029566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ce->event, ce, 0, 0, 0)); 2039566063dSJacob Faibussowitsch PetscCall(PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r * Nf])); 2049566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ce->event, ce, 0, 0, 0)); 205f2ed2dc7SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 206f2ed2dc7SMatthew G. Knepley PetscSection s, fs; 207f2ed2dc7SMatthew G. Knepley PetscInt lsize; 208f2ed2dc7SMatthew G. Knepley 209f2ed2dc7SMatthew G. Knepley /* Could use DMGetOutputDM() to add in Dirichlet dofs */ 2109566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm[r], &s)); 2119566063dSJacob Faibussowitsch PetscCall(PetscSectionGetField(s, f, &fs)); 2129566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstrainedStorageSize(fs, &lsize)); 213462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&lsize, &ce->dofs[r * Nf + f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)ts))); 2149566063dSJacob Faibussowitsch PetscCall(PetscLogEventSetDof(ce->event, f, ce->dofs[r * Nf + f])); 2159566063dSJacob Faibussowitsch PetscCall(PetscLogEventSetError(ce->event, f, ce->errors[r * Nf + f])); 216f2ed2dc7SMatthew G. Knepley } 217f2ed2dc7SMatthew G. Knepley /* Monitor */ 2189566063dSJacob Faibussowitsch PetscCall(PetscConvEstMonitorDefault(ce, r)); 219f2ed2dc7SMatthew G. Knepley if (!r) { 220f2ed2dc7SMatthew G. Knepley /* PCReset() does not wipe out the level structure */ 221f2ed2dc7SMatthew G. Knepley SNES snes; 222f2ed2dc7SMatthew G. Knepley KSP ksp; 223f2ed2dc7SMatthew G. Knepley PC pc; 224f2ed2dc7SMatthew G. Knepley 2259566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 2269566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 2279566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 2289566063dSJacob Faibussowitsch PetscCall(PCMGGetLevels(pc, &oldnlev)); 229f2ed2dc7SMatthew G. Knepley } 230f2ed2dc7SMatthew G. Knepley /* Cleanup */ 2319566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 232f2ed2dc7SMatthew G. Knepley } 233ae11cfe1SMatthew G. Knepley PetscCall(DMTSGetIFunctionLocal(dm[r - 1], &ifunc, &fctx)); 234ae11cfe1SMatthew G. Knepley PetscCall(DMTSGetIJacobianLocal(dm[r - 1], &ijac, &jctx)); 235ae11cfe1SMatthew G. Knepley PetscCall(DMTSGetRHSFunctionLocal(dm[r - 1], &rhsfunc, &rctx)); 23648a46eb9SPierre Jolivet for (r = 1; r <= Nr; ++r) PetscCall(DMDestroy(&dm[r])); 237f2ed2dc7SMatthew G. Knepley /* Fit convergence rate */ 2389566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nr + 1, &x, Nr + 1, &y)); 239f2ed2dc7SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 240f2ed2dc7SMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 2417809adefSMatthew G. Knepley x[r] = PetscLog10Real(ce->dofs[r * Nf + f]); 242f2ed2dc7SMatthew G. Knepley y[r] = PetscLog10Real(ce->errors[r * Nf + f]); 243f2ed2dc7SMatthew G. Knepley } 2449566063dSJacob Faibussowitsch PetscCall(PetscLinearRegression(Nr + 1, x, y, &slope, &intercept)); 245f2ed2dc7SMatthew G. Knepley /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */ 246f2ed2dc7SMatthew G. Knepley alpha[f] = -slope * dim; 247f2ed2dc7SMatthew G. Knepley } 2489566063dSJacob Faibussowitsch PetscCall(PetscFree2(x, y)); 2499566063dSJacob Faibussowitsch PetscCall(PetscFree(dm)); 250f2ed2dc7SMatthew G. Knepley /* Restore solver */ 2519566063dSJacob Faibussowitsch PetscCall(TSReset(ts)); 252f2ed2dc7SMatthew G. Knepley { 253f2ed2dc7SMatthew G. Knepley /* PCReset() does not wipe out the level structure */ 254f2ed2dc7SMatthew G. Knepley SNES snes; 255f2ed2dc7SMatthew G. Knepley KSP ksp; 256f2ed2dc7SMatthew G. Knepley PC pc; 257f2ed2dc7SMatthew G. Knepley 2589566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 2599566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 2609566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 2619566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels(pc, oldnlev, NULL)); 2629566063dSJacob Faibussowitsch PetscCall(DMSetRefineLevel(ce->idm, oldlevel)); /* The damn DMCoarsen() calls in PCMG can reset this */ 263f2ed2dc7SMatthew G. Knepley } 2649566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, ce->idm)); 2659566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(ce->idm, DMPlexTSComputeBoundary, ctx)); 26648a46eb9SPierre Jolivet if (ifunc) PetscCall(DMTSSetIFunctionLocal(ce->idm, ifunc, fctx)); 26748a46eb9SPierre Jolivet if (ijac) PetscCall(DMTSSetIJacobianLocal(ce->idm, ijac, jctx)); 26848a46eb9SPierre Jolivet if (rhsfunc) PetscCall(DMTSSetRHSFunctionLocal(ce->idm, rhsfunc, rctx)); 2699566063dSJacob Faibussowitsch PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_ITERATING)); 2709566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 2719566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 2729566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 2739566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, uInitial)); 2749566063dSJacob Faibussowitsch PetscCall(PetscConvEstComputeInitialGuess(ce, 0, NULL, uInitial)); 275a7e8c07fSStefano Zampini PetscCall(VecDestroy(&uInitial)); 2763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 277f2ed2dc7SMatthew G. Knepley } 278f2ed2dc7SMatthew G. Knepley 279d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstUseTS(PetscConvEst ce, PetscBool checkTemporal) 280d71ae5a4SJacob Faibussowitsch { 281900f6b5bSMatthew G. Knepley PetscFunctionBegin; 282900f6b5bSMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 283900f6b5bSMatthew G. Knepley ce->ops->setsolver = PetscConvEstSetTS_Private; 284900f6b5bSMatthew G. Knepley ce->ops->initguess = PetscConvEstInitGuessTS_Private; 285900f6b5bSMatthew G. Knepley ce->ops->computeerror = PetscConvEstComputeErrorTS_Private; 286f2ed2dc7SMatthew G. Knepley if (checkTemporal) { 287f2ed2dc7SMatthew G. Knepley ce->ops->getconvrate = PetscConvEstGetConvRateTS_Temporal_Private; 288f2ed2dc7SMatthew G. Knepley } else { 289f2ed2dc7SMatthew G. Knepley ce->ops->getconvrate = PetscConvEstGetConvRateTS_Spatial_Private; 290f2ed2dc7SMatthew G. Knepley } 2913ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 292900f6b5bSMatthew G. Knepley } 293