Actual source code: adaptcfl.c
petsc-3.8.4 2018-03-24
1: #include <petsc/private/tsimpl.h>
3: static PetscErrorCode TSAdaptChoose_CFL(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte,PetscReal *wltea,PetscReal *wlter)
4: {
5: PetscErrorCode ierr;
6: PetscReal hcfl,cfltimestep,ccfl;
7: PetscInt ncandidates;
8: const PetscReal *ccflarray;
11: TSGetCFLTime(ts,&cfltimestep);
12: TSAdaptCandidatesGet(adapt,&ncandidates,NULL,NULL,&ccflarray,NULL);
13: ccfl = (ncandidates > 0) ? ccflarray[0] : 1.0;
15: if (!adapt->always_accept) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Step rejection not implemented. The CFL implementation is incomplete/unusable");
17: /* Determine whether the step is accepted of rejected */
18: *accept = PETSC_TRUE;
19: if (h > cfltimestep * ccfl) {
20: if (adapt->always_accept) {
21: 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);
22: } else {
23: 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);
24: *accept = PETSC_FALSE;
25: }
26: }
28: /* The optimal new step based purely on CFL constraint for this step. */
29: hcfl = adapt->safety * cfltimestep * ccfl;
30: if (hcfl < adapt->dt_min) {
31: 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)adapt->safety,(double)adapt->dt_min,(double)ccfl);
32: }
34: *next_sc = 0;
35: *next_h = PetscClipInterval(hcfl,adapt->dt_min,adapt->dt_max);
36: *wlte = -1; /* Weighted local truncation error was not evaluated */
37: *wltea = -1; /* Weighted absolute local truncation error was not evaluated */
38: *wlter = -1; /* Weighted relative local truncation error was not evaluated */
39: return(0);
40: }
42: /*MC
43: TSADAPTCFL - CFL adaptive controller for time stepping
45: Level: intermediate
47: .seealso: TS, TSAdapt, TSGetAdapt()
48: M*/
49: PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt adapt)
50: {
52: adapt->ops->choose = TSAdaptChoose_CFL;
53: return(0);
54: }