Actual source code: adaptcfl.c
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: PetscReal hcfl, cfltimestep, ccfl;
6: PetscInt ncandidates;
7: const PetscReal *ccflarray;
9: PetscFunctionBegin;
10: PetscCall(TSGetCFLTime(ts, &cfltimestep));
11: PetscCall(TSAdaptCandidatesGet(adapt, &ncandidates, NULL, NULL, &ccflarray, NULL));
12: ccfl = (ncandidates > 0) ? ccflarray[0] : 1.0;
14: PetscCheck(adapt->always_accept, PetscObjectComm((PetscObject)adapt), PETSC_ERR_SUP, "Step rejection not implemented. The CFL implementation is incomplete/unusable");
16: /* Determine whether the step is accepted of rejected */
17: *accept = PETSC_TRUE;
18: if (h > cfltimestep * ccfl) {
19: if (adapt->always_accept) {
20: PetscCall(PetscInfo(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));
21: } else {
22: PetscCall(PetscInfo(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));
23: *accept = PETSC_FALSE;
24: }
25: }
27: /* The optimal new step based purely on CFL constraint for this step. */
28: hcfl = adapt->safety * cfltimestep * ccfl;
29: if (hcfl < adapt->dt_min) {
30: PetscCall(PetscInfo(adapt, "Cannot satisfy CFL constraint %g (with %g safety) at minimum time step %g with method coefficient %g, proceeding anyway\n", (double)cfltimestep, (double)adapt->safety, (double)adapt->dt_min, (double)ccfl));
31: }
33: *next_sc = 0;
34: *next_h = PetscClipInterval(hcfl, adapt->dt_min, adapt->dt_max);
35: *wlte = -1; /* Weighted local truncation error was not evaluated */
36: *wltea = -1; /* Weighted absolute local truncation error was not evaluated */
37: *wlter = -1; /* Weighted relative local truncation error was not evaluated */
38: PetscFunctionReturn(PETSC_SUCCESS);
39: }
41: /*MC
42: TSADAPTCFL - CFL adaptive controller for time stepping
44: Level: intermediate
46: .seealso: [](ch_ts), `TS`, `TSAdapt`, `TSGetAdapt()`, `TSAdaptType`
47: M*/
48: PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt adapt)
49: {
50: PetscFunctionBegin;
51: adapt->ops->choose = TSAdaptChoose_CFL;
52: PetscFunctionReturn(PETSC_SUCCESS);
53: }