xref: /petsc/src/ts/adapt/impls/glee/adaptglee.c (revision e91eccc2778f456fc991f5a9b142a3a67738bfd5)
1ae2316d0SEmil Constantinescu #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/
2ae2316d0SEmil Constantinescu 
3ae2316d0SEmil Constantinescu typedef struct {
41917a363SLisandro Dalcin   Vec E,Y;
5726095e4SEmil Constantinescu } TSAdapt_GLEE;
6ae2316d0SEmil Constantinescu 
75e4ed32fSEmil 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)
8ae2316d0SEmil Constantinescu {
9657c1e31SEmil Constantinescu   TSAdapt_GLEE  *glee = (TSAdapt_GLEE*)adapt->data;
100a01e1b2SEmil Constantinescu   TSType         time_scheme;      /* Type of time-integration scheme        */
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;
15df3bac28SEmil Constantinescu   PetscBool      bGTEMethod=PETSC_FALSE;
16ae2316d0SEmil Constantinescu 
17ae2316d0SEmil Constantinescu   PetscFunctionBegin;
18df3bac28SEmil Constantinescu 
195e4ed32fSEmil Constantinescu   *next_sc = 0; /* Reuse the same order scheme */
201917a363SLisandro Dalcin   safety = adapt->safety;
210a01e1b2SEmil Constantinescu   ierr = TSGetType(ts,&time_scheme);CHKERRQ(ierr);
22df3bac28SEmil Constantinescu   if (!strcmp(time_scheme,TSGLEE)) bGTEMethod=PETSC_TRUE;
23df3bac28SEmil Constantinescu   order = adapt->candidates.order[0];
24df3bac28SEmil Constantinescu 
25df3bac28SEmil Constantinescu   if (bGTEMethod){/* the method is of GLEE type */
268a175baeSEmil Constantinescu     ierr = TSGetSolution(ts,&X);CHKERRQ(ierr);
278a175baeSEmil Constantinescu     if (!glee->E) {ierr = VecDuplicate(X,&glee->E);CHKERRQ(ierr);}
288a175baeSEmil Constantinescu     E    = glee->E;
298a175baeSEmil Constantinescu     ierr = TSGetTimeError(ts,0,&E);CHKERRQ(ierr);
30df3bac28SEmil Constantinescu     /* this should be called with Y (the solution at the beginning of the step)*/
318a175baeSEmil Constantinescu     ierr = TSErrorWeightedENorm(ts,E,X,X,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr);
320a01e1b2SEmil Constantinescu   } else {
33df3bac28SEmil Constantinescu     /* the method is NOT of GLEE type; use the stantard basic augmented by separate atol and rtol */
34ae2316d0SEmil Constantinescu     ierr = TSGetSolution(ts,&X);CHKERRQ(ierr);
35657c1e31SEmil Constantinescu     if (!glee->Y) {ierr = VecDuplicate(X,&glee->Y);CHKERRQ(ierr);}
36657c1e31SEmil Constantinescu     Y     = glee->Y;
37ae2316d0SEmil Constantinescu     ierr  = TSEvaluateStep(ts,order-1,Y,NULL);CHKERRQ(ierr);
387453f775SEmil Constantinescu     ierr  = TSErrorWeightedNorm(ts,X,Y,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr);
390a01e1b2SEmil Constantinescu   }
405e4ed32fSEmil Constantinescu 
415e4ed32fSEmil Constantinescu   if (enorm < 0) {
425e4ed32fSEmil Constantinescu     *accept  = PETSC_TRUE;
435e4ed32fSEmil Constantinescu     *next_h  = h;            /* Reuse the old step */
44df3bac28SEmil Constantinescu     *wlte    = -1;           /* Weighted error was not evaluated */
45df3bac28SEmil Constantinescu     *wltea   = -1;           /* Weighted absolute error was not evaluated */
46df3bac28SEmil Constantinescu     *wlter   = -1;           /* Weighted relative error was not evaluated */
475e4ed32fSEmil Constantinescu     PetscFunctionReturn(0);
485e4ed32fSEmil Constantinescu   }
495e4ed32fSEmil Constantinescu 
50df3bac28SEmil Constantinescu   if (enorm > 1. || enorma > 1. || enormr > 1.) {
511917a363SLisandro Dalcin     if (!*accept) safety *= adapt->reject_safety; /* The last attempt also failed, shorten more aggressively */
52ae2316d0SEmil Constantinescu     if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) {
53df3bac28SEmil 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);
54ae2316d0SEmil Constantinescu       *accept = PETSC_TRUE;
55bf997491SLisandro Dalcin     } else if (adapt->always_accept) {
56df3bac28SEmil 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);
57ae2316d0SEmil Constantinescu       *accept = PETSC_TRUE;
58ae2316d0SEmil Constantinescu     } else {
59df3bac28SEmil 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);
60ae2316d0SEmil Constantinescu       *accept = PETSC_FALSE;
61ae2316d0SEmil Constantinescu     }
62ae2316d0SEmil Constantinescu   } else {
63df3bac28SEmil 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);
64ae2316d0SEmil Constantinescu     *accept = PETSC_TRUE;
65ae2316d0SEmil Constantinescu   }
66ae2316d0SEmil Constantinescu 
67df3bac28SEmil Constantinescu   if (bGTEMethod){
68df3bac28SEmil Constantinescu     /* The optimal new step based on the current global truncation error. */
69df3bac28SEmil Constantinescu     if (enorm > 0) {
70df3bac28SEmil Constantinescu       /* factor based on the absolute tolerance */
71df3bac28SEmil Constantinescu       hfac_ltea = safety * PetscPowReal(1./enorma,((PetscReal)1)/order);
72df3bac28SEmil Constantinescu       /* factor based on the relative tolerance */
73df3bac28SEmil Constantinescu       hfac_lter = safety * PetscPowReal(1./enormr,((PetscReal)1)/order);
74df3bac28SEmil Constantinescu       /* pick the minimum time step among the relative and absolute tolerances */
75df3bac28SEmil Constantinescu       hfac_lte  = PetscMin(hfac_ltea,hfac_lter);
76df3bac28SEmil Constantinescu     } else {
77ae2316d0SEmil Constantinescu       hfac_lte = safety * PETSC_INFINITY;
78df3bac28SEmil Constantinescu     }
791917a363SLisandro Dalcin     h_lte = h * PetscClipInterval(hfac_lte,adapt->clip[0],adapt->clip[1]);
80ae2316d0SEmil Constantinescu     *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max);
81df3bac28SEmil Constantinescu   } else {
82df3bac28SEmil Constantinescu     /* The optimal new step based purely on local truncation error for this step. */
83df3bac28SEmil Constantinescu     if (enorm > 0) {
84df3bac28SEmil Constantinescu       /* factor based on the absolute tolerance */
85df3bac28SEmil Constantinescu       hfac_ltea = safety * PetscPowReal(enorma,((PetscReal)-1)/order);
86df3bac28SEmil Constantinescu       /* factor based on the relative tolerance */
87df3bac28SEmil Constantinescu       hfac_lter = safety * PetscPowReal(enormr,((PetscReal)-1)/order);
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 {
91df3bac28SEmil Constantinescu       hfac_lte = safety * PETSC_INFINITY;
92df3bac28SEmil Constantinescu     }
931917a363SLisandro Dalcin     h_lte = h * PetscClipInterval(hfac_lte,adapt->clip[0],adapt->clip[1]);
94df3bac28SEmil Constantinescu     *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max);
95df3bac28SEmil Constantinescu   }
96ae2316d0SEmil Constantinescu   *wlte   = enorm;
975e4ed32fSEmil Constantinescu   *wltea  = enorma;
985e4ed32fSEmil Constantinescu   *wlter  = enormr;
99ae2316d0SEmil Constantinescu   PetscFunctionReturn(0);
100ae2316d0SEmil Constantinescu }
101ae2316d0SEmil Constantinescu 
102726095e4SEmil Constantinescu static PetscErrorCode TSAdaptReset_GLEE(TSAdapt adapt)
103ae2316d0SEmil Constantinescu {
104657c1e31SEmil Constantinescu   TSAdapt_GLEE  *glee = (TSAdapt_GLEE*)adapt->data;
105ae2316d0SEmil Constantinescu   PetscErrorCode ierr;
106ae2316d0SEmil Constantinescu 
107ae2316d0SEmil Constantinescu   PetscFunctionBegin;
108657c1e31SEmil Constantinescu   ierr = VecDestroy(&glee->Y);CHKERRQ(ierr);
1091917a363SLisandro Dalcin   ierr = VecDestroy(&glee->E);CHKERRQ(ierr);
110ae2316d0SEmil Constantinescu   PetscFunctionReturn(0);
111ae2316d0SEmil Constantinescu }
112ae2316d0SEmil Constantinescu 
113726095e4SEmil Constantinescu static PetscErrorCode TSAdaptDestroy_GLEE(TSAdapt adapt)
114ae2316d0SEmil Constantinescu {
115ae2316d0SEmil Constantinescu   PetscErrorCode ierr;
116ae2316d0SEmil Constantinescu 
117ae2316d0SEmil Constantinescu   PetscFunctionBegin;
118726095e4SEmil Constantinescu   ierr = TSAdaptReset_GLEE(adapt);CHKERRQ(ierr);
119ae2316d0SEmil Constantinescu   ierr = PetscFree(adapt->data);CHKERRQ(ierr);
120ae2316d0SEmil Constantinescu   PetscFunctionReturn(0);
121ae2316d0SEmil Constantinescu }
122ae2316d0SEmil Constantinescu 
123ae2316d0SEmil Constantinescu /*MC
124657c1e31SEmil Constantinescu    TSADAPTGLEE - GLEE adaptive controller for time stepping
125ae2316d0SEmil Constantinescu 
126ae2316d0SEmil Constantinescu    Level: intermediate
127ae2316d0SEmil Constantinescu 
128*e91eccc2SStefano Zampini .seealso: TS, TSAdapt, TSGetAdapt()
129ae2316d0SEmil Constantinescu M*/
130726095e4SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt adapt)
131ae2316d0SEmil Constantinescu {
132ae2316d0SEmil Constantinescu   PetscErrorCode ierr;
1331917a363SLisandro Dalcin   TSAdapt_GLEE  *glee;
134ae2316d0SEmil Constantinescu 
135ae2316d0SEmil Constantinescu   PetscFunctionBegin;
1361917a363SLisandro Dalcin   ierr = PetscNewLog(adapt,&glee);CHKERRQ(ierr);
1371917a363SLisandro Dalcin   adapt->data         = (void*)glee;
138726095e4SEmil Constantinescu   adapt->ops->choose  = TSAdaptChoose_GLEE;
1391917a363SLisandro Dalcin   adapt->ops->reset   = TSAdaptReset_GLEE;
140726095e4SEmil Constantinescu   adapt->ops->destroy = TSAdaptDestroy_GLEE;
141ae2316d0SEmil Constantinescu   PetscFunctionReturn(0);
142ae2316d0SEmil Constantinescu }
143