Actual source code: adaptbasic.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 clip[2];            /* admissible decrease/increase factors */
  6:   PetscReal safety;             /* safety factor relative to target error */
  7:   PetscReal reject_safety;      /* extra safety factor if the last step was rejected */
  8:   Vec       Y;
  9: } TSAdapt_Basic;

 13: static PetscErrorCode TSAdaptChoose_Basic(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte)
 14: {
 15:   TSAdapt_Basic     *basic = (TSAdapt_Basic*)adapt->data;
 16:   PetscErrorCode    ierr;
 17:   Vec               X,Y;
 18:   PetscReal         enorm,hfac_lte,h_lte,safety;
 19:   PetscInt          order,stepno;

 22:   TSGetTimeStepNumber(ts,&stepno);
 23:   TSGetSolution(ts,&X);
 24:   if (!basic->Y) {VecDuplicate(X,&basic->Y);}
 25:   Y = basic->Y;
 26:   order = adapt->candidates.order[0];
 27:   TSEvaluateStep(ts,order-1,Y,PETSC_NULL);

 29:   safety = basic->safety;
 30:   TSErrorNormWRMS(ts,Y,&enorm);
 31:   if (enorm > 1.) {
 32:     if (!*accept) safety *= basic->reject_safety; /* The last attempt also failed, shorten more aggressively */
 33:     if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) {
 34:       PetscInfo2(adapt,"Estimated scaled local truncation error %G, accepting because step size %G is at minimum\n",enorm,h);
 35:       *accept = PETSC_TRUE;
 36:     } else if (basic->always_accept) {
 37:       PetscInfo2(adapt,"Estimated scaled local truncation error %G, accepting step of size %G because always_accept is set\n",enorm,h);
 38:       *accept = PETSC_TRUE;
 39:     } else {
 40:       PetscInfo2(adapt,"Estimated scaled local truncation error %G, rejecting step of size %G\n",enorm,h);
 41:       *accept = PETSC_FALSE;
 42:     }
 43:   } else {
 44:     PetscInfo2(adapt,"Estimated scaled local truncation error %G, accepting step of size %G\n",enorm,h);
 45:     *accept = PETSC_TRUE;
 46:   }

 48:   /* The optimal new step based purely on local truncation error for this step. */
 49:   hfac_lte = safety * PetscRealPart(PetscPowScalar((PetscScalar)enorm,(PetscReal)(-1./order)));
 50:   h_lte = h * PetscClipInterval(hfac_lte,basic->clip[0],basic->clip[1]);

 52:   *next_sc = 0;
 53:   *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max);
 54:   *wlte = enorm;
 55:   return(0);
 56: }

 60: static PetscErrorCode TSAdaptDestroy_Basic(TSAdapt adapt)
 61: {
 62:   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;

 66:   VecDestroy(&basic->Y);
 67:   PetscFree(adapt->data);
 68:   return(0);
 69: }

 73: static PetscErrorCode TSAdaptSetFromOptions_Basic(TSAdapt adapt)
 74: {
 75:   TSAdapt_Basic  *basic = (TSAdapt_Basic*)adapt->data;
 77:   PetscInt       two;
 78:   PetscBool      set;

 81:   PetscOptionsHead("Basic adaptive controller options");
 82:   two = 2;
 83:   PetscOptionsRealArray("-ts_adapt_basic_clip","Admissible decrease/increase in step size","",basic->clip,&two,&set);
 84:   if (set && (two != 2 || basic->clip[0] > basic->clip[1]))
 85:     SETERRQ(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Must give exactly two values to -ts_adapt_basic_clip");
 86:   PetscOptionsReal("-ts_adapt_basic_safety","Safety factor relative to target error","",basic->safety,&basic->safety,PETSC_NULL);
 87:   PetscOptionsReal("-ts_adapt_basic_reject_safety","Extra safety factor to apply if the last step was rejected","",basic->reject_safety,&basic->reject_safety,PETSC_NULL);
 88:   PetscOptionsBool("-ts_adapt_basic_always_accept","Always accept the step regardless of whether local truncation error meets goal","",basic->always_accept,&basic->always_accept,PETSC_NULL);
 89:   PetscOptionsTail();
 90:   return(0);
 91: }

 93: EXTERN_C_BEGIN
 96: /*MC
 97:    TSADAPTBASIC - Basic adaptive controller for time stepping

 99:    Level: intermediate

101: .seealso: TS, TSAdapt, TSSetAdapt()
102: M*/
103: PetscErrorCode TSAdaptCreate_Basic(TSAdapt adapt)
104: {
106:   TSAdapt_Basic *a;

109:   PetscNewLog(adapt,TSAdapt_Basic,&a);
110:   adapt->data = (void*)a;
111:   adapt->ops->choose         = TSAdaptChoose_Basic;
112:   adapt->ops->setfromoptions = TSAdaptSetFromOptions_Basic;
113:   adapt->ops->destroy        = TSAdaptDestroy_Basic;

115:   a->clip[0]       = 0.1;
116:   a->clip[1]       = 10.;
117:   a->safety        = 0.9;
118:   a->reject_safety = 0.5;
119:   a->always_accept = PETSC_FALSE;
120:   return(0);
121: }
122: EXTERN_C_END