xref: /petsc/src/ts/adapt/impls/cfl/adaptcfl.c (revision 9566063d113dddea24716c546802770db7481bc0)
1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
28d59e960SJed Brown 
35e4ed32fSEmil Constantinescu static PetscErrorCode TSAdaptChoose_CFL(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte,PetscReal *wltea,PetscReal *wlter)
48d59e960SJed Brown {
51566a47fSLisandro Dalcin   PetscReal       hcfl,cfltimestep,ccfl;
61566a47fSLisandro Dalcin   PetscInt        ncandidates;
71566a47fSLisandro Dalcin   const PetscReal *ccflarray;
88d59e960SJed Brown 
98d59e960SJed Brown   PetscFunctionBegin;
10*9566063dSJacob Faibussowitsch   PetscCall(TSGetCFLTime(ts,&cfltimestep));
11*9566063dSJacob 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) {
20*9566063dSJacob 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 {
22*9566063dSJacob 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) {
30*9566063dSJacob Faibussowitsch     PetscCall(PetscInfo(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));
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 */
388d59e960SJed Brown   PetscFunctionReturn(0);
398d59e960SJed Brown }
408d59e960SJed Brown 
418d59e960SJed Brown /*MC
428d59e960SJed Brown    TSADAPTCFL - CFL adaptive controller for time stepping
438d59e960SJed Brown 
448d59e960SJed Brown    Level: intermediate
458d59e960SJed Brown 
46e91eccc2SStefano Zampini .seealso: TS, TSAdapt, TSGetAdapt()
478d59e960SJed Brown M*/
488cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt adapt)
498d59e960SJed Brown {
508d59e960SJed Brown   PetscFunctionBegin;
518d59e960SJed Brown   adapt->ops->choose = TSAdaptChoose_CFL;
528d59e960SJed Brown   PetscFunctionReturn(0);
538d59e960SJed Brown }
54