xref: /petsc/src/ts/adapt/impls/cfl/adaptcfl.c (revision 1566a47f5f7ddc4d41547989035e2da0a3965616)
1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
28d59e960SJed Brown 
38d59e960SJed Brown typedef struct {
4*1566a47fSLisandro Dalcin   PetscReal safety;         /* safety factor relative to target CFL constraint */
58d59e960SJed Brown   PetscBool always_accept;
68d59e960SJed Brown } TSAdapt_CFL;
78d59e960SJed Brown 
88d59e960SJed Brown #undef __FUNCT__
98d59e960SJed Brown #define __FUNCT__ "TSAdaptChoose_CFL"
100b99f514SJed Brown static PetscErrorCode TSAdaptChoose_CFL(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte)
118d59e960SJed Brown {
128d59e960SJed Brown   TSAdapt_CFL     *cfl = (TSAdapt_CFL*)adapt->data;
138d59e960SJed Brown   PetscErrorCode  ierr;
14*1566a47fSLisandro Dalcin   PetscReal       hcfl,cfltimestep,ccfl;
15*1566a47fSLisandro Dalcin   PetscInt        ncandidates;
16*1566a47fSLisandro Dalcin   const PetscReal *ccflarray;
178d59e960SJed Brown 
188d59e960SJed Brown   PetscFunctionBegin;
19*1566a47fSLisandro Dalcin   ierr = TSGetCFLTime(ts,&cfltimestep);CHKERRQ(ierr);
20*1566a47fSLisandro Dalcin   ierr = TSAdaptCandidatesGet(adapt,&ncandidates,NULL,NULL,&ccflarray,NULL);CHKERRQ(ierr);
21*1566a47fSLisandro Dalcin   ccfl = (ncandidates > 0) ? ccflarray[0] : 1.0;
228d59e960SJed Brown 
23*1566a47fSLisandro Dalcin   /* Determine whether the step is accepted of rejected */
24*1566a47fSLisandro Dalcin   *accept = PETSC_TRUE;
25*1566a47fSLisandro Dalcin   if (h > cfltimestep * ccfl) {
268d59e960SJed Brown     if (cfl->always_accept) {
27*1566a47fSLisandro 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);
288d59e960SJed Brown     } else {
29*1566a47fSLisandro 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);
308d59e960SJed Brown       *accept  = PETSC_FALSE;
318d59e960SJed Brown     }
328d59e960SJed Brown   }
338d59e960SJed Brown 
34*1566a47fSLisandro Dalcin   /* The optimal new step based purely on CFL constraint for this step. */
35*1566a47fSLisandro Dalcin   hcfl = cfl->safety * cfltimestep * ccfl;
36*1566a47fSLisandro Dalcin   if (hcfl < adapt->dt_min) {
37*1566a47fSLisandro 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);
38*1566a47fSLisandro Dalcin   }
39*1566a47fSLisandro Dalcin 
408d59e960SJed Brown   *next_sc = 0;
418d59e960SJed Brown   *next_h  = PetscClipInterval(hcfl,adapt->dt_min,adapt->dt_max);
420b99f514SJed Brown   *wlte    = -1;  /* Weighted local truncation error was not evaluated */
438d59e960SJed Brown   PetscFunctionReturn(0);
448d59e960SJed Brown }
458d59e960SJed Brown 
468d59e960SJed Brown #undef __FUNCT__
478d59e960SJed Brown #define __FUNCT__ "TSAdaptDestroy_CFL"
488d59e960SJed Brown static PetscErrorCode TSAdaptDestroy_CFL(TSAdapt adapt)
498d59e960SJed Brown {
508d59e960SJed Brown   PetscErrorCode ierr;
518d59e960SJed Brown 
528d59e960SJed Brown   PetscFunctionBegin;
538d59e960SJed Brown   ierr = PetscFree(adapt->data);CHKERRQ(ierr);
548d59e960SJed Brown   PetscFunctionReturn(0);
558d59e960SJed Brown }
568d59e960SJed Brown 
578d59e960SJed Brown #undef __FUNCT__
588d59e960SJed Brown #define __FUNCT__ "TSAdaptSetFromOptions_CFL"
594416b707SBarry Smith static PetscErrorCode TSAdaptSetFromOptions_CFL(PetscOptionItems *PetscOptionsObject,TSAdapt adapt)
608d59e960SJed Brown {
618d59e960SJed Brown   TSAdapt_CFL    *cfl = (TSAdapt_CFL*)adapt->data;
628d59e960SJed Brown   PetscErrorCode ierr;
638d59e960SJed Brown 
648d59e960SJed Brown   PetscFunctionBegin;
65e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"CFL adaptive controller options");CHKERRQ(ierr);
66*1566a47fSLisandro Dalcin   ierr = PetscOptionsReal("-ts_adapt_cfl_safety","Safety factor relative to target CFL constraint","",cfl->safety,&cfl->safety,NULL);CHKERRQ(ierr);
67*1566a47fSLisandro Dalcin   ierr = PetscOptionsBool("-ts_adapt_cfl_always_accept","Always accept the step regardless of whether CFL constraint meets goal","",cfl->always_accept,&cfl->always_accept,NULL);CHKERRQ(ierr);
688d59e960SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
698d59e960SJed Brown   PetscFunctionReturn(0);
708d59e960SJed Brown }
718d59e960SJed Brown 
728d59e960SJed Brown #undef __FUNCT__
738d59e960SJed Brown #define __FUNCT__ "TSAdaptCreate_CFL"
748d59e960SJed Brown /*MC
758d59e960SJed Brown    TSADAPTCFL - CFL adaptive controller for time stepping
768d59e960SJed Brown 
778d59e960SJed Brown    Level: intermediate
788d59e960SJed Brown 
798d59e960SJed Brown .seealso: TS, TSAdapt, TSSetAdapt()
808d59e960SJed Brown M*/
818cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt adapt)
828d59e960SJed Brown {
838d59e960SJed Brown   PetscErrorCode ierr;
848d59e960SJed Brown   TSAdapt_CFL    *a;
858d59e960SJed Brown 
868d59e960SJed Brown   PetscFunctionBegin;
87b00a9115SJed Brown   ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr);
888d59e960SJed Brown   adapt->data = (void*)a;
89*1566a47fSLisandro Dalcin 
908d59e960SJed Brown   adapt->ops->choose         = TSAdaptChoose_CFL;
918d59e960SJed Brown   adapt->ops->setfromoptions = TSAdaptSetFromOptions_CFL;
928d59e960SJed Brown   adapt->ops->destroy        = TSAdaptDestroy_CFL;
938d59e960SJed Brown 
948d59e960SJed Brown   a->safety        = 0.9;
958d59e960SJed Brown   a->always_accept = PETSC_FALSE;
968d59e960SJed Brown   PetscFunctionReturn(0);
978d59e960SJed Brown }
98