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