xref: /petsc/src/ts/adapt/interface/tsadapt.c (revision 1c3436cf04c11f6e65a053bdffa658cacfbd163a) !
1 
2 #include <private/tsimpl.h> /*I  "petscts.h" I*/
3 
4 static PetscFList TSAdaptList;
5 static PetscBool  TSAdaptPackageInitialized;
6 static PetscBool  TSAdaptRegisterAllCalled;
7 static PetscClassId TSADAPT_CLASSID;
8 
9 EXTERN_C_BEGIN
10 PetscErrorCode  TSAdaptCreate_Basic(TSAdapt);
11 EXTERN_C_END
12 
13 #undef __FUNCT__
14 #define __FUNCT__ "TSAdaptRegister"
15 /*@C
16    TSAdaptRegister - see TSAdaptRegisterDynamic()
17 
18    Level: advanced
19 @*/
20 PetscErrorCode  TSAdaptRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(TSAdapt))
21 {
22   PetscErrorCode ierr;
23   char           fullname[PETSC_MAX_PATH_LEN];
24 
25   PetscFunctionBegin;
26   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
27   ierr = PetscFListAdd(&TSAdaptList,sname,fullname,(void(*)(void))function);CHKERRQ(ierr);
28   PetscFunctionReturn(0);
29 }
30 
31 #undef __FUNCT__
32 #define __FUNCT__ "TSAdaptRegisterAll"
33 /*@C
34   TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt
35 
36   Not Collective
37 
38   Level: advanced
39 
40 .keywords: TSAdapt, register, all
41 
42 .seealso: TSAdaptRegisterDestroy()
43 @*/
44 PetscErrorCode  TSAdaptRegisterAll(const char path[])
45 {
46   PetscErrorCode ierr;
47 
48   PetscFunctionBegin;
49   ierr = TSAdaptRegisterDynamic(TSADAPTBASIC,path,"TSAdaptCreate_Basic",TSAdaptCreate_Basic);CHKERRQ(ierr);
50   PetscFunctionReturn(0);
51 }
52 
53 #undef __FUNCT__
54 #define __FUNCT__ "TSAdaptFinalizePackage"
55 /*@C
56   TSFinalizePackage - This function destroys everything in the TS package. It is
57   called from PetscFinalize().
58 
59   Level: developer
60 
61 .keywords: Petsc, destroy, package
62 .seealso: PetscFinalize()
63 @*/
64 PetscErrorCode  TSAdaptFinalizePackage(void)
65 {
66   PetscFunctionBegin;
67   TSAdaptPackageInitialized = PETSC_FALSE;
68   TSAdaptRegisterAllCalled  = PETSC_FALSE;
69   TSAdaptList               = PETSC_NULL;
70   PetscFunctionReturn(0);
71 }
72 
73 #undef __FUNCT__
74 #define __FUNCT__ "TSAdaptInitializePackage"
75 /*@C
76   TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is
77   called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to
78   TSCreate_GL() when using static libraries.
79 
80   Input Parameter:
81   path - The dynamic library path, or PETSC_NULL
82 
83   Level: developer
84 
85 .keywords: TSAdapt, initialize, package
86 .seealso: PetscInitialize()
87 @*/
88 PetscErrorCode  TSAdaptInitializePackage(const char path[])
89 {
90   PetscErrorCode ierr;
91 
92   PetscFunctionBegin;
93   if (TSAdaptPackageInitialized) PetscFunctionReturn(0);
94   TSAdaptPackageInitialized = PETSC_TRUE;
95   ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr);
96   ierr = TSAdaptRegisterAll(path);CHKERRQ(ierr);
97   ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr);
98   PetscFunctionReturn(0);
99 }
100 
101 #undef __FUNCT__
102 #define __FUNCT__ "TSAdaptRegisterDestroy"
103 /*@C
104    TSAdaptRegisterDestroy - Frees the list of adaptivity schemes that were registered by TSAdaptRegister()/TSAdaptRegisterDynamic().
105 
106    Not Collective
107 
108    Level: advanced
109 
110 .keywords: TSAdapt, register, destroy
111 .seealso: TSAdaptRegister(), TSAdaptRegisterAll(), TSAdaptRegisterDynamic()
112 @*/
113 PetscErrorCode  TSAdaptRegisterDestroy(void)
114 {
115   PetscErrorCode ierr;
116 
117   PetscFunctionBegin;
118   ierr = PetscFListDestroy(&TSAdaptList);CHKERRQ(ierr);
119   TSAdaptRegisterAllCalled = PETSC_FALSE;
120   PetscFunctionReturn(0);
121 }
122 
123 
124 #undef __FUNCT__
125 #define __FUNCT__ "TSAdaptSetType"
126 PetscErrorCode  TSAdaptSetType(TSAdapt adapt,const TSAdaptType type)
127 {
128   PetscErrorCode ierr,(*r)(TSAdapt);
129 
130   PetscFunctionBegin;
131   ierr = PetscFListFind(TSAdaptList,((PetscObject)adapt)->comm,type,PETSC_TRUE,(void(**)(void))&r);CHKERRQ(ierr);
132   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type);
133   if (((PetscObject)adapt)->type_name) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);}
134   ierr = (*r)(adapt);CHKERRQ(ierr);
135   ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr);
136   PetscFunctionReturn(0);
137 }
138 
139 #undef __FUNCT__
140 #define __FUNCT__ "TSAdaptSetOptionsPrefix"
141 PetscErrorCode  TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[])
142 {
143   PetscErrorCode ierr;
144 
145   PetscFunctionBegin;
146   ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr);
147   PetscFunctionReturn(0);
148 }
149 
150 #undef __FUNCT__
151 #define __FUNCT__ "TSAdaptView"
152 PetscErrorCode  TSAdaptView(TSAdapt adapt,PetscViewer viewer)
153 {
154   PetscErrorCode ierr;
155   PetscBool      iascii;
156 
157   PetscFunctionBegin;
158   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
159   if (iascii) {
160     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer,"TSAdapt Object");CHKERRQ(ierr);
161     ierr = PetscViewerASCIIPrintf(viewer,"number of candidates %D\n",adapt->candidates.n);CHKERRQ(ierr);
162     if (adapt->ops->view) {
163       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
164       ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
165       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
166     }
167   } else {
168     SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Viewer type %s not supported",((PetscObject)viewer)->type_name);
169   }
170   PetscFunctionReturn(0);
171 }
172 
173 #undef __FUNCT__
174 #define __FUNCT__ "TSAdaptDestroy"
175 PetscErrorCode  TSAdaptDestroy(TSAdapt *adapt)
176 {
177   PetscErrorCode ierr;
178 
179   PetscFunctionBegin;
180   if (!*adapt) PetscFunctionReturn(0);
181   PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1);
182   if (--((PetscObject)(*adapt))->refct > 0) {*adapt = 0; PetscFunctionReturn(0);}
183   if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);}
184   ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr);
185   ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr);
186   PetscFunctionReturn(0);
187 }
188 
189 #undef __FUNCT__
190 #define __FUNCT__ "TSAdaptSetMonitor"
191 /*@
192    TSAdaptSetMonitor - Monitor the choices made by the adaptive controller
193 
194    Collective on TSAdapt
195 
196    Input Arguments:
197 +  adapt - adaptive controller context
198 -  flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable
199 
200    Level: intermediate
201 
202 .seealso: TSAdaptChoose()
203 @*/
204 PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg)
205 {
206   PetscErrorCode ierr;
207 
208   PetscFunctionBegin;
209   if (flg) {
210     if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(((PetscObject)adapt)->comm,"stdout",&adapt->monitor);CHKERRQ(ierr);}
211   } else {
212     ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr);
213   }
214   PetscFunctionReturn(0);
215 }
216 
217 #undef __FUNCT__
218 #define __FUNCT__ "TSAdaptSetStepLimits"
219 /*@
220    TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller
221 
222    Logically Collective
223 
224    Input Arguments:
225 +  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
226 .  hmin - minimum time step
227 -  hmax - maximum time step
228 
229    Options Database Keys:
230 +  -ts_adapt_dt_min - minimum time step
231 -  -ts_adapt_dt_max - maximum time step
232 
233    Level: intermediate
234 
235 .seealso: TSAdapt
236 @*/
237 PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax)
238 {
239 
240   PetscFunctionBegin;
241   if (hmin != PETSC_DECIDE) adapt->dt_min = hmin;
242   if (hmax != PETSC_DECIDE) adapt->dt_max = hmax;
243   PetscFunctionReturn(0);
244 }
245 
246 #undef __FUNCT__
247 #define __FUNCT__ "TSAdaptSetFromOptions"
248 /*@
249    TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options.
250 
251    Collective on TSAdapt
252 
253    Input Parameter:
254 .  adapt - the TSAdapt context
255 
256    Options Database Keys:
257 .  -ts_adapt_type <type> - basic
258 
259    Level: advanced
260 
261    Notes:
262    This function is automatically called by TSSetFromOptions()
263 
264 .keywords: TS, TSGetAdapt(), TSAdaptSetType()
265 
266 .seealso: TSGetType()
267 @*/
268 PetscErrorCode  TSAdaptSetFromOptions(TSAdapt adapt)
269 {
270   PetscErrorCode ierr;
271   char           type[256] = TSADAPTBASIC;
272   PetscBool      set,flg;
273 
274   PetscFunctionBegin;
275   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
276   * function can only be called from inside TSSetFromOptions_GL()  */
277   ierr = PetscOptionsHead("TS Adaptivity options");CHKERRQ(ierr);
278   ierr = PetscOptionsList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList,
279                           ((PetscObject)adapt)->type_name?((PetscObject)adapt)->type_name:type,type,sizeof type,&flg);CHKERRQ(ierr);
280   if (flg || !((PetscObject)adapt)->type_name) {
281     ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr);
282   }
283   if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(adapt);CHKERRQ(ierr);}
284   ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",adapt->dt_min,&adapt->dt_min,PETSC_NULL);CHKERRQ(ierr);
285   ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",adapt->dt_max,&adapt->dt_max,PETSC_NULL);CHKERRQ(ierr);
286   ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
287   if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);}
288   ierr = PetscOptionsTail();CHKERRQ(ierr);
289   PetscFunctionReturn(0);
290 }
291 
292 #undef __FUNCT__
293 #define __FUNCT__ "TSAdaptCandidatesClear"
294 /*@
295    TSAdaptCandidatesClear - clear any previously set candidate schemes
296 
297    Logically Collective
298 
299    Input Argument:
300 .  adapt - adaptive controller
301 
302    Level: developer
303 
304 .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose()
305 @*/
306 PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
307 {
308   PetscErrorCode ierr;
309 
310   PetscFunctionBegin;
311   ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr);
312   PetscFunctionReturn(0);
313 }
314 
315 #undef __FUNCT__
316 #define __FUNCT__ "TSAdaptCandidateAdd"
317 /*@C
318    TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from
319 
320    Logically Collective
321 
322    Input Arguments:
323 +  adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate()
324 .  name - name of the candidate scheme to add
325 .  order - order of the candidate scheme
326 .  stageorder - stage order of the candidate scheme
327 .  leadingerror - leading error coefficient of the candidate scheme
328 .  cost - relative measure of the amount of work required for the candidate scheme
329 -  inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme
330 
331    Note:
332    This routine is not available in Fortran.
333 
334    Level: developer
335 
336 .seealso: TSAdaptCandidatesClear(), TSAdaptChoose()
337 @*/
338 PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal leadingerror,PetscReal cost,PetscBool inuse)
339 {
340   PetscInt c;
341 
342   PetscFunctionBegin;
343   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
344   if (order < 1) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order);
345   if (inuse) {
346     if (adapt->candidates.inuse_set) SETERRQ(((PetscObject)adapt)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()");
347     adapt->candidates.inuse_set = PETSC_TRUE;
348   }
349   /* first slot if this is the current scheme, otherwise the next available slot */
350   c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;
351   adapt->candidates.name[c]         = name;
352   adapt->candidates.order[c]        = order;
353   adapt->candidates.stageorder[c]   = stageorder;
354   adapt->candidates.leadingerror[c] = leadingerror;
355   adapt->candidates.cost[c]         = cost;
356   adapt->candidates.n++;
357   PetscFunctionReturn(0);
358 }
359 
360 #undef __FUNCT__
361 #define __FUNCT__ "TSAdaptChoose"
362 /*@C
363    TSAdaptChoose - choose which method and step size to use for the next step
364 
365    Logically Collective
366 
367    Input Arguments:
368 +  adapt - adaptive contoller
369 -  h - current step size
370 
371    Output Arguments:
372 +  next_sc - scheme to use for the next step
373 .  next_h - step size to use for the next step
374 -  accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size
375 
376    Note:
377    The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is
378    being retried after an initial rejection.
379 
380    Level: developer
381 
382 .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd()
383 @*/
384 PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept)
385 {
386   PetscErrorCode ierr;
387 
388   PetscFunctionBegin;
389   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
390   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
391   PetscValidIntPointer(next_sc,4);
392   PetscValidPointer(next_h,5);
393   PetscValidIntPointer(accept,6);
394   if (adapt->candidates.n < 1) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_WRONGSTATE,"%D candidates have been registered",adapt->candidates.n);
395   if (!adapt->candidates.inuse_set) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_WRONGSTATE,"The current in-use scheme is not among the %D candidates",adapt->candidates.n);
396   ierr = (*adapt->ops->choose)(adapt,ts,h,next_sc,next_h,accept);CHKERRQ(ierr);
397   if (*next_sc < 0 || adapt->candidates.n <= *next_sc) SETERRQ2(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Chosen scheme %D not in valid range 0..%D",*next_sc,adapt->candidates.n-1);
398   if (!(*next_h > 0.)) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %G must be positive",*next_h);
399 
400   if (adapt->monitor) {
401     ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
402     ierr = PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt '%s': step %s  family='%s' scheme=%D:'%s' dt=%G\n",((PetscObject)adapt)->type_name,accept?"accepted":"rejected",((PetscObject)ts)->type_name,*next_sc,adapt->candidates.name[*next_sc],*next_h);CHKERRQ(ierr);
403     ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
404   }
405   PetscFunctionReturn(0);
406 }
407 
408 #undef __FUNCT__
409 #define __FUNCT__ "TSAdaptCreate"
410 /*@
411   TSAdaptCreate - create an adaptive controller context for time stepping
412 
413   Collective on MPI_Comm
414 
415   Input Parameter:
416 . comm - The communicator
417 
418   Output Parameter:
419 . adapt - new TSAdapt object
420 
421   Level: developer
422 
423   Notes:
424   TSAdapt creation is handled by TS, so users should not need to call this function.
425 
426 .keywords: TSAdapt, create
427 .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy()
428 @*/
429 PetscErrorCode  TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt)
430 {
431   PetscErrorCode ierr;
432   TSAdapt adapt;
433 
434   PetscFunctionBegin;
435   *inadapt = 0;
436   ierr = PetscHeaderCreate(adapt,_p_TSAdapt,struct _TSAdaptOps,TSADAPT_CLASSID,0,"TSAdapt","General Linear adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr);
437 
438   adapt->dt_min = 1e-20;
439   adapt->dt_max = 1e50;
440 
441   *inadapt = adapt;
442   PetscFunctionReturn(0);
443 }
444