Actual source code: adaptbasic.c

  1: #include <petsc/private/tsimpl.h>
  2: #include <petscdm.h>

  4: static PetscErrorCode TSAdaptChoose_Basic(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte,PetscReal *wltea,PetscReal *wlter)
  5: {
  6:   Vec            Y;
  7:   DM             dm;
  8:   PetscInt       order  = PETSC_DECIDE;
  9:   PetscReal      enorm  = -1;
 10:   PetscReal      enorma,enormr;
 11:   PetscReal      safety = adapt->safety;
 12:   PetscReal      hfac_lte,h_lte;

 14:   *next_sc = 0;   /* Reuse the same order scheme */
 15:   *wltea   = -1;  /* Weighted absolute local truncation error is not used */
 16:   *wlter   = -1;  /* Weighted relative local truncation error is not used */

 18:   if (ts->ops->evaluatewlte) {
 19:     TSEvaluateWLTE(ts,adapt->wnormtype,&order,&enorm);
 21:   } else if (ts->ops->evaluatestep) {
 24:     order = adapt->candidates.order[0];
 25:     TSGetDM(ts,&dm);
 26:     DMGetGlobalVector(dm,&Y);
 27:     TSEvaluateStep(ts,order-1,Y,NULL);
 28:     TSErrorWeightedNorm(ts,ts->vec_sol,Y,adapt->wnormtype,&enorm,&enorma,&enormr);
 29:     DMRestoreGlobalVector(dm,&Y);
 30:   }

 32:   if (enorm < 0) {
 33:     *accept  = PETSC_TRUE;
 34:     *next_h  = h;            /* Reuse the old step */
 35:     *wlte    = -1;           /* Weighted local truncation error was not evaluated */
 36:     return 0;
 37:   }

 39:   /* Determine whether the step is accepted of rejected */
 40:   if (enorm > 1) {
 41:     if (!*accept) safety *= adapt->reject_safety; /* The last attempt also failed, shorten more aggressively */
 42:     if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) {
 43:       PetscInfo(adapt,"Estimated scaled local truncation error %g, accepting because step size %g is at minimum\n",(double)enorm,(double)h);
 44:       *accept = PETSC_TRUE;
 45:     } else if (adapt->always_accept) {
 46:       PetscInfo(adapt,"Estimated scaled local truncation error %g, accepting step of size %g because always_accept is set\n",(double)enorm,(double)h);
 47:       *accept = PETSC_TRUE;
 48:     } else {
 49:       PetscInfo(adapt,"Estimated scaled local truncation error %g, rejecting step of size %g\n",(double)enorm,(double)h);
 50:       *accept = PETSC_FALSE;
 51:     }
 52:   } else {
 53:     PetscInfo(adapt,"Estimated scaled local truncation error %g, accepting step of size %g\n",(double)enorm,(double)h);
 54:     *accept = PETSC_TRUE;
 55:   }

 57:   /* The optimal new step based purely on local truncation error for this step. */
 58:   if (enorm > 0)
 59:     hfac_lte = safety * PetscPowReal(enorm,((PetscReal)-1)/order);
 60:   else
 61:     hfac_lte = safety * PETSC_INFINITY;
 62:   if (adapt->timestepjustdecreased) {
 63:     hfac_lte = PetscMin(hfac_lte,1.0);
 64:     adapt->timestepjustdecreased--;
 65:   }
 66:   h_lte = h * PetscClipInterval(hfac_lte,adapt->clip[0],adapt->clip[1]);

 68:   *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max);
 69:   *wlte   = enorm;
 70:   return 0;
 71: }

 73: /*MC
 74:    TSADAPTBASIC - Basic adaptive controller for time stepping

 76:    Level: intermediate

 78: .seealso: TS, TSAdapt, TSGetAdapt()
 79: M*/
 80: PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt adapt)
 81: {
 82:   adapt->ops->choose = TSAdaptChoose_Basic;
 83:   return 0;
 84: }