Actual source code: adaptcfl.c

petsc-3.3-p7 2013-05-11
  1: #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/

  3: typedef struct {
  4:   PetscBool always_accept;
  5:   PetscReal safety;             /* safety factor relative to target error */
  6: } TSAdapt_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;

 20:   TSGetTimeStepNumber(ts,&stepno);
 21:   TSGetCFLTime(ts,&cfltime);
 22:   TSAdaptCandidatesGet(adapt,&ncandidates,&order,PETSC_NULL,&ccfl,PETSC_NULL);

 24:   hcfl = cfl->safety * cfltime * ccfl[0];
 25:   if (hcfl < adapt->dt_min) {
 26:     PetscInfo4(adapt,"Cannot satisfy CFL constraint %G (with %G safety) at minimum time step %G with method coefficient %G, proceding anyway\n",cfltime,cfl->safety,adapt->dt_min,ccfl[0]);
 27:   }

 29:   if (h > cfltime * ccfl[0]) {
 30:     if (cfl->always_accept) {
 31:       PetscInfo3(adapt,"Step length %G with scheme of CFL coefficient %G did not satisfy user-provided CFL constraint %G, proceeding anyway\n",h,ccfl[0],cfltime);
 32:     } else {
 33:       PetscInfo3(adapt,"Step length %G with scheme of CFL coefficient %G did not satisfy user-provided CFL constraint %G, step REJECTED\n",h,ccfl[0],cfltime);
 34:       *next_sc = 0;
 35:       *next_h = PetscClipInterval(hcfl,adapt->dt_min,adapt->dt_max);
 36:       *accept = PETSC_FALSE;
 37:     }
 38:   }

 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:   return(0);
 45: }

 49: static PetscErrorCode TSAdaptDestroy_CFL(TSAdapt adapt)
 50: {

 54:   PetscFree(adapt->data);
 55:   return(0);
 56: }

 60: static PetscErrorCode TSAdaptSetFromOptions_CFL(TSAdapt adapt)
 61: {
 62:   TSAdapt_CFL  *cfl = (TSAdapt_CFL*)adapt->data;

 66:   PetscOptionsHead("CFL adaptive controller options");
 67:   PetscOptionsReal("-ts_adapt_cfl_safety","Safety factor relative to target error","",cfl->safety,&cfl->safety,PETSC_NULL);
 68:   PetscOptionsBool("-ts_adapt_cfl_always_accept","Always accept the step regardless of whether local truncation error meets goal","",cfl->always_accept,&cfl->always_accept,PETSC_NULL);
 69:   if (!cfl->always_accept) SETERRQ(((PetscObject)adapt)->comm,PETSC_ERR_SUP,"step rejection not implemented yet");
 70:   PetscOptionsTail();
 71:   return(0);
 72: }

 74: EXTERN_C_BEGIN
 77: /*MC
 78:    TSADAPTCFL - CFL adaptive controller for time stepping

 80:    Level: intermediate

 82: .seealso: TS, TSAdapt, TSSetAdapt()
 83: M*/
 84: PetscErrorCode TSAdaptCreate_CFL(TSAdapt adapt)
 85: {
 87:   TSAdapt_CFL *a;

 90:   PetscNewLog(adapt,TSAdapt_CFL,&a);
 91:   adapt->data = (void*)a;
 92:   adapt->ops->choose         = TSAdaptChoose_CFL;
 93:   adapt->ops->setfromoptions = TSAdaptSetFromOptions_CFL;
 94:   adapt->ops->destroy        = TSAdaptDestroy_CFL;

 96:   a->safety        = 0.9;
 97:   a->always_accept = PETSC_FALSE;
 98:   return(0);
 99: }
100: EXTERN_C_END