1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2*86e171c2SStefano Zampini #include <petscdm.h> 384df9cb4SJed Brown 45e4ed32fSEmil Constantinescu static PetscErrorCode TSAdaptChoose_Basic(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte,PetscReal *wltea,PetscReal *wlter) 584df9cb4SJed Brown { 6*86e171c2SStefano Zampini Vec Y; 7*86e171c2SStefano Zampini DM dm; 8b1f5bccaSLisandro Dalcin PetscInt order = PETSC_DECIDE; 9b1f5bccaSLisandro Dalcin PetscReal enorm = -1; 107453f775SEmil Constantinescu PetscReal enorma,enormr; 111917a363SLisandro Dalcin PetscReal safety = adapt->safety; 12b1f5bccaSLisandro Dalcin PetscReal hfac_lte,h_lte; 131c3436cfSJed Brown PetscErrorCode ierr; 1484df9cb4SJed Brown 1584df9cb4SJed Brown PetscFunctionBegin; 161566a47fSLisandro Dalcin *next_sc = 0; /* Reuse the same order scheme */ 17bf997491SLisandro Dalcin *wltea = -1; /* Weighted absolute local truncation error is not used */ 18bf997491SLisandro Dalcin *wlter = -1; /* Weighted relative local truncation error is not used */ 191566a47fSLisandro Dalcin 201ebf93c6SLisandro Dalcin if (ts->ops->evaluatewlte) { 211566a47fSLisandro Dalcin ierr = TSEvaluateWLTE(ts,adapt->wnormtype,&order,&enorm);CHKERRQ(ierr); 221566a47fSLisandro Dalcin if (enorm >= 0 && order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed error order %D must be positive",order); 231566a47fSLisandro Dalcin } else if (ts->ops->evaluatestep) { 241566a47fSLisandro Dalcin if (adapt->candidates.n < 1) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"No candidate has been registered"); 251566a47fSLisandro Dalcin if (!adapt->candidates.inuse_set) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"The current in-use scheme is not among the %D candidates",adapt->candidates.n); 261566a47fSLisandro Dalcin order = adapt->candidates.order[0]; 27*86e171c2SStefano Zampini ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 28*86e171c2SStefano Zampini ierr = DMGetGlobalVector(dm,&Y);CHKERRQ(ierr); 29*86e171c2SStefano Zampini ierr = TSEvaluateStep(ts,order-1,Y,NULL);CHKERRQ(ierr); 30*86e171c2SStefano Zampini ierr = TSErrorWeightedNorm(ts,ts->vec_sol,Y,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr); 31*86e171c2SStefano Zampini ierr = DMRestoreGlobalVector(dm,&Y);CHKERRQ(ierr); 321ebf93c6SLisandro Dalcin } 331c3436cfSJed Brown 341566a47fSLisandro Dalcin if (enorm < 0) { 351566a47fSLisandro Dalcin *accept = PETSC_TRUE; 361566a47fSLisandro Dalcin *next_h = h; /* Reuse the old step */ 371566a47fSLisandro Dalcin *wlte = -1; /* Weighted local truncation error was not evaluated */ 381566a47fSLisandro Dalcin PetscFunctionReturn(0); 391566a47fSLisandro Dalcin } 401566a47fSLisandro Dalcin 411566a47fSLisandro Dalcin /* Determine whether the step is accepted of rejected */ 421566a47fSLisandro Dalcin if (enorm > 1) { 431917a363SLisandro Dalcin if (!*accept) safety *= adapt->reject_safety; /* The last attempt also failed, shorten more aggressively */ 44fd94acc0SJed Brown if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) { 4557622a8eSBarry Smith ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting because step size %g is at minimum\n",(double)enorm,(double)h);CHKERRQ(ierr); 46fd94acc0SJed Brown *accept = PETSC_TRUE; 47bf997491SLisandro Dalcin } else if (adapt->always_accept) { 4857622a8eSBarry Smith ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g because always_accept is set\n",(double)enorm,(double)h);CHKERRQ(ierr); 497e1ba4dcSJed Brown *accept = PETSC_TRUE; 501c3436cfSJed Brown } else { 5157622a8eSBarry Smith ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, rejecting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr); 52fd94acc0SJed Brown *accept = PETSC_FALSE; 53fd94acc0SJed Brown } 54fd94acc0SJed Brown } else { 5557622a8eSBarry Smith ierr = PetscInfo2(adapt,"Estimated scaled local truncation error %g, accepting step of size %g\n",(double)enorm,(double)h);CHKERRQ(ierr); 56fd94acc0SJed Brown *accept = PETSC_TRUE; 571c3436cfSJed Brown } 581c3436cfSJed Brown 591c3436cfSJed Brown /* The optimal new step based purely on local truncation error for this step. */ 601566a47fSLisandro Dalcin if (enorm > 0) 611566a47fSLisandro Dalcin hfac_lte = safety * PetscPowReal(enorm,((PetscReal)-1)/order); 621566a47fSLisandro Dalcin else 63bf0c4ff7SBarry Smith hfac_lte = safety * PETSC_INFINITY; 64de50f1caSBarry Smith if (adapt->timestepjustdecreased){ 65a191cbb8SBarry Smith hfac_lte = PetscMin(hfac_lte,1.0); 66de50f1caSBarry Smith adapt->timestepjustdecreased--; 67a191cbb8SBarry Smith } 681917a363SLisandro Dalcin h_lte = h * PetscClipInterval(hfac_lte,adapt->clip[0],adapt->clip[1]); 691c3436cfSJed Brown 701c3436cfSJed Brown *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max); 710b99f514SJed Brown *wlte = enorm; 7284df9cb4SJed Brown PetscFunctionReturn(0); 7384df9cb4SJed Brown } 7484df9cb4SJed Brown 7584df9cb4SJed Brown /*MC 7684df9cb4SJed Brown TSADAPTBASIC - Basic adaptive controller for time stepping 7784df9cb4SJed Brown 7884df9cb4SJed Brown Level: intermediate 7984df9cb4SJed Brown 80e91eccc2SStefano Zampini .seealso: TS, TSAdapt, TSGetAdapt() 8184df9cb4SJed Brown M*/ 828cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt adapt) 8384df9cb4SJed Brown { 8484df9cb4SJed Brown PetscFunctionBegin; 851c3436cfSJed Brown adapt->ops->choose = TSAdaptChoose_Basic; 861b9b13dfSLisandro Dalcin PetscFunctionReturn(0); 871b9b13dfSLisandro Dalcin } 88