Actual source code: adaptglee.c

  1: #include <petsc/private/tsimpl.h>
  2: #include <petscdm.h>

  4: typedef struct {
  5:   Vec Y;
  6: } TSAdapt_GLEE;

  8: static PetscErrorCode TSAdaptChoose_GLEE(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept, PetscReal *wlte, PetscReal *wltea, PetscReal *wlter)
  9: {
 10:   TSAdapt_GLEE *glee = (TSAdapt_GLEE *)adapt->data;
 11:   Vec           X, Y, E;
 12:   PetscReal     enorm, enorma, enormr, hfac_lte, hfac_ltea, hfac_lter, h_lte, safety;
 13:   PetscInt      order;
 14:   PetscBool     bGTEMethod;

 16:   PetscFunctionBegin;
 17:   *next_sc = 0; /* Reuse the same order scheme */
 18:   safety   = adapt->safety;
 19:   PetscCall(PetscObjectTypeCompare((PetscObject)ts, TSGLEE, &bGTEMethod));
 20:   order = adapt->candidates.order[0];

 22:   if (bGTEMethod) { /* the method is of GLEE type */
 23:     DM dm;

 25:     PetscCall(TSGetSolution(ts, &X));
 26:     if (!glee->Y && adapt->glee_use_local) PetscCall(VecDuplicate(X, &glee->Y)); /*create vector to store previous step global error*/
 27:     PetscCall(TSGetDM(ts, &dm));
 28:     PetscCall(DMGetGlobalVector(dm, &E));
 29:     PetscCall(TSGetTimeError(ts, 0, &E));

 31:     if (adapt->glee_use_local) PetscCall(VecAXPY(E, -1.0, glee->Y)); /* local error = current error - previous step error */

 33:     /* this should be called with the solution at the beginning of the step too*/
 34:     PetscCall(TSErrorWeightedENorm(ts, E, X, X, adapt->wnormtype, &enorm, &enorma, &enormr));
 35:     PetscCall(DMRestoreGlobalVector(dm, &E));
 36:   } else {
 37:     /* the method is NOT of GLEE type; use the stantard basic augmented by separate atol and rtol */
 38:     PetscCall(TSGetSolution(ts, &X));
 39:     if (!glee->Y) PetscCall(VecDuplicate(X, &glee->Y));
 40:     Y = glee->Y;
 41:     PetscCall(TSEvaluateStep(ts, order - 1, Y, NULL));
 42:     PetscCall(TSErrorWeightedNorm(ts, X, Y, adapt->wnormtype, &enorm, &enorma, &enormr));
 43:   }

 45:   if (enorm < 0) {
 46:     *accept = PETSC_TRUE;
 47:     *next_h = h;  /* Reuse the old step */
 48:     *wlte   = -1; /* Weighted error was not evaluated */
 49:     *wltea  = -1; /* Weighted absolute error was not evaluated */
 50:     *wlter  = -1; /* Weighted relative error was not evaluated */
 51:     PetscFunctionReturn(PETSC_SUCCESS);
 52:   }

 54:   if (enorm > 1. || enorma > 1. || enormr > 1.) {
 55:     if (!*accept) safety *= adapt->reject_safety; /* The last attempt also failed, shorten more aggressively */
 56:     if (h < (1 + PETSC_SQRT_MACHINE_EPSILON) * adapt->dt_min) {
 57:       PetscCall(PetscInfo(adapt, "Estimated scaled truncation error [combined, absolute, relative]] [%g, %g, %g], accepting because step size %g is at minimum\n", (double)enorm, (double)enorma, (double)enormr, (double)h));
 58:       *accept = PETSC_TRUE;
 59:     } else if (adapt->always_accept) {
 60:       PetscCall(PetscInfo(adapt, "Estimated scaled truncation error [combined, absolute, relative]] [%g, %g, %g], accepting step of size %g because always_accept is set\n", (double)enorm, (double)enorma, (double)enormr, (double)h));
 61:       *accept = PETSC_TRUE;
 62:     } else {
 63:       PetscCall(PetscInfo(adapt, "Estimated scaled truncation error [combined, absolute, relative]] [%g, %g, %g], rejecting step of size %g\n", (double)enorm, (double)enorma, (double)enormr, (double)h));
 64:       *accept = PETSC_FALSE;
 65:     }
 66:   } else {
 67:     PetscCall(PetscInfo(adapt, "Estimated scaled truncation error [combined, absolute, relative] [%g, %g, %g], accepting step of size %g\n", (double)enorm, (double)enorma, (double)enormr, (double)h));
 68:     *accept = PETSC_TRUE;
 69:   }

 71:   if (bGTEMethod) {
 72:     if (*accept == PETSC_TRUE && adapt->glee_use_local) {
 73:       /* If step is accepted, then overwrite previous step error with the current error to be used on the next step */
 74:       /* WARNING: if the adapters are composable, then the accept test will not be reliable*/
 75:       PetscCall(TSGetTimeError(ts, 0, &glee->Y));
 76:     }

 78:     /* The optimal new step based on the current global truncation error. */
 79:     if (enorm > 0) {
 80:       /* factor based on the absolute tolerance */
 81:       hfac_ltea = safety * PetscPowReal(1. / enorma, ((PetscReal)1) / (order + 1));
 82:       /* factor based on the relative tolerance */
 83:       hfac_lter = safety * PetscPowReal(1. / enormr, ((PetscReal)1) / (order + 1));
 84:       /* pick the minimum time step among the relative and absolute tolerances */
 85:       hfac_lte = PetscMin(hfac_ltea, hfac_lter);
 86:     } else {
 87:       hfac_lte = safety * PETSC_INFINITY;
 88:     }
 89:     h_lte   = h * PetscClipInterval(hfac_lte, adapt->clip[0], adapt->clip[1]);
 90:     *next_h = PetscClipInterval(h_lte, adapt->dt_min, adapt->dt_max);
 91:   } else {
 92:     /* The optimal new step based purely on local truncation error for this step. */
 93:     if (enorm > 0) {
 94:       /* factor based on the absolute tolerance */
 95:       hfac_ltea = safety * PetscPowReal(enorma, ((PetscReal)-1) / order);
 96:       /* factor based on the relative tolerance */
 97:       hfac_lter = safety * PetscPowReal(enormr, ((PetscReal)-1) / order);
 98:       /* pick the minimum time step among the relative and absolute tolerances */
 99:       hfac_lte = PetscMin(hfac_ltea, hfac_lter);
100:     } else {
101:       hfac_lte = safety * PETSC_INFINITY;
102:     }
103:     h_lte   = h * PetscClipInterval(hfac_lte, adapt->clip[0], adapt->clip[1]);
104:     *next_h = PetscClipInterval(h_lte, adapt->dt_min, adapt->dt_max);
105:   }
106:   *wlte  = enorm;
107:   *wltea = enorma;
108:   *wlter = enormr;
109:   PetscFunctionReturn(PETSC_SUCCESS);
110: }

