xref: /petsc/src/ts/adapt/impls/glee/adaptglee.c (revision 1c167fc2a0c466620c6784b5ffb486042a891a3a)
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);
27*1c167fc2SEmil Constantinescu     if (!glee->Y && adapt->glee_use_local) {
28*1c167fc2SEmil Constantinescu       ierr = VecDuplicate(X,&glee->Y);CHKERRQ(ierr);/*create vector to store previous step global error*/
29*1c167fc2SEmil Constantinescu       ierr = VecZeroEntries(glee->Y);CHKERRQ(ierr); /*set error to zero on the first step - may not work if error is not zero initially*/
30*1c167fc2SEmil Constantinescu     }
318a175baeSEmil Constantinescu     if (!glee->E) {ierr = VecDuplicate(X,&glee->E);CHKERRQ(ierr);}
328a175baeSEmil Constantinescu     E    = glee->E;
338a175baeSEmil Constantinescu     ierr = TSGetTimeError(ts,0,&E);CHKERRQ(ierr);
34*1c167fc2SEmil Constantinescu 
35*1c167fc2SEmil Constantinescu     if (adapt->glee_use_local) {ierr = VecAXPY(E,-1.0,glee->Y);CHKERRQ(ierr);} /* local error = current error - previous step error */
36*1c167fc2SEmil Constantinescu 
37*1c167fc2SEmil 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);
390a01e1b2SEmil Constantinescu   } else {
40df3bac28SEmil Constantinescu     /* the method is NOT of GLEE type; use the stantard basic augmented by separate atol and rtol */
41ae2316d0SEmil Constantinescu     ierr = TSGetSolution(ts,&X);CHKERRQ(ierr);
42657c1e31SEmil Constantinescu     if (!glee->Y) {ierr = VecDuplicate(X,&glee->Y);CHKERRQ(ierr);}
43657c1e31SEmil Constantinescu     Y     = glee->Y;
44ae2316d0SEmil Constantinescu     ierr  = TSEvaluateStep(ts,order-1,Y,NULL);CHKERRQ(ierr);
457453f775SEmil Constantinescu     ierr  = TSErrorWeightedNorm(ts,X,Y,adapt->wnormtype,&enorm,&enorma,&enormr);CHKERRQ(ierr);
460a01e1b2SEmil Constantinescu   }
475e4ed32fSEmil Constantinescu 
485e4ed32fSEmil Constantinescu   if (enorm < 0) {
495e4ed32fSEmil Constantinescu     *accept  = PETSC_TRUE;
505e4ed32fSEmil Constantinescu     *next_h  = h;            /* Reuse the old step */
51df3bac28SEmil Constantinescu     *wlte    = -1;           /* Weighted error was not evaluated */
52df3bac28SEmil Constantinescu     *wltea   = -1;           /* Weighted absolute error was not evaluated */
53df3bac28SEmil Constantinescu     *wlter   = -1;           /* Weighted relative error was not evaluated */
545e4ed32fSEmil Constantinescu     PetscFunctionReturn(0);
555e4ed32fSEmil Constantinescu   }
565e4ed32fSEmil Constantinescu 
57df3bac28SEmil Constantinescu   if (enorm > 1. || enorma > 1. || enormr > 1.) {
581917a363SLisandro Dalcin     if (!*accept) safety *= adapt->reject_safety; /* The last attempt also failed, shorten more aggressively */
59ae2316d0SEmil Constantinescu     if (h < (1 + PETSC_SQRT_MACHINE_EPSILON)*adapt->dt_min) {
60df3bac28SEmil 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);
61ae2316d0SEmil Constantinescu       *accept = PETSC_TRUE;
62bf997491SLisandro Dalcin     } else if (adapt->always_accept) {
63df3bac28SEmil 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);
64ae2316d0SEmil Constantinescu       *accept = PETSC_TRUE;
65ae2316d0SEmil Constantinescu     } else {
66df3bac28SEmil 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);
67ae2316d0SEmil Constantinescu       *accept = PETSC_FALSE;
68ae2316d0SEmil Constantinescu     }
69ae2316d0SEmil Constantinescu   } else {
70df3bac28SEmil 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);
71ae2316d0SEmil Constantinescu     *accept = PETSC_TRUE;
72ae2316d0SEmil Constantinescu   }
73ae2316d0SEmil Constantinescu 
74df3bac28SEmil Constantinescu   if (bGTEMethod){
75*1c167fc2SEmil Constantinescu     if (*accept == PETSC_TRUE && adapt->glee_use_local) {
76*1c167fc2SEmil Constantinescu       /* If step is accepted, then overwrite previous step error with the current error to be used on the next step */
77*1c167fc2SEmil Constantinescu       /* WARNING: if the adapters are composable, then the accept test will not be reliable*/
78*1c167fc2SEmil Constantinescu       ierr = TSGetTimeError(ts,0,&(glee->Y));CHKERRQ(ierr);
79*1c167fc2SEmil Constantinescu     }
80*1c167fc2SEmil Constantinescu 
81df3bac28SEmil Constantinescu     /* The optimal new step based on the current global truncation error. */
82df3bac28SEmil Constantinescu     if (enorm > 0) {
83df3bac28SEmil Constantinescu       /* factor based on the absolute tolerance */
84*1c167fc2SEmil Constantinescu       hfac_ltea = safety * PetscPowReal(1./enorma,((PetscReal)1)/(order+1));
85df3bac28SEmil Constantinescu       /* factor based on the relative tolerance */
86*1c167fc2SEmil Constantinescu       hfac_lter = safety * PetscPowReal(1./enormr,((PetscReal)1)/(order+1));
87df3bac28SEmil Constantinescu       /* pick the minimum time step among the relative and absolute tolerances */
88df3bac28SEmil Constantinescu       hfac_lte  = PetscMin(hfac_ltea,hfac_lter);
89df3bac28SEmil Constantinescu     } else {
90ae2316d0SEmil Constantinescu       hfac_lte = safety * PETSC_INFINITY;
91df3bac28SEmil Constantinescu     }
921917a363SLisandro Dalcin     h_lte = h * PetscClipInterval(hfac_lte,adapt->clip[0],adapt->clip[1]);
93ae2316d0SEmil Constantinescu     *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max);
94df3bac28SEmil Constantinescu   } else {
95df3bac28SEmil Constantinescu     /* The optimal new step based purely on local truncation error for this step. */
96df3bac28SEmil Constantinescu     if (enorm > 0) {
97df3bac28SEmil Constantinescu       /* factor based on the absolute tolerance */
98df3bac28SEmil Constantinescu       hfac_ltea = safety * PetscPowReal(enorma,((PetscReal)-1)/order);
99df3bac28SEmil Constantinescu       /* factor based on the relative tolerance */
100df3bac28SEmil Constantinescu       hfac_lter = safety * PetscPowReal(enormr,((PetscReal)-1)/order);
101df3bac28SEmil Constantinescu       /* pick the minimum time step among the relative and absolute tolerances */
102df3bac28SEmil Constantinescu       hfac_lte  = PetscMin(hfac_ltea,hfac_lter);
103df3bac28SEmil Constantinescu     } else {
104df3bac28SEmil Constantinescu       hfac_lte = safety * PETSC_INFINITY;
105df3bac28SEmil Constantinescu     }
1061917a363SLisandro Dalcin     h_lte = h * PetscClipInterval(hfac_lte,adapt->clip[0],adapt->clip[1]);
107df3bac28SEmil Constantinescu     *next_h = PetscClipInterval(h_lte,adapt->dt_min,adapt->dt_max);
108df3bac28SEmil Constantinescu   }
109ae2316d0SEmil Constantinescu   *wlte   = enorm;
1105e4ed32fSEmil Constantinescu   *wltea  = enorma;
1115e4ed32fSEmil Constantinescu   *wlter  = enormr;
112ae2316d0SEmil Constantinescu   PetscFunctionReturn(0);
113ae2316d0SEmil Constantinescu }
114ae2316d0SEmil Constantinescu 
115726095e4SEmil Constantinescu static PetscErrorCode TSAdaptReset_GLEE(TSAdapt adapt)
116ae2316d0SEmil Constantinescu {
117657c1e31SEmil Constantinescu   TSAdapt_GLEE  *glee = (TSAdapt_GLEE*)adapt->data;
118ae2316d0SEmil Constantinescu   PetscErrorCode ierr;
119ae2316d0SEmil Constantinescu 
120ae2316d0SEmil Constantinescu   PetscFunctionBegin;
121657c1e31SEmil Constantinescu   ierr = VecDestroy(&glee->Y);CHKERRQ(ierr);
1221917a363SLisandro Dalcin   ierr = VecDestroy(&glee->E);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