xref: /petsc/src/ts/adapt/interface/tsadapt.c (revision d0288e4feab9efbf18cfa81e5796fd5fb5e8e9b7)
184df9cb4SJed Brown 
2af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I  "petscts.h" I*/
384df9cb4SJed Brown 
41b9b13dfSLisandro Dalcin PetscClassId TSADAPT_CLASSID;
51b9b13dfSLisandro Dalcin 
6140e18c1SBarry Smith static PetscFunctionList TSAdaptList;
784df9cb4SJed Brown static PetscBool         TSAdaptPackageInitialized;
884df9cb4SJed Brown static PetscBool         TSAdaptRegisterAllCalled;
984df9cb4SJed Brown 
10726095e4SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt);
118cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_None(TSAdapt);
121566a47fSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt);
138cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt);
1484df9cb4SJed Brown 
1584df9cb4SJed Brown /*@C
161c84c290SBarry Smith    TSAdaptRegister -  adds a TSAdapt implementation
171c84c290SBarry Smith 
181c84c290SBarry Smith    Not Collective
191c84c290SBarry Smith 
201c84c290SBarry Smith    Input Parameters:
211c84c290SBarry Smith +  name_scheme - name of user-defined adaptivity scheme
221c84c290SBarry Smith -  routine_create - routine to create method context
231c84c290SBarry Smith 
241c84c290SBarry Smith    Notes:
251c84c290SBarry Smith    TSAdaptRegister() may be called multiple times to add several user-defined families.
261c84c290SBarry Smith 
271c84c290SBarry Smith    Sample usage:
281c84c290SBarry Smith .vb
29bdf89e91SBarry Smith    TSAdaptRegister("my_scheme",MySchemeCreate);
301c84c290SBarry Smith .ve
311c84c290SBarry Smith 
321c84c290SBarry Smith    Then, your scheme can be chosen with the procedural interface via
331c84c290SBarry Smith $     TSAdaptSetType(ts,"my_scheme")
341c84c290SBarry Smith    or at runtime via the option
351c84c290SBarry Smith $     -ts_adapt_type my_scheme
3684df9cb4SJed Brown 
3784df9cb4SJed Brown    Level: advanced
381c84c290SBarry Smith 
391c84c290SBarry Smith .keywords: TSAdapt, register
401c84c290SBarry Smith 
411c84c290SBarry Smith .seealso: TSAdaptRegisterAll()
4284df9cb4SJed Brown @*/
43bdf89e91SBarry Smith PetscErrorCode  TSAdaptRegister(const char sname[],PetscErrorCode (*function)(TSAdapt))
4484df9cb4SJed Brown {
4584df9cb4SJed Brown   PetscErrorCode ierr;
4684df9cb4SJed Brown 
4784df9cb4SJed Brown   PetscFunctionBegin;
48a240a19fSJed Brown   ierr = PetscFunctionListAdd(&TSAdaptList,sname,function);CHKERRQ(ierr);
4984df9cb4SJed Brown   PetscFunctionReturn(0);
5084df9cb4SJed Brown }
5184df9cb4SJed Brown 
5284df9cb4SJed Brown /*@C
5384df9cb4SJed Brown   TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt
5484df9cb4SJed Brown 
5584df9cb4SJed Brown   Not Collective
5684df9cb4SJed Brown 
5784df9cb4SJed Brown   Level: advanced
5884df9cb4SJed Brown 
5984df9cb4SJed Brown .keywords: TSAdapt, register, all
6084df9cb4SJed Brown 
6184df9cb4SJed Brown .seealso: TSAdaptRegisterDestroy()
6284df9cb4SJed Brown @*/
63607a6623SBarry Smith PetscErrorCode  TSAdaptRegisterAll(void)
6484df9cb4SJed Brown {
6584df9cb4SJed Brown   PetscErrorCode ierr;
6684df9cb4SJed Brown 
6784df9cb4SJed Brown   PetscFunctionBegin;
680f51fdf8SToby Isaac   if (TSAdaptRegisterAllCalled) PetscFunctionReturn(0);
690f51fdf8SToby Isaac   TSAdaptRegisterAllCalled = PETSC_TRUE;
70726095e4SEmil Constantinescu   ierr = TSAdaptRegister(TSADAPTGLEE ,TSAdaptCreate_GLEE);CHKERRQ(ierr);
71bdf89e91SBarry Smith   ierr = TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None);CHKERRQ(ierr);
721566a47fSLisandro Dalcin   ierr = TSAdaptRegister(TSADAPTBASIC,TSAdaptCreate_Basic);CHKERRQ(ierr);
73bdf89e91SBarry Smith   ierr = TSAdaptRegister(TSADAPTCFL,  TSAdaptCreate_CFL);CHKERRQ(ierr);
7484df9cb4SJed Brown   PetscFunctionReturn(0);
7584df9cb4SJed Brown }
7684df9cb4SJed Brown 
7784df9cb4SJed Brown /*@C
78560360afSLisandro Dalcin   TSAdaptFinalizePackage - This function destroys everything in the TS package. It is
7984df9cb4SJed Brown   called from PetscFinalize().
8084df9cb4SJed Brown 
8184df9cb4SJed Brown   Level: developer
8284df9cb4SJed Brown 
8384df9cb4SJed Brown .keywords: Petsc, destroy, package
8484df9cb4SJed Brown .seealso: PetscFinalize()
8584df9cb4SJed Brown @*/
8684df9cb4SJed Brown PetscErrorCode  TSAdaptFinalizePackage(void)
8784df9cb4SJed Brown {
8837e93019SBarry Smith   PetscErrorCode ierr;
8937e93019SBarry Smith 
9084df9cb4SJed Brown   PetscFunctionBegin;
9137e93019SBarry Smith   ierr = PetscFunctionListDestroy(&TSAdaptList);CHKERRQ(ierr);
9284df9cb4SJed Brown   TSAdaptPackageInitialized = PETSC_FALSE;
9384df9cb4SJed Brown   TSAdaptRegisterAllCalled  = PETSC_FALSE;
9484df9cb4SJed Brown   PetscFunctionReturn(0);
9584df9cb4SJed Brown }
9684df9cb4SJed Brown 
9784df9cb4SJed Brown /*@C
9884df9cb4SJed Brown   TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is
9984df9cb4SJed Brown   called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to
10026d28e4eSEmil Constantinescu   TSCreate_GLLE() when using static libraries.
10184df9cb4SJed Brown 
10284df9cb4SJed Brown   Level: developer
10384df9cb4SJed Brown 
10484df9cb4SJed Brown .keywords: TSAdapt, initialize, package
10584df9cb4SJed Brown .seealso: PetscInitialize()
10684df9cb4SJed Brown @*/
107607a6623SBarry Smith PetscErrorCode  TSAdaptInitializePackage(void)
10884df9cb4SJed Brown {
10984df9cb4SJed Brown   PetscErrorCode ierr;
11084df9cb4SJed Brown 
11184df9cb4SJed Brown   PetscFunctionBegin;
11284df9cb4SJed Brown   if (TSAdaptPackageInitialized) PetscFunctionReturn(0);
11384df9cb4SJed Brown   TSAdaptPackageInitialized = PETSC_TRUE;
11484df9cb4SJed Brown   ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr);
115607a6623SBarry Smith   ierr = TSAdaptRegisterAll();CHKERRQ(ierr);
11684df9cb4SJed Brown   ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr);
11784df9cb4SJed Brown   PetscFunctionReturn(0);
11884df9cb4SJed Brown }
11984df9cb4SJed Brown 
1207eef6055SBarry Smith /*@C
1217eef6055SBarry Smith   TSAdaptSetType - sets the approach used for the error adapter, currently there is only TSADAPTBASIC and TSADAPTNONE
1227eef6055SBarry Smith 
1237eef6055SBarry Smith   Logicially Collective on TSAdapt
1247eef6055SBarry Smith 
1257eef6055SBarry Smith   Input Parameter:
126*d0288e4fSLisandro Dalcin + adapt - the TS adapter, most likely obtained with TSGetAdapt()
1277eef6055SBarry Smith - type - either  TSADAPTBASIC or TSADAPTNONE
1287eef6055SBarry Smith 
1297eef6055SBarry Smith   Options Database:
1307eef6055SBarry Smith .  -ts_adapt_type basic or none - to setting the adapter type
1317eef6055SBarry Smith 
1327eef6055SBarry Smith   Level: intermediate
1337eef6055SBarry Smith 
1347eef6055SBarry Smith .keywords: TSAdapt, create
1357eef6055SBarry Smith 
1367eef6055SBarry Smith .seealso: TSGetAdapt(), TSAdaptDestroy(), TSAdaptType, TSAdaptGetType()
1377eef6055SBarry Smith @*/
13819fd82e9SBarry Smith PetscErrorCode  TSAdaptSetType(TSAdapt adapt,TSAdaptType type)
13984df9cb4SJed Brown {
140ef749922SLisandro Dalcin   PetscBool      match;
14184df9cb4SJed Brown   PetscErrorCode ierr,(*r)(TSAdapt);
14284df9cb4SJed Brown 
14384df9cb4SJed Brown   PetscFunctionBegin;
1444782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
145ef749922SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)adapt,type,&match);CHKERRQ(ierr);
146ef749922SLisandro Dalcin   if (match) PetscFunctionReturn(0);
1471c9cd337SJed Brown   ierr = PetscFunctionListFind(TSAdaptList,type,&r);CHKERRQ(ierr);
14884df9cb4SJed Brown   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type);
149ef749922SLisandro Dalcin   if (adapt->ops->destroy) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);}
1504cefc2ffSBarry Smith   ierr = PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));CHKERRQ(ierr);
15184df9cb4SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr);
15268ae941aSLisandro Dalcin   ierr = (*r)(adapt);CHKERRQ(ierr);
15384df9cb4SJed Brown   PetscFunctionReturn(0);
15484df9cb4SJed Brown }
15584df9cb4SJed Brown 
156*d0288e4fSLisandro Dalcin /*@C
157*d0288e4fSLisandro Dalcin   TSAdaptGetType - gets the TS adapter method type (as a string).
158*d0288e4fSLisandro Dalcin 
159*d0288e4fSLisandro Dalcin   Not Collective
160*d0288e4fSLisandro Dalcin 
161*d0288e4fSLisandro Dalcin   Input Parameter:
162*d0288e4fSLisandro Dalcin . adapt - The TS adapter, most likely obtained with TSGetAdapt()
163*d0288e4fSLisandro Dalcin 
164*d0288e4fSLisandro Dalcin   Output Parameter:
165*d0288e4fSLisandro Dalcin . type - The name of TS adapter method
166*d0288e4fSLisandro Dalcin 
167*d0288e4fSLisandro Dalcin   Level: intermediate
168*d0288e4fSLisandro Dalcin 
169*d0288e4fSLisandro Dalcin .keywords: TSAdapt, get, type
170*d0288e4fSLisandro Dalcin .seealso TSAdaptSetType()
171*d0288e4fSLisandro Dalcin @*/
172*d0288e4fSLisandro Dalcin PetscErrorCode TSAdaptGetType(TSAdapt adapt,TSAdaptType *type)
173*d0288e4fSLisandro Dalcin {
174*d0288e4fSLisandro Dalcin   PetscFunctionBegin;
175*d0288e4fSLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
176*d0288e4fSLisandro Dalcin   PetscValidPointer(type,2);
177*d0288e4fSLisandro Dalcin   *type = ((PetscObject)adapt)->type_name;
178*d0288e4fSLisandro Dalcin   PetscFunctionReturn(0);
179*d0288e4fSLisandro Dalcin }
180*d0288e4fSLisandro Dalcin 
18184df9cb4SJed Brown PetscErrorCode  TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[])
18284df9cb4SJed Brown {
18384df9cb4SJed Brown   PetscErrorCode ierr;
18484df9cb4SJed Brown 
18584df9cb4SJed Brown   PetscFunctionBegin;
1864782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
18784df9cb4SJed Brown   ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr);
18884df9cb4SJed Brown   PetscFunctionReturn(0);
18984df9cb4SJed Brown }
19084df9cb4SJed Brown 
191ad6bc421SBarry Smith /*@C
192ad6bc421SBarry Smith   TSAdaptLoad - Loads a TSAdapt that has been stored in binary  with TSAdaptView().
193ad6bc421SBarry Smith 
194ad6bc421SBarry Smith   Collective on PetscViewer
195ad6bc421SBarry Smith 
196ad6bc421SBarry Smith   Input Parameters:
197ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or
198ad6bc421SBarry Smith            some related function before a call to TSAdaptLoad().
199ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
200ad6bc421SBarry Smith            HDF5 file viewer, obtained from PetscViewerHDF5Open()
201ad6bc421SBarry Smith 
202ad6bc421SBarry Smith    Level: intermediate
203ad6bc421SBarry Smith 
204ad6bc421SBarry Smith   Notes:
205ad6bc421SBarry Smith    The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored.
206ad6bc421SBarry Smith 
207ad6bc421SBarry Smith   Notes for advanced users:
208ad6bc421SBarry Smith   Most users should not need to know the details of the binary storage
209ad6bc421SBarry Smith   format, since TSAdaptLoad() and TSAdaptView() completely hide these details.
210ad6bc421SBarry Smith   But for anyone who's interested, the standard binary matrix storage
211ad6bc421SBarry Smith   format is
212ad6bc421SBarry Smith .vb
213ad6bc421SBarry Smith      has not yet been determined
214ad6bc421SBarry Smith .ve
215ad6bc421SBarry Smith 
216ad6bc421SBarry Smith .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad()
217ad6bc421SBarry Smith @*/
2184782b174SLisandro Dalcin PetscErrorCode  TSAdaptLoad(TSAdapt adapt,PetscViewer viewer)
219ad6bc421SBarry Smith {
220ad6bc421SBarry Smith   PetscErrorCode ierr;
221ad6bc421SBarry Smith   PetscBool      isbinary;
222ad6bc421SBarry Smith   char           type[256];
223ad6bc421SBarry Smith 
224ad6bc421SBarry Smith   PetscFunctionBegin;
2254782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
226ad6bc421SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
227ad6bc421SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
228ad6bc421SBarry Smith   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
229ad6bc421SBarry Smith 
230060da220SMatthew G. Knepley   ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr);
2314782b174SLisandro Dalcin   ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr);
2324782b174SLisandro Dalcin   if (adapt->ops->load) {
2334782b174SLisandro Dalcin     ierr = (*adapt->ops->load)(adapt,viewer);CHKERRQ(ierr);
234ad6bc421SBarry Smith   }
235ad6bc421SBarry Smith   PetscFunctionReturn(0);
236ad6bc421SBarry Smith }
237ad6bc421SBarry Smith 
23884df9cb4SJed Brown PetscErrorCode  TSAdaptView(TSAdapt adapt,PetscViewer viewer)
23984df9cb4SJed Brown {
24084df9cb4SJed Brown   PetscErrorCode ierr;
241ad6bc421SBarry Smith   PetscBool      iascii,isbinary;
24284df9cb4SJed Brown 
24384df9cb4SJed Brown   PetscFunctionBegin;
2444782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
2454782b174SLisandro Dalcin   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer);CHKERRQ(ierr);}
2464782b174SLisandro Dalcin   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
2474782b174SLisandro Dalcin   PetscCheckSameComm(adapt,1,viewer,2);
248251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
249ad6bc421SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
25084df9cb4SJed Brown   if (iascii) {
251dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr);
25284df9cb4SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  number of candidates %D\n",adapt->candidates.n);CHKERRQ(ierr);
253bf997491SLisandro Dalcin     if (adapt->always_accept) {ierr = PetscViewerASCIIPrintf(viewer,"  always accepting steps\n");CHKERRQ(ierr);}
254bc0b8257SSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  maximum allowed timestep %g\n",(double)adapt->dt_max);CHKERRQ(ierr);
255bc0b8257SSatish Balay     ierr = PetscViewerASCIIPrintf(viewer,"  minimum allowed timestep %g\n",(double)adapt->dt_min);CHKERRQ(ierr);
25684df9cb4SJed Brown     if (adapt->ops->view) {
25784df9cb4SJed Brown       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
25884df9cb4SJed Brown       ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
25984df9cb4SJed Brown       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
26084df9cb4SJed Brown     }
261ad6bc421SBarry Smith   } else if (isbinary) {
262ad6bc421SBarry Smith     char type[256];
263ad6bc421SBarry Smith 
264ad6bc421SBarry Smith     /* need to save FILE_CLASS_ID for adapt class */
265ad6bc421SBarry Smith     ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr);
266ad6bc421SBarry Smith     ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr);
267bbd56ea5SKarl Rupp   } else if (adapt->ops->view) {
268f2c2a1b9SBarry Smith     ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
269f2c2a1b9SBarry Smith   }
27084df9cb4SJed Brown   PetscFunctionReturn(0);
27184df9cb4SJed Brown }
27284df9cb4SJed Brown 
273099af0a3SLisandro Dalcin /*@
274099af0a3SLisandro Dalcin    TSAdaptReset - Resets a TSAdapt context.
275099af0a3SLisandro Dalcin 
276099af0a3SLisandro Dalcin    Collective on TS
277099af0a3SLisandro Dalcin 
278099af0a3SLisandro Dalcin    Input Parameter:
279099af0a3SLisandro Dalcin .  adapt - the TSAdapt context obtained from TSAdaptCreate()
280099af0a3SLisandro Dalcin 
281099af0a3SLisandro Dalcin    Level: developer
282099af0a3SLisandro Dalcin 
283099af0a3SLisandro Dalcin .seealso: TSAdaptCreate(), TSAdaptDestroy()
284099af0a3SLisandro Dalcin @*/
285099af0a3SLisandro Dalcin PetscErrorCode  TSAdaptReset(TSAdapt adapt)
286099af0a3SLisandro Dalcin {
287099af0a3SLisandro Dalcin   PetscErrorCode ierr;
288099af0a3SLisandro Dalcin 
289099af0a3SLisandro Dalcin   PetscFunctionBegin;
290099af0a3SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
291099af0a3SLisandro Dalcin   if (adapt->ops->reset) {ierr = (*adapt->ops->reset)(adapt);CHKERRQ(ierr);}
292099af0a3SLisandro Dalcin   PetscFunctionReturn(0);
293099af0a3SLisandro Dalcin }
294099af0a3SLisandro Dalcin 
29584df9cb4SJed Brown PetscErrorCode  TSAdaptDestroy(TSAdapt *adapt)
29684df9cb4SJed Brown {
29784df9cb4SJed Brown   PetscErrorCode ierr;
29884df9cb4SJed Brown 
29984df9cb4SJed Brown   PetscFunctionBegin;
30084df9cb4SJed Brown   if (!*adapt) PetscFunctionReturn(0);
30184df9cb4SJed Brown   PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1);
3024782b174SLisandro Dalcin   if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);}
303099af0a3SLisandro Dalcin 
304099af0a3SLisandro Dalcin   ierr = TSAdaptReset(*adapt);CHKERRQ(ierr);
305099af0a3SLisandro Dalcin 
30684df9cb4SJed Brown   if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);}
3071c3436cfSJed Brown   ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr);
30884df9cb4SJed Brown   ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr);
30984df9cb4SJed Brown   PetscFunctionReturn(0);
31084df9cb4SJed Brown }
31184df9cb4SJed Brown 
3121c3436cfSJed Brown /*@
3131c3436cfSJed Brown    TSAdaptSetMonitor - Monitor the choices made by the adaptive controller
3141c3436cfSJed Brown 
3151c3436cfSJed Brown    Collective on TSAdapt
3161c3436cfSJed Brown 
3171c3436cfSJed Brown    Input Arguments:
3181c3436cfSJed Brown +  adapt - adaptive controller context
3191c3436cfSJed Brown -  flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable
3201c3436cfSJed Brown 
321bf997491SLisandro Dalcin    Options Database Keys:
322bf997491SLisandro Dalcin .  -ts_adapt_monitor
323bf997491SLisandro Dalcin 
3241c3436cfSJed Brown    Level: intermediate
3251c3436cfSJed Brown 
3261c3436cfSJed Brown .seealso: TSAdaptChoose()
3271c3436cfSJed Brown @*/
3281c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg)
3291c3436cfSJed Brown {
3301c3436cfSJed Brown   PetscErrorCode ierr;
3311c3436cfSJed Brown 
3321c3436cfSJed Brown   PetscFunctionBegin;
3334782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
3344782b174SLisandro Dalcin   PetscValidLogicalCollectiveBool(adapt,flg,2);
3351c3436cfSJed Brown   if (flg) {
336ce94432eSBarry Smith     if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);}
3371c3436cfSJed Brown   } else {
3381c3436cfSJed Brown     ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr);
3391c3436cfSJed Brown   }
3401c3436cfSJed Brown   PetscFunctionReturn(0);
3411c3436cfSJed Brown }
3421c3436cfSJed Brown 
3430873808bSJed Brown /*@C
344bf997491SLisandro Dalcin    TSAdaptSetCheckStage - Set a callback to check convergence for a stage
3450873808bSJed Brown 
3460873808bSJed Brown    Logically collective on TSAdapt
3470873808bSJed Brown 
3480873808bSJed Brown    Input Arguments:
3490873808bSJed Brown +  adapt - adaptive controller context
3500873808bSJed Brown -  func - stage check function
3510873808bSJed Brown 
3520873808bSJed Brown    Arguments of func:
3530873808bSJed Brown $  PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept)
3540873808bSJed Brown 
3550873808bSJed Brown +  adapt - adaptive controller context
3560873808bSJed Brown .  ts - time stepping context
3570873808bSJed Brown -  accept - pending choice of whether to accept, can be modified by this routine
3580873808bSJed Brown 
3590873808bSJed Brown    Level: advanced
3600873808bSJed Brown 
3610873808bSJed Brown .seealso: TSAdaptChoose()
3620873808bSJed Brown @*/
363b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscReal,Vec,PetscBool*))
3640873808bSJed Brown {
3650873808bSJed Brown 
3660873808bSJed Brown   PetscFunctionBegin;
3670873808bSJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
36868ae941aSLisandro Dalcin   adapt->checkstage = func;
3690873808bSJed Brown   PetscFunctionReturn(0);
3700873808bSJed Brown }
3710873808bSJed Brown 
3721c3436cfSJed Brown /*@
373bf997491SLisandro Dalcin    TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of
374bf997491SLisandro Dalcin    any error or stability condition not meeting the prescribed goal.
375bf997491SLisandro Dalcin 
376bf997491SLisandro Dalcin 
377bf997491SLisandro Dalcin    Input Arguments:
378bf997491SLisandro Dalcin +  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
379bf997491SLisandro Dalcin -  flag - whether to always accept steps
380bf997491SLisandro Dalcin 
381bf997491SLisandro Dalcin    Options Database Keys:
382bf997491SLisandro Dalcin .  -ts_adapt_always_accept
383bf997491SLisandro Dalcin 
384bf997491SLisandro Dalcin    Level: intermediate
385bf997491SLisandro Dalcin 
386bf997491SLisandro Dalcin .seealso: TSAdapt, TSAdaptChoose()
387bf997491SLisandro Dalcin @*/
388bf997491SLisandro Dalcin PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt,PetscBool flag)
389bf997491SLisandro Dalcin {
390bf997491SLisandro Dalcin   PetscFunctionBegin;
391bf997491SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
392bf997491SLisandro Dalcin   PetscValidLogicalCollectiveBool(adapt,flag,2);
393bf997491SLisandro Dalcin   adapt->always_accept = flag;
394bf997491SLisandro Dalcin   PetscFunctionReturn(0);
395bf997491SLisandro Dalcin }
396bf997491SLisandro Dalcin 
397bf997491SLisandro Dalcin /*@
3981c3436cfSJed Brown    TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller
3991c3436cfSJed Brown 
4001c3436cfSJed Brown    Logically Collective
4011c3436cfSJed Brown 
4021c3436cfSJed Brown    Input Arguments:
403552698daSJed Brown +  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
4041c3436cfSJed Brown .  hmin - minimum time step
4051c3436cfSJed Brown -  hmax - maximum time step
4061c3436cfSJed Brown 
4071c3436cfSJed Brown    Options Database Keys:
4081c3436cfSJed Brown +  -ts_adapt_dt_min - minimum time step
4091c3436cfSJed Brown -  -ts_adapt_dt_max - maximum time step
4101c3436cfSJed Brown 
4111c3436cfSJed Brown    Level: intermediate
4121c3436cfSJed Brown 
413bf997491SLisandro Dalcin .seealso: TSAdapt, TSAdaptChoose()
4141c3436cfSJed Brown @*/
4151c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax)
4161c3436cfSJed Brown {
4171c3436cfSJed Brown 
4181c3436cfSJed Brown   PetscFunctionBegin;
4194782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
420b1f5bccaSLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt,hmin,2);
421b1f5bccaSLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt,hmax,3);
422b1f5bccaSLisandro Dalcin   if (hmin != PETSC_DEFAULT && hmin < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmin);
423b1f5bccaSLisandro Dalcin   if (hmax != PETSC_DEFAULT && hmax < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmax);
424b1f5bccaSLisandro Dalcin   if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin;
425b1f5bccaSLisandro Dalcin   if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax;
426b1f5bccaSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
427b1f5bccaSLisandro Dalcin   hmin = adapt->dt_min;
428b1f5bccaSLisandro Dalcin   hmax = adapt->dt_max;
429b1f5bccaSLisandro Dalcin   if (hmax <= hmin) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Maximum time step %g must geather than minimum time step %g",(double)hmax,(double)hmin);
430b1f5bccaSLisandro Dalcin #endif
4311c3436cfSJed Brown   PetscFunctionReturn(0);
4321c3436cfSJed Brown }
4331c3436cfSJed Brown 
434e55864a3SBarry Smith /*
43584df9cb4SJed Brown    TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options.
43684df9cb4SJed Brown 
43784df9cb4SJed Brown    Collective on TSAdapt
43884df9cb4SJed Brown 
43984df9cb4SJed Brown    Input Parameter:
44084df9cb4SJed Brown .  adapt - the TSAdapt context
44184df9cb4SJed Brown 
44284df9cb4SJed Brown    Options Database Keys:
44323de1d84SBarry Smith +  -ts_adapt_type <type> - basic
444bf997491SLisandro Dalcin .  -ts_adapt_always_accept - always accept steps regardless of error/stability goals
44523de1d84SBarry Smith .  -ts_adapt_dt_min <min> - minimum timestep to use
44623de1d84SBarry Smith .  -ts_adapt_dt_max <max> - maximum timestep to use
44723de1d84SBarry Smith .  -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails
44823de1d84SBarry Smith -  -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates
44984df9cb4SJed Brown 
45084df9cb4SJed Brown    Level: advanced
45184df9cb4SJed Brown 
45284df9cb4SJed Brown    Notes:
45384df9cb4SJed Brown    This function is automatically called by TSSetFromOptions()
45484df9cb4SJed Brown 
45523de1d84SBarry Smith .keywords: TS, TSGetAdapt(), TSAdaptSetType(), TSAdaptSetStepLimits()
45684df9cb4SJed Brown 
45784df9cb4SJed Brown .seealso: TSGetType()
458e55864a3SBarry Smith */
4594416b707SBarry Smith PetscErrorCode  TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt)
46084df9cb4SJed Brown {
46184df9cb4SJed Brown   PetscErrorCode ierr;
46284df9cb4SJed Brown   char           type[256] = TSADAPTBASIC;
463bf997491SLisandro Dalcin   PetscReal      hmin,hmax;
4641c3436cfSJed Brown   PetscBool      set,flg;
46584df9cb4SJed Brown 
46684df9cb4SJed Brown   PetscFunctionBegin;
4674782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
46884df9cb4SJed Brown   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
4691566a47fSLisandro Dalcin    * function can only be called from inside TSSetFromOptions()  */
470e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");CHKERRQ(ierr);
47123de1d84SBarry Smith   ierr = PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList,((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);CHKERRQ(ierr);
47284df9cb4SJed Brown   if (flg || !((PetscObject)adapt)->type_name) {
47384df9cb4SJed Brown     ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr);
47484df9cb4SJed Brown   }
475bf997491SLisandro Dalcin   ierr = PetscOptionsBool("-ts_adapt_always_accept","Always accept the step","TSAdaptSetAlwaysAccept",adapt->always_accept,&flg,&set);CHKERRQ(ierr);
476bf997491SLisandro Dalcin   if (set) {ierr = TSAdaptSetAlwaysAccept(adapt,flg);CHKERRQ(ierr);}
477bf997491SLisandro Dalcin   ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",hmin=adapt->dt_min,&hmin,&set);CHKERRQ(ierr);
478bf997491SLisandro Dalcin   ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",hmax=adapt->dt_max,&hmax,&flg);CHKERRQ(ierr);
479bf997491SLisandro Dalcin   if (set || flg) {ierr = TSAdaptSetStepLimits(adapt,hmin,hmax);CHKERRQ(ierr);}
4800298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_adapt_scale_solve_failed","Scale step by this factor if solve fails","",adapt->scale_solve_failed,&adapt->scale_solve_failed,NULL);CHKERRQ(ierr);
4811c3436cfSJed Brown   ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
482bf997491SLisandro Dalcin   if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);}
4837619abb3SShri   ierr = PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);CHKERRQ(ierr);
4847619abb3SShri   if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported");
485e55864a3SBarry Smith   if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);}
48684df9cb4SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
48784df9cb4SJed Brown   PetscFunctionReturn(0);
48884df9cb4SJed Brown }
48984df9cb4SJed Brown 
49084df9cb4SJed Brown /*@
49184df9cb4SJed Brown    TSAdaptCandidatesClear - clear any previously set candidate schemes
49284df9cb4SJed Brown 
49384df9cb4SJed Brown    Logically Collective
49484df9cb4SJed Brown 
49584df9cb4SJed Brown    Input Argument:
49684df9cb4SJed Brown .  adapt - adaptive controller
49784df9cb4SJed Brown 
49884df9cb4SJed Brown    Level: developer
49984df9cb4SJed Brown 
50084df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose()
50184df9cb4SJed Brown @*/
50284df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
50384df9cb4SJed Brown {
50484df9cb4SJed Brown   PetscErrorCode ierr;
50584df9cb4SJed Brown 
50684df9cb4SJed Brown   PetscFunctionBegin;
5074782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
50884df9cb4SJed Brown   ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr);
50984df9cb4SJed Brown   PetscFunctionReturn(0);
51084df9cb4SJed Brown }
51184df9cb4SJed Brown 
51284df9cb4SJed Brown /*@C
51384df9cb4SJed Brown    TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from
51484df9cb4SJed Brown 
51584df9cb4SJed Brown    Logically Collective
51684df9cb4SJed Brown 
51784df9cb4SJed Brown    Input Arguments:
518552698daSJed Brown +  adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate()
51984df9cb4SJed Brown .  name - name of the candidate scheme to add
52084df9cb4SJed Brown .  order - order of the candidate scheme
52184df9cb4SJed Brown .  stageorder - stage order of the candidate scheme
5228d59e960SJed Brown .  ccfl - stability coefficient relative to explicit Euler, used for CFL constraints
52384df9cb4SJed Brown .  cost - relative measure of the amount of work required for the candidate scheme
52484df9cb4SJed Brown -  inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme
52584df9cb4SJed Brown 
52684df9cb4SJed Brown    Note:
52784df9cb4SJed Brown    This routine is not available in Fortran.
52884df9cb4SJed Brown 
52984df9cb4SJed Brown    Level: developer
53084df9cb4SJed Brown 
53184df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose()
53284df9cb4SJed Brown @*/
5338d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse)
53484df9cb4SJed Brown {
53584df9cb4SJed Brown   PetscInt c;
53684df9cb4SJed Brown 
53784df9cb4SJed Brown   PetscFunctionBegin;
53884df9cb4SJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
539ce94432eSBarry Smith   if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order);
54084df9cb4SJed Brown   if (inuse) {
541ce94432eSBarry Smith     if (adapt->candidates.inuse_set) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()");
54284df9cb4SJed Brown     adapt->candidates.inuse_set = PETSC_TRUE;
54384df9cb4SJed Brown   }
5441c3436cfSJed Brown   /* first slot if this is the current scheme, otherwise the next available slot */
5451c3436cfSJed Brown   c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;
546bbd56ea5SKarl Rupp 
54784df9cb4SJed Brown   adapt->candidates.name[c]       = name;
54884df9cb4SJed Brown   adapt->candidates.order[c]      = order;
54984df9cb4SJed Brown   adapt->candidates.stageorder[c] = stageorder;
5508d59e960SJed Brown   adapt->candidates.ccfl[c]       = ccfl;
55184df9cb4SJed Brown   adapt->candidates.cost[c]       = cost;
55284df9cb4SJed Brown   adapt->candidates.n++;
55384df9cb4SJed Brown   PetscFunctionReturn(0);
55484df9cb4SJed Brown }
55584df9cb4SJed Brown 
5568d59e960SJed Brown /*@C
5578d59e960SJed Brown    TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost
5588d59e960SJed Brown 
5598d59e960SJed Brown    Not Collective
5608d59e960SJed Brown 
5618d59e960SJed Brown    Input Arguments:
5628d59e960SJed Brown .  adapt - time step adaptivity context
5638d59e960SJed Brown 
5648d59e960SJed Brown    Output Arguments:
5658d59e960SJed Brown +  n - number of candidate schemes, always at least 1
5668d59e960SJed Brown .  order - the order of each candidate scheme
5678d59e960SJed Brown .  stageorder - the stage order of each candidate scheme
5688d59e960SJed Brown .  ccfl - the CFL coefficient of each scheme
5698d59e960SJed Brown -  cost - the relative cost of each scheme
5708d59e960SJed Brown 
5718d59e960SJed Brown    Level: developer
5728d59e960SJed Brown 
5738d59e960SJed Brown    Note:
5748d59e960SJed Brown    The current scheme is always returned in the first slot
5758d59e960SJed Brown 
5768d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose()
5778d59e960SJed Brown @*/
5788d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost)
5798d59e960SJed Brown {
5808d59e960SJed Brown   PetscFunctionBegin;
5818d59e960SJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
5828d59e960SJed Brown   if (n) *n = adapt->candidates.n;
5838d59e960SJed Brown   if (order) *order = adapt->candidates.order;
5848d59e960SJed Brown   if (stageorder) *stageorder = adapt->candidates.stageorder;
5858d59e960SJed Brown   if (ccfl) *ccfl = adapt->candidates.ccfl;
5868d59e960SJed Brown   if (cost) *cost = adapt->candidates.cost;
5878d59e960SJed Brown   PetscFunctionReturn(0);
5888d59e960SJed Brown }
5898d59e960SJed Brown 
59084df9cb4SJed Brown /*@C
59184df9cb4SJed Brown    TSAdaptChoose - choose which method and step size to use for the next step
59284df9cb4SJed Brown 
59384df9cb4SJed Brown    Logically Collective
59484df9cb4SJed Brown 
59584df9cb4SJed Brown    Input Arguments:
59684df9cb4SJed Brown +  adapt - adaptive contoller
59784df9cb4SJed Brown -  h - current step size
59884df9cb4SJed Brown 
59984df9cb4SJed Brown    Output Arguments:
6001566a47fSLisandro Dalcin +  next_sc - optional, scheme to use for the next step
60184df9cb4SJed Brown .  next_h - step size to use for the next step
60284df9cb4SJed Brown -  accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size
60384df9cb4SJed Brown 
6041c3436cfSJed Brown    Note:
6051c3436cfSJed Brown    The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is
6061c3436cfSJed Brown    being retried after an initial rejection.
6071c3436cfSJed Brown 
60884df9cb4SJed Brown    Level: developer
60984df9cb4SJed Brown 
61084df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd()
61184df9cb4SJed Brown @*/
61284df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept)
61384df9cb4SJed Brown {
61484df9cb4SJed Brown   PetscErrorCode ierr;
6151566a47fSLisandro Dalcin   PetscInt       ncandidates = adapt->candidates.n;
6161566a47fSLisandro Dalcin   PetscInt       scheme = 0;
6170b99f514SJed Brown   PetscReal      wlte = -1.0;
6185e4ed32fSEmil Constantinescu   PetscReal      wltea = -1.0;
6195e4ed32fSEmil Constantinescu   PetscReal      wlter = -1.0;
62084df9cb4SJed Brown 
62184df9cb4SJed Brown   PetscFunctionBegin;
62284df9cb4SJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
62384df9cb4SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
6241566a47fSLisandro Dalcin   if (next_sc) PetscValidIntPointer(next_sc,4);
62584df9cb4SJed Brown   PetscValidPointer(next_h,5);
62684df9cb4SJed Brown   PetscValidIntPointer(accept,6);
6271566a47fSLisandro Dalcin   if (next_sc) *next_sc = 0;
6281566a47fSLisandro Dalcin 
6291566a47fSLisandro Dalcin   /* Do not mess with adaptivity while handling events*/
6301566a47fSLisandro Dalcin   if (ts->event && ts->event->status != TSEVENT_NONE) {
6311566a47fSLisandro Dalcin     *next_h = h;
6321566a47fSLisandro Dalcin     *accept = PETSC_TRUE;
6331566a47fSLisandro Dalcin     PetscFunctionReturn(0);
6341566a47fSLisandro Dalcin   }
6351566a47fSLisandro Dalcin 
6365e4ed32fSEmil Constantinescu   ierr = (*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte,&wltea,&wlter);CHKERRQ(ierr);
6371566a47fSLisandro Dalcin   if (scheme < 0 || (ncandidates > 0 && scheme >= ncandidates)) SETERRQ2(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Chosen scheme %D not in valid range 0..%D",scheme,ncandidates-1);
6381566a47fSLisandro Dalcin   if (*next_h < 0) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h);
6391566a47fSLisandro Dalcin   if (next_sc) *next_sc = scheme;
6401566a47fSLisandro Dalcin 
64149354f04SShri Abhyankar   if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
64249354f04SShri Abhyankar     /* Reduce time step if it overshoots max time */
6431566a47fSLisandro Dalcin     if (ts->ptime + ts->time_step + *next_h >= ts->max_time) {
6441566a47fSLisandro Dalcin       PetscReal next_dt = ts->max_time - (ts->ptime + ts->time_step);
6458709b12bSShri Abhyankar       if (next_dt > PETSC_SMALL) *next_h = next_dt;
64649354f04SShri Abhyankar       else ts->reason = TS_CONVERGED_TIME;
64749354f04SShri Abhyankar     }
64849354f04SShri Abhyankar   }
6491c3436cfSJed Brown 
6501c3436cfSJed Brown   if (adapt->monitor) {
6511566a47fSLisandro Dalcin     const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : "";
6521c3436cfSJed Brown     ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
6530b99f514SJed Brown     if (wlte < 0) {
6540dea241dSBarry Smith       ierr = PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt %s %s %D:%s step %3D %s t=%-11g+%10.3e dt=%-10.3e\n",((PetscObject)adapt)->type_name,((PetscObject)ts)->type_name,scheme,sc_name,ts->steps,*accept ? "accepted" : "rejected",(double)ts->ptime,(double)h,(double)*next_h);CHKERRQ(ierr);
6550b99f514SJed Brown     } else {
6560dea241dSBarry Smith       ierr = PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt %s %s %D:%s step %3D %s t=%-11g+%10.3e dt=%-10.3e wlte=%5.3g  wltea=%5.3g wlter=%5.3g\n",((PetscObject)adapt)->type_name,((PetscObject)ts)->type_name,scheme,sc_name,ts->steps,*accept ? "accepted" : "rejected",(double)ts->ptime,(double)h,(double)*next_h,(double)wlte,(double)wltea,(double)wlter);CHKERRQ(ierr);
6570b99f514SJed Brown     }
6581c3436cfSJed Brown     ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
6591c3436cfSJed Brown   }
66084df9cb4SJed Brown   PetscFunctionReturn(0);
66184df9cb4SJed Brown }
66284df9cb4SJed Brown 
66397335746SJed Brown /*@
66497335746SJed Brown    TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails)
66597335746SJed Brown 
66697335746SJed Brown    Collective
66797335746SJed Brown 
66897335746SJed Brown    Input Arguments:
66997335746SJed Brown +  adapt - adaptive controller context
670b295832fSPierre Barbier de Reuille .  ts - time stepper
671b295832fSPierre Barbier de Reuille .  t - Current simulation time
672b295832fSPierre Barbier de Reuille -  Y - Current solution vector
67397335746SJed Brown 
67497335746SJed Brown    Output Arguments:
67597335746SJed Brown .  accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject
67697335746SJed Brown 
67797335746SJed Brown    Level: developer
67897335746SJed Brown 
67997335746SJed Brown .seealso:
68097335746SJed Brown @*/
681b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept)
68297335746SJed Brown {
68397335746SJed Brown   PetscErrorCode      ierr;
6841566a47fSLisandro Dalcin   SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;
68597335746SJed Brown 
68697335746SJed Brown   PetscFunctionBegin;
6874782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
6884782b174SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
6894782b174SLisandro Dalcin   PetscValidIntPointer(accept,3);
6901566a47fSLisandro Dalcin 
6911566a47fSLisandro Dalcin   if (ts->snes) {ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);}
69297335746SJed Brown   if (snesreason < 0) {
69397335746SJed Brown     *accept = PETSC_FALSE;
6946de24e2aSJed Brown     if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) {
69597335746SJed Brown       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
69648c19aefSLisandro Dalcin       ierr = PetscInfo2(ts,"Step=%D, nonlinear solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr);
69797335746SJed Brown       if (adapt->monitor) {
69897335746SJed Brown         ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
6990dea241dSBarry Smith         ierr = PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt %s step %3D stage rejected t=%-11g+%10.3e, nonlinear solve failures %D greater than current TS allowed\n",((PetscObject)adapt)->type_name,ts->steps,(double)ts->ptime,(double)ts->time_step,ts->num_snes_failures);CHKERRQ(ierr);
70097335746SJed Brown         ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
70197335746SJed Brown       }
702cb9d8021SPierre Barbier de Reuille     }
703cb9d8021SPierre Barbier de Reuille   } else {
7041566a47fSLisandro Dalcin     *accept = PETSC_TRUE;
705b295832fSPierre Barbier de Reuille     ierr = TSFunctionDomainError(ts,t,Y,accept);CHKERRQ(ierr);
706cb9d8021SPierre Barbier de Reuille     if(*accept && adapt->checkstage) {
707b295832fSPierre Barbier de Reuille       ierr = (*adapt->checkstage)(adapt,ts,t,Y,accept);CHKERRQ(ierr);
708cb9d8021SPierre Barbier de Reuille     }
709cb9d8021SPierre Barbier de Reuille   }
710cb9d8021SPierre Barbier de Reuille 
7111566a47fSLisandro Dalcin   if(!(*accept) && !ts->reason) {
7121566a47fSLisandro Dalcin     PetscReal dt,new_dt;
7131566a47fSLisandro Dalcin     ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
714cb9d8021SPierre Barbier de Reuille     new_dt = dt * adapt->scale_solve_failed;
71597335746SJed Brown     ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr);
71697335746SJed Brown     if (adapt->monitor) {
71797335746SJed Brown       ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
7180dea241dSBarry Smith       ierr = PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt %s step %3D stage rejected t=%-11g+%10.3e retrying with dt=%-10.3e\n",((PetscObject)adapt)->type_name,ts->steps,(double)ts->ptime,(double)dt,(double)new_dt);CHKERRQ(ierr);
71997335746SJed Brown       ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
72097335746SJed Brown     }
72197335746SJed Brown   }
72297335746SJed Brown   PetscFunctionReturn(0);
72397335746SJed Brown }
72497335746SJed Brown 
72584df9cb4SJed Brown /*@
72684df9cb4SJed Brown   TSAdaptCreate - create an adaptive controller context for time stepping
72784df9cb4SJed Brown 
72884df9cb4SJed Brown   Collective on MPI_Comm
72984df9cb4SJed Brown 
73084df9cb4SJed Brown   Input Parameter:
73184df9cb4SJed Brown . comm - The communicator
73284df9cb4SJed Brown 
73384df9cb4SJed Brown   Output Parameter:
73484df9cb4SJed Brown . adapt - new TSAdapt object
73584df9cb4SJed Brown 
73684df9cb4SJed Brown   Level: developer
73784df9cb4SJed Brown 
73884df9cb4SJed Brown   Notes:
73984df9cb4SJed Brown   TSAdapt creation is handled by TS, so users should not need to call this function.
74084df9cb4SJed Brown 
74184df9cb4SJed Brown .keywords: TSAdapt, create
742552698daSJed Brown .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy()
74384df9cb4SJed Brown @*/
74484df9cb4SJed Brown PetscErrorCode  TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt)
74584df9cb4SJed Brown {
74684df9cb4SJed Brown   PetscErrorCode ierr;
74784df9cb4SJed Brown   TSAdapt        adapt;
74884df9cb4SJed Brown 
74984df9cb4SJed Brown   PetscFunctionBegin;
7503b3bcf4cSLisandro Dalcin   PetscValidPointer(inadapt,1);
7513b3bcf4cSLisandro Dalcin   *inadapt = NULL;
7523b3bcf4cSLisandro Dalcin   ierr = TSAdaptInitializePackage();CHKERRQ(ierr);
7533b3bcf4cSLisandro Dalcin 
75473107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr);
7551c3436cfSJed Brown 
756bf997491SLisandro Dalcin   adapt->always_accept      = PETSC_FALSE;
7571c3436cfSJed Brown   adapt->dt_min             = 1e-20;
758bc0b8257SSatish Balay #if defined(PETSC_USE_REAL_SINGLE)
759bc0b8257SSatish Balay   adapt->dt_max             = 1e35;
760bc0b8257SSatish Balay #else
7611c3436cfSJed Brown   adapt->dt_max             = 1e50;
762bc0b8257SSatish Balay #endif
76397335746SJed Brown   adapt->scale_solve_failed = 0.25;
7647619abb3SShri   adapt->wnormtype          = NORM_2;
7657eef6055SBarry Smith   ierr = TSAdaptSetType(adapt,TSADAPTBASIC);CHKERRQ(ierr);
7661c3436cfSJed Brown 
76784df9cb4SJed Brown   *inadapt = adapt;
76884df9cb4SJed Brown   PetscFunctionReturn(0);
76984df9cb4SJed Brown }
770