xref: /petsc/src/ts/adapt/interface/tsadapt.c (revision 560360af1e0197128c3ad6271bbfa2e76ad345d4)
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 
108cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_None(TSAdapt);
111566a47fSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt);
128cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt);
1384df9cb4SJed Brown 
1484df9cb4SJed Brown #undef __FUNCT__
1584df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegister"
1684df9cb4SJed Brown /*@C
171c84c290SBarry Smith    TSAdaptRegister -  adds a TSAdapt implementation
181c84c290SBarry Smith 
191c84c290SBarry Smith    Not Collective
201c84c290SBarry Smith 
211c84c290SBarry Smith    Input Parameters:
221c84c290SBarry Smith +  name_scheme - name of user-defined adaptivity scheme
231c84c290SBarry Smith -  routine_create - routine to create method context
241c84c290SBarry Smith 
251c84c290SBarry Smith    Notes:
261c84c290SBarry Smith    TSAdaptRegister() may be called multiple times to add several user-defined families.
271c84c290SBarry Smith 
281c84c290SBarry Smith    Sample usage:
291c84c290SBarry Smith .vb
30bdf89e91SBarry Smith    TSAdaptRegister("my_scheme",MySchemeCreate);
311c84c290SBarry Smith .ve
321c84c290SBarry Smith 
331c84c290SBarry Smith    Then, your scheme can be chosen with the procedural interface via
341c84c290SBarry Smith $     TSAdaptSetType(ts,"my_scheme")
351c84c290SBarry Smith    or at runtime via the option
361c84c290SBarry Smith $     -ts_adapt_type my_scheme
3784df9cb4SJed Brown 
3884df9cb4SJed Brown    Level: advanced
391c84c290SBarry Smith 
401c84c290SBarry Smith .keywords: TSAdapt, register
411c84c290SBarry Smith 
421c84c290SBarry Smith .seealso: TSAdaptRegisterAll()
4384df9cb4SJed Brown @*/
44bdf89e91SBarry Smith PetscErrorCode  TSAdaptRegister(const char sname[],PetscErrorCode (*function)(TSAdapt))
4584df9cb4SJed Brown {
4684df9cb4SJed Brown   PetscErrorCode ierr;
4784df9cb4SJed Brown 
4884df9cb4SJed Brown   PetscFunctionBegin;
49a240a19fSJed Brown   ierr = PetscFunctionListAdd(&TSAdaptList,sname,function);CHKERRQ(ierr);
5084df9cb4SJed Brown   PetscFunctionReturn(0);
5184df9cb4SJed Brown }
5284df9cb4SJed Brown 
5384df9cb4SJed Brown #undef __FUNCT__
5484df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegisterAll"
5584df9cb4SJed Brown /*@C
5684df9cb4SJed Brown   TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt
5784df9cb4SJed Brown 
5884df9cb4SJed Brown   Not Collective
5984df9cb4SJed Brown 
6084df9cb4SJed Brown   Level: advanced
6184df9cb4SJed Brown 
6284df9cb4SJed Brown .keywords: TSAdapt, register, all
6384df9cb4SJed Brown 
6484df9cb4SJed Brown .seealso: TSAdaptRegisterDestroy()
6584df9cb4SJed Brown @*/
66607a6623SBarry Smith PetscErrorCode  TSAdaptRegisterAll(void)
6784df9cb4SJed Brown {
6884df9cb4SJed Brown   PetscErrorCode ierr;
6984df9cb4SJed Brown 
7084df9cb4SJed Brown   PetscFunctionBegin;
710f51fdf8SToby Isaac   if (TSAdaptRegisterAllCalled) PetscFunctionReturn(0);
720f51fdf8SToby Isaac   TSAdaptRegisterAllCalled = PETSC_TRUE;
73bdf89e91SBarry Smith   ierr = TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None);CHKERRQ(ierr);
741566a47fSLisandro Dalcin   ierr = TSAdaptRegister(TSADAPTBASIC,TSAdaptCreate_Basic);CHKERRQ(ierr);
75bdf89e91SBarry Smith   ierr = TSAdaptRegister(TSADAPTCFL,  TSAdaptCreate_CFL);CHKERRQ(ierr);
7684df9cb4SJed Brown   PetscFunctionReturn(0);
7784df9cb4SJed Brown }
7884df9cb4SJed Brown 
7984df9cb4SJed Brown #undef __FUNCT__
8084df9cb4SJed Brown #define __FUNCT__ "TSAdaptFinalizePackage"
8184df9cb4SJed Brown /*@C
82*560360afSLisandro Dalcin   TSAdaptFinalizePackage - This function destroys everything in the TS package. It is
8384df9cb4SJed Brown   called from PetscFinalize().
8484df9cb4SJed Brown 
8584df9cb4SJed Brown   Level: developer
8684df9cb4SJed Brown 
8784df9cb4SJed Brown .keywords: Petsc, destroy, package
8884df9cb4SJed Brown .seealso: PetscFinalize()
8984df9cb4SJed Brown @*/
9084df9cb4SJed Brown PetscErrorCode  TSAdaptFinalizePackage(void)
9184df9cb4SJed Brown {
9237e93019SBarry Smith   PetscErrorCode ierr;
9337e93019SBarry Smith 
9484df9cb4SJed Brown   PetscFunctionBegin;
9537e93019SBarry Smith   ierr = PetscFunctionListDestroy(&TSAdaptList);CHKERRQ(ierr);
9684df9cb4SJed Brown   TSAdaptPackageInitialized = PETSC_FALSE;
9784df9cb4SJed Brown   TSAdaptRegisterAllCalled  = PETSC_FALSE;
9884df9cb4SJed Brown   PetscFunctionReturn(0);
9984df9cb4SJed Brown }
10084df9cb4SJed Brown 
10184df9cb4SJed Brown #undef __FUNCT__
10284df9cb4SJed Brown #define __FUNCT__ "TSAdaptInitializePackage"
10384df9cb4SJed Brown /*@C
10484df9cb4SJed Brown   TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is
10584df9cb4SJed Brown   called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to
10684df9cb4SJed Brown   TSCreate_GL() when using static libraries.
10784df9cb4SJed Brown 
10884df9cb4SJed Brown   Level: developer
10984df9cb4SJed Brown 
11084df9cb4SJed Brown .keywords: TSAdapt, initialize, package
11184df9cb4SJed Brown .seealso: PetscInitialize()
11284df9cb4SJed Brown @*/
113607a6623SBarry Smith PetscErrorCode  TSAdaptInitializePackage(void)
11484df9cb4SJed Brown {
11584df9cb4SJed Brown   PetscErrorCode ierr;
11684df9cb4SJed Brown 
11784df9cb4SJed Brown   PetscFunctionBegin;
11884df9cb4SJed Brown   if (TSAdaptPackageInitialized) PetscFunctionReturn(0);
11984df9cb4SJed Brown   TSAdaptPackageInitialized = PETSC_TRUE;
12084df9cb4SJed Brown   ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr);
121607a6623SBarry Smith   ierr = TSAdaptRegisterAll();CHKERRQ(ierr);
12284df9cb4SJed Brown   ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr);
12384df9cb4SJed Brown   PetscFunctionReturn(0);
12484df9cb4SJed Brown }
12584df9cb4SJed Brown 
12684df9cb4SJed Brown #undef __FUNCT__
12784df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetType"
12819fd82e9SBarry Smith PetscErrorCode  TSAdaptSetType(TSAdapt adapt,TSAdaptType type)
12984df9cb4SJed Brown {
130ef749922SLisandro Dalcin   PetscBool      match;
13184df9cb4SJed Brown   PetscErrorCode ierr,(*r)(TSAdapt);
13284df9cb4SJed Brown 
13384df9cb4SJed Brown   PetscFunctionBegin;
1344782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
135ef749922SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)adapt,type,&match);CHKERRQ(ierr);
136ef749922SLisandro Dalcin   if (match) PetscFunctionReturn(0);
1371c9cd337SJed Brown   ierr = PetscFunctionListFind(TSAdaptList,type,&r);CHKERRQ(ierr);
13884df9cb4SJed Brown   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type);
139ef749922SLisandro Dalcin   if (adapt->ops->destroy) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);}
1404cefc2ffSBarry Smith   ierr = PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));CHKERRQ(ierr);
14184df9cb4SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr);
14268ae941aSLisandro Dalcin   ierr = (*r)(adapt);CHKERRQ(ierr);
14384df9cb4SJed Brown   PetscFunctionReturn(0);
14484df9cb4SJed Brown }
14584df9cb4SJed Brown 
14684df9cb4SJed Brown #undef __FUNCT__
14784df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetOptionsPrefix"
14884df9cb4SJed Brown PetscErrorCode  TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[])
14984df9cb4SJed Brown {
15084df9cb4SJed Brown   PetscErrorCode ierr;
15184df9cb4SJed Brown 
15284df9cb4SJed Brown   PetscFunctionBegin;
1534782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
15484df9cb4SJed Brown   ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr);
15584df9cb4SJed Brown   PetscFunctionReturn(0);
15684df9cb4SJed Brown }
15784df9cb4SJed Brown 
15884df9cb4SJed Brown #undef __FUNCT__
159ad6bc421SBarry Smith #define __FUNCT__ "TSAdaptLoad"
160ad6bc421SBarry Smith /*@C
161ad6bc421SBarry Smith   TSAdaptLoad - Loads a TSAdapt that has been stored in binary  with TSAdaptView().
162ad6bc421SBarry Smith 
163ad6bc421SBarry Smith   Collective on PetscViewer
164ad6bc421SBarry Smith 
165ad6bc421SBarry Smith   Input Parameters:
166ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or
167ad6bc421SBarry Smith            some related function before a call to TSAdaptLoad().
168ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
169ad6bc421SBarry Smith            HDF5 file viewer, obtained from PetscViewerHDF5Open()
170ad6bc421SBarry Smith 
171ad6bc421SBarry Smith    Level: intermediate
172ad6bc421SBarry Smith 
173ad6bc421SBarry Smith   Notes:
174ad6bc421SBarry Smith    The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored.
175ad6bc421SBarry Smith 
176ad6bc421SBarry Smith   Notes for advanced users:
177ad6bc421SBarry Smith   Most users should not need to know the details of the binary storage
178ad6bc421SBarry Smith   format, since TSAdaptLoad() and TSAdaptView() completely hide these details.
179ad6bc421SBarry Smith   But for anyone who's interested, the standard binary matrix storage
180ad6bc421SBarry Smith   format is
181ad6bc421SBarry Smith .vb
182ad6bc421SBarry Smith      has not yet been determined
183ad6bc421SBarry Smith .ve
184ad6bc421SBarry Smith 
185ad6bc421SBarry Smith .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad()
186ad6bc421SBarry Smith @*/
1874782b174SLisandro Dalcin PetscErrorCode  TSAdaptLoad(TSAdapt adapt,PetscViewer viewer)
188ad6bc421SBarry Smith {
189ad6bc421SBarry Smith   PetscErrorCode ierr;
190ad6bc421SBarry Smith   PetscBool      isbinary;
191ad6bc421SBarry Smith   char           type[256];
192ad6bc421SBarry Smith 
193ad6bc421SBarry Smith   PetscFunctionBegin;
1944782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
195ad6bc421SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
196ad6bc421SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
197ad6bc421SBarry Smith   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
198ad6bc421SBarry Smith 
199060da220SMatthew G. Knepley   ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr);
2004782b174SLisandro Dalcin   ierr = TSAdaptSetType(adapt, type);CHKERRQ(ierr);
2014782b174SLisandro Dalcin   if (adapt->ops->load) {
2024782b174SLisandro Dalcin     ierr = (*adapt->ops->load)(adapt,viewer);CHKERRQ(ierr);
203ad6bc421SBarry Smith   }
204ad6bc421SBarry Smith   PetscFunctionReturn(0);
205ad6bc421SBarry Smith }
206ad6bc421SBarry Smith 
207ad6bc421SBarry Smith #undef __FUNCT__
20884df9cb4SJed Brown #define __FUNCT__ "TSAdaptView"
20984df9cb4SJed Brown PetscErrorCode  TSAdaptView(TSAdapt adapt,PetscViewer viewer)
21084df9cb4SJed Brown {
21184df9cb4SJed Brown   PetscErrorCode ierr;
212ad6bc421SBarry Smith   PetscBool      iascii,isbinary;
21384df9cb4SJed Brown 
21484df9cb4SJed Brown   PetscFunctionBegin;
2154782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
2164782b174SLisandro Dalcin   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer);CHKERRQ(ierr);}
2174782b174SLisandro Dalcin   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
2184782b174SLisandro Dalcin   PetscCheckSameComm(adapt,1,viewer,2);
219251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
220ad6bc421SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
22184df9cb4SJed Brown   if (iascii) {
222dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr);
22384df9cb4SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  number of candidates %D\n",adapt->candidates.n);CHKERRQ(ierr);
22484df9cb4SJed Brown     if (adapt->ops->view) {
22584df9cb4SJed Brown       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
22684df9cb4SJed Brown       ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
22784df9cb4SJed Brown       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
22884df9cb4SJed Brown     }
229ad6bc421SBarry Smith   } else if (isbinary) {
230ad6bc421SBarry Smith     char type[256];
231ad6bc421SBarry Smith 
232ad6bc421SBarry Smith     /* need to save FILE_CLASS_ID for adapt class */
233ad6bc421SBarry Smith     ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr);
234ad6bc421SBarry Smith     ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr);
235bbd56ea5SKarl Rupp   } else if (adapt->ops->view) {
236f2c2a1b9SBarry Smith     ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
237f2c2a1b9SBarry Smith   }
23884df9cb4SJed Brown   PetscFunctionReturn(0);
23984df9cb4SJed Brown }
24084df9cb4SJed Brown 
24184df9cb4SJed Brown #undef __FUNCT__
242099af0a3SLisandro Dalcin #define __FUNCT__ "TSAdaptReset"
243099af0a3SLisandro Dalcin /*@
244099af0a3SLisandro Dalcin    TSAdaptReset - Resets a TSAdapt context.
245099af0a3SLisandro Dalcin 
246099af0a3SLisandro Dalcin    Collective on TS
247099af0a3SLisandro Dalcin 
248099af0a3SLisandro Dalcin    Input Parameter:
249099af0a3SLisandro Dalcin .  adapt - the TSAdapt context obtained from TSAdaptCreate()
250099af0a3SLisandro Dalcin 
251099af0a3SLisandro Dalcin    Level: developer
252099af0a3SLisandro Dalcin 
253099af0a3SLisandro Dalcin .seealso: TSAdaptCreate(), TSAdaptDestroy()
254099af0a3SLisandro Dalcin @*/
255099af0a3SLisandro Dalcin PetscErrorCode  TSAdaptReset(TSAdapt adapt)
256099af0a3SLisandro Dalcin {
257099af0a3SLisandro Dalcin   PetscErrorCode ierr;
258099af0a3SLisandro Dalcin 
259099af0a3SLisandro Dalcin   PetscFunctionBegin;
260099af0a3SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
261099af0a3SLisandro Dalcin   if (adapt->ops->reset) {ierr = (*adapt->ops->reset)(adapt);CHKERRQ(ierr);}
262099af0a3SLisandro Dalcin   PetscFunctionReturn(0);
263099af0a3SLisandro Dalcin }
264099af0a3SLisandro Dalcin 
265099af0a3SLisandro Dalcin #undef __FUNCT__
26684df9cb4SJed Brown #define __FUNCT__ "TSAdaptDestroy"
26784df9cb4SJed Brown PetscErrorCode  TSAdaptDestroy(TSAdapt *adapt)
26884df9cb4SJed Brown {
26984df9cb4SJed Brown   PetscErrorCode ierr;
27084df9cb4SJed Brown 
27184df9cb4SJed Brown   PetscFunctionBegin;
27284df9cb4SJed Brown   if (!*adapt) PetscFunctionReturn(0);
27384df9cb4SJed Brown   PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1);
2744782b174SLisandro Dalcin   if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);}
275099af0a3SLisandro Dalcin 
276099af0a3SLisandro Dalcin   ierr = TSAdaptReset(*adapt);CHKERRQ(ierr);
277099af0a3SLisandro Dalcin 
27884df9cb4SJed Brown   if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);}
2791c3436cfSJed Brown   ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr);
28084df9cb4SJed Brown   ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr);
28184df9cb4SJed Brown   PetscFunctionReturn(0);
28284df9cb4SJed Brown }
28384df9cb4SJed Brown 
28484df9cb4SJed Brown #undef __FUNCT__
2851c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetMonitor"
2861c3436cfSJed Brown /*@
2871c3436cfSJed Brown    TSAdaptSetMonitor - Monitor the choices made by the adaptive controller
2881c3436cfSJed Brown 
2891c3436cfSJed Brown    Collective on TSAdapt
2901c3436cfSJed Brown 
2911c3436cfSJed Brown    Input Arguments:
2921c3436cfSJed Brown +  adapt - adaptive controller context
2931c3436cfSJed Brown -  flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable
2941c3436cfSJed Brown 
2951c3436cfSJed Brown    Level: intermediate
2961c3436cfSJed Brown 
2971c3436cfSJed Brown .seealso: TSAdaptChoose()
2981c3436cfSJed Brown @*/
2991c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg)
3001c3436cfSJed Brown {
3011c3436cfSJed Brown   PetscErrorCode ierr;
3021c3436cfSJed Brown 
3031c3436cfSJed Brown   PetscFunctionBegin;
3044782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
3054782b174SLisandro Dalcin   PetscValidLogicalCollectiveBool(adapt,flg,2);
3061c3436cfSJed Brown   if (flg) {
307ce94432eSBarry Smith     if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);}
3081c3436cfSJed Brown   } else {
3091c3436cfSJed Brown     ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr);
3101c3436cfSJed Brown   }
3111c3436cfSJed Brown   PetscFunctionReturn(0);
3121c3436cfSJed Brown }
3131c3436cfSJed Brown 
3141c3436cfSJed Brown #undef __FUNCT__
3150873808bSJed Brown #define __FUNCT__ "TSAdaptSetCheckStage"
3160873808bSJed Brown /*@C
3170873808bSJed Brown    TSAdaptSetCheckStage - set a callback to check convergence for a stage
3180873808bSJed Brown 
3190873808bSJed Brown    Logically collective on TSAdapt
3200873808bSJed Brown 
3210873808bSJed Brown    Input Arguments:
3220873808bSJed Brown +  adapt - adaptive controller context
3230873808bSJed Brown -  func - stage check function
3240873808bSJed Brown 
3250873808bSJed Brown    Arguments of func:
3260873808bSJed Brown $  PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept)
3270873808bSJed Brown 
3280873808bSJed Brown +  adapt - adaptive controller context
3290873808bSJed Brown .  ts - time stepping context
3300873808bSJed Brown -  accept - pending choice of whether to accept, can be modified by this routine
3310873808bSJed Brown 
3320873808bSJed Brown    Level: advanced
3330873808bSJed Brown 
3340873808bSJed Brown .seealso: TSAdaptChoose()
3350873808bSJed Brown @*/
336b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscReal,Vec,PetscBool*))
3370873808bSJed Brown {
3380873808bSJed Brown 
3390873808bSJed Brown   PetscFunctionBegin;
3400873808bSJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
34168ae941aSLisandro Dalcin   adapt->checkstage = func;
3420873808bSJed Brown   PetscFunctionReturn(0);
3430873808bSJed Brown }
3440873808bSJed Brown 
3450873808bSJed Brown #undef __FUNCT__
3461c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetStepLimits"
3471c3436cfSJed Brown /*@
3481c3436cfSJed Brown    TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller
3491c3436cfSJed Brown 
3501c3436cfSJed Brown    Logically Collective
3511c3436cfSJed Brown 
3521c3436cfSJed Brown    Input Arguments:
353552698daSJed Brown +  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
3541c3436cfSJed Brown .  hmin - minimum time step
3551c3436cfSJed Brown -  hmax - maximum time step
3561c3436cfSJed Brown 
3571c3436cfSJed Brown    Options Database Keys:
3581c3436cfSJed Brown +  -ts_adapt_dt_min - minimum time step
3591c3436cfSJed Brown -  -ts_adapt_dt_max - maximum time step
3601c3436cfSJed Brown 
3611c3436cfSJed Brown    Level: intermediate
3621c3436cfSJed Brown 
3631c3436cfSJed Brown .seealso: TSAdapt
3641c3436cfSJed Brown @*/
3651c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax)
3661c3436cfSJed Brown {
3671c3436cfSJed Brown 
3681c3436cfSJed Brown   PetscFunctionBegin;
3694782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
3701c3436cfSJed Brown   if (hmin != PETSC_DECIDE) adapt->dt_min = hmin;
3711c3436cfSJed Brown   if (hmax != PETSC_DECIDE) adapt->dt_max = hmax;
3721c3436cfSJed Brown   PetscFunctionReturn(0);
3731c3436cfSJed Brown }
3741c3436cfSJed Brown 
3751c3436cfSJed Brown #undef __FUNCT__
37684df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetFromOptions"
377e55864a3SBarry Smith /*
37884df9cb4SJed Brown    TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options.
37984df9cb4SJed Brown 
38084df9cb4SJed Brown    Collective on TSAdapt
38184df9cb4SJed Brown 
38284df9cb4SJed Brown    Input Parameter:
38384df9cb4SJed Brown .  adapt - the TSAdapt context
38484df9cb4SJed Brown 
38584df9cb4SJed Brown    Options Database Keys:
38684df9cb4SJed Brown .  -ts_adapt_type <type> - basic
38784df9cb4SJed Brown 
38884df9cb4SJed Brown    Level: advanced
38984df9cb4SJed Brown 
39084df9cb4SJed Brown    Notes:
39184df9cb4SJed Brown    This function is automatically called by TSSetFromOptions()
39284df9cb4SJed Brown 
393552698daSJed Brown .keywords: TS, TSGetAdapt(), TSAdaptSetType()
39484df9cb4SJed Brown 
39584df9cb4SJed Brown .seealso: TSGetType()
396e55864a3SBarry Smith */
3974416b707SBarry Smith PetscErrorCode  TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt)
39884df9cb4SJed Brown {
39984df9cb4SJed Brown   PetscErrorCode ierr;
40084df9cb4SJed Brown   char           type[256] = TSADAPTBASIC;
4011c3436cfSJed Brown   PetscBool      set,flg;
40284df9cb4SJed Brown 
40384df9cb4SJed Brown   PetscFunctionBegin;
4044782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
40584df9cb4SJed Brown   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
4061566a47fSLisandro Dalcin    * function can only be called from inside TSSetFromOptions()  */
407e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");CHKERRQ(ierr);
408a264d7a6SBarry Smith   ierr = PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList,
4098caf3d72SBarry Smith                           ((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);CHKERRQ(ierr);
41084df9cb4SJed Brown   if (flg || !((PetscObject)adapt)->type_name) {
41184df9cb4SJed Brown     ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr);
41284df9cb4SJed Brown   }
4130298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",adapt->dt_min,&adapt->dt_min,NULL);CHKERRQ(ierr);
4140298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",adapt->dt_max,&adapt->dt_max,NULL);CHKERRQ(ierr);
4150298fd71SBarry 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);
4161c3436cfSJed Brown   ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
4177619abb3SShri   ierr = PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);CHKERRQ(ierr);
4187619abb3SShri   if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported");
4191c3436cfSJed Brown   if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);}
420e55864a3SBarry Smith   if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);}
42184df9cb4SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
42284df9cb4SJed Brown   PetscFunctionReturn(0);
42384df9cb4SJed Brown }
42484df9cb4SJed Brown 
42584df9cb4SJed Brown #undef __FUNCT__
42684df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidatesClear"
42784df9cb4SJed Brown /*@
42884df9cb4SJed Brown    TSAdaptCandidatesClear - clear any previously set candidate schemes
42984df9cb4SJed Brown 
43084df9cb4SJed Brown    Logically Collective
43184df9cb4SJed Brown 
43284df9cb4SJed Brown    Input Argument:
43384df9cb4SJed Brown .  adapt - adaptive controller
43484df9cb4SJed Brown 
43584df9cb4SJed Brown    Level: developer
43684df9cb4SJed Brown 
43784df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose()
43884df9cb4SJed Brown @*/
43984df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
44084df9cb4SJed Brown {
44184df9cb4SJed Brown   PetscErrorCode ierr;
44284df9cb4SJed Brown 
44384df9cb4SJed Brown   PetscFunctionBegin;
4444782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
44584df9cb4SJed Brown   ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr);
44684df9cb4SJed Brown   PetscFunctionReturn(0);
44784df9cb4SJed Brown }
44884df9cb4SJed Brown 
44984df9cb4SJed Brown #undef __FUNCT__
45084df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidateAdd"
45184df9cb4SJed Brown /*@C
45284df9cb4SJed Brown    TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from
45384df9cb4SJed Brown 
45484df9cb4SJed Brown    Logically Collective
45584df9cb4SJed Brown 
45684df9cb4SJed Brown    Input Arguments:
457552698daSJed Brown +  adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate()
45884df9cb4SJed Brown .  name - name of the candidate scheme to add
45984df9cb4SJed Brown .  order - order of the candidate scheme
46084df9cb4SJed Brown .  stageorder - stage order of the candidate scheme
4618d59e960SJed Brown .  ccfl - stability coefficient relative to explicit Euler, used for CFL constraints
46284df9cb4SJed Brown .  cost - relative measure of the amount of work required for the candidate scheme
46384df9cb4SJed Brown -  inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme
46484df9cb4SJed Brown 
46584df9cb4SJed Brown    Note:
46684df9cb4SJed Brown    This routine is not available in Fortran.
46784df9cb4SJed Brown 
46884df9cb4SJed Brown    Level: developer
46984df9cb4SJed Brown 
47084df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose()
47184df9cb4SJed Brown @*/
4728d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse)
47384df9cb4SJed Brown {
47484df9cb4SJed Brown   PetscInt c;
47584df9cb4SJed Brown 
47684df9cb4SJed Brown   PetscFunctionBegin;
47784df9cb4SJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
478ce94432eSBarry Smith   if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order);
47984df9cb4SJed Brown   if (inuse) {
480ce94432eSBarry 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()");
48184df9cb4SJed Brown     adapt->candidates.inuse_set = PETSC_TRUE;
48284df9cb4SJed Brown   }
4831c3436cfSJed Brown   /* first slot if this is the current scheme, otherwise the next available slot */
4841c3436cfSJed Brown   c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;
485bbd56ea5SKarl Rupp 
48684df9cb4SJed Brown   adapt->candidates.name[c]       = name;
48784df9cb4SJed Brown   adapt->candidates.order[c]      = order;
48884df9cb4SJed Brown   adapt->candidates.stageorder[c] = stageorder;
4898d59e960SJed Brown   adapt->candidates.ccfl[c]       = ccfl;
49084df9cb4SJed Brown   adapt->candidates.cost[c]       = cost;
49184df9cb4SJed Brown   adapt->candidates.n++;
49284df9cb4SJed Brown   PetscFunctionReturn(0);
49384df9cb4SJed Brown }
49484df9cb4SJed Brown 
49584df9cb4SJed Brown #undef __FUNCT__
4968d59e960SJed Brown #define __FUNCT__ "TSAdaptCandidatesGet"
4978d59e960SJed Brown /*@C
4988d59e960SJed Brown    TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost
4998d59e960SJed Brown 
5008d59e960SJed Brown    Not Collective
5018d59e960SJed Brown 
5028d59e960SJed Brown    Input Arguments:
5038d59e960SJed Brown .  adapt - time step adaptivity context
5048d59e960SJed Brown 
5058d59e960SJed Brown    Output Arguments:
5068d59e960SJed Brown +  n - number of candidate schemes, always at least 1
5078d59e960SJed Brown .  order - the order of each candidate scheme
5088d59e960SJed Brown .  stageorder - the stage order of each candidate scheme
5098d59e960SJed Brown .  ccfl - the CFL coefficient of each scheme
5108d59e960SJed Brown -  cost - the relative cost of each scheme
5118d59e960SJed Brown 
5128d59e960SJed Brown    Level: developer
5138d59e960SJed Brown 
5148d59e960SJed Brown    Note:
5158d59e960SJed Brown    The current scheme is always returned in the first slot
5168d59e960SJed Brown 
5178d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose()
5188d59e960SJed Brown @*/
5198d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost)
5208d59e960SJed Brown {
5218d59e960SJed Brown   PetscFunctionBegin;
5228d59e960SJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
5238d59e960SJed Brown   if (n) *n = adapt->candidates.n;
5248d59e960SJed Brown   if (order) *order = adapt->candidates.order;
5258d59e960SJed Brown   if (stageorder) *stageorder = adapt->candidates.stageorder;
5268d59e960SJed Brown   if (ccfl) *ccfl = adapt->candidates.ccfl;
5278d59e960SJed Brown   if (cost) *cost = adapt->candidates.cost;
5288d59e960SJed Brown   PetscFunctionReturn(0);
5298d59e960SJed Brown }
5308d59e960SJed Brown 
5318d59e960SJed Brown #undef __FUNCT__
53284df9cb4SJed Brown #define __FUNCT__ "TSAdaptChoose"
53384df9cb4SJed Brown /*@C
53484df9cb4SJed Brown    TSAdaptChoose - choose which method and step size to use for the next step
53584df9cb4SJed Brown 
53684df9cb4SJed Brown    Logically Collective
53784df9cb4SJed Brown 
53884df9cb4SJed Brown    Input Arguments:
53984df9cb4SJed Brown +  adapt - adaptive contoller
54084df9cb4SJed Brown -  h - current step size
54184df9cb4SJed Brown 
54284df9cb4SJed Brown    Output Arguments:
5431566a47fSLisandro Dalcin +  next_sc - optional, scheme to use for the next step
54484df9cb4SJed Brown .  next_h - step size to use for the next step
54584df9cb4SJed Brown -  accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size
54684df9cb4SJed Brown 
5471c3436cfSJed Brown    Note:
5481c3436cfSJed Brown    The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is
5491c3436cfSJed Brown    being retried after an initial rejection.
5501c3436cfSJed Brown 
55184df9cb4SJed Brown    Level: developer
55284df9cb4SJed Brown 
55384df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd()
55484df9cb4SJed Brown @*/
55584df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept)
55684df9cb4SJed Brown {
55784df9cb4SJed Brown   PetscErrorCode ierr;
5581566a47fSLisandro Dalcin   PetscInt       ncandidates = adapt->candidates.n;
5591566a47fSLisandro Dalcin   PetscInt       scheme = 0;
5600b99f514SJed Brown   PetscReal      wlte = -1.0;
56184df9cb4SJed Brown 
56284df9cb4SJed Brown   PetscFunctionBegin;
56384df9cb4SJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
56484df9cb4SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
5651566a47fSLisandro Dalcin   if (next_sc) PetscValidIntPointer(next_sc,4);
56684df9cb4SJed Brown   PetscValidPointer(next_h,5);
56784df9cb4SJed Brown   PetscValidIntPointer(accept,6);
5681566a47fSLisandro Dalcin   if (next_sc) *next_sc = 0;
5691566a47fSLisandro Dalcin 
5701566a47fSLisandro Dalcin   /* Do not mess with adaptivity while handling events*/
5711566a47fSLisandro Dalcin   if (ts->event && ts->event->status != TSEVENT_NONE) {
5721566a47fSLisandro Dalcin     *next_h = h;
5731566a47fSLisandro Dalcin     *accept = PETSC_TRUE;
5741566a47fSLisandro Dalcin     PetscFunctionReturn(0);
5751566a47fSLisandro Dalcin   }
5761566a47fSLisandro Dalcin 
5771566a47fSLisandro Dalcin   ierr = (*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte);CHKERRQ(ierr);
5781566a47fSLisandro 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);
5791566a47fSLisandro Dalcin   if (*next_h < 0) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h);
5801566a47fSLisandro Dalcin   if (next_sc) *next_sc = scheme;
5811566a47fSLisandro Dalcin 
58249354f04SShri Abhyankar   if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
58349354f04SShri Abhyankar     /* Reduce time step if it overshoots max time */
5841566a47fSLisandro Dalcin     if (ts->ptime + ts->time_step + *next_h >= ts->max_time) {
5851566a47fSLisandro Dalcin       PetscReal next_dt = ts->max_time - (ts->ptime + ts->time_step);
5868709b12bSShri Abhyankar       if (next_dt > PETSC_SMALL) *next_h = next_dt;
58749354f04SShri Abhyankar       else ts->reason = TS_CONVERGED_TIME;
58849354f04SShri Abhyankar     }
58949354f04SShri Abhyankar   }
5901c3436cfSJed Brown 
5911c3436cfSJed Brown   if (adapt->monitor) {
5921566a47fSLisandro Dalcin     const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : "";
5931c3436cfSJed Brown     ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
5940b99f514SJed Brown     if (wlte < 0) {
5951566a47fSLisandro Dalcin       ierr = PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt '%s': step %3D %s t=%-11g+%10.3e family='%s' scheme=%D:'%s' dt=%-10.3e\n",((PetscObject)adapt)->type_name,ts->steps,*accept ? "accepted" : "rejected",(double)ts->ptime,(double)h,((PetscObject)ts)->type_name,scheme,sc_name,(double)*next_h);CHKERRQ(ierr);
5960b99f514SJed Brown     } else {
5971566a47fSLisandro Dalcin       ierr = PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt '%s': step %3D %s t=%-11g+%10.3e wlte=%5.3g family='%s' scheme=%D:'%s' dt=%-10.3e\n",((PetscObject)adapt)->type_name,ts->steps,*accept ? "accepted" : "rejected",(double)ts->ptime,(double)h,(double)wlte,((PetscObject)ts)->type_name,scheme,sc_name,(double)*next_h);CHKERRQ(ierr);
5980b99f514SJed Brown     }
5991c3436cfSJed Brown     ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
6001c3436cfSJed Brown   }
60184df9cb4SJed Brown   PetscFunctionReturn(0);
60284df9cb4SJed Brown }
60384df9cb4SJed Brown 
60484df9cb4SJed Brown #undef __FUNCT__
60597335746SJed Brown #define __FUNCT__ "TSAdaptCheckStage"
60697335746SJed Brown /*@
60797335746SJed Brown    TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails)
60897335746SJed Brown 
60997335746SJed Brown    Collective
61097335746SJed Brown 
61197335746SJed Brown    Input Arguments:
61297335746SJed Brown +  adapt - adaptive controller context
613b295832fSPierre Barbier de Reuille .  ts - time stepper
614b295832fSPierre Barbier de Reuille .  t - Current simulation time
615b295832fSPierre Barbier de Reuille -  Y - Current solution vector
61697335746SJed Brown 
61797335746SJed Brown    Output Arguments:
61897335746SJed Brown .  accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject
61997335746SJed Brown 
62097335746SJed Brown    Level: developer
62197335746SJed Brown 
62297335746SJed Brown .seealso:
62397335746SJed Brown @*/
624b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept)
62597335746SJed Brown {
62697335746SJed Brown   PetscErrorCode      ierr;
6271566a47fSLisandro Dalcin   SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;
62897335746SJed Brown 
62997335746SJed Brown   PetscFunctionBegin;
6304782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
6314782b174SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
6324782b174SLisandro Dalcin   PetscValidIntPointer(accept,3);
6331566a47fSLisandro Dalcin 
6341566a47fSLisandro Dalcin   if (ts->snes) {ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);}
63597335746SJed Brown   if (snesreason < 0) {
63697335746SJed Brown     *accept = PETSC_FALSE;
6376de24e2aSJed Brown     if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) {
63897335746SJed Brown       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
63948c19aefSLisandro 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);
64097335746SJed Brown       if (adapt->monitor) {
64197335746SJed Brown         ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
6421566a47fSLisandro Dalcin         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);
64397335746SJed Brown         ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
64497335746SJed Brown       }
645cb9d8021SPierre Barbier de Reuille     }
646cb9d8021SPierre Barbier de Reuille   } else {
6471566a47fSLisandro Dalcin     *accept = PETSC_TRUE;
648b295832fSPierre Barbier de Reuille     ierr = TSFunctionDomainError(ts,t,Y,accept);CHKERRQ(ierr);
649cb9d8021SPierre Barbier de Reuille     if(*accept && adapt->checkstage) {
650b295832fSPierre Barbier de Reuille       ierr = (*adapt->checkstage)(adapt,ts,t,Y,accept);CHKERRQ(ierr);
651cb9d8021SPierre Barbier de Reuille     }
652cb9d8021SPierre Barbier de Reuille   }
653cb9d8021SPierre Barbier de Reuille 
6541566a47fSLisandro Dalcin   if(!(*accept) && !ts->reason) {
6551566a47fSLisandro Dalcin     PetscReal dt,new_dt;
6561566a47fSLisandro Dalcin     ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
657cb9d8021SPierre Barbier de Reuille     new_dt = dt * adapt->scale_solve_failed;
65897335746SJed Brown     ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr);
65997335746SJed Brown     if (adapt->monitor) {
66097335746SJed Brown       ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
66197335746SJed Brown       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);
66297335746SJed Brown       ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
66397335746SJed Brown     }
66497335746SJed Brown   }
66597335746SJed Brown   PetscFunctionReturn(0);
66697335746SJed Brown }
66797335746SJed Brown 
66897335746SJed Brown #undef __FUNCT__
66984df9cb4SJed Brown #define __FUNCT__ "TSAdaptCreate"
67084df9cb4SJed Brown /*@
67184df9cb4SJed Brown   TSAdaptCreate - create an adaptive controller context for time stepping
67284df9cb4SJed Brown 
67384df9cb4SJed Brown   Collective on MPI_Comm
67484df9cb4SJed Brown 
67584df9cb4SJed Brown   Input Parameter:
67684df9cb4SJed Brown . comm - The communicator
67784df9cb4SJed Brown 
67884df9cb4SJed Brown   Output Parameter:
67984df9cb4SJed Brown . adapt - new TSAdapt object
68084df9cb4SJed Brown 
68184df9cb4SJed Brown   Level: developer
68284df9cb4SJed Brown 
68384df9cb4SJed Brown   Notes:
68484df9cb4SJed Brown   TSAdapt creation is handled by TS, so users should not need to call this function.
68584df9cb4SJed Brown 
68684df9cb4SJed Brown .keywords: TSAdapt, create
687552698daSJed Brown .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy()
68884df9cb4SJed Brown @*/
68984df9cb4SJed Brown PetscErrorCode  TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt)
69084df9cb4SJed Brown {
69184df9cb4SJed Brown   PetscErrorCode ierr;
69284df9cb4SJed Brown   TSAdapt        adapt;
69384df9cb4SJed Brown 
69484df9cb4SJed Brown   PetscFunctionBegin;
6953b3bcf4cSLisandro Dalcin   PetscValidPointer(inadapt,1);
6963b3bcf4cSLisandro Dalcin   *inadapt = NULL;
6973b3bcf4cSLisandro Dalcin   ierr = TSAdaptInitializePackage();CHKERRQ(ierr);
6983b3bcf4cSLisandro Dalcin 
69973107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr);
7001c3436cfSJed Brown 
7011c3436cfSJed Brown   adapt->dt_min             = 1e-20;
7021c3436cfSJed Brown   adapt->dt_max             = 1e50;
70397335746SJed Brown   adapt->scale_solve_failed = 0.25;
7047619abb3SShri   adapt->wnormtype          = NORM_2;
7051c3436cfSJed Brown 
70684df9cb4SJed Brown   *inadapt = adapt;
70784df9cb4SJed Brown   PetscFunctionReturn(0);
70884df9cb4SJed Brown }
709