xref: /petsc/src/ts/utils/tsconvest.c (revision 9566063d113dddea24716c546802770db7481bc0)
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 
7900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstSetTS_Private(PetscConvEst ce, PetscObject solver)
8900f6b5bSMatthew G. Knepley {
9900f6b5bSMatthew G. Knepley   PetscClassId   id;
10900f6b5bSMatthew G. Knepley 
11900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
12*9566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetClassId(ce->solver, &id));
133c633725SBarry Smith   PetscCheck(id == TS_CLASSID,PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "Solver was not a TS");
14*9566063dSJacob Faibussowitsch   PetscCall(TSGetDM((TS) ce->solver, &ce->idm));
15900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
16900f6b5bSMatthew G. Knepley }
17900f6b5bSMatthew G. Knepley 
18900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstInitGuessTS_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u)
19900f6b5bSMatthew G. Knepley {
20900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
21*9566063dSJacob Faibussowitsch   PetscCall(TSComputeInitialCondition((TS) ce->solver, u));
22900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
23900f6b5bSMatthew G. Knepley }
24900f6b5bSMatthew G. Knepley 
25900f6b5bSMatthew G. Knepley static PetscErrorCode PetscConvEstComputeErrorTS_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
26900f6b5bSMatthew G. Knepley {
27f2cacb80SMatthew G. Knepley   TS               ts = (TS) ce->solver;
28f2cacb80SMatthew G. Knepley   PetscErrorCode (*exactError)(TS, Vec, Vec);
29900f6b5bSMatthew G. Knepley 
30900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
31*9566063dSJacob Faibussowitsch   PetscCall(TSGetComputeExactError(ts, &exactError));
32f2cacb80SMatthew G. Knepley   if (exactError) {
33f2cacb80SMatthew G. Knepley     Vec      e;
34f2cacb80SMatthew G. Knepley     PetscInt f;
35f2cacb80SMatthew G. Knepley 
36*9566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(u, &e));
37*9566063dSJacob Faibussowitsch     PetscCall(TSComputeExactError(ts, u, e));
38*9566063dSJacob Faibussowitsch     PetscCall(VecNorm(e, NORM_2, errors));
39f2cacb80SMatthew G. Knepley     for (f = 1; f < ce->Nf; ++f) errors[f] = errors[0];
40*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&e));
41f2cacb80SMatthew G. Knepley   } else {
42f2cacb80SMatthew G. Knepley     PetscReal t;
43f2cacb80SMatthew G. Knepley 
44*9566063dSJacob Faibussowitsch     PetscCall(TSGetSolveTime(ts, &t));
45*9566063dSJacob Faibussowitsch     PetscCall(DMComputeL2FieldDiff(dm, t, ce->exactSol, ce->ctxs, u, errors));
46f2cacb80SMatthew G. Knepley   }
47900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
48900f6b5bSMatthew G. Knepley }
49900f6b5bSMatthew G. Knepley 
50f2ed2dc7SMatthew G. Knepley static PetscErrorCode PetscConvEstGetConvRateTS_Temporal_Private(PetscConvEst ce, PetscReal alpha[])
51900f6b5bSMatthew G. Knepley {
52900f6b5bSMatthew G. Knepley   TS             ts = (TS) ce->solver;
53900f6b5bSMatthew G. Knepley   Vec            u;
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;
58*9566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &u));
59*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(Nr+1, &dt));
60*9566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt[0]));
61*9566063dSJacob Faibussowitsch   PetscCall(TSGetMaxSteps(ts, &oNs));
62900f6b5bSMatthew G. Knepley   Ns   = oNs;
63900f6b5bSMatthew G. Knepley   for (r = 0; r <= Nr; ++r) {
64900f6b5bSMatthew G. Knepley     if (r > 0) {
652e61be88SMatthew G. Knepley       dt[r] = dt[r-1]/ce->r;
662e61be88SMatthew G. Knepley       Ns    = PetscCeilReal(Ns*ce->r);
67900f6b5bSMatthew G. Knepley     }
68*9566063dSJacob Faibussowitsch     PetscCall(TSSetTime(ts, 0.0));
69*9566063dSJacob Faibussowitsch     PetscCall(TSSetStepNumber(ts, 0));
70*9566063dSJacob Faibussowitsch     PetscCall(TSSetTimeStep(ts, dt[r]));
71*9566063dSJacob Faibussowitsch     PetscCall(TSSetMaxSteps(ts, Ns));
72*9566063dSJacob Faibussowitsch     PetscCall(PetscConvEstComputeInitialGuess(ce, r, NULL, u));
73*9566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, u));
74*9566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(ce->event, ce, 0, 0, 0));
75*9566063dSJacob Faibussowitsch     PetscCall(PetscConvEstComputeError(ce, r, ce->idm, u, &ce->errors[r*Nf]));
76*9566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(ce->event, ce, 0, 0, 0));
77f2ed2dc7SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
789327cc11SMatthew G. Knepley       ce->dofs[r*Nf+f] = 1.0/dt[r];
79*9566063dSJacob Faibussowitsch       PetscCall(PetscLogEventSetDof(ce->event, f, ce->dofs[r*Nf+f]));
80*9566063dSJacob Faibussowitsch       PetscCall(PetscLogEventSetError(ce->event, f, ce->errors[r*Nf+f]));
81900f6b5bSMatthew G. Knepley     }
82f2cacb80SMatthew G. Knepley     /* Monitor */
83*9566063dSJacob Faibussowitsch     PetscCall(PetscConvEstMonitorDefault(ce, r));
84900f6b5bSMatthew G. Knepley   }
85900f6b5bSMatthew G. Knepley   /* Fit convergence rate */
86900f6b5bSMatthew G. Knepley   if (Nr) {
87*9566063dSJacob Faibussowitsch     PetscCall(PetscMalloc2(Nr+1, &x, Nr+1, &y));
88900f6b5bSMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
89900f6b5bSMatthew G. Knepley       for (r = 0; r <= Nr; ++r) {
90f2cacb80SMatthew G. Knepley         x[r] = PetscLog10Real(dt[r]);
91900f6b5bSMatthew G. Knepley         y[r] = PetscLog10Real(ce->errors[r*Nf+f]);
92900f6b5bSMatthew G. Knepley       }
93*9566063dSJacob Faibussowitsch       PetscCall(PetscLinearRegression(Nr+1, x, y, &slope, &intercept));
94900f6b5bSMatthew G. Knepley       /* Since lg err = s lg dt + b */
95900f6b5bSMatthew G. Knepley       alpha[f] = slope;
96900f6b5bSMatthew G. Knepley     }
97*9566063dSJacob Faibussowitsch     PetscCall(PetscFree2(x, y));
98900f6b5bSMatthew G. Knepley   }
99900f6b5bSMatthew G. Knepley   /* Reset solver */
100*9566063dSJacob Faibussowitsch   PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_ITERATING));
101*9566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
102*9566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts, 0));
103*9566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt[0]));
104*9566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, oNs));
105*9566063dSJacob Faibussowitsch   PetscCall(PetscConvEstComputeInitialGuess(ce, 0, NULL, u));
106*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(dt));
107900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
108900f6b5bSMatthew G. Knepley }
109900f6b5bSMatthew G. Knepley 
110f2ed2dc7SMatthew G. Knepley static PetscErrorCode PetscConvEstGetConvRateTS_Spatial_Private(PetscConvEst ce, PetscReal alpha[])
111f2ed2dc7SMatthew G. Knepley {
112f2ed2dc7SMatthew G. Knepley   TS             ts = (TS) ce->solver;
113f2ed2dc7SMatthew G. Knepley   Vec            uInitial;
114f2ed2dc7SMatthew G. Knepley   DM            *dm;
115f2ed2dc7SMatthew G. Knepley   PetscObject    disc;
116f2ed2dc7SMatthew G. Knepley   PetscReal     *x, *y, slope, intercept;
1177809adefSMatthew G. Knepley   PetscInt       Nr = ce->Nr, r, Nf = ce->Nf, f, dim, oldlevel, oldnlev;
118f2ed2dc7SMatthew G. Knepley   void          *ctx;
119f2ed2dc7SMatthew G. Knepley 
120f2ed2dc7SMatthew G. Knepley   PetscFunctionBegin;
1213c633725SBarry Smith   PetscCheck(ce->r == 2.0,PetscObjectComm((PetscObject) ce), PETSC_ERR_SUP, "Only refinement factor 2 is currently supported (not %g)", (double) ce->r);
122*9566063dSJacob Faibussowitsch   PetscCall(DMGetDimension(ce->idm, &dim));
123*9566063dSJacob Faibussowitsch   PetscCall(DMGetApplicationContext(ce->idm, &ctx));
124*9566063dSJacob Faibussowitsch   PetscCall(DMPlexSetRefinementUniform(ce->idm, PETSC_TRUE));
125*9566063dSJacob Faibussowitsch   PetscCall(DMGetRefineLevel(ce->idm, &oldlevel));
126*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1((Nr+1), &dm));
127*9566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &uInitial));
128f2ed2dc7SMatthew G. Knepley   /* Loop over meshes */
129f2ed2dc7SMatthew G. Knepley   dm[0] = ce->idm;
130f2ed2dc7SMatthew G. Knepley   for (r = 0; r <= Nr; ++r) {
131f2ed2dc7SMatthew G. Knepley     Vec           u;
132f2ed2dc7SMatthew G. Knepley #if defined(PETSC_USE_LOG)
133f2ed2dc7SMatthew G. Knepley     PetscLogStage stage;
134f2ed2dc7SMatthew G. Knepley #endif
135f2ed2dc7SMatthew G. Knepley     char          stageName[PETSC_MAX_PATH_LEN];
136f2ed2dc7SMatthew G. Knepley     const char   *dmname, *uname;
137f2ed2dc7SMatthew G. Knepley 
138*9566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(stageName, PETSC_MAX_PATH_LEN-1, "ConvEst Refinement Level %D", r));
139608e5a7aSMatthew G. Knepley #if defined(PETSC_USE_LOG)
140*9566063dSJacob Faibussowitsch     PetscCall(PetscLogStageGetId(stageName, &stage));
141*9566063dSJacob Faibussowitsch     if (stage < 0) PetscCall(PetscLogStageRegister(stageName, &stage));
142608e5a7aSMatthew G. Knepley #endif
143*9566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePush(stage));
144f2ed2dc7SMatthew G. Knepley     if (r > 0) {
145b2df8587SMatthew G. Knepley       if (!ce->noRefine) {
146*9566063dSJacob Faibussowitsch         PetscCall(DMRefine(dm[r-1], MPI_COMM_NULL, &dm[r]));
147*9566063dSJacob Faibussowitsch         PetscCall(DMSetCoarseDM(dm[r], dm[r-1]));
148b2df8587SMatthew G. Knepley       } else {
149b2df8587SMatthew G. Knepley         DM cdm, rcdm;
150b2df8587SMatthew G. Knepley 
151*9566063dSJacob Faibussowitsch         PetscCall(DMClone(dm[r-1], &dm[r]));
152*9566063dSJacob Faibussowitsch         PetscCall(DMCopyDisc(dm[r-1], dm[r]));
153*9566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinateDM(dm[r-1], &cdm));
154*9566063dSJacob Faibussowitsch         PetscCall(DMGetCoordinateDM(dm[r],   &rcdm));
155*9566063dSJacob Faibussowitsch         PetscCall(DMCopyDisc(cdm, rcdm));
156b2df8587SMatthew G. Knepley       }
157*9566063dSJacob Faibussowitsch       PetscCall(DMCopyTransform(ce->idm, dm[r]));
158*9566063dSJacob Faibussowitsch       PetscCall(PetscObjectGetName((PetscObject) dm[r-1], &dmname));
159*9566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetName((PetscObject) dm[r], dmname));
160f2ed2dc7SMatthew G. Knepley       for (f = 0; f <= Nf; ++f) {
1618cda7954SMatthew G. Knepley         PetscErrorCode (*nspconstr)(DM, PetscInt, PetscInt, MatNullSpace *);
162f2ed2dc7SMatthew G. Knepley 
163*9566063dSJacob Faibussowitsch         PetscCall(DMGetNullSpaceConstructor(dm[r-1], f, &nspconstr));
164*9566063dSJacob Faibussowitsch         PetscCall(DMSetNullSpaceConstructor(dm[r],   f,  nspconstr));
165f2ed2dc7SMatthew G. Knepley       }
166f2ed2dc7SMatthew G. Knepley     }
167*9566063dSJacob Faibussowitsch     PetscCall(DMViewFromOptions(dm[r], NULL, "-conv_dm_view"));
168f2ed2dc7SMatthew G. Knepley     /* Create solution */
169*9566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(dm[r], &u));
170*9566063dSJacob Faibussowitsch     PetscCall(DMGetField(dm[r], 0, NULL, &disc));
171*9566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetName(disc, &uname));
172*9566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetName((PetscObject) u, uname));
173f2ed2dc7SMatthew G. Knepley     /* Setup solver */
174*9566063dSJacob Faibussowitsch     PetscCall(TSReset(ts));
175*9566063dSJacob Faibussowitsch     PetscCall(TSSetDM(ts, dm[r]));
176*9566063dSJacob Faibussowitsch     PetscCall(DMTSSetBoundaryLocal(dm[r], DMPlexTSComputeBoundary, ctx));
177*9566063dSJacob Faibussowitsch     PetscCall(DMTSSetIFunctionLocal(dm[r], DMPlexTSComputeIFunctionFEM, ctx));
178*9566063dSJacob Faibussowitsch     PetscCall(DMTSSetIJacobianLocal(dm[r], DMPlexTSComputeIJacobianFEM, ctx));
179*9566063dSJacob Faibussowitsch     PetscCall(TSSetTime(ts, 0.0));
180*9566063dSJacob Faibussowitsch     PetscCall(TSSetStepNumber(ts, 0));
181*9566063dSJacob Faibussowitsch     PetscCall(TSSetFromOptions(ts));
182f2ed2dc7SMatthew G. Knepley     /* Create initial guess */
183*9566063dSJacob Faibussowitsch     PetscCall(PetscConvEstComputeInitialGuess(ce, r, dm[r], u));
184*9566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts, u));
185*9566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(ce->event, ce, 0, 0, 0));
186*9566063dSJacob Faibussowitsch     PetscCall(PetscConvEstComputeError(ce, r, dm[r], u, &ce->errors[r*Nf]));
187*9566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(ce->event, ce, 0, 0, 0));
188f2ed2dc7SMatthew G. Knepley     for (f = 0; f < Nf; ++f) {
189f2ed2dc7SMatthew G. Knepley       PetscSection s, fs;
190f2ed2dc7SMatthew G. Knepley       PetscInt     lsize;
191f2ed2dc7SMatthew G. Knepley 
192f2ed2dc7SMatthew G. Knepley       /* Could use DMGetOutputDM() to add in Dirichlet dofs */
193*9566063dSJacob Faibussowitsch       PetscCall(DMGetLocalSection(dm[r], &s));
194*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetField(s, f, &fs));
195*9566063dSJacob Faibussowitsch       PetscCall(PetscSectionGetConstrainedStorageSize(fs, &lsize));
196*9566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Allreduce(&lsize, &ce->dofs[r*Nf+f], 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) ts)));
197*9566063dSJacob Faibussowitsch       PetscCall(PetscLogEventSetDof(ce->event, f, ce->dofs[r*Nf+f]));
198*9566063dSJacob Faibussowitsch       PetscCall(PetscLogEventSetError(ce->event, f, ce->errors[r*Nf+f]));
199f2ed2dc7SMatthew G. Knepley     }
200f2ed2dc7SMatthew G. Knepley     /* Monitor */
201*9566063dSJacob Faibussowitsch     PetscCall(PetscConvEstMonitorDefault(ce, r));
202f2ed2dc7SMatthew G. Knepley     if (!r) {
203f2ed2dc7SMatthew G. Knepley       /* PCReset() does not wipe out the level structure */
204f2ed2dc7SMatthew G. Knepley       SNES snes;
205f2ed2dc7SMatthew G. Knepley       KSP  ksp;
206f2ed2dc7SMatthew G. Knepley       PC   pc;
207f2ed2dc7SMatthew G. Knepley 
208*9566063dSJacob Faibussowitsch       PetscCall(TSGetSNES(ts, &snes));
209*9566063dSJacob Faibussowitsch       PetscCall(SNESGetKSP(snes, &ksp));
210*9566063dSJacob Faibussowitsch       PetscCall(KSPGetPC(ksp, &pc));
211*9566063dSJacob Faibussowitsch       PetscCall(PCMGGetLevels(pc, &oldnlev));
212f2ed2dc7SMatthew G. Knepley     }
213f2ed2dc7SMatthew G. Knepley     /* Cleanup */
214*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&u));
215*9566063dSJacob Faibussowitsch     PetscCall(PetscLogStagePop());
216f2ed2dc7SMatthew G. Knepley   }
217f2ed2dc7SMatthew G. Knepley   for (r = 1; r <= Nr; ++r) {
218*9566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&dm[r]));
219f2ed2dc7SMatthew G. Knepley   }
220f2ed2dc7SMatthew G. Knepley   /* Fit convergence rate */
221*9566063dSJacob Faibussowitsch   PetscCall(PetscMalloc2(Nr+1, &x, Nr+1, &y));
222f2ed2dc7SMatthew G. Knepley   for (f = 0; f < Nf; ++f) {
223f2ed2dc7SMatthew G. Knepley     for (r = 0; r <= Nr; ++r) {
2247809adefSMatthew G. Knepley       x[r] = PetscLog10Real(ce->dofs[r*Nf+f]);
225f2ed2dc7SMatthew G. Knepley       y[r] = PetscLog10Real(ce->errors[r*Nf+f]);
226f2ed2dc7SMatthew G. Knepley     }
227*9566063dSJacob Faibussowitsch     PetscCall(PetscLinearRegression(Nr+1, x, y, &slope, &intercept));
228f2ed2dc7SMatthew G. Knepley     /* Since h^{-dim} = N, lg err = s lg N + b = -s dim lg h + b */
229f2ed2dc7SMatthew G. Knepley     alpha[f] = -slope * dim;
230f2ed2dc7SMatthew G. Knepley   }
231*9566063dSJacob Faibussowitsch   PetscCall(PetscFree2(x, y));
232*9566063dSJacob Faibussowitsch   PetscCall(PetscFree(dm));
233f2ed2dc7SMatthew G. Knepley   /* Restore solver */
234*9566063dSJacob Faibussowitsch   PetscCall(TSReset(ts));
235f2ed2dc7SMatthew G. Knepley   {
236f2ed2dc7SMatthew G. Knepley     /* PCReset() does not wipe out the level structure */
237f2ed2dc7SMatthew G. Knepley     SNES snes;
238f2ed2dc7SMatthew G. Knepley     KSP  ksp;
239f2ed2dc7SMatthew G. Knepley     PC   pc;
240f2ed2dc7SMatthew G. Knepley 
241*9566063dSJacob Faibussowitsch     PetscCall(TSGetSNES(ts, &snes));
242*9566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(snes, &ksp));
243*9566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
244*9566063dSJacob Faibussowitsch     PetscCall(PCMGSetLevels(pc, oldnlev, NULL));
245*9566063dSJacob Faibussowitsch     PetscCall(DMSetRefineLevel(ce->idm, oldlevel)); /* The damn DMCoarsen() calls in PCMG can reset this */
246f2ed2dc7SMatthew G. Knepley   }
247*9566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, ce->idm));
248*9566063dSJacob Faibussowitsch   PetscCall(DMTSSetBoundaryLocal(ce->idm, DMPlexTSComputeBoundary, ctx));
249*9566063dSJacob Faibussowitsch   PetscCall(DMTSSetIFunctionLocal(ce->idm, DMPlexTSComputeIFunctionFEM, ctx));
250*9566063dSJacob Faibussowitsch   PetscCall(DMTSSetIJacobianLocal(ce->idm, DMPlexTSComputeIJacobianFEM, ctx));
251*9566063dSJacob Faibussowitsch   PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_ITERATING));
252*9566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
253*9566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts, 0));
254*9566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
255*9566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, uInitial));
256*9566063dSJacob Faibussowitsch   PetscCall(PetscConvEstComputeInitialGuess(ce, 0, NULL, uInitial));
257f2ed2dc7SMatthew G. Knepley   PetscFunctionReturn(0);
258f2ed2dc7SMatthew G. Knepley }
259f2ed2dc7SMatthew G. Knepley 
260f2ed2dc7SMatthew G. Knepley PetscErrorCode PetscConvEstUseTS(PetscConvEst ce, PetscBool checkTemporal)
261900f6b5bSMatthew G. Knepley {
262900f6b5bSMatthew G. Knepley   PetscFunctionBegin;
263900f6b5bSMatthew G. Knepley   PetscValidHeaderSpecific(ce, PETSC_OBJECT_CLASSID, 1);
264900f6b5bSMatthew G. Knepley   ce->ops->setsolver     = PetscConvEstSetTS_Private;
265900f6b5bSMatthew G. Knepley   ce->ops->initguess     = PetscConvEstInitGuessTS_Private;
266900f6b5bSMatthew G. Knepley   ce->ops->computeerror  = PetscConvEstComputeErrorTS_Private;
267f2ed2dc7SMatthew G. Knepley   if (checkTemporal) {
268f2ed2dc7SMatthew G. Knepley     ce->ops->getconvrate = PetscConvEstGetConvRateTS_Temporal_Private;
269f2ed2dc7SMatthew G. Knepley   } else {
270f2ed2dc7SMatthew G. Knepley     ce->ops->getconvrate = PetscConvEstGetConvRateTS_Spatial_Private;
271f2ed2dc7SMatthew G. Knepley   }
272900f6b5bSMatthew G. Knepley   PetscFunctionReturn(0);
273900f6b5bSMatthew G. Knepley }
274