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