1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
28d59e960SJed Brown
TSAdaptChoose_CFL(TSAdapt adapt,TS ts,PetscReal h,PetscInt * next_sc,PetscReal * next_h,PetscBool * accept,PetscReal * wlte,PetscReal * wltea,PetscReal * wlter)3d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSAdaptChoose_CFL(TSAdapt adapt, TS ts, PetscReal h, PetscInt *next_sc, PetscReal *next_h, PetscBool *accept, PetscReal *wlte, PetscReal *wltea, PetscReal *wlter)
4d71ae5a4SJacob Faibussowitsch {
51566a47fSLisandro Dalcin PetscReal hcfl, cfltimestep, ccfl;
61566a47fSLisandro Dalcin PetscInt ncandidates;
71566a47fSLisandro Dalcin const PetscReal *ccflarray;
88d59e960SJed Brown
98d59e960SJed Brown PetscFunctionBegin;
109566063dSJacob Faibussowitsch PetscCall(TSGetCFLTime(ts, &cfltimestep));
119566063dSJacob Faibussowitsch PetscCall(TSAdaptCandidatesGet(adapt, &ncandidates, NULL, NULL, &ccflarray, NULL));
121566a47fSLisandro Dalcin ccfl = (ncandidates > 0) ? ccflarray[0] : 1.0;
138d59e960SJed Brown
143c633725SBarry Smith PetscCheck(adapt->always_accept, PetscObjectComm((PetscObject)adapt), PETSC_ERR_SUP, "Step rejection not implemented. The CFL implementation is incomplete/unusable");
151917a363SLisandro Dalcin
161566a47fSLisandro Dalcin /* Determine whether the step is accepted of rejected */
171566a47fSLisandro Dalcin *accept = PETSC_TRUE;
181566a47fSLisandro Dalcin if (h > cfltimestep * ccfl) {
19bf997491SLisandro Dalcin if (adapt->always_accept) {
209566063dSJacob Faibussowitsch 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));
218d59e960SJed Brown } else {
229566063dSJacob Faibussowitsch 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));
238d59e960SJed Brown *accept = PETSC_FALSE;
248d59e960SJed Brown }
258d59e960SJed Brown }
268d59e960SJed Brown
271566a47fSLisandro Dalcin /* The optimal new step based purely on CFL constraint for this step. */
281917a363SLisandro Dalcin hcfl = adapt->safety * cfltimestep * ccfl;
291566a47fSLisandro Dalcin if (hcfl < adapt->dt_min) {
30aaa8cc7dSPierre Jolivet 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));
311566a47fSLisandro Dalcin }
321566a47fSLisandro Dalcin
338d59e960SJed Brown *next_sc = 0;
348d59e960SJed Brown *next_h = PetscClipInterval(hcfl, adapt->dt_min, adapt->dt_max);
350b99f514SJed Brown *wlte = -1; /* Weighted local truncation error was not evaluated */
365e4ed32fSEmil Constantinescu *wltea = -1; /* Weighted absolute local truncation error was not evaluated */
375e4ed32fSEmil Constantinescu *wlter = -1; /* Weighted relative local truncation error was not evaluated */
383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
398d59e960SJed Brown }
408d59e960SJed Brown
418d59e960SJed Brown /*MC
428d59e960SJed Brown TSADAPTCFL - CFL adaptive controller for time stepping
438d59e960SJed Brown
448d59e960SJed Brown Level: intermediate
458d59e960SJed Brown
46*efa39862SBarry Smith .seealso: [](ch_ts), [](sec_ts_error_control), `TS`, `TSAdapt`, `TSGetAdapt()`, `TSAdaptType`
478d59e960SJed Brown M*/
TSAdaptCreate_CFL(TSAdapt adapt)48d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt adapt)
49d71ae5a4SJacob Faibussowitsch {
508d59e960SJed Brown PetscFunctionBegin;
518d59e960SJed Brown adapt->ops->choose = TSAdaptChoose_CFL;
523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
538d59e960SJed Brown }
54