xref: /petsc/src/ts/impls/implicit/glle/glleadapt.c (revision 5b6bfdb9644f185dbf5e5a09b808ec241507e1e7)
1 
2 #include <../src/ts/impls/implicit/glle/glle.h> /*I  "petscts.h" I*/
3 
4 static PetscFunctionList TSGLLEAdaptList;
5 static PetscBool         TSGLLEAdaptPackageInitialized;
6 static PetscBool         TSGLLEAdaptRegisterAllCalled;
7 static PetscClassId      TSGLLEADAPT_CLASSID;
8 
9 struct _TSGLLEAdaptOps {
10   PetscErrorCode (*choose)(TSGLLEAdapt,PetscInt,const PetscInt[],const PetscReal[],const PetscReal[],PetscInt,PetscReal,PetscReal,PetscInt*,PetscReal*,PetscBool*);
11   PetscErrorCode (*destroy)(TSGLLEAdapt);
12   PetscErrorCode (*view)(TSGLLEAdapt,PetscViewer);
13   PetscErrorCode (*setfromoptions)(PetscOptionItems*,TSGLLEAdapt);
14 };
15 
16 struct _p_TSGLLEAdapt {
17   PETSCHEADER(struct _TSGLLEAdaptOps);
18   void *data;
19 };
20 
21 PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate_None(TSGLLEAdapt);
22 PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate_Size(TSGLLEAdapt);
23 PETSC_EXTERN PetscErrorCode TSGLLEAdaptCreate_Both(TSGLLEAdapt);
24 
25 /*@C
26    TSGLLEAdaptRegister -  adds a TSGLLEAdapt implementation
27 
28    Not Collective
29 
30    Input Parameters:
31 +  name_scheme - name of user-defined adaptivity scheme
32 -  routine_create - routine to create method context
33 
34    Notes:
35    TSGLLEAdaptRegister() may be called multiple times to add several user-defined families.
36 
37    Sample usage:
38 .vb
39    TSGLLEAdaptRegister("my_scheme",MySchemeCreate);
40 .ve
41 
42    Then, your scheme can be chosen with the procedural interface via
43 $     TSGLLEAdaptSetType(ts,"my_scheme")
44    or at runtime via the option
45 $     -ts_adapt_type my_scheme
46 
47    Level: advanced
48 
49 .keywords: TSGLLEAdapt, register
50 
51 .seealso: TSGLLEAdaptRegisterAll()
52 @*/
53 PetscErrorCode  TSGLLEAdaptRegister(const char sname[],PetscErrorCode (*function)(TSGLLEAdapt))
54 {
55   PetscErrorCode ierr;
56 
57   PetscFunctionBegin;
58   ierr = PetscFunctionListAdd(&TSGLLEAdaptList,sname,function);CHKERRQ(ierr);
59   PetscFunctionReturn(0);
60 }
61 
62 /*@C
63   TSGLLEAdaptRegisterAll - Registers all of the adaptivity schemes in TSGLLEAdapt
64 
65   Not Collective
66 
67   Level: advanced
68 
69 .keywords: TSGLLEAdapt, register, all
70 
71 .seealso: TSGLLEAdaptRegisterDestroy()
72 @*/
73 PetscErrorCode  TSGLLEAdaptRegisterAll(void)
74 {
75   PetscErrorCode ierr;
76 
77   PetscFunctionBegin;
78   if (TSGLLEAdaptRegisterAllCalled) PetscFunctionReturn(0);
79   TSGLLEAdaptRegisterAllCalled = PETSC_TRUE;
80   ierr = TSGLLEAdaptRegister(TSGLLEADAPT_NONE,TSGLLEAdaptCreate_None);CHKERRQ(ierr);
81   ierr = TSGLLEAdaptRegister(TSGLLEADAPT_SIZE,TSGLLEAdaptCreate_Size);CHKERRQ(ierr);
82   ierr = TSGLLEAdaptRegister(TSGLLEADAPT_BOTH,TSGLLEAdaptCreate_Both);CHKERRQ(ierr);
83   PetscFunctionReturn(0);
84 }
85 
86 /*@C
87   TSGLLEFinalizePackage - This function destroys everything in the TSGLLE package. It is
88   called from PetscFinalize().
89 
90   Level: developer
91 
92 .keywords: Petsc, destroy, package
93 .seealso: PetscFinalize()
94 @*/
95 PetscErrorCode  TSGLLEAdaptFinalizePackage(void)
96 {
97   PetscErrorCode ierr;
98 
99   PetscFunctionBegin;
100   ierr = PetscFunctionListDestroy(&TSGLLEAdaptList);CHKERRQ(ierr);
101   TSGLLEAdaptPackageInitialized = PETSC_FALSE;
102   TSGLLEAdaptRegisterAllCalled  = PETSC_FALSE;
103   PetscFunctionReturn(0);
104 }
105 
106 /*@C
107   TSGLLEAdaptInitializePackage - This function initializes everything in the TSGLLEAdapt package. It is
108   called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to
109   TSCreate_GLLE() when using static libraries.
110 
111   Level: developer
112 
113 .keywords: TSGLLEAdapt, initialize, package
114 .seealso: PetscInitialize()
115 @*/
116 PetscErrorCode  TSGLLEAdaptInitializePackage(void)
117 {
118   PetscErrorCode ierr;
119 
120   PetscFunctionBegin;
121   if (TSGLLEAdaptPackageInitialized) PetscFunctionReturn(0);
122   TSGLLEAdaptPackageInitialized = PETSC_TRUE;
123   ierr = PetscClassIdRegister("TSGLLEAdapt",&TSGLLEADAPT_CLASSID);CHKERRQ(ierr);
124   ierr = TSGLLEAdaptRegisterAll();CHKERRQ(ierr);
125   ierr = PetscRegisterFinalize(TSGLLEAdaptFinalizePackage);CHKERRQ(ierr);
126   PetscFunctionReturn(0);
127 }
128 
129 PetscErrorCode  TSGLLEAdaptSetType(TSGLLEAdapt adapt,TSGLLEAdaptType type)
130 {
131   PetscErrorCode ierr,(*r)(TSGLLEAdapt);
132 
133   PetscFunctionBegin;
134   ierr = PetscFunctionListFind(TSGLLEAdaptList,type,&r);CHKERRQ(ierr);
135   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSGLLEAdapt type \"%s\" given",type);
136   if (((PetscObject)adapt)->type_name) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);}
137   ierr = (*r)(adapt);CHKERRQ(ierr);
138   ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr);
139   PetscFunctionReturn(0);
140 }
141 
142 PetscErrorCode  TSGLLEAdaptSetOptionsPrefix(TSGLLEAdapt adapt,const char prefix[])
143 {
144   PetscErrorCode ierr;
145 
146   PetscFunctionBegin;
147   ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr);
148   PetscFunctionReturn(0);
149 }
150 
151 PetscErrorCode  TSGLLEAdaptView(TSGLLEAdapt adapt,PetscViewer viewer)
152 {
153   PetscErrorCode ierr;
154   PetscBool      iascii;
155 
156   PetscFunctionBegin;
157   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
158   if (iascii) {
159     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr);
160     if (adapt->ops->view) {
161       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
162       ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
163       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
164     }
165   }
166   PetscFunctionReturn(0);
167 }
168 
169 PetscErrorCode  TSGLLEAdaptDestroy(TSGLLEAdapt *adapt)
170 {
171   PetscErrorCode ierr;
172 
173   PetscFunctionBegin;
174   if (!*adapt) PetscFunctionReturn(0);
175   PetscValidHeaderSpecific(*adapt,TSGLLEADAPT_CLASSID,1);
176   if (--((PetscObject)(*adapt))->refct > 0) {*adapt = 0; PetscFunctionReturn(0);}
177   if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);}
178   ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr);
179   PetscFunctionReturn(0);
180 }
181 
182 PetscErrorCode  TSGLLEAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSGLLEAdapt adapt)
183 {
184   PetscErrorCode ierr;
185   char           type[256] = TSGLLEADAPT_BOTH;
186   PetscBool      flg;
187 
188   PetscFunctionBegin;
189   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TSGLLE, but currently this
190   * function can only be called from inside TSSetFromOptions_GLLE()  */
191   ierr = PetscOptionsHead(PetscOptionsObject,"TSGLLE Adaptivity options");CHKERRQ(ierr);
192   ierr = PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSGLLEAdaptSetType",TSGLLEAdaptList,
193                           ((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);CHKERRQ(ierr);
194   if (flg || !((PetscObject)adapt)->type_name) {
195     ierr = TSGLLEAdaptSetType(adapt,type);CHKERRQ(ierr);
196   }
197   if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);}
198   ierr = PetscOptionsTail();CHKERRQ(ierr);
199   PetscFunctionReturn(0);
200 }
201 
202 PetscErrorCode  TSGLLEAdaptChoose(TSGLLEAdapt adapt,PetscInt n,const PetscInt orders[],const PetscReal errors[],const PetscReal cost[],PetscInt cur,PetscReal h,PetscReal tleft,PetscInt *next_sc,PetscReal *next_h,PetscBool  *finish)
203 {
204   PetscErrorCode ierr;
205 
206   PetscFunctionBegin;
207   PetscValidHeaderSpecific(adapt,TSGLLEADAPT_CLASSID,1);
208   PetscValidIntPointer(orders,3);
209   PetscValidPointer(errors,4);
210   PetscValidPointer(cost,5);
211   PetscValidIntPointer(next_sc,9);
212   PetscValidPointer(next_h,10);
213   PetscValidIntPointer(finish,11);
214   ierr = (*adapt->ops->choose)(adapt,n,orders,errors,cost,cur,h,tleft,next_sc,next_h,finish);CHKERRQ(ierr);
215   PetscFunctionReturn(0);
216 }
217 
218 PetscErrorCode  TSGLLEAdaptCreate(MPI_Comm comm,TSGLLEAdapt *inadapt)
219 {
220   PetscErrorCode ierr;
221   TSGLLEAdapt      adapt;
222 
223   PetscFunctionBegin;
224   *inadapt = NULL;
225   ierr     = PetscHeaderCreate(adapt,TSGLLEADAPT_CLASSID,"TSGLLEAdapt","General Linear adaptivity","TS",comm,TSGLLEAdaptDestroy,TSGLLEAdaptView);CHKERRQ(ierr);
226   *inadapt = adapt;
227   PetscFunctionReturn(0);
228 }
229 
230 
231 /*
232 *  Implementations
233 */
234 
235 static PetscErrorCode TSGLLEAdaptDestroy_JustFree(TSGLLEAdapt adapt)
236 {
237   PetscErrorCode ierr;
238 
239   PetscFunctionBegin;
240   ierr = PetscFree(adapt->data);CHKERRQ(ierr);
241   PetscFunctionReturn(0);
242 }
243 
244 /* -------------------------------- None ----------------------------------- */
245 typedef struct {
246   PetscInt  scheme;
247   PetscReal h;
248 } TSGLLEAdapt_None;
249 
250 static PetscErrorCode TSGLLEAdaptChoose_None(TSGLLEAdapt adapt,PetscInt n,const PetscInt orders[],const PetscReal errors[],const PetscReal cost[],PetscInt cur,PetscReal h,PetscReal tleft,PetscInt *next_sc,PetscReal *next_h,PetscBool  *finish)
251 {
252 
253   PetscFunctionBegin;
254   *next_sc = cur;
255   *next_h  = h;
256   if (*next_h > tleft) {
257     *finish = PETSC_TRUE;
258     *next_h = tleft;
259   } else *finish = PETSC_FALSE;
260   PetscFunctionReturn(0);
261 }
262 
263 PetscErrorCode  TSGLLEAdaptCreate_None(TSGLLEAdapt adapt)
264 {
265   PetscErrorCode ierr;
266   TSGLLEAdapt_None *a;
267 
268   PetscFunctionBegin;
269   ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr);
270   adapt->data         = (void*)a;
271   adapt->ops->choose  = TSGLLEAdaptChoose_None;
272   adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree;
273   PetscFunctionReturn(0);
274 }
275 
276 /* -------------------------------- Size ----------------------------------- */
277 typedef struct {
278   PetscReal desired_h;
279 } TSGLLEAdapt_Size;
280 
281 
282 static PetscErrorCode TSGLLEAdaptChoose_Size(TSGLLEAdapt adapt,PetscInt n,const PetscInt orders[],const PetscReal errors[],const PetscReal cost[],PetscInt cur,PetscReal h,PetscReal tleft,PetscInt *next_sc,PetscReal *next_h,PetscBool  *finish)
283 {
284   TSGLLEAdapt_Size *sz = (TSGLLEAdapt_Size*)adapt->data;
285   PetscReal      dec = 0.2,inc = 5.0,safe = 0.9,optimal,last_desired_h;
286 
287   PetscFunctionBegin;
288   *next_sc = cur;
289   optimal  = PetscPowReal((PetscReal)errors[cur],(PetscReal)-1./(safe*orders[cur]));
290   /* Step sizes oscillate when there is no smoothing.  Here we use a geometric mean of the current step size and the
291   * one that would have been taken (without smoothing) on the last step. */
292   last_desired_h = sz->desired_h;
293   sz->desired_h  = h*PetscMax(dec,PetscMin(inc,optimal)); /* Trim to [dec,inc] */
294 
295   /* Normally only happens on the first step */
296   if (last_desired_h > 1e-14) *next_h = PetscSqrtReal(last_desired_h * sz->desired_h);
297   else *next_h = sz->desired_h;
298 
299   if (*next_h > tleft) {
300     *finish = PETSC_TRUE;
301     *next_h = tleft;
302   } else *finish = PETSC_FALSE;
303   PetscFunctionReturn(0);
304 }
305 
306 PetscErrorCode  TSGLLEAdaptCreate_Size(TSGLLEAdapt adapt)
307 {
308   PetscErrorCode ierr;
309   TSGLLEAdapt_Size *a;
310 
311   PetscFunctionBegin;
312   ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr);
313   adapt->data         = (void*)a;
314   adapt->ops->choose  = TSGLLEAdaptChoose_Size;
315   adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree;
316   PetscFunctionReturn(0);
317 }
318 
319 /* -------------------------------- Both ----------------------------------- */
320 typedef struct {
321   PetscInt  count_at_order;
322   PetscReal desired_h;
323 } TSGLLEAdapt_Both;
324 
325 
326 static PetscErrorCode TSGLLEAdaptChoose_Both(TSGLLEAdapt adapt,PetscInt n,const PetscInt orders[],const PetscReal errors[],const PetscReal cost[],PetscInt cur,PetscReal h,PetscReal tleft,PetscInt *next_sc,PetscReal *next_h,PetscBool  *finish)
327 {
328   TSGLLEAdapt_Both *both = (TSGLLEAdapt_Both*)adapt->data;
329   PetscErrorCode ierr;
330   PetscReal      dec = 0.2,inc = 5.0,safe = 0.9;
331   struct {PetscInt id; PetscReal h,eff;} best={-1,0,0},trial={-1,0,0},current={-1,0,0};
332   PetscInt i;
333 
334   PetscFunctionBegin;
335   for (i=0; i<n; i++) {
336     PetscReal optimal;
337     trial.id  = i;
338     optimal   = PetscPowReal((PetscReal)errors[i],(PetscReal)-1./(safe*orders[i]));
339     trial.h   = h*optimal;
340     trial.eff = trial.h/cost[i];
341     if (trial.eff > best.eff) {ierr = PetscMemcpy(&best,&trial,sizeof(trial));CHKERRQ(ierr);}
342     if (i == cur) {ierr = PetscMemcpy(&current,&trial,sizeof(trial));CHKERRQ(ierr);}
343   }
344   /* Only switch orders if the scheme offers significant benefits over the current one.
345   When the scheme is not changing, only change step size if it offers significant benefits. */
346   if (best.eff < 1.2*current.eff || both->count_at_order < orders[cur]+2) {
347     PetscReal last_desired_h;
348     *next_sc        = current.id;
349     last_desired_h  = both->desired_h;
350     both->desired_h = PetscMax(h*dec,PetscMin(h*inc,current.h));
351     *next_h         = (both->count_at_order > 0)
352                       ? PetscSqrtReal(last_desired_h * both->desired_h)
353                       : both->desired_h;
354     both->count_at_order++;
355   } else {
356     PetscReal rat = cost[best.id]/cost[cur];
357     *next_sc = best.id;
358     *next_h  = PetscMax(h*rat*dec,PetscMin(h*rat*inc,best.h));
359     both->count_at_order = 0;
360     both->desired_h      = best.h;
361   }
362 
363   if (*next_h > tleft) {
364     *finish = PETSC_TRUE;
365     *next_h = tleft;
366   } else *finish = PETSC_FALSE;
367   PetscFunctionReturn(0);
368 }
369 
370 PetscErrorCode TSGLLEAdaptCreate_Both(TSGLLEAdapt adapt)
371 {
372   PetscErrorCode ierr;
373   TSGLLEAdapt_Both *a;
374 
375   PetscFunctionBegin;
376   ierr = PetscNewLog(adapt,&a);CHKERRQ(ierr);
377   adapt->data         = (void*)a;
378   adapt->ops->choose  = TSGLLEAdaptChoose_Both;
379   adapt->ops->destroy = TSGLLEAdaptDestroy_JustFree;
380   PetscFunctionReturn(0);
381 }
382