xref: /petsc/src/ts/adapt/impls/cfl/adaptcfl.c (revision 8d59e96071807dd8d22c4ae77f4b4dc641017dcf)
1*8d59e960SJed Brown #include <private/tsimpl.h> /*I "petscts.h" I*/
2*8d59e960SJed Brown 
3*8d59e960SJed Brown typedef struct {
4*8d59e960SJed Brown   PetscBool always_accept;
5*8d59e960SJed Brown   PetscReal safety;             /* safety factor relative to target error */
6*8d59e960SJed Brown } TSAdapt_CFL;
7*8d59e960SJed Brown 
8*8d59e960SJed Brown #undef __FUNCT__
9*8d59e960SJed Brown #define __FUNCT__ "TSAdaptChoose_CFL"
10*8d59e960SJed Brown static PetscErrorCode TSAdaptChoose_CFL(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept)
11*8d59e960SJed Brown {
12*8d59e960SJed Brown   TSAdapt_CFL     *cfl = (TSAdapt_CFL*)adapt->data;
13*8d59e960SJed Brown   PetscErrorCode  ierr;
14*8d59e960SJed Brown   PetscReal       hcfl,cfltime;
15*8d59e960SJed Brown   PetscInt        stepno,ncandidates;
16*8d59e960SJed Brown   const PetscInt  *order;
17*8d59e960SJed Brown   const PetscReal *ccfl;
18*8d59e960SJed Brown 
19*8d59e960SJed Brown   PetscFunctionBegin;
20*8d59e960SJed Brown   ierr = TSGetTimeStepNumber(ts,&stepno);CHKERRQ(ierr);
21*8d59e960SJed Brown   ierr = TSGetCFLTime(ts,&cfltime);CHKERRQ(ierr);
22*8d59e960SJed Brown   ierr = TSAdaptCandidatesGet(adapt,&ncandidates,&order,PETSC_NULL,&ccfl,PETSC_NULL);CHKERRQ(ierr);
23*8d59e960SJed Brown 
24*8d59e960SJed Brown   hcfl = cfl->safety * cfltime * ccfl[0];
25*8d59e960SJed Brown   if (hcfl < adapt->dt_min) {
26*8d59e960SJed Brown     ierr = PetscInfo4(adapt,"Cannot satisfy CFL constraint %G (with %G safety) at minimum time step %G with method coefficient %G, proceding anyway\n",cfltime,cfl->safety,adapt->dt_min,ccfl[0]);CHKERRQ(ierr);
27*8d59e960SJed Brown   }
28*8d59e960SJed Brown 
29*8d59e960SJed Brown   if (h > cfltime * ccfl[0]) {
30*8d59e960SJed Brown     if (cfl->always_accept) {
31*8d59e960SJed Brown       ierr = PetscInfo3(adapt,"Step length %G with scheme of CFL coefficient %G did not satisfy user-provided CFL constraint %G, proceeding anyway\n",h,ccfl[0],cfltime);CHKERRQ(ierr);
32*8d59e960SJed Brown     } else {
33*8d59e960SJed Brown       ierr = PetscInfo3(adapt,"Step length %G with scheme of CFL coefficient %G did not satisfy user-provided CFL constraint %G, step REJECTED\n",h,ccfl[0],cfltime);CHKERRQ(ierr);
34*8d59e960SJed Brown       *next_sc = 0;
35*8d59e960SJed Brown       *next_h = PetscClipInterval(hcfl,adapt->dt_min,adapt->dt_max);
36*8d59e960SJed Brown       *accept = PETSC_FALSE;
37*8d59e960SJed Brown     }
38*8d59e960SJed Brown   }
39*8d59e960SJed Brown 
40*8d59e960SJed Brown   *next_sc = 0;
41*8d59e960SJed Brown   *next_h = PetscClipInterval(hcfl,adapt->dt_min,adapt->dt_max);
42*8d59e960SJed Brown   *accept = PETSC_TRUE;
43*8d59e960SJed Brown   PetscFunctionReturn(0);
44*8d59e960SJed Brown }
45*8d59e960SJed Brown 
46*8d59e960SJed Brown #undef __FUNCT__
47*8d59e960SJed Brown #define __FUNCT__ "TSAdaptDestroy_CFL"
48*8d59e960SJed Brown static PetscErrorCode TSAdaptDestroy_CFL(TSAdapt adapt)
49*8d59e960SJed Brown {
50*8d59e960SJed Brown   PetscErrorCode ierr;
51*8d59e960SJed Brown 
52*8d59e960SJed Brown   PetscFunctionBegin;
53*8d59e960SJed Brown   ierr = PetscFree(adapt->data);CHKERRQ(ierr);
54*8d59e960SJed Brown   PetscFunctionReturn(0);
55*8d59e960SJed Brown }
56*8d59e960SJed Brown 
57*8d59e960SJed Brown #undef __FUNCT__
58*8d59e960SJed Brown #define __FUNCT__ "TSAdaptSetFromOptions_CFL"
59*8d59e960SJed Brown static PetscErrorCode TSAdaptSetFromOptions_CFL(TSAdapt adapt)
60*8d59e960SJed Brown {
61*8d59e960SJed Brown   TSAdapt_CFL  *cfl = (TSAdapt_CFL*)adapt->data;
62*8d59e960SJed Brown   PetscErrorCode ierr;
63*8d59e960SJed Brown 
64*8d59e960SJed Brown   PetscFunctionBegin;
65*8d59e960SJed Brown   ierr = PetscOptionsHead("CFL adaptive controller options");CHKERRQ(ierr);
66*8d59e960SJed Brown   ierr = PetscOptionsReal("-ts_adapt_cfl_safety","Safety factor relative to target error","",cfl->safety,&cfl->safety,PETSC_NULL);CHKERRQ(ierr);
67*8d59e960SJed Brown   ierr = PetscOptionsBool("-ts_adapt_cfl_always_accept","Always accept the step regardless of whether local truncation error meets goal","",cfl->always_accept,&cfl->always_accept,PETSC_NULL);CHKERRQ(ierr);
68*8d59e960SJed Brown   if (!cfl->always_accept) SETERRQ(((PetscObject)adapt)->comm,PETSC_ERR_SUP,"step rejection not implemented yet");
69*8d59e960SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
70*8d59e960SJed Brown   PetscFunctionReturn(0);
71*8d59e960SJed Brown }
72*8d59e960SJed Brown 
73*8d59e960SJed Brown EXTERN_C_BEGIN
74*8d59e960SJed Brown #undef __FUNCT__
75*8d59e960SJed Brown #define __FUNCT__ "TSAdaptCreate_CFL"
76*8d59e960SJed Brown /*MC
77*8d59e960SJed Brown    TSADAPTCFL - CFL adaptive controller for time stepping
78*8d59e960SJed Brown 
79*8d59e960SJed Brown    Level: intermediate
80*8d59e960SJed Brown 
81*8d59e960SJed Brown .seealso: TS, TSAdapt, TSSetAdapt()
82*8d59e960SJed Brown M*/
83*8d59e960SJed Brown PetscErrorCode TSAdaptCreate_CFL(TSAdapt adapt)
84*8d59e960SJed Brown {
85*8d59e960SJed Brown   PetscErrorCode ierr;
86*8d59e960SJed Brown   TSAdapt_CFL *a;
87*8d59e960SJed Brown 
88*8d59e960SJed Brown   PetscFunctionBegin;
89*8d59e960SJed Brown   ierr = PetscNewLog(adapt,TSAdapt_CFL,&a);CHKERRQ(ierr);
90*8d59e960SJed Brown   adapt->data = (void*)a;
91*8d59e960SJed Brown   adapt->ops->choose         = TSAdaptChoose_CFL;
92*8d59e960SJed Brown   adapt->ops->setfromoptions = TSAdaptSetFromOptions_CFL;
93*8d59e960SJed Brown   adapt->ops->destroy        = TSAdaptDestroy_CFL;
94*8d59e960SJed Brown 
95*8d59e960SJed Brown   a->safety        = 0.9;
96*8d59e960SJed Brown   a->always_accept = PETSC_FALSE;
97*8d59e960SJed Brown   PetscFunctionReturn(0);
98*8d59e960SJed Brown }
99*8d59e960SJed Brown EXTERN_C_END
100