112: static PetscErrorCode TSAdaptReset_GLEE(TSAdapt adapt)
113: {
114:   TSAdapt_GLEE *glee = (TSAdapt_GLEE *)adapt->data;

116:   PetscFunctionBegin;
117:   PetscCall(VecDestroy(&glee->Y));
118:   PetscFunctionReturn(PETSC_SUCCESS);
119: }

121: static PetscErrorCode TSAdaptDestroy_GLEE(TSAdapt adapt)
122: {
123:   PetscFunctionBegin;
124:   PetscCall(TSAdaptReset_GLEE(adapt));
125:   PetscCall(PetscFree(adapt->data));
126:   PetscFunctionReturn(PETSC_SUCCESS);
127: }

129: /*MC
130:    TSADAPTGLEE - GLEE adaptive controller for time stepping

132:    Level: intermediate

134: .seealso: [](ch_ts), [](sec_ts_error_control), `TS`, `TSAdapt`, `TSGetAdapt()`, `TSAdaptType`
135: M*/
136: PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt adapt)
137: {
138:   TSAdapt_GLEE *glee;

140:   PetscFunctionBegin;
141:   PetscCall(PetscNew(&glee));
142:   adapt->data         = (void *)glee;
143:   adapt->ops->choose  = TSAdaptChoose_GLEE;
144:   adapt->ops->reset   = TSAdaptReset_GLEE;
145:   adapt->ops->destroy = TSAdaptDestroy_GLEE;
146:   PetscFunctionReturn(PETSC_SUCCESS);
147: }