xref: /petsc/src/ts/adapt/impls/basic/adaptbasic.c (revision 1917a363a0108e5a01ea00f477c9a5d39d8d3ccc)
1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
284df9cb4SJed Brown 
384df9cb4SJed Brown typedef struct {
41c3436cfSJed Brown   Vec Y;
584df9cb4SJed Brown } TSAdapt_Basic;
684df9cb4SJed Brown 
75e4ed32fSEmil Constantinescu static PetscErrorCode TSAdaptChoose_Basic(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte,PetscReal *wltea,PetscReal *wlter)
884df9cb4SJed Brown {
91c3436cfSJed Brown   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;
10b1f5bccaSLisandro Dalcin   PetscInt       order  = PETSC_DECIDE;
11b1f5bccaSLisandro Dalcin   PetscReal      enorm  = -1;
127453f775SEmil Constantinescu   PetscReal      enorma,enormr;
13*1917a363SLisandro Dalcin   PetscReal      safety = adapt->safety;
14b1f5bccaSLisandro Dalcin   PetscReal      hfac_lte,h_lte;
151c3436cfSJed Brown   PetscErrorCode ierr;
1684df9cb4SJed Brown 
1784df9cb4SJed Brown   PetscFunctionBegin;
181566a47fSLisandro Dalcin   *next_sc = 0;   /* Reuse the same order scheme */
19bf997491SLisandro Dalcin   *wltea   = -1;  /* Weighted absolute local truncation error is not used */
20bf997491SLisandro Dalcin   *wlter   = -1;  /* Weighted relative local truncation error is not used */
211566a47fSLisandro Dalcin 
221ebf93c6SLisandro Dalcin   if (ts->ops->evaluatewlte) {
231566a47fSLisandro Dalcin     ierr = TSEvaluateWLTE(ts,adapt->wnormtype,&order,&enorm);CHKERRQ(ierr);
241566a47fSLisandro Dalcin     if (enorm >= 0 && order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed error order %D must be positive",order);
251566a47fSLisandro Dalcin   } else if (ts->ops->evaluatestep) {
261566a47fSLisandro Dalcin     if (adapt->candidates.n < 1) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"No candidate has been registered");
271566a47fSLisandro Dalcin     if (!adapt->candidates.inuse_set) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"The current in-use scheme is not among the %D candidates",adapt->candidates.n);
281566a47fSLisandro Dalcin     if (!basic->Y) {ierr = VecDuplicate(ts->vec_sol,&basic->Y);CHKERRQ(ierr);}
291566a47fSLisandro Dalcin     order = adapt->candidates.order[0];
301ebf93c6SLisandro Dalcin     ierr = TSEvaluateStep(ts,order-1,basic->Y,NULL);CHKERRQ(ierr);
317453f775SEmil Constantinescu     ierr = TSErrorWeightedNorm(ts,ts->vec_sol,basic->Y,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr);
321ebf93c6SLisandro Dalcin   }
331c3436cfSJed Brown 
341566a47fSLisandro Dalcin   if (enorm < 0) {
351566a47fSLisandro Dalcin     *accept  = PETSC_TRUE;
361566a47fSLisandro Dalcin     *next_h  = h;            /* Reuse the old step */
371566a47fSLisandro Dalcin     *wlte    = -1;           /* Weighted local truncation error was not evaluated */
381566a47fSLisandro Dalcin     PetscFunctionReturn(0);
391566a47fSLisandro Dalcin   }
401566a47fSLisandro Dalcin 
411566a47fSLisandro Dalcin   /* Determine whether the step is accepted of rejected */
421566a47fSLisandro Dalcin   if (enorm > 1) {
43*1917a363SLisandro Dalcin     if (!*accept) safety *= adapt->reject_safety; /* The last attempt also failed, shorten more aggressively */
44fd94acc0SJed Brown     if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) {
4557622a8eSBarry Smith       ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting because step size %g is at minimum\n",(double)enorm,(double)h);CHKERRQ(ierr);
46fd94acc0SJed Brown       *accept = PETSC_TRUE;
47bf997491SLisandro Dalcin     } else if (adapt->always_accept) {
4857622a8eSBarry Smith       ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g because always_accept is set\n",(double)enorm,(double)h);CHKERRQ(ierr);
497e1ba4dcSJed Brown       *accept = PETSC_TRUE;
501c3436cfSJed Brown     } else {
5157622a8eSBarry Smith       ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, rejecting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr);
52fd94acc0SJed Brown       *accept = PETSC_FALSE;
53fd94acc0SJed Brown     }
54fd94acc0SJed Brown   } else {
5557622a8eSBarry Smith     ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr);
56fd94acc0SJed Brown     *accept = PETSC_TRUE;
571c3436cfSJed Brown   }
581c3436cfSJed Brown 
591c3436cfSJed Brown   /* The optimal new step based purely on local truncation error for this step. */
601566a47fSLisandro Dalcin   if (enorm > 0)
611566a47fSLisandro Dalcin     hfac_lte = safety * PetscPowReal(enorm,((PetscReal)-1)/order);
621566a47fSLisandro Dalcin   else
63bf0c4ff7SBarry Smith     hfac_lte = safety * PETSC_INFINITY;
64*1917a363SLisandro 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 
71099af0a3SLisandro Dalcin static PetscErrorCode TSAdaptReset_Basic(TSAdapt adapt)
7284df9cb4SJed Brown {
731c3436cfSJed Brown   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;
7484df9cb4SJed Brown   PetscErrorCode ierr;
7584df9cb4SJed Brown 
7684df9cb4SJed Brown   PetscFunctionBegin;
771c3436cfSJed Brown   ierr = VecDestroy(&basic->Y);CHKERRQ(ierr);
78099af0a3SLisandro Dalcin   PetscFunctionReturn(0);
79099af0a3SLisandro Dalcin }
80099af0a3SLisandro Dalcin 
81099af0a3SLisandro Dalcin static PetscErrorCode TSAdaptDestroy_Basic(TSAdapt adapt)
82099af0a3SLisandro Dalcin {
83099af0a3SLisandro Dalcin   PetscErrorCode ierr;
84099af0a3SLisandro Dalcin 
85099af0a3SLisandro Dalcin   PetscFunctionBegin;
86099af0a3SLisandro Dalcin   ierr = TSAdaptReset_Basic(adapt);CHKERRQ(ierr);
8784df9cb4SJed Brown   ierr = PetscFree(adapt->data);CHKERRQ(ierr);
8884df9cb4SJed Brown   PetscFunctionReturn(0);
8984df9cb4SJed Brown }
9084df9cb4SJed Brown 
9184df9cb4SJed Brown /*MC
9284df9cb4SJed Brown    TSADAPTBASIC - Basic adaptive controller for time stepping
9384df9cb4SJed Brown 
9484df9cb4SJed Brown    Level: intermediate
9584df9cb4SJed Brown 
9684df9cb4SJed Brown .seealso: TS, TSAdapt, TSSetAdapt()
9784df9cb4SJed Brown M*/
988cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt adapt)
9984df9cb4SJed Brown {
10084df9cb4SJed Brown   PetscErrorCode ierr;
101*1917a363SLisandro Dalcin   TSAdapt_Basic  *basic;
10284df9cb4SJed Brown 
10384df9cb4SJed Brown   PetscFunctionBegin;
104*1917a363SLisandro Dalcin   ierr= PetscNewLog(adapt,&basic);CHKERRQ(ierr);
105*1917a363SLisandro Dalcin   adapt->data         = (void*)basic;
1061c3436cfSJed Brown   adapt->ops->choose  = TSAdaptChoose_Basic;
107*1917a363SLisandro Dalcin   adapt->ops->reset   = TSAdaptReset_Basic;
10884df9cb4SJed Brown   adapt->ops->destroy = TSAdaptDestroy_Basic;
1091b9b13dfSLisandro Dalcin   PetscFunctionReturn(0);
1101b9b13dfSLisandro Dalcin }
111