xref: /petsc/src/ts/adapt/impls/cfl/adaptcfl.c (revision e6e75211d226c622f451867f53ce5d558649ff4f)
1 #include <petsc/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,PetscReal *wlte)
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,NULL,&ccfl,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",(double)cfltime,(double)cfl->safety,(double)adapt->dt_min,(double)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",(double)h,(double)ccfl[0],(double)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",(double)h,(double)ccfl[0],(double)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   *wlte    = -1;                /* Weighted local truncation error was not evaluated */
44   PetscFunctionReturn(0);
45 }
46 
47 #undef __FUNCT__
48 #define __FUNCT__ "TSAdaptDestroy_CFL"
49 static PetscErrorCode TSAdaptDestroy_CFL(TSAdapt adapt)
50 {
51   PetscErrorCode ierr;
52 
53   PetscFunctionBegin;
54   ierr = PetscFree(adapt->data);CHKERRQ(ierr);
55   PetscFunctionReturn(0);
56 }
57 
58 #undef __FUNCT__
59 #define __FUNCT__ "TSAdaptSetFromOptions_CFL"
60 static PetscErrorCode TSAdaptSetFromOptions_CFL(PetscOptions *PetscOptionsObject,TSAdapt adapt)
61 {
62   TSAdapt_CFL    *cfl = (TSAdapt_CFL*)adapt->data;
63   PetscErrorCode ierr;
64 
65   PetscFunctionBegin;
66   ierr = PetscOptionsHead(PetscOptionsObject,"CFL adaptive controller options");CHKERRQ(ierr);
67   ierr = PetscOptionsReal("-ts_adapt_cfl_safety","Safety factor relative to target error","",cfl->safety,&cfl->safety,NULL);CHKERRQ(ierr);
68   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,NULL);CHKERRQ(ierr);
69   if (!cfl->always_accept) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"step rejection not implemented yet");
70   ierr = PetscOptionsTail();CHKERRQ(ierr);
71   PetscFunctionReturn(0);
72 }
73 
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 PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt adapt)
84 {
85   PetscErrorCode ierr;
86   TSAdapt_CFL    *a;
87 
88   PetscFunctionBegin;
89   ierr                       = PetscNewLog(adapt,&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