xref: /petsc/src/ts/utils/tsconvest.c (revision a7e8c07f09f75bc4f67bdfe093b39d3875753d0d)
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 *);
127ae11cfe1SMatthew G. Knepley   void *fctx, *jctx, *rctx;
128f2ed2dc7SMatthew G. Knepley   void *ctx;
129f2ed2dc7SMatthew G. Knepley 
130f2ed2dc7SMatthew G. Knepley   PetscFunctionBegin;
1313c633725SBarry Smith   PetscCheck(ce->r == 2.0, PetscObjectComm((PetscObject)ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double)ce->r);
1329566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(ce->idm, &dim));
1339566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(ce->idm, &ctx));
1349566063dSJacob Faibussowitsch   PetscCall(DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE));
1359566063dSJacob Faibussowitsch   PetscCall(DMGetRefineLevel(ce->idm, &oldlevel));
1369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((Nr + 1), &dm));
1379566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &uInitial));
138*a7e8c07fSStefano Zampini   PetscCall(PetscObjectReference((PetscObject)uInitial));
139*a7e8c07fSStefano Zampini 
140f2ed2dc7SMatthew G. Knepley   /* Loop over meshes */
141f2ed2dc7SMatthew G. Knepley   dm[0] = ce->idm;
142f2ed2dc7SMatthew G. Knepley   for (r = 0; r <= Nr; ++r) {
143f2ed2dc7SMatthew G. Knepley     Vec u;
144f2ed2dc7SMatthew G. Knepley #if defined(PETSC_USE_LOG)
145f2ed2dc7SMatthew G. Knepley     PetscLogStage stage;
146f2ed2dc7SMatthew G. Knepley #endif
147f2ed2dc7SMatthew G. Knepley     char        stageName[PETSC_MAX_PATH_LEN];
148f2ed2dc7SMatthew G. Knepley     const char *dmname, *uname;
149f2ed2dc7SMatthew G. Knepley 
15063a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN - 1, "ConvEst Refinement Level %" PetscInt_FMT, r));
151608e5a7aSMatthew G. Knepley #if defined(PETSC_USE_LOG)
1529566063dSJacob Faibussowitsch     PetscCall(PetscLogStageGetId(stageName, &stage));
1539566063dSJacob Faibussowitsch     if (stage < 0) PetscCall(PetscLogStageRegister(stageName, &stage));
154608e5a7aSMatthew G. Knepley #endif
1559566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage));
156f2ed2dc7SMatthew G. Knepley     if (r > 0) {
157b2df8587SMatthew G. Knepley       if (!ce->noRefine) {
1589566063dSJacob Faibussowitsch         PetscCall(DMRefine(dm[r - 1], MPI_COMM_NULL, &dm[r]));
1599566063dSJacob Faibussowitsch         PetscCall(DMSetCoarseDM(dm[r], dm[r - 1]));
160b2df8587SMatthew G. Knepley       } else {
161b2df8587SMatthew G. Knepley         DM cdm, rcdm;
162b2df8587SMatthew G. Knepley 
1639566063dSJacob Faibussowitsch         PetscCall(DMClone(dm[r - 1], &dm[r]));
1649566063dSJacob Faibussowitsch         PetscCall(DMCopyDisc(dm[r - 1], dm[r]));
1659566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinateDM(dm[r - 1], &cdm));
1669566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinateDM(dm[r], &rcdm));
1679566063dSJacob Faibussowitsch         PetscCall(DMCopyDisc(cdm, rcdm));
168b2df8587SMatthew G. Knepley       }
1699566063dSJacob Faibussowitsch       PetscCall(DMCopyTransform(ce->idm, dm[r]));
1709566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject)dm[r - 1], &dmname));
1719566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject)dm[r], dmname));
172f2ed2dc7SMatthew G. Knepley       for (f = 0; f <= Nf; ++f) {
1738cda7954SMatthew G. Knepley         PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *);
174f2ed2dc7SMatthew G. Knepley 
1759566063dSJacob Faibussowitsch         PetscCall(DMGetNullSpaceConstructor(dm[r - 1], f, &nspconstr));
1769566063dSJacob Faibussowitsch         PetscCall(DMSetNullSpaceConstructor(dm[r], f, nspconstr));
177f2ed2dc7SMatthew G. Knepley       }
178f2ed2dc7SMatthew G. Knepley     }
1799566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(dm[r], NULL, "-conv_dm_view"));
180f2ed2dc7SMatthew G. Knepley     /* Create solution */
1819566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dm[r], &u));
1829566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm[r], 0, NULL, &disc));
1839566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName(disc, &uname));
1849566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject)u, uname));
185f2ed2dc7SMatthew G. Knepley     /* Setup solver */
1869566063dSJacob Faibussowitsch     PetscCall(TSReset(ts));
1879566063dSJacob Faibussowitsch     PetscCall(TSSetDM(ts, dm[r]));
1889566063dSJacob Faibussowitsch     PetscCall(DMTSSetBoundaryLocal(dm[r], DMPlexTSComputeBoundary, ctx));
189ae11cfe1SMatthew G. Knepley     if (r > 0) {
190ae11cfe1SMatthew G. Knepley       PetscCall(DMTSGetIFunctionLocal(dm[r - 1], &ifunc, &fctx));
191ae11cfe1SMatthew G. Knepley       PetscCall(DMTSGetIJacobianLocal(dm[r - 1], &ijac, &jctx));
192ae11cfe1SMatthew G. Knepley       PetscCall(DMTSGetRHSFunctionLocal(dm[r - 1], &rhsfunc, &rctx));
19348a46eb9SPierre Jolivet       if (ifunc) PetscCall(DMTSSetIFunctionLocal(dm[r], ifunc, fctx));
19448a46eb9SPierre Jolivet       if (ijac) PetscCall(DMTSSetIJacobianLocal(dm[r], ijac, jctx));
19548a46eb9SPierre Jolivet       if (rhsfunc) PetscCall(DMTSSetRHSFunctionLocal(dm[r], rhsfunc, rctx));
196ae11cfe1SMatthew G. Knepley     }
1979566063dSJacob Faibussowitsch     PetscCall(TSSetTime(ts, 0.0));
1989566063dSJacob Faibussowitsch     PetscCall(TSSetStepNumber(ts, 0));
1999566063dSJacob Faibussowitsch     PetscCall(TSSetFromOptions(ts));
20049979724SMatthew G. Knepley     PetscCall(TSSetSolution(ts, u));
20149979724SMatthew G. Knepley     PetscCall(VecDestroy(&u));
202f2ed2dc7SMatthew G. Knepley     /* Create initial guess */
20349979724SMatthew G. Knepley     PetscCall(TSGetSolution(ts, &u));
2049566063dSJacob Faibussowitsch     PetscCall(PetscConvEstComputeInitialGuess(ce, r, dm[r], u));
20549979724SMatthew G. Knepley     PetscCall(TSSolve(ts, NULL));
20649979724SMatthew G. Knepley     PetscCall(TSGetSolution(ts, &u));
2079566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(ce->event, ce, 0, 0, 0));
2089566063dSJacob Faibussowitsch     PetscCall(PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r * Nf]));
2099566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(ce->event, ce, 0, 0, 0));
210f2ed2dc7SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
211f2ed2dc7SMatthew G. Knepley       PetscSection s, fs;
212f2ed2dc7SMatthew G. Knepley       PetscInt     lsize;
213f2ed2dc7SMatthew G. Knepley 
214f2ed2dc7SMatthew G. Knepley       /* Could use DMGetOutputDM() to add in Dirichlet dofs */
2159566063dSJacob Faibussowitsch       PetscCall(DMGetLocalSection(dm[r], &s));
2169566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetField(s, f, &fs));
2179566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstrainedStorageSize(fs, &lsize));
218712fec58SPierre Jolivet       PetscCall(MPIU_Allreduce(&lsize, &ce->dofs[r * Nf + f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)ts)));
2199566063dSJacob Faibussowitsch       PetscCall(PetscLogEventSetDof(ce->event, f, ce->dofs[r * Nf + f]));
2209566063dSJacob Faibussowitsch       PetscCall(PetscLogEventSetError(ce->event, f, ce->errors[r * Nf + f]));
221f2ed2dc7SMatthew G. Knepley     }
222f2ed2dc7SMatthew G. Knepley     /* Monitor */
2239566063dSJacob Faibussowitsch     PetscCall(PetscConvEstMonitorDefault(ce, r));
224f2ed2dc7SMatthew G. Knepley     if (!r) {
225f2ed2dc7SMatthew G. Knepley       /* PCReset() does not wipe out the level structure */
226f2ed2dc7SMatthew G. Knepley       SNES snes;
227f2ed2dc7SMatthew G. Knepley       KSP  ksp;
228f2ed2dc7SMatthew G. Knepley       PC   pc;
229f2ed2dc7SMatthew G. Knepley 
2309566063dSJacob Faibussowitsch       PetscCall(TSGetSNES(ts, &snes));
2319566063dSJacob Faibussowitsch       PetscCall(SNESGetKSP(snes, &ksp));
2329566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(ksp, &pc));
2339566063dSJacob Faibussowitsch       PetscCall(PCMGGetLevels(pc, &oldnlev));
234f2ed2dc7SMatthew G. Knepley     }
235f2ed2dc7SMatthew G. Knepley     /* Cleanup */
2369566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
237f2ed2dc7SMatthew G. Knepley   }
238ae11cfe1SMatthew G. Knepley   PetscCall(DMTSGetIFunctionLocal(dm[r - 1], &ifunc, &fctx));
239ae11cfe1SMatthew G. Knepley   PetscCall(DMTSGetIJacobianLocal(dm[r - 1], &ijac, &jctx));
240ae11cfe1SMatthew G. Knepley   PetscCall(DMTSGetRHSFunctionLocal(dm[r - 1], &rhsfunc, &rctx));
24148a46eb9SPierre Jolivet   for (r = 1; r <= Nr; ++r) PetscCall(DMDestroy(&dm[r]));
242f2ed2dc7SMatthew G. Knepley   /* Fit convergence rate */
2439566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nr + 1, &x, Nr + 1, &y));
244f2ed2dc7SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
245f2ed2dc7SMatthew G. Knepley     for (r = 0; r <= Nr; ++r) {
2467809adefSMatthew G. Knepley       x[r] = PetscLog10Real(ce->dofs[r * Nf + f]);
247f2ed2dc7SMatthew G. Knepley       y[r] = PetscLog10Real(ce->errors[r * Nf + f]);
248f2ed2dc7SMatthew G. Knepley     }
2499566063dSJacob Faibussowitsch     PetscCall(PetscLinearRegression(Nr + 1, x, y, &slope, &intercept));
250f2ed2dc7SMatthew G. Knepley     /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */
251f2ed2dc7SMatthew G. Knepley     alpha[f] = -slope * dim;
252f2ed2dc7SMatthew G. Knepley   }
2539566063dSJacob Faibussowitsch   PetscCall(PetscFree2(x, y));
2549566063dSJacob Faibussowitsch   PetscCall(PetscFree(dm));
255f2ed2dc7SMatthew G. Knepley   /* Restore solver */
2569566063dSJacob Faibussowitsch   PetscCall(TSReset(ts));
257f2ed2dc7SMatthew G. Knepley   {
258f2ed2dc7SMatthew G. Knepley     /* PCReset() does not wipe out the level structure */
259f2ed2dc7SMatthew G. Knepley     SNES snes;
260f2ed2dc7SMatthew G. Knepley     KSP  ksp;
261f2ed2dc7SMatthew G. Knepley     PC   pc;
262f2ed2dc7SMatthew G. Knepley 
2639566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts, &snes));
2649566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes, &ksp));
2659566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
2669566063dSJacob Faibussowitsch     PetscCall(PCMGSetLevels(pc, oldnlev, NULL));
2679566063dSJacob Faibussowitsch     PetscCall(DMSetRefineLevel(ce->idm, oldlevel)); /* The damn DMCoarsen() calls in PCMG can reset this */
268f2ed2dc7SMatthew G. Knepley   }
2699566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, ce->idm));
2709566063dSJacob Faibussowitsch   PetscCall(DMTSSetBoundaryLocal(ce->idm, DMPlexTSComputeBoundary, ctx));
27148a46eb9SPierre Jolivet   if (ifunc) PetscCall(DMTSSetIFunctionLocal(ce->idm, ifunc, fctx));
27248a46eb9SPierre Jolivet   if (ijac) PetscCall(DMTSSetIJacobianLocal(ce->idm, ijac, jctx));
27348a46eb9SPierre Jolivet   if (rhsfunc) PetscCall(DMTSSetRHSFunctionLocal(ce->idm, rhsfunc, rctx));
2749566063dSJacob Faibussowitsch   PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_ITERATING));
2759566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
2769566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts, 0));
2779566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
2789566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, uInitial));
2799566063dSJacob Faibussowitsch   PetscCall(PetscConvEstComputeInitialGuess(ce, 0, NULL, uInitial));
280*a7e8c07fSStefano Zampini   PetscCall(VecDestroy(&uInitial));
2813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
282f2ed2dc7SMatthew G. Knepley }
283f2ed2dc7SMatthew G. Knepley 
284d71ae5a4SJacob Faibussowitsch PetscErrorCode PetscConvEstUseTS(PetscConvEst ce, PetscBool checkTemporal)
285d71ae5a4SJacob Faibussowitsch {
286900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
287900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
288900f6b5bSMatthew G. Knepley   ce->ops->setsolver    = PetscConvEstSetTS_Private;
289900f6b5bSMatthew G. Knepley   ce->ops->initguess    = PetscConvEstInitGuessTS_Private;
290900f6b5bSMatthew G. Knepley   ce->ops->computeerror = PetscConvEstComputeErrorTS_Private;
291f2ed2dc7SMatthew G. Knepley   if (checkTemporal) {
292f2ed2dc7SMatthew G. Knepley     ce->ops->getconvrate = PetscConvEstGetConvRateTS_Temporal_Private;
293f2ed2dc7SMatthew G. Knepley   } else {
294f2ed2dc7SMatthew G. Knepley     ce->ops->getconvrate = PetscConvEstGetConvRateTS_Spatial_Private;
295f2ed2dc7SMatthew G. Knepley   }
2963ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
297900f6b5bSMatthew G. Knepley }
298