1ae2316d0SEmil Constantinescu #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2*86e171c2SStefano Zampini #include <petscdm.h> 3ae2316d0SEmil Constantinescu 4ae2316d0SEmil Constantinescu typedef struct { 5*86e171c2SStefano Zampini Vec Y; 6726095e4SEmil Constantinescu } TSAdapt_GLEE; 7ae2316d0SEmil Constantinescu 85e4ed32fSEmil Constantinescu static PetscErrorCode TSAdaptChoose_GLEE(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept,PetscReal *wlte,PetscReal *wltea,PetscReal *wlter) 9ae2316d0SEmil Constantinescu { 10657c1e31SEmil Constantinescu TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 11ae2316d0SEmil Constantinescu PetscErrorCode ierr; 128a175baeSEmil Constantinescu Vec X,Y,E; 13df3bac28SEmil Constantinescu PetscReal enorm,enorma,enormr,hfac_lte,hfac_ltea,hfac_lter,h_lte,safety; 1480275a0aSLisandro Dalcin PetscInt order; 15*86e171c2SStefano Zampini PetscBool bGTEMethod; 16ae2316d0SEmil Constantinescu 17ae2316d0SEmil Constantinescu PetscFunctionBegin; 185e4ed32fSEmil Constantinescu *next_sc = 0; /* Reuse the same order scheme */ 191917a363SLisandro Dalcin safety = adapt->safety; 20*86e171c2SStefano Zampini ierr = PetscObjectTypeCompare((PetscObject)ts,TSGLEE,&bGTEMethod);CHKERRQ(ierr); 21df3bac28SEmil Constantinescu order = adapt->candidates.order[0]; 22df3bac28SEmil Constantinescu 23df3bac28SEmil Constantinescu if (bGTEMethod){/* the method is of GLEE type */ 24*86e171c2SStefano Zampini DM dm; 25*86e171c2SStefano Zampini 268a175baeSEmil Constantinescu ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 271c167fc2SEmil Constantinescu if (!glee->Y && adapt->glee_use_local) { 281c167fc2SEmil Constantinescu ierr = VecDuplicate(X,&glee->Y);CHKERRQ(ierr);/*create vector to store previous step global error*/ 291c167fc2SEmil Constantinescu ierr = VecZeroEntries(glee->Y);CHKERRQ(ierr); /*set error to zero on the first step - may not work if error is not zero initially*/ 301c167fc2SEmil Constantinescu } 31*86e171c2SStefano Zampini ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 32*86e171c2SStefano Zampini ierr = DMGetGlobalVector(dm,&E);CHKERRQ(ierr); 338a175baeSEmil Constantinescu ierr = TSGetTimeError(ts,0,&E);CHKERRQ(ierr); 341c167fc2SEmil Constantinescu 351c167fc2SEmil Constantinescu if (adapt->glee_use_local) {ierr = VecAXPY(E,-1.0,glee->Y);CHKERRQ(ierr);} /* local error = current error - previous step error */ 361c167fc2SEmil Constantinescu 371c167fc2SEmil Constantinescu /* this should be called with the solution at the beginning of the step too*/ 388a175baeSEmil Constantinescu ierr = TSErrorWeightedENorm(ts,E,X,X,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr); 39*86e171c2SStefano Zampini ierr = DMRestoreGlobalVector(dm,&E);CHKERRQ(ierr); 400a01e1b2SEmil Constantinescu } else { 41df3bac28SEmil Constantinescu /* the method is NOT of GLEE type; use the stantard basic augmented by separate atol and rtol */ 42ae2316d0SEmil Constantinescu ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 43657c1e31SEmil Constantinescu if (!glee->Y) {ierr = VecDuplicate(X,&glee->Y);CHKERRQ(ierr);} 44657c1e31SEmil Constantinescu Y = glee->Y; 45ae2316d0SEmil Constantinescu ierr = TSEvaluateStep(ts,order-1,Y,NULL);CHKERRQ(ierr); 467453f775SEmil Constantinescu ierr = TSErrorWeightedNorm(ts,X,Y,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr); 470a01e1b2SEmil Constantinescu } 485e4ed32fSEmil Constantinescu 495e4ed32fSEmil Constantinescu if (enorm < 0) { 505e4ed32fSEmil Constantinescu *accept = PETSC_TRUE; 515e4ed32fSEmil Constantinescu *next_h = h; /* Reuse the old step */ 52df3bac28SEmil Constantinescu *wlte = -1; /* Weighted error was not evaluated */ 53df3bac28SEmil Constantinescu *wltea = -1; /* Weighted absolute error was not evaluated */ 54df3bac28SEmil Constantinescu *wlter = -1; /* Weighted relative error was not evaluated */ 555e4ed32fSEmil Constantinescu PetscFunctionReturn(0); 565e4ed32fSEmil Constantinescu } 575e4ed32fSEmil Constantinescu 58df3bac28SEmil Constantinescu if (enorm > 1. || enorma > 1. || enormr > 1.) { 591917a363SLisandro Dalcin if (!*accept) safety *= adapt->reject_safety; /* The last attempt also failed, shorten more aggressively */ 60ae2316d0SEmil Constantinescu if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) { 61df3bac28SEmil Constantinescu ierr = PetscInfo4(adapt,"Estimated scaled truncation error [combined, absolute, relative]] [%g, %g, %g], accepting because step size %g is at minimum\n",(double)enorm,(double)enorma,(double)enormr,(double)h);CHKERRQ(ierr); 62ae2316d0SEmil Constantinescu *accept = PETSC_TRUE; 63bf997491SLisandro Dalcin } else if (adapt->always_accept) { 64df3bac28SEmil Constantinescu ierr = PetscInfo4(adapt,"Estimated scaled truncation error [combined, absolute, relative]] [%g, %g, %g], accepting step of size %g because always_accept is set\n",(double)enorm,(double)enorma,(double)enormr,(double)h);CHKERRQ(ierr); 65ae2316d0SEmil Constantinescu *accept = PETSC_TRUE; 66ae2316d0SEmil Constantinescu } else { 67df3bac28SEmil Constantinescu ierr = PetscInfo4(adapt,"Estimated scaled truncation error [combined, absolute, relative]] [%g, %g, %g], rejecting step of size %g\n",(double)enorm,(double)enorma,(double)enormr,(double)h);CHKERRQ(ierr); 68ae2316d0SEmil Constantinescu *accept = PETSC_FALSE; 69ae2316d0SEmil Constantinescu } 70ae2316d0SEmil Constantinescu } else { 71df3bac28SEmil Constantinescu ierr = PetscInfo4(adapt,"Estimated scaled truncation error [combined, absolute, relative] [%g, %g, %g], accepting step of size %g\n",(double)enorm,(double)enorma,(double)enormr,(double)h);CHKERRQ(ierr); 72ae2316d0SEmil Constantinescu *accept = PETSC_TRUE; 73ae2316d0SEmil Constantinescu } 74ae2316d0SEmil Constantinescu 75df3bac28SEmil Constantinescu if (bGTEMethod){ 761c167fc2SEmil Constantinescu if (*accept == PETSC_TRUE && adapt->glee_use_local) { 771c167fc2SEmil Constantinescu /* If step is accepted, then overwrite previous step error with the current error to be used on the next step */ 781c167fc2SEmil Constantinescu /* WARNING: if the adapters are composable, then the accept test will not be reliable*/ 79*86e171c2SStefano Zampini ierr = TSGetTimeError(ts,0,&glee->Y);CHKERRQ(ierr); 801c167fc2SEmil Constantinescu } 811c167fc2SEmil Constantinescu 82df3bac28SEmil Constantinescu /* The optimal new step based on the current global truncation error. */ 83df3bac28SEmil Constantinescu if (enorm > 0) { 84df3bac28SEmil Constantinescu /* factor based on the absolute tolerance */ 851c167fc2SEmil Constantinescu hfac_ltea = safety * PetscPowReal(1./enorma,((PetscReal)1)/(order+1)); 86df3bac28SEmil Constantinescu /* factor based on the relative tolerance */ 871c167fc2SEmil Constantinescu hfac_lter = safety * PetscPowReal(1./enormr,((PetscReal)1)/(order+1)); 88df3bac28SEmil Constantinescu /* pick the minimum time step among the relative and absolute tolerances */ 89df3bac28SEmil Constantinescu hfac_lte = PetscMin(hfac_ltea,hfac_lter); 90df3bac28SEmil Constantinescu } else { 91ae2316d0SEmil Constantinescu hfac_lte = safety * PETSC_INFINITY; 92df3bac28SEmil Constantinescu } 931917a363SLisandro Dalcin h_lte = h * PetscClipInterval(hfac_lte,adapt->clip[0],adapt->clip[1]); 94ae2316d0SEmil Constantinescu *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max); 95df3bac28SEmil Constantinescu } else { 96df3bac28SEmil Constantinescu /* The optimal new step based purely on local truncation error for this step. */ 97df3bac28SEmil Constantinescu if (enorm > 0) { 98df3bac28SEmil Constantinescu /* factor based on the absolute tolerance */ 99df3bac28SEmil Constantinescu hfac_ltea = safety * PetscPowReal(enorma,((PetscReal)-1)/order); 100df3bac28SEmil Constantinescu /* factor based on the relative tolerance */ 101df3bac28SEmil Constantinescu hfac_lter = safety * PetscPowReal(enormr,((PetscReal)-1)/order); 102df3bac28SEmil Constantinescu /* pick the minimum time step among the relative and absolute tolerances */ 103df3bac28SEmil Constantinescu hfac_lte = PetscMin(hfac_ltea,hfac_lter); 104df3bac28SEmil Constantinescu } else { 105df3bac28SEmil Constantinescu hfac_lte = safety * PETSC_INFINITY; 106df3bac28SEmil Constantinescu } 1071917a363SLisandro Dalcin h_lte = h * PetscClipInterval(hfac_lte,adapt->clip[0],adapt->clip[1]); 108df3bac28SEmil Constantinescu *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max); 109df3bac28SEmil Constantinescu } 110ae2316d0SEmil Constantinescu *wlte = enorm; 1115e4ed32fSEmil Constantinescu *wltea = enorma; 1125e4ed32fSEmil Constantinescu *wlter = enormr; 113ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 114ae2316d0SEmil Constantinescu } 115ae2316d0SEmil Constantinescu 116726095e4SEmil Constantinescu static PetscErrorCode TSAdaptReset_GLEE(TSAdapt adapt) 117ae2316d0SEmil Constantinescu { 118657c1e31SEmil Constantinescu TSAdapt_GLEE *glee = (TSAdapt_GLEE*)adapt->data; 119ae2316d0SEmil Constantinescu PetscErrorCode ierr; 120ae2316d0SEmil Constantinescu 121ae2316d0SEmil Constantinescu PetscFunctionBegin; 122657c1e31SEmil Constantinescu ierr = VecDestroy(&glee->Y);CHKERRQ(ierr); 123ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 124ae2316d0SEmil Constantinescu } 125ae2316d0SEmil Constantinescu 126726095e4SEmil Constantinescu static PetscErrorCode TSAdaptDestroy_GLEE(TSAdapt adapt) 127ae2316d0SEmil Constantinescu { 128ae2316d0SEmil Constantinescu PetscErrorCode ierr; 129ae2316d0SEmil Constantinescu 130ae2316d0SEmil Constantinescu PetscFunctionBegin; 131726095e4SEmil Constantinescu ierr = TSAdaptReset_GLEE(adapt);CHKERRQ(ierr); 132ae2316d0SEmil Constantinescu ierr = PetscFree(adapt->data);CHKERRQ(ierr); 133ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 134ae2316d0SEmil Constantinescu } 135ae2316d0SEmil Constantinescu 136ae2316d0SEmil Constantinescu /*MC 137657c1e31SEmil Constantinescu TSADAPTGLEE - GLEE adaptive controller for time stepping 138ae2316d0SEmil Constantinescu 139ae2316d0SEmil Constantinescu Level: intermediate 140ae2316d0SEmil Constantinescu 141e91eccc2SStefano Zampini .seealso: TS, TSAdapt, TSGetAdapt() 142ae2316d0SEmil Constantinescu M*/ 143726095e4SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt adapt) 144ae2316d0SEmil Constantinescu { 145ae2316d0SEmil Constantinescu PetscErrorCode ierr; 1461917a363SLisandro Dalcin TSAdapt_GLEE *glee; 147ae2316d0SEmil Constantinescu 148ae2316d0SEmil Constantinescu PetscFunctionBegin; 1491917a363SLisandro Dalcin ierr = PetscNewLog(adapt,&glee);CHKERRQ(ierr); 1501917a363SLisandro Dalcin adapt->data = (void*)glee; 151726095e4SEmil Constantinescu adapt->ops->choose = TSAdaptChoose_GLEE; 1521917a363SLisandro Dalcin adapt->ops->reset = TSAdaptReset_GLEE; 153726095e4SEmil Constantinescu adapt->ops->destroy = TSAdaptDestroy_GLEE; 154ae2316d0SEmil Constantinescu PetscFunctionReturn(0); 155ae2316d0SEmil Constantinescu } 156