xref: /petsc/src/ts/adapt/impls/basic/adaptbasic.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
286e171c2SStefano Zampini #include <petscdm.h>
384df9cb4SJed Brown 
4*9371c9d4SSatish Balay static PetscErrorCode TSAdaptChoose_Basic(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept, PetscReal *wlte, PetscReal *wltea, PetscReal *wlter) {
586e171c2SStefano Zampini   Vec       Y;
686e171c2SStefano Zampini   DM        dm;
7b1f5bccaSLisandro Dalcin   PetscInt  order = PETSC_DECIDE;
8b1f5bccaSLisandro Dalcin   PetscReal enorm = -1;
97453f775SEmil Constantinescu   PetscReal enorma, enormr;
101917a363SLisandro Dalcin   PetscReal safety = adapt->safety;
11b1f5bccaSLisandro Dalcin   PetscReal hfac_lte, h_lte;
1284df9cb4SJed Brown 
1384df9cb4SJed Brown   PetscFunctionBegin;
141566a47fSLisandro Dalcin   *next_sc = 0;  /* Reuse the same order scheme */
15bf997491SLisandro Dalcin   *wltea   = -1; /* Weighted absolute local truncation error is not used */
16bf997491SLisandro Dalcin   *wlter   = -1; /* Weighted relative local truncation error is not used */
171566a47fSLisandro Dalcin 
181ebf93c6SLisandro Dalcin   if (ts->ops->evaluatewlte) {
199566063dSJacob Faibussowitsch     PetscCall(TSEvaluateWLTE(ts, adapt->wnormtype, &order, &enorm));
20cad9d221SBarry Smith     PetscCheck(enorm < 0 || order >= 1, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_OUTOFRANGE, "Computed error order %" PetscInt_FMT " must be positive", order);
211566a47fSLisandro Dalcin   } else if (ts->ops->evaluatestep) {
2208401ef6SPierre Jolivet     PetscCheck(adapt->candidates.n >= 1, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONGSTATE, "No candidate has been registered");
2363a3b9bcSJacob Faibussowitsch     PetscCheck(adapt->candidates.inuse_set, PetscObjectComm((PetscObject)adapt), PETSC_ERR_ARG_WRONGSTATE, "The current in-use scheme is not among the %" PetscInt_FMT " candidates", adapt->candidates.n);
241566a47fSLisandro Dalcin     order = adapt->candidates.order[0];
259566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts, &dm));
269566063dSJacob Faibussowitsch     PetscCall(DMGetGlobalVector(dm, &Y));
279566063dSJacob Faibussowitsch     PetscCall(TSEvaluateStep(ts, order - 1, Y, NULL));
289566063dSJacob Faibussowitsch     PetscCall(TSErrorWeightedNorm(ts, ts->vec_sol, Y, adapt->wnormtype, &enorm, &enorma, &enormr));
299566063dSJacob Faibussowitsch     PetscCall(DMRestoreGlobalVector(dm, &Y));
301ebf93c6SLisandro Dalcin   }
311c3436cfSJed Brown 
321566a47fSLisandro Dalcin   if (enorm < 0) {
331566a47fSLisandro Dalcin     *accept = PETSC_TRUE;
341566a47fSLisandro Dalcin     *next_h = h;  /* Reuse the old step */
351566a47fSLisandro Dalcin     *wlte   = -1; /* Weighted local truncation error was not evaluated */
361566a47fSLisandro Dalcin     PetscFunctionReturn(0);
371566a47fSLisandro Dalcin   }
381566a47fSLisandro Dalcin 
391566a47fSLisandro Dalcin   /* Determine whether the step is accepted of rejected */
401566a47fSLisandro Dalcin   if (enorm > 1) {
411917a363SLisandro Dalcin     if (!*accept) safety *= adapt->reject_safety; /* The last attempt also failed, shorten more aggressively */
42fd94acc0SJed Brown     if (h < (1 + PETSC_SQRT_MACHINE_EPSILON) * adapt->dt_min) {
439566063dSJacob Faibussowitsch       PetscCall(PetscInfo(adapt, "Estimated scaled local truncation error %g, accepting because step size %g is at minimum\n", (double)enorm, (double)h));
44fd94acc0SJed Brown       *accept = PETSC_TRUE;
45bf997491SLisandro Dalcin     } else if (adapt->always_accept) {
469566063dSJacob Faibussowitsch       PetscCall(PetscInfo(adapt, "Estimated scaled local truncation error %g, accepting step of size %g because always_accept is set\n", (double)enorm, (double)h));
477e1ba4dcSJed Brown       *accept = PETSC_TRUE;
481c3436cfSJed Brown     } else {
499566063dSJacob Faibussowitsch       PetscCall(PetscInfo(adapt, "Estimated scaled local truncation error %g, rejecting step of size %g\n", (double)enorm, (double)h));
50fd94acc0SJed Brown       *accept = PETSC_FALSE;
51fd94acc0SJed Brown     }
52fd94acc0SJed Brown   } else {
539566063dSJacob Faibussowitsch     PetscCall(PetscInfo(adapt, "Estimated scaled local truncation error %g, accepting step of size %g\n", (double)enorm, (double)h));
54fd94acc0SJed Brown     *accept = PETSC_TRUE;
551c3436cfSJed Brown   }
561c3436cfSJed Brown 
571c3436cfSJed Brown   /* The optimal new step based purely on local truncation error for this step. */
58*9371c9d4SSatish Balay   if (enorm > 0) hfac_lte = safety * PetscPowReal(enorm, ((PetscReal)-1) / order);
59*9371c9d4SSatish Balay   else hfac_lte = safety * PETSC_INFINITY;
60de50f1caSBarry Smith   if (adapt->timestepjustdecreased) {
61a191cbb8SBarry Smith     hfac_lte = PetscMin(hfac_lte, 1.0);
62de50f1caSBarry Smith     adapt->timestepjustdecreased--;
63a191cbb8SBarry Smith   }
641917a363SLisandro Dalcin   h_lte = h * PetscClipInterval(hfac_lte, adapt->clip[0], adapt->clip[1]);
651c3436cfSJed Brown 
661c3436cfSJed Brown   *next_h = PetscClipInterval(h_lte, adapt->dt_min, adapt->dt_max);
670b99f514SJed Brown   *wlte   = enorm;
6884df9cb4SJed Brown   PetscFunctionReturn(0);
6984df9cb4SJed Brown }
7084df9cb4SJed Brown 
7184df9cb4SJed Brown /*MC
7284df9cb4SJed Brown    TSADAPTBASIC - Basic adaptive controller for time stepping
7384df9cb4SJed Brown 
7484df9cb4SJed Brown    Level: intermediate
7584df9cb4SJed Brown 
76db781477SPatrick Sanan .seealso: `TS`, `TSAdapt`, `TSGetAdapt()`
7784df9cb4SJed Brown M*/
78*9371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt adapt) {
7984df9cb4SJed Brown   PetscFunctionBegin;
801c3436cfSJed Brown   adapt->ops->choose = TSAdaptChoose_Basic;
811b9b13dfSLisandro Dalcin   PetscFunctionReturn(0);
821b9b13dfSLisandro Dalcin }
83