xref: /petsc/src/ts/adapt/impls/cfl/adaptcfl.c (revision bf997491f620a56e3f945815f638946a12caa405)
1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
28d59e960SJed Brown 
38d59e960SJed Brown typedef struct {
41566a47fSLisandro Dalcin   PetscReal safety;         /* safety factor relative to target CFL constraint */
58d59e960SJed Brown } TSAdapt_CFL;
68d59e960SJed Brown 
75e4ed32fSEmil Constantinescu static PetscErrorCode TSAdaptChoose_CFL(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte,PetscReal *wltea,PetscReal *wlter)
88d59e960SJed Brown {
98d59e960SJed Brown   TSAdapt_CFL     *cfl = (TSAdapt_CFL*)adapt->data;
108d59e960SJed Brown   PetscErrorCode  ierr;
111566a47fSLisandro Dalcin   PetscReal       hcfl,cfltimestep,ccfl;
121566a47fSLisandro Dalcin   PetscInt        ncandidates;
131566a47fSLisandro Dalcin   const PetscReal *ccflarray;
148d59e960SJed Brown 
158d59e960SJed Brown   PetscFunctionBegin;
161566a47fSLisandro Dalcin   ierr = TSGetCFLTime(ts,&cfltimestep);CHKERRQ(ierr);
171566a47fSLisandro Dalcin   ierr = TSAdaptCandidatesGet(adapt,&ncandidates,NULL,NULL,&ccflarray,NULL);CHKERRQ(ierr);
181566a47fSLisandro Dalcin   ccfl = (ncandidates > 0) ? ccflarray[0] : 1.0;
198d59e960SJed Brown 
201566a47fSLisandro Dalcin   /* Determine whether the step is accepted of rejected */
211566a47fSLisandro Dalcin   *accept = PETSC_TRUE;
221566a47fSLisandro Dalcin   if (h > cfltimestep * ccfl) {
23*bf997491SLisandro Dalcin     if (adapt->always_accept) {
241566a47fSLisandro Dalcin       ierr = PetscInfo3(adapt,"Step length %g with scheme of CFL coefficient %g did not satisfy user-provided CFL constraint %g, proceeding anyway\n",(double)h,(double)ccfl,(double)cfltimestep);CHKERRQ(ierr);
258d59e960SJed Brown     } else {
261566a47fSLisandro Dalcin       ierr = PetscInfo3(adapt,"Step length %g with scheme of CFL coefficient %g did not satisfy user-provided CFL constraint %g, step REJECTED\n",(double)h,(double)ccfl,(double)cfltimestep);CHKERRQ(ierr);
278d59e960SJed Brown       *accept  = PETSC_FALSE;
288d59e960SJed Brown     }
298d59e960SJed Brown   }
308d59e960SJed Brown 
311566a47fSLisandro Dalcin   /* The optimal new step based purely on CFL constraint for this step. */
321566a47fSLisandro Dalcin   hcfl = cfl->safety * cfltimestep * ccfl;
331566a47fSLisandro Dalcin   if (hcfl < adapt->dt_min) {
341566a47fSLisandro Dalcin     ierr = PetscInfo4(adapt,"Cannot satisfy CFL constraint %g (with %g safety) at minimum time step %g with method coefficient %g, proceding anyway\n",(double)cfltimestep,(double)cfl->safety,(double)adapt->dt_min,(double)ccfl);CHKERRQ(ierr);
351566a47fSLisandro Dalcin   }
361566a47fSLisandro Dalcin 
378d59e960SJed Brown   *next_sc = 0;
388d59e960SJed Brown   *next_h  = PetscClipInterval(hcfl,adapt->dt_min,adapt->dt_max);
390b99f514SJed Brown   *wlte    = -1;   /* Weighted local truncation error was not evaluated */
405e4ed32fSEmil Constantinescu   *wltea   = -1;   /* Weighted absolute local truncation error was not evaluated */
415e4ed32fSEmil Constantinescu   *wlter   = -1;   /* Weighted relative local truncation error was not evaluated */
428d59e960SJed Brown   PetscFunctionReturn(0);
438d59e960SJed Brown }
448d59e960SJed Brown 
458d59e960SJed Brown static PetscErrorCode TSAdaptDestroy_CFL(TSAdapt adapt)
468d59e960SJed Brown {
478d59e960SJed Brown   PetscErrorCode ierr;
488d59e960SJed Brown 
498d59e960SJed Brown   PetscFunctionBegin;
508d59e960SJed Brown   ierr = PetscFree(adapt->data);CHKERRQ(ierr);
518d59e960SJed Brown   PetscFunctionReturn(0);
528d59e960SJed Brown }
538d59e960SJed Brown 
544416b707SBarry Smith static PetscErrorCode TSAdaptSetFromOptions_CFL(PetscOptionItems *PetscOptionsObject,TSAdapt adapt)
558d59e960SJed Brown {
568d59e960SJed Brown   TSAdapt_CFL    *cfl = (TSAdapt_CFL*)adapt->data;
578d59e960SJed Brown   PetscErrorCode ierr;
588d59e960SJed Brown 
598d59e960SJed Brown   PetscFunctionBegin;
60e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"CFL adaptive controller options");CHKERRQ(ierr);
611566a47fSLisandro Dalcin   ierr = PetscOptionsReal("-ts_adapt_cfl_safety","Safety factor relative to target CFL constraint","",cfl->safety,&cfl->safety,NULL);CHKERRQ(ierr);
62c524d1a0SJed Brown   /* The TS implementations do not currently communicate CFL information to the controller.  There is a placeholder, but
63c524d1a0SJed Brown    * we do not believe it to provide sufficiently rich information.  That means the CFL adaptor is incomplete and
64c524d1a0SJed Brown    * unusable.  Do not delete the guard below unless you have finished the implementation. */
65*bf997491SLisandro Dalcin   if (!adapt->always_accept) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Step rejection not implemented. The CFL implementation is incomplete/unusable");
668d59e960SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
678d59e960SJed Brown   PetscFunctionReturn(0);
688d59e960SJed Brown }
698d59e960SJed Brown 
708d59e960SJed Brown /*MC
718d59e960SJed Brown    TSADAPTCFL - CFL adaptive controller for time stepping
728d59e960SJed Brown 
738d59e960SJed Brown    Level: intermediate
748d59e960SJed Brown 
758d59e960SJed Brown .seealso: TS, TSAdapt, TSSetAdapt()
768d59e960SJed Brown M*/
778cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt adapt)
788d59e960SJed Brown {
798d59e960SJed Brown   PetscErrorCode ierr;
808d59e960SJed Brown   TSAdapt_CFL    *a;
818d59e960SJed Brown 
828d59e960SJed Brown   PetscFunctionBegin;
83b00a9115SJed Brown   ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr);
848d59e960SJed Brown   adapt->data = (void*)a;
851566a47fSLisandro Dalcin 
868d59e960SJed Brown   adapt->ops->choose         = TSAdaptChoose_CFL;
878d59e960SJed Brown   adapt->ops->setfromoptions = TSAdaptSetFromOptions_CFL;
888d59e960SJed Brown   adapt->ops->destroy        = TSAdaptDestroy_CFL;
898d59e960SJed Brown 
908d59e960SJed Brown   a->safety = 0.9;
918d59e960SJed Brown   PetscFunctionReturn(0);
928d59e960SJed Brown }
93