1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 286e171c2SStefano 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 { 686e171c2SStefano Zampini Vec Y; 786e171c2SStefano 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; 1384df9cb4SJed Brown 1484df9cb4SJed Brown PetscFunctionBegin; 151566a47fSLisandro Dalcin *next_sc = 0; /* Reuse the same order scheme */ 16bf997491SLisandro Dalcin *wltea = -1; /* Weighted absolute local truncation error is not used */ 17bf997491SLisandro Dalcin *wlter = -1; /* Weighted relative local truncation error is not used */ 181566a47fSLisandro Dalcin 191ebf93c6SLisandro Dalcin if (ts->ops->evaluatewlte) { 20*9566063dSJacob Faibussowitsch PetscCall(TSEvaluateWLTE(ts,adapt->wnormtype,&order,&enorm)); 212c71b3e2SJacob Faibussowitsch PetscCheckFalse(enorm >= 0 && order < 1,PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed error order %D must be positive",order); 221566a47fSLisandro Dalcin } else if (ts->ops->evaluatestep) { 232c71b3e2SJacob Faibussowitsch PetscCheckFalse(adapt->candidates.n < 1,PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"No candidate has been registered"); 243c633725SBarry Smith PetscCheck(adapt->candidates.inuse_set,PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"The current in-use scheme is not among the %D candidates",adapt->candidates.n); 251566a47fSLisandro Dalcin order = adapt->candidates.order[0]; 26*9566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 27*9566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm,&Y)); 28*9566063dSJacob Faibussowitsch PetscCall(TSEvaluateStep(ts,order-1,Y,NULL)); 29*9566063dSJacob Faibussowitsch PetscCall(TSErrorWeightedNorm(ts,ts->vec_sol,Y,adapt->wnormtype,&enorm,&enorma,&enormr)); 30*9566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm,&Y)); 311ebf93c6SLisandro Dalcin } 321c3436cfSJed Brown 331566a47fSLisandro Dalcin if (enorm < 0) { 341566a47fSLisandro Dalcin *accept = PETSC_TRUE; 351566a47fSLisandro Dalcin *next_h = h; /* Reuse the old step */ 361566a47fSLisandro Dalcin *wlte = -1; /* Weighted local truncation error was not evaluated */ 371566a47fSLisandro Dalcin PetscFunctionReturn(0); 381566a47fSLisandro Dalcin } 391566a47fSLisandro Dalcin 401566a47fSLisandro Dalcin /* Determine whether the step is accepted of rejected */ 411566a47fSLisandro Dalcin if (enorm > 1) { 421917a363SLisandro Dalcin if (!*accept) safety *= adapt->reject_safety; /* The last attempt also failed, shorten more aggressively */ 43fd94acc0SJed Brown if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) { 44*9566063dSJacob Faibussowitsch PetscCall(PetscInfo(adapt,"Estimated scaled local truncation error %g, accepting because step size %g is at minimum\n",(double)enorm,(double)h)); 45fd94acc0SJed Brown *accept = PETSC_TRUE; 46bf997491SLisandro Dalcin } else if (adapt->always_accept) { 47*9566063dSJacob Faibussowitsch PetscCall(PetscInfo(adapt,"Estimated scaled local truncation error %g, accepting step of size %g because always_accept is set\n",(double)enorm,(double)h)); 487e1ba4dcSJed Brown *accept = PETSC_TRUE; 491c3436cfSJed Brown } else { 50*9566063dSJacob Faibussowitsch PetscCall(PetscInfo(adapt,"Estimated scaled local truncation error %g, rejecting step of size %g\n",(double)enorm,(double)h)); 51fd94acc0SJed Brown *accept = PETSC_FALSE; 52fd94acc0SJed Brown } 53fd94acc0SJed Brown } else { 54*9566063dSJacob Faibussowitsch PetscCall(PetscInfo(adapt,"Estimated scaled local truncation error %g, accepting step of size %g\n",(double)enorm,(double)h)); 55fd94acc0SJed Brown *accept = PETSC_TRUE; 561c3436cfSJed Brown } 571c3436cfSJed Brown 581c3436cfSJed Brown /* The optimal new step based purely on local truncation error for this step. */ 591566a47fSLisandro Dalcin if (enorm > 0) 601566a47fSLisandro Dalcin hfac_lte = safety * PetscPowReal(enorm,((PetscReal)-1)/order); 611566a47fSLisandro Dalcin else 62bf0c4ff7SBarry Smith hfac_lte = safety * PETSC_INFINITY; 63de50f1caSBarry Smith if (adapt->timestepjustdecreased) { 64a191cbb8SBarry Smith hfac_lte = PetscMin(hfac_lte,1.0); 65de50f1caSBarry Smith adapt->timestepjustdecreased--; 66a191cbb8SBarry Smith } 671917a363SLisandro Dalcin h_lte = h * PetscClipInterval(hfac_lte,adapt->clip[0],adapt->clip[1]); 681c3436cfSJed Brown 691c3436cfSJed Brown *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max); 700b99f514SJed Brown *wlte = enorm; 7184df9cb4SJed Brown PetscFunctionReturn(0); 7284df9cb4SJed Brown } 7384df9cb4SJed Brown 7484df9cb4SJed Brown /*MC 7584df9cb4SJed Brown TSADAPTBASIC - Basic adaptive controller for time stepping 7684df9cb4SJed Brown 7784df9cb4SJed Brown Level: intermediate 7884df9cb4SJed Brown 79e91eccc2SStefano Zampini .seealso: TS, TSAdapt, TSGetAdapt() 8084df9cb4SJed Brown M*/ 818cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt adapt) 8284df9cb4SJed Brown { 8384df9cb4SJed Brown PetscFunctionBegin; 841c3436cfSJed Brown adapt->ops->choose = TSAdaptChoose_Basic; 851b9b13dfSLisandro Dalcin PetscFunctionReturn(0); 861b9b13dfSLisandro Dalcin } 87