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 79371c9d4SSatish Balay static PetscErrorCode PetscConvEstSetTS_Private(PetscConvEst ce, PetscObject solver) { 8900f6b5bSMatthew G. Knepley PetscClassId id; 9900f6b5bSMatthew G. Knepley 10900f6b5bSMatthew G. Knepley PetscFunctionBegin; 119566063dSJacob Faibussowitsch PetscCall(PetscObjectGetClassId(ce->solver, &id)); 123c633725SBarry Smith PetscCheck(id == TS_CLASSID, PetscObjectComm((PetscObject)ce), PETSC_ERR_ARG_WRONG, "Solver was not a TS"); 139566063dSJacob Faibussowitsch PetscCall(TSGetDM((TS)ce->solver, &ce->idm)); 14900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 15900f6b5bSMatthew G. Knepley } 16900f6b5bSMatthew G. Knepley 179371c9d4SSatish Balay static PetscErrorCode PetscConvEstInitGuessTS_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u) { 18900f6b5bSMatthew G. Knepley PetscFunctionBegin; 199566063dSJacob Faibussowitsch PetscCall(TSComputeInitialCondition((TS)ce->solver, u)); 20900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 21900f6b5bSMatthew G. Knepley } 22900f6b5bSMatthew G. Knepley 239371c9d4SSatish Balay static PetscErrorCode PetscConvEstComputeErrorTS_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[]) { 24f2cacb80SMatthew G. Knepley TS ts = (TS)ce->solver; 25f2cacb80SMatthew G. Knepley PetscErrorCode (*exactError)(TS, Vec, Vec); 26900f6b5bSMatthew G. Knepley 27900f6b5bSMatthew G. Knepley PetscFunctionBegin; 289566063dSJacob Faibussowitsch PetscCall(TSGetComputeExactError(ts, &exactError)); 29f2cacb80SMatthew G. Knepley if (exactError) { 30f2cacb80SMatthew G. Knepley Vec e; 31f2cacb80SMatthew G. Knepley PetscInt f; 32f2cacb80SMatthew G. Knepley 339566063dSJacob Faibussowitsch PetscCall(VecDuplicate(u, &e)); 349566063dSJacob Faibussowitsch PetscCall(TSComputeExactError(ts, u, e)); 359566063dSJacob Faibussowitsch PetscCall(VecNorm(e, NORM_2, errors)); 36f2cacb80SMatthew G. Knepley for (f = 1; f < ce->Nf; ++f) errors[f] = errors[0]; 379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&e)); 38f2cacb80SMatthew G. Knepley } else { 39f2cacb80SMatthew G. Knepley PetscReal t; 40f2cacb80SMatthew G. Knepley 419566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &t)); 429566063dSJacob Faibussowitsch PetscCall(DMComputeL2FieldDiff(dm, t, ce->exactSol, ce->ctxs, u, errors)); 43f2cacb80SMatthew G. Knepley } 44900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 45900f6b5bSMatthew G. Knepley } 46900f6b5bSMatthew G. Knepley 479371c9d4SSatish Balay static PetscErrorCode PetscConvEstGetConvRateTS_Temporal_Private(PetscConvEst ce, PetscReal alpha[]) { 48900f6b5bSMatthew G. Knepley TS ts = (TS)ce->solver; 4949979724SMatthew G. Knepley Vec u, u0; 50900f6b5bSMatthew G. Knepley PetscReal *dt, *x, *y, slope, intercept; 51900f6b5bSMatthew G. Knepley PetscInt Ns, oNs, Nf = ce->Nf, f, Nr = ce->Nr, r; 52900f6b5bSMatthew G. Knepley 53900f6b5bSMatthew G. Knepley PetscFunctionBegin; 549566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(Nr + 1, &dt)); 559566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt[0])); 569566063dSJacob Faibussowitsch PetscCall(TSGetMaxSteps(ts, &oNs)); 5749979724SMatthew G. Knepley PetscCall(TSGetSolution(ts, &u0)); 5849979724SMatthew G. Knepley PetscCall(PetscObjectReference((PetscObject)u0)); 59900f6b5bSMatthew G. Knepley Ns = oNs; 60900f6b5bSMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 61900f6b5bSMatthew G. Knepley if (r > 0) { 622e61be88SMatthew G. Knepley dt[r] = dt[r - 1] / ce->r; 632e61be88SMatthew G. Knepley Ns = PetscCeilReal(Ns * ce->r); 64900f6b5bSMatthew G. Knepley } 659566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 669566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 679566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt[r])); 689566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, Ns)); 6949979724SMatthew G. Knepley PetscCall(TSGetSolution(ts, &u)); 709566063dSJacob Faibussowitsch PetscCall(PetscConvEstComputeInitialGuess(ce, r, NULL, u)); 7149979724SMatthew G. Knepley PetscCall(TSSolve(ts, NULL)); 7249979724SMatthew G. Knepley PetscCall(TSGetSolution(ts, &u)); 739566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ce->event, ce, 0, 0, 0)); 749566063dSJacob Faibussowitsch PetscCall(PetscConvEstComputeError(ce, r, ce->idm, u, &ce->errors[r * Nf])); 759566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ce->event, ce, 0, 0, 0)); 76f2ed2dc7SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 779327cc11SMatthew G. Knepley ce->dofs[r * Nf + f] = 1.0 / dt[r]; 789566063dSJacob Faibussowitsch PetscCall(PetscLogEventSetDof(ce->event, f, ce->dofs[r * Nf + f])); 799566063dSJacob Faibussowitsch PetscCall(PetscLogEventSetError(ce->event, f, ce->errors[r * Nf + f])); 80900f6b5bSMatthew G. Knepley } 81f2cacb80SMatthew G. Knepley /* Monitor */ 829566063dSJacob Faibussowitsch PetscCall(PetscConvEstMonitorDefault(ce, r)); 83900f6b5bSMatthew G. Knepley } 84900f6b5bSMatthew G. Knepley /* Fit convergence rate */ 85900f6b5bSMatthew G. Knepley if (Nr) { 869566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nr + 1, &x, Nr + 1, &y)); 87900f6b5bSMatthew G. Knepley for (f = 0; f < Nf; ++f) { 88900f6b5bSMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 89f2cacb80SMatthew G. Knepley x[r] = PetscLog10Real(dt[r]); 90900f6b5bSMatthew G. Knepley y[r] = PetscLog10Real(ce->errors[r * Nf + f]); 91900f6b5bSMatthew G. Knepley } 929566063dSJacob Faibussowitsch PetscCall(PetscLinearRegression(Nr + 1, x, y, &slope, &intercept)); 93900f6b5bSMatthew G. Knepley /* Since lg err = s lg dt + b */ 94900f6b5bSMatthew G. Knepley alpha[f] = slope; 95900f6b5bSMatthew G. Knepley } 969566063dSJacob Faibussowitsch PetscCall(PetscFree2(x, y)); 97900f6b5bSMatthew G. Knepley } 98900f6b5bSMatthew G. Knepley /* Reset solver */ 9949979724SMatthew G. Knepley PetscCall(TSReset(ts)); 1009566063dSJacob Faibussowitsch PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_ITERATING)); 1019566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 1029566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 1039566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, dt[0])); 1049566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts, oNs)); 10549979724SMatthew G. Knepley PetscCall(TSSetSolution(ts, u0)); 10649979724SMatthew G. Knepley PetscCall(PetscConvEstComputeInitialGuess(ce, 0, NULL, u0)); 10749979724SMatthew G. Knepley PetscCall(VecDestroy(&u0)); 1089566063dSJacob Faibussowitsch PetscCall(PetscFree(dt)); 109900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 110900f6b5bSMatthew G. Knepley } 111900f6b5bSMatthew G. Knepley 1129371c9d4SSatish Balay static PetscErrorCode PetscConvEstGetConvRateTS_Spatial_Private(PetscConvEst ce, PetscReal alpha[]) { 113f2ed2dc7SMatthew G. Knepley TS ts = (TS)ce->solver; 114f2ed2dc7SMatthew G. Knepley Vec uInitial; 115f2ed2dc7SMatthew G. Knepley DM *dm; 116f2ed2dc7SMatthew G. Knepley PetscObject disc; 117f2ed2dc7SMatthew G. Knepley PetscReal *x, *y, slope, intercept; 1187809adefSMatthew G. Knepley PetscInt Nr = ce->Nr, r, Nf = ce->Nf, f, dim, oldlevel, oldnlev; 119ae11cfe1SMatthew G. Knepley PetscErrorCode (*ifunc)(DM, PetscReal, Vec, Vec, Vec, void *); 120ae11cfe1SMatthew G. Knepley PetscErrorCode (*ijac)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 121ae11cfe1SMatthew G. Knepley PetscErrorCode (*rhsfunc)(DM, PetscReal, Vec, Vec, void *); 122ae11cfe1SMatthew G. Knepley void *fctx, *jctx, *rctx; 123f2ed2dc7SMatthew G. Knepley void *ctx; 124f2ed2dc7SMatthew G. Knepley 125f2ed2dc7SMatthew G. Knepley PetscFunctionBegin; 1263c633725SBarry Smith PetscCheck(ce->r == 2.0, PetscObjectComm((PetscObject)ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double)ce->r); 1279566063dSJacob Faibussowitsch PetscCall(DMGetDimension(ce->idm, &dim)); 1289566063dSJacob Faibussowitsch PetscCall(DMGetApplicationContext(ce->idm, &ctx)); 1299566063dSJacob Faibussowitsch PetscCall(DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE)); 1309566063dSJacob Faibussowitsch PetscCall(DMGetRefineLevel(ce->idm, &oldlevel)); 1319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1((Nr + 1), &dm)); 1329566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &uInitial)); 133f2ed2dc7SMatthew G. Knepley /* Loop over meshes */ 134f2ed2dc7SMatthew G. Knepley dm[0] = ce->idm; 135f2ed2dc7SMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 136f2ed2dc7SMatthew G. Knepley Vec u; 137f2ed2dc7SMatthew G. Knepley #if defined(PETSC_USE_LOG) 138f2ed2dc7SMatthew G. Knepley PetscLogStage stage; 139f2ed2dc7SMatthew G. Knepley #endif 140f2ed2dc7SMatthew G. Knepley char stageName[PETSC_MAX_PATH_LEN]; 141f2ed2dc7SMatthew G. Knepley const char *dmname, *uname; 142f2ed2dc7SMatthew G. Knepley 14363a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN - 1, "ConvEst Refinement Level %" PetscInt_FMT, r)); 144608e5a7aSMatthew G. Knepley #if defined(PETSC_USE_LOG) 1459566063dSJacob Faibussowitsch PetscCall(PetscLogStageGetId(stageName, &stage)); 1469566063dSJacob Faibussowitsch if (stage < 0) PetscCall(PetscLogStageRegister(stageName, &stage)); 147608e5a7aSMatthew G. Knepley #endif 1489566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(stage)); 149f2ed2dc7SMatthew G. Knepley if (r > 0) { 150b2df8587SMatthew G. Knepley if (!ce->noRefine) { 1519566063dSJacob Faibussowitsch PetscCall(DMRefine(dm[r - 1], MPI_COMM_NULL, &dm[r])); 1529566063dSJacob Faibussowitsch PetscCall(DMSetCoarseDM(dm[r], dm[r - 1])); 153b2df8587SMatthew G. Knepley } else { 154b2df8587SMatthew G. Knepley DM cdm, rcdm; 155b2df8587SMatthew G. Knepley 1569566063dSJacob Faibussowitsch PetscCall(DMClone(dm[r - 1], &dm[r])); 1579566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(dm[r - 1], dm[r])); 1589566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm[r - 1], &cdm)); 1599566063dSJacob Faibussowitsch PetscCall(DMGetCoordinateDM(dm[r], &rcdm)); 1609566063dSJacob Faibussowitsch PetscCall(DMCopyDisc(cdm, rcdm)); 161b2df8587SMatthew G. Knepley } 1629566063dSJacob Faibussowitsch PetscCall(DMCopyTransform(ce->idm, dm[r])); 1639566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName((PetscObject)dm[r - 1], &dmname)); 1649566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)dm[r], dmname)); 165f2ed2dc7SMatthew G. Knepley for (f = 0; f <= Nf; ++f) { 1668cda7954SMatthew G. Knepley PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *); 167f2ed2dc7SMatthew G. Knepley 1689566063dSJacob Faibussowitsch PetscCall(DMGetNullSpaceConstructor(dm[r - 1], f, &nspconstr)); 1699566063dSJacob Faibussowitsch PetscCall(DMSetNullSpaceConstructor(dm[r], f, nspconstr)); 170f2ed2dc7SMatthew G. Knepley } 171f2ed2dc7SMatthew G. Knepley } 1729566063dSJacob Faibussowitsch PetscCall(DMViewFromOptions(dm[r], NULL, "-conv_dm_view")); 173f2ed2dc7SMatthew G. Knepley /* Create solution */ 1749566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dm[r], &u)); 1759566063dSJacob Faibussowitsch PetscCall(DMGetField(dm[r], 0, NULL, &disc)); 1769566063dSJacob Faibussowitsch PetscCall(PetscObjectGetName(disc, &uname)); 1779566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)u, uname)); 178f2ed2dc7SMatthew G. Knepley /* Setup solver */ 1799566063dSJacob Faibussowitsch PetscCall(TSReset(ts)); 1809566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, dm[r])); 1819566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(dm[r], DMPlexTSComputeBoundary, ctx)); 182ae11cfe1SMatthew G. Knepley if (r > 0) { 183ae11cfe1SMatthew G. Knepley PetscCall(DMTSGetIFunctionLocal(dm[r - 1], &ifunc, &fctx)); 184ae11cfe1SMatthew G. Knepley PetscCall(DMTSGetIJacobianLocal(dm[r - 1], &ijac, &jctx)); 185ae11cfe1SMatthew G. Knepley PetscCall(DMTSGetRHSFunctionLocal(dm[r - 1], &rhsfunc, &rctx)); 186*48a46eb9SPierre Jolivet if (ifunc) PetscCall(DMTSSetIFunctionLocal(dm[r], ifunc, fctx)); 187*48a46eb9SPierre Jolivet if (ijac) PetscCall(DMTSSetIJacobianLocal(dm[r], ijac, jctx)); 188*48a46eb9SPierre Jolivet if (rhsfunc) PetscCall(DMTSSetRHSFunctionLocal(dm[r], rhsfunc, rctx)); 189ae11cfe1SMatthew G. Knepley } 1909566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 1919566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 1929566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 19349979724SMatthew G. Knepley PetscCall(TSSetSolution(ts, u)); 19449979724SMatthew G. Knepley PetscCall(VecDestroy(&u)); 195f2ed2dc7SMatthew G. Knepley /* Create initial guess */ 19649979724SMatthew G. Knepley PetscCall(TSGetSolution(ts, &u)); 1979566063dSJacob Faibussowitsch PetscCall(PetscConvEstComputeInitialGuess(ce, r, dm[r], u)); 19849979724SMatthew G. Knepley PetscCall(TSSolve(ts, NULL)); 19949979724SMatthew G. Knepley PetscCall(TSGetSolution(ts, &u)); 2009566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(ce->event, ce, 0, 0, 0)); 2019566063dSJacob Faibussowitsch PetscCall(PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r * Nf])); 2029566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(ce->event, ce, 0, 0, 0)); 203f2ed2dc7SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 204f2ed2dc7SMatthew G. Knepley PetscSection s, fs; 205f2ed2dc7SMatthew G. Knepley PetscInt lsize; 206f2ed2dc7SMatthew G. Knepley 207f2ed2dc7SMatthew G. Knepley /* Could use DMGetOutputDM() to add in Dirichlet dofs */ 2089566063dSJacob Faibussowitsch PetscCall(DMGetLocalSection(dm[r], &s)); 2099566063dSJacob Faibussowitsch PetscCall(PetscSectionGetField(s, f, &fs)); 2109566063dSJacob Faibussowitsch PetscCall(PetscSectionGetConstrainedStorageSize(fs, &lsize)); 2119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Allreduce(&lsize, &ce->dofs[r * Nf + f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)ts))); 2129566063dSJacob Faibussowitsch PetscCall(PetscLogEventSetDof(ce->event, f, ce->dofs[r * Nf + f])); 2139566063dSJacob Faibussowitsch PetscCall(PetscLogEventSetError(ce->event, f, ce->errors[r * Nf + f])); 214f2ed2dc7SMatthew G. Knepley } 215f2ed2dc7SMatthew G. Knepley /* Monitor */ 2169566063dSJacob Faibussowitsch PetscCall(PetscConvEstMonitorDefault(ce, r)); 217f2ed2dc7SMatthew G. Knepley if (!r) { 218f2ed2dc7SMatthew G. Knepley /* PCReset() does not wipe out the level structure */ 219f2ed2dc7SMatthew G. Knepley SNES snes; 220f2ed2dc7SMatthew G. Knepley KSP ksp; 221f2ed2dc7SMatthew G. Knepley PC pc; 222f2ed2dc7SMatthew G. Knepley 2239566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 2249566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 2259566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 2269566063dSJacob Faibussowitsch PetscCall(PCMGGetLevels(pc, &oldnlev)); 227f2ed2dc7SMatthew G. Knepley } 228f2ed2dc7SMatthew G. Knepley /* Cleanup */ 2299566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 230f2ed2dc7SMatthew G. Knepley } 231ae11cfe1SMatthew G. Knepley PetscCall(DMTSGetIFunctionLocal(dm[r - 1], &ifunc, &fctx)); 232ae11cfe1SMatthew G. Knepley PetscCall(DMTSGetIJacobianLocal(dm[r - 1], &ijac, &jctx)); 233ae11cfe1SMatthew G. Knepley PetscCall(DMTSGetRHSFunctionLocal(dm[r - 1], &rhsfunc, &rctx)); 234*48a46eb9SPierre Jolivet for (r = 1; r <= Nr; ++r) PetscCall(DMDestroy(&dm[r])); 235f2ed2dc7SMatthew G. Knepley /* Fit convergence rate */ 2369566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(Nr + 1, &x, Nr + 1, &y)); 237f2ed2dc7SMatthew G. Knepley for (f = 0; f < Nf; ++f) { 238f2ed2dc7SMatthew G. Knepley for (r = 0; r <= Nr; ++r) { 2397809adefSMatthew G. Knepley x[r] = PetscLog10Real(ce->dofs[r * Nf + f]); 240f2ed2dc7SMatthew G. Knepley y[r] = PetscLog10Real(ce->errors[r * Nf + f]); 241f2ed2dc7SMatthew G. Knepley } 2429566063dSJacob Faibussowitsch PetscCall(PetscLinearRegression(Nr + 1, x, y, &slope, &intercept)); 243f2ed2dc7SMatthew G. Knepley /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */ 244f2ed2dc7SMatthew G. Knepley alpha[f] = -slope * dim; 245f2ed2dc7SMatthew G. Knepley } 2469566063dSJacob Faibussowitsch PetscCall(PetscFree2(x, y)); 2479566063dSJacob Faibussowitsch PetscCall(PetscFree(dm)); 248f2ed2dc7SMatthew G. Knepley /* Restore solver */ 2499566063dSJacob Faibussowitsch PetscCall(TSReset(ts)); 250f2ed2dc7SMatthew G. Knepley { 251f2ed2dc7SMatthew G. Knepley /* PCReset() does not wipe out the level structure */ 252f2ed2dc7SMatthew G. Knepley SNES snes; 253f2ed2dc7SMatthew G. Knepley KSP ksp; 254f2ed2dc7SMatthew G. Knepley PC pc; 255f2ed2dc7SMatthew G. Knepley 2569566063dSJacob Faibussowitsch PetscCall(TSGetSNES(ts, &snes)); 2579566063dSJacob Faibussowitsch PetscCall(SNESGetKSP(snes, &ksp)); 2589566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 2599566063dSJacob Faibussowitsch PetscCall(PCMGSetLevels(pc, oldnlev, NULL)); 2609566063dSJacob Faibussowitsch PetscCall(DMSetRefineLevel(ce->idm, oldlevel)); /* The damn DMCoarsen() calls in PCMG can reset this */ 261f2ed2dc7SMatthew G. Knepley } 2629566063dSJacob Faibussowitsch PetscCall(TSSetDM(ts, ce->idm)); 2639566063dSJacob Faibussowitsch PetscCall(DMTSSetBoundaryLocal(ce->idm, DMPlexTSComputeBoundary, ctx)); 264*48a46eb9SPierre Jolivet if (ifunc) PetscCall(DMTSSetIFunctionLocal(ce->idm, ifunc, fctx)); 265*48a46eb9SPierre Jolivet if (ijac) PetscCall(DMTSSetIJacobianLocal(ce->idm, ijac, jctx)); 266*48a46eb9SPierre Jolivet if (rhsfunc) PetscCall(DMTSSetRHSFunctionLocal(ce->idm, rhsfunc, rctx)); 2679566063dSJacob Faibussowitsch PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_ITERATING)); 2689566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 2699566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 2709566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 2719566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, uInitial)); 2729566063dSJacob Faibussowitsch PetscCall(PetscConvEstComputeInitialGuess(ce, 0, NULL, uInitial)); 273f2ed2dc7SMatthew G. Knepley PetscFunctionReturn(0); 274f2ed2dc7SMatthew G. Knepley } 275f2ed2dc7SMatthew G. Knepley 2769371c9d4SSatish Balay PetscErrorCode PetscConvEstUseTS(PetscConvEst ce, PetscBool checkTemporal) { 277900f6b5bSMatthew G. Knepley PetscFunctionBegin; 278900f6b5bSMatthew G. Knepley PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1); 279900f6b5bSMatthew G. Knepley ce->ops->setsolver = PetscConvEstSetTS_Private; 280900f6b5bSMatthew G. Knepley ce->ops->initguess = PetscConvEstInitGuessTS_Private; 281900f6b5bSMatthew G. Knepley ce->ops->computeerror = PetscConvEstComputeErrorTS_Private; 282f2ed2dc7SMatthew G. Knepley if (checkTemporal) { 283f2ed2dc7SMatthew G. Knepley ce->ops->getconvrate = PetscConvEstGetConvRateTS_Temporal_Private; 284f2ed2dc7SMatthew G. Knepley } else { 285f2ed2dc7SMatthew G. Knepley ce->ops->getconvrate = PetscConvEstGetConvRateTS_Spatial_Private; 286f2ed2dc7SMatthew G. Knepley } 287900f6b5bSMatthew G. Knepley PetscFunctionReturn(0); 288900f6b5bSMatthew G. Knepley } 289