Actual source code: tsconvest.c
petsc-3.13.3 2020-07-01
1: #include <petscconvest.h>
2: #include <petscts.h>
4: #include <petsc/private/petscconvestimpl.h>
6: static PetscErrorCode PetscConvEstSetTS_Private(PetscConvEst ce, PetscObject solver)
7: {
8: PetscClassId id;
12: PetscObjectGetClassId(ce->solver, &id);
13: if (id != TS_CLASSID) SETERRQ(PetscObjectComm((PetscObject) ce), PETSC_ERR_ARG_WRONG, "Solver was not a TS");
14: TSGetDM((TS) ce->solver, &ce->idm);
15: return(0);
16: }
18: static PetscErrorCode PetscConvEstInitGuessTS_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u)
19: {
23: TSComputeInitialCondition((TS) ce->solver, u);
24: return(0);
25: }
27: static PetscErrorCode PetscConvEstComputeErrorTS_Private(PetscConvEst ce, PetscInt r, DM dm, Vec u, PetscReal errors[])
28: {
29: Vec e;
33: VecDuplicate(u, &e);
34: TSComputeExactError((TS) ce->solver, u, e);
35: VecNorm(e, NORM_2, errors);
36: VecDestroy(&e);
37: return(0);
38: }
40: static PetscErrorCode PetscConvEstGetConvRateTS_Private(PetscConvEst ce, PetscReal alpha[])
41: {
42: TS ts = (TS) ce->solver;
43: Vec u;
44: PetscReal *dt, *x, *y, slope, intercept;
45: PetscInt Ns, oNs, Nf = ce->Nf, f, Nr = ce->Nr, r;
49: TSGetSolution(ts, &u);
50: PetscMalloc1(Nr+1, &dt);
51: TSGetTimeStep(ts, &dt[0]);
52: TSGetMaxSteps(ts, &oNs);
53: Ns = oNs;
54: for (r = 0; r <= Nr; ++r) {
55: if (r > 0) {
56: dt[r] = dt[r-1]/ce->r;
57: Ns = PetscCeilReal(Ns*ce->r);
58: }
59: TSSetTime(ts, 0.0);
60: TSSetStepNumber(ts, 0);
61: TSSetTimeStep(ts, dt[r]);
62: TSSetMaxSteps(ts, Ns);
63: PetscConvEstComputeInitialGuess(ce, r, NULL, u);
64: TSSolve(ts, u);
65: PetscLogEventBegin(ce->event, ce, 0, 0, 0);
66: PetscConvEstComputeError(ce, r, NULL, u, &ce->errors[r*Nf]);
67: PetscLogEventEnd(ce->event, ce, 0, 0, 0);
68: for (f = 0; f < ce->Nf; ++f) {
69: PetscLogEventSetDof(ce->event, f, 1.0/dt[r]);
70: PetscLogEventSetError(ce->event, f, ce->errors[r*Nf+f]);
71: }
72: }
73: /* Fit convergence rate */
74: if (Nr) {
75: PetscMalloc2(Nr+1, &x, Nr+1, &y);
76: for (f = 0; f < Nf; ++f) {
77: for (r = 0; r <= Nr; ++r) {
78: x[r] = PetscLog10Real(dt[r*Nf+f]);
79: y[r] = PetscLog10Real(ce->errors[r*Nf+f]);
80: }
81: PetscLinearRegression(Nr+1, x, y, &slope, &intercept);
82: /* Since lg err = s lg dt + b */
83: alpha[f] = slope;
84: }
85: PetscFree2(x, y);
86: }
87: /* Reset solver */
88: TSSetConvergedReason(ts, TS_CONVERGED_ITERATING);
89: TSSetTime(ts, 0.0);
90: TSSetStepNumber(ts, 0);
91: TSSetTimeStep(ts, dt[0]);
92: TSSetMaxSteps(ts, oNs);
93: PetscConvEstComputeInitialGuess(ce, 0, NULL, u);
94: PetscFree(dt);
95: return(0);
96: }
98: PetscErrorCode PetscConvEstUseTS(PetscConvEst ce)
99: {
102: ce->ops->setsolver = PetscConvEstSetTS_Private;
103: ce->ops->initguess = PetscConvEstInitGuessTS_Private;
104: ce->ops->computeerror = PetscConvEstComputeErrorTS_Private;
105: ce->ops->getconvrate = PetscConvEstGetConvRateTS_Private;
106: return(0);
107: }