xref: /petsc/src/ts/adapt/interface/tsadapt.c (revision 7eef605549633c35e83b43f04cbaa60e08860c91)
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
82560360afSLisandro 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"
128*7eef6055SBarry Smith /*@C
129*7eef6055SBarry Smith   TSAdaptSetType - sets the approach used for the error adapter, currently there is only TSADAPTBASIC and TSADAPTNONE
130*7eef6055SBarry Smith 
131*7eef6055SBarry Smith   Logicially Collective on TSAdapt
132*7eef6055SBarry Smith 
133*7eef6055SBarry Smith   Input Parameter:
134*7eef6055SBarry Smith + adapt - the TS error adapter, most likely obtained with TSGetAdapt()
135*7eef6055SBarry Smith - type - either  TSADAPTBASIC or TSADAPTNONE
136*7eef6055SBarry Smith 
137*7eef6055SBarry Smith   Options Database:
138*7eef6055SBarry Smith .  -ts_adapt_type basic or none - to setting the adapter type
139*7eef6055SBarry Smith 
140*7eef6055SBarry Smith   Level: intermediate
141*7eef6055SBarry Smith 
142*7eef6055SBarry Smith .keywords: TSAdapt, create
143*7eef6055SBarry Smith 
144*7eef6055SBarry Smith .seealso: TSGetAdapt(), TSAdaptDestroy(), TSAdaptType, TSAdaptGetType()
145*7eef6055SBarry Smith @*/
14619fd82e9SBarry Smith PetscErrorCode  TSAdaptSetType(TSAdapt adapt,TSAdaptType type)
14784df9cb4SJed Brown {
148ef749922SLisandro Dalcin   PetscBool      match;
14984df9cb4SJed Brown   PetscErrorCode ierr,(*r)(TSAdapt);
15084df9cb4SJed Brown 
15184df9cb4SJed Brown   PetscFunctionBegin;
1524782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
153ef749922SLisandro Dalcin   ierr = PetscObjectTypeCompare((PetscObject)adapt,type,&match);CHKERRQ(ierr);
154ef749922SLisandro Dalcin   if (match) PetscFunctionReturn(0);
1551c9cd337SJed Brown   ierr = PetscFunctionListFind(TSAdaptList,type,&r);CHKERRQ(ierr);
15684df9cb4SJed Brown   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type);
157ef749922SLisandro Dalcin   if (adapt->ops->destroy) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);}
1584cefc2ffSBarry Smith   ierr = PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));CHKERRQ(ierr);
15984df9cb4SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr);
16068ae941aSLisandro Dalcin   ierr = (*r)(adapt);CHKERRQ(ierr);
16184df9cb4SJed Brown   PetscFunctionReturn(0);
16284df9cb4SJed Brown }
16384df9cb4SJed Brown 
16484df9cb4SJed Brown #undef __FUNCT__
16584df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetOptionsPrefix"
16684df9cb4SJed Brown PetscErrorCode  TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[])
16784df9cb4SJed Brown {
16884df9cb4SJed Brown   PetscErrorCode ierr;
16984df9cb4SJed Brown 
17084df9cb4SJed Brown   PetscFunctionBegin;
1714782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
17284df9cb4SJed Brown   ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr);
17384df9cb4SJed Brown   PetscFunctionReturn(0);
17484df9cb4SJed Brown }
17584df9cb4SJed Brown 
17684df9cb4SJed Brown #undef __FUNCT__
177ad6bc421SBarry Smith #define __FUNCT__ "TSAdaptLoad"
178ad6bc421SBarry Smith /*@C
179ad6bc421SBarry Smith   TSAdaptLoad - Loads a TSAdapt that has been stored in binary  with TSAdaptView().
180ad6bc421SBarry Smith 
181ad6bc421SBarry Smith   Collective on PetscViewer
182ad6bc421SBarry Smith 
183ad6bc421SBarry Smith   Input Parameters:
184ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or
185ad6bc421SBarry Smith            some related function before a call to TSAdaptLoad().
186ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
187ad6bc421SBarry Smith            HDF5 file viewer, obtained from PetscViewerHDF5Open()
188ad6bc421SBarry Smith 
189ad6bc421SBarry Smith    Level: intermediate
190ad6bc421SBarry Smith 
191ad6bc421SBarry Smith   Notes:
192ad6bc421SBarry Smith    The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored.
193ad6bc421SBarry Smith 
194ad6bc421SBarry Smith   Notes for advanced users:
195ad6bc421SBarry Smith   Most users should not need to know the details of the binary storage
196ad6bc421SBarry Smith   format, since TSAdaptLoad() and TSAdaptView() completely hide these details.
197ad6bc421SBarry Smith   But for anyone who's interested, the standard binary matrix storage
198ad6bc421SBarry Smith   format is
199ad6bc421SBarry Smith .vb
200ad6bc421SBarry Smith      has not yet been determined
201ad6bc421SBarry Smith .ve
202ad6bc421SBarry Smith 
203ad6bc421SBarry Smith .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad()
204ad6bc421SBarry Smith @*/
2054782b174SLisandro Dalcin PetscErrorCode  TSAdaptLoad(TSAdapt adapt,PetscViewer viewer)
206ad6bc421SBarry Smith {
207ad6bc421SBarry Smith   PetscErrorCode ierr;
208ad6bc421SBarry Smith   PetscBool      isbinary;
209ad6bc421SBarry Smith   char           type[256];
210ad6bc421SBarry Smith 
211ad6bc421SBarry Smith   PetscFunctionBegin;
2124782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
213ad6bc421SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
214ad6bc421SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
215ad6bc421SBarry Smith   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
216ad6bc421SBarry Smith 
217060da220SMatthew G. Knepley   ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr);
2184782b174SLisandro Dalcin   ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr);
2194782b174SLisandro Dalcin   if (adapt->ops->load) {
2204782b174SLisandro Dalcin     ierr = (*adapt->ops->load)(adapt,viewer);CHKERRQ(ierr);
221ad6bc421SBarry Smith   }
222ad6bc421SBarry Smith   PetscFunctionReturn(0);
223ad6bc421SBarry Smith }
224ad6bc421SBarry Smith 
225ad6bc421SBarry Smith #undef __FUNCT__
22684df9cb4SJed Brown #define __FUNCT__ "TSAdaptView"
22784df9cb4SJed Brown PetscErrorCode  TSAdaptView(TSAdapt adapt,PetscViewer viewer)
22884df9cb4SJed Brown {
22984df9cb4SJed Brown   PetscErrorCode ierr;
230ad6bc421SBarry Smith   PetscBool      iascii,isbinary;
23184df9cb4SJed Brown 
23284df9cb4SJed Brown   PetscFunctionBegin;
2334782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
2344782b174SLisandro Dalcin   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer);CHKERRQ(ierr);}
2354782b174SLisandro Dalcin   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
2364782b174SLisandro Dalcin   PetscCheckSameComm(adapt,1,viewer,2);
237251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
238ad6bc421SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
23984df9cb4SJed Brown   if (iascii) {
240dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr);
24184df9cb4SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  number of candidates %D\n",adapt->candidates.n);CHKERRQ(ierr);
24284df9cb4SJed Brown     if (adapt->ops->view) {
24384df9cb4SJed Brown       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
24484df9cb4SJed Brown       ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
24584df9cb4SJed Brown       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
24684df9cb4SJed Brown     }
247ad6bc421SBarry Smith   } else if (isbinary) {
248ad6bc421SBarry Smith     char type[256];
249ad6bc421SBarry Smith 
250ad6bc421SBarry Smith     /* need to save FILE_CLASS_ID for adapt class */
251ad6bc421SBarry Smith     ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr);
252ad6bc421SBarry Smith     ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr);
253bbd56ea5SKarl Rupp   } else if (adapt->ops->view) {
254f2c2a1b9SBarry Smith     ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
255f2c2a1b9SBarry Smith   }
25684df9cb4SJed Brown   PetscFunctionReturn(0);
25784df9cb4SJed Brown }
25884df9cb4SJed Brown 
25984df9cb4SJed Brown #undef __FUNCT__
260099af0a3SLisandro Dalcin #define __FUNCT__ "TSAdaptReset"
261099af0a3SLisandro Dalcin /*@
262099af0a3SLisandro Dalcin    TSAdaptReset - Resets a TSAdapt context.
263099af0a3SLisandro Dalcin 
264099af0a3SLisandro Dalcin    Collective on TS
265099af0a3SLisandro Dalcin 
266099af0a3SLisandro Dalcin    Input Parameter:
267099af0a3SLisandro Dalcin .  adapt - the TSAdapt context obtained from TSAdaptCreate()
268099af0a3SLisandro Dalcin 
269099af0a3SLisandro Dalcin    Level: developer
270099af0a3SLisandro Dalcin 
271099af0a3SLisandro Dalcin .seealso: TSAdaptCreate(), TSAdaptDestroy()
272099af0a3SLisandro Dalcin @*/
273099af0a3SLisandro Dalcin PetscErrorCode  TSAdaptReset(TSAdapt adapt)
274099af0a3SLisandro Dalcin {
275099af0a3SLisandro Dalcin   PetscErrorCode ierr;
276099af0a3SLisandro Dalcin 
277099af0a3SLisandro Dalcin   PetscFunctionBegin;
278099af0a3SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
279099af0a3SLisandro Dalcin   if (adapt->ops->reset) {ierr = (*adapt->ops->reset)(adapt);CHKERRQ(ierr);}
280099af0a3SLisandro Dalcin   PetscFunctionReturn(0);
281099af0a3SLisandro Dalcin }
282099af0a3SLisandro Dalcin 
283099af0a3SLisandro Dalcin #undef __FUNCT__
28484df9cb4SJed Brown #define __FUNCT__ "TSAdaptDestroy"
28584df9cb4SJed Brown PetscErrorCode  TSAdaptDestroy(TSAdapt *adapt)
28684df9cb4SJed Brown {
28784df9cb4SJed Brown   PetscErrorCode ierr;
28884df9cb4SJed Brown 
28984df9cb4SJed Brown   PetscFunctionBegin;
29084df9cb4SJed Brown   if (!*adapt) PetscFunctionReturn(0);
29184df9cb4SJed Brown   PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1);
2924782b174SLisandro Dalcin   if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);}
293099af0a3SLisandro Dalcin 
294099af0a3SLisandro Dalcin   ierr = TSAdaptReset(*adapt);CHKERRQ(ierr);
295099af0a3SLisandro Dalcin 
29684df9cb4SJed Brown   if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);}
2971c3436cfSJed Brown   ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr);
29884df9cb4SJed Brown   ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr);
29984df9cb4SJed Brown   PetscFunctionReturn(0);
30084df9cb4SJed Brown }
30184df9cb4SJed Brown 
30284df9cb4SJed Brown #undef __FUNCT__
3031c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetMonitor"
3041c3436cfSJed Brown /*@
3051c3436cfSJed Brown    TSAdaptSetMonitor - Monitor the choices made by the adaptive controller
3061c3436cfSJed Brown 
3071c3436cfSJed Brown    Collective on TSAdapt
3081c3436cfSJed Brown 
3091c3436cfSJed Brown    Input Arguments:
3101c3436cfSJed Brown +  adapt - adaptive controller context
3111c3436cfSJed Brown -  flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable
3121c3436cfSJed Brown 
3131c3436cfSJed Brown    Level: intermediate
3141c3436cfSJed Brown 
3151c3436cfSJed Brown .seealso: TSAdaptChoose()
3161c3436cfSJed Brown @*/
3171c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg)
3181c3436cfSJed Brown {
3191c3436cfSJed Brown   PetscErrorCode ierr;
3201c3436cfSJed Brown 
3211c3436cfSJed Brown   PetscFunctionBegin;
3224782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
3234782b174SLisandro Dalcin   PetscValidLogicalCollectiveBool(adapt,flg,2);
3241c3436cfSJed Brown   if (flg) {
325ce94432eSBarry Smith     if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);}
3261c3436cfSJed Brown   } else {
3271c3436cfSJed Brown     ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr);
3281c3436cfSJed Brown   }
3291c3436cfSJed Brown   PetscFunctionReturn(0);
3301c3436cfSJed Brown }
3311c3436cfSJed Brown 
3321c3436cfSJed Brown #undef __FUNCT__
3330873808bSJed Brown #define __FUNCT__ "TSAdaptSetCheckStage"
3340873808bSJed Brown /*@C
3350873808bSJed Brown    TSAdaptSetCheckStage - set a callback to check convergence for a stage
3360873808bSJed Brown 
3370873808bSJed Brown    Logically collective on TSAdapt
3380873808bSJed Brown 
3390873808bSJed Brown    Input Arguments:
3400873808bSJed Brown +  adapt - adaptive controller context
3410873808bSJed Brown -  func - stage check function
3420873808bSJed Brown 
3430873808bSJed Brown    Arguments of func:
3440873808bSJed Brown $  PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept)
3450873808bSJed Brown 
3460873808bSJed Brown +  adapt - adaptive controller context
3470873808bSJed Brown .  ts - time stepping context
3480873808bSJed Brown -  accept - pending choice of whether to accept, can be modified by this routine
3490873808bSJed Brown 
3500873808bSJed Brown    Level: advanced
3510873808bSJed Brown 
3520873808bSJed Brown .seealso: TSAdaptChoose()
3530873808bSJed Brown @*/
354b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscReal,Vec,PetscBool*))
3550873808bSJed Brown {
3560873808bSJed Brown 
3570873808bSJed Brown   PetscFunctionBegin;
3580873808bSJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
35968ae941aSLisandro Dalcin   adapt->checkstage = func;
3600873808bSJed Brown   PetscFunctionReturn(0);
3610873808bSJed Brown }
3620873808bSJed Brown 
3630873808bSJed Brown #undef __FUNCT__
3641c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetStepLimits"
3651c3436cfSJed Brown /*@
3661c3436cfSJed Brown    TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller
3671c3436cfSJed Brown 
3681c3436cfSJed Brown    Logically Collective
3691c3436cfSJed Brown 
3701c3436cfSJed Brown    Input Arguments:
371552698daSJed Brown +  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
3721c3436cfSJed Brown .  hmin - minimum time step
3731c3436cfSJed Brown -  hmax - maximum time step
3741c3436cfSJed Brown 
3751c3436cfSJed Brown    Options Database Keys:
3761c3436cfSJed Brown +  -ts_adapt_dt_min - minimum time step
3771c3436cfSJed Brown -  -ts_adapt_dt_max - maximum time step
3781c3436cfSJed Brown 
3791c3436cfSJed Brown    Level: intermediate
3801c3436cfSJed Brown 
3811c3436cfSJed Brown .seealso: TSAdapt
3821c3436cfSJed Brown @*/
3831c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax)
3841c3436cfSJed Brown {
3851c3436cfSJed Brown 
3861c3436cfSJed Brown   PetscFunctionBegin;
3874782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
388b1f5bccaSLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt,hmin,2);
389b1f5bccaSLisandro Dalcin   PetscValidLogicalCollectiveReal(adapt,hmax,3);
390b1f5bccaSLisandro 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);
391b1f5bccaSLisandro 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);
392b1f5bccaSLisandro Dalcin   if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin;
393b1f5bccaSLisandro Dalcin   if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax;
394b1f5bccaSLisandro Dalcin #if defined(PETSC_USE_DEBUG)
395b1f5bccaSLisandro Dalcin   hmin = adapt->dt_min;
396b1f5bccaSLisandro Dalcin   hmax = adapt->dt_max;
397b1f5bccaSLisandro 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);
398b1f5bccaSLisandro Dalcin #endif
3991c3436cfSJed Brown   PetscFunctionReturn(0);
4001c3436cfSJed Brown }
4011c3436cfSJed Brown 
4021c3436cfSJed Brown #undef __FUNCT__
40384df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetFromOptions"
404e55864a3SBarry Smith /*
40584df9cb4SJed Brown    TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options.
40684df9cb4SJed Brown 
40784df9cb4SJed Brown    Collective on TSAdapt
40884df9cb4SJed Brown 
40984df9cb4SJed Brown    Input Parameter:
41084df9cb4SJed Brown .  adapt - the TSAdapt context
41184df9cb4SJed Brown 
41284df9cb4SJed Brown    Options Database Keys:
41384df9cb4SJed Brown .  -ts_adapt_type <type> - basic
41484df9cb4SJed Brown 
41584df9cb4SJed Brown    Level: advanced
41684df9cb4SJed Brown 
41784df9cb4SJed Brown    Notes:
41884df9cb4SJed Brown    This function is automatically called by TSSetFromOptions()
41984df9cb4SJed Brown 
420552698daSJed Brown .keywords: TS, TSGetAdapt(), TSAdaptSetType()
42184df9cb4SJed Brown 
42284df9cb4SJed Brown .seealso: TSGetType()
423e55864a3SBarry Smith */
4244416b707SBarry Smith PetscErrorCode  TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt)
42584df9cb4SJed Brown {
42684df9cb4SJed Brown   PetscErrorCode ierr;
42784df9cb4SJed Brown   char           type[256] = TSADAPTBASIC;
4281c3436cfSJed Brown   PetscBool      set,flg;
42984df9cb4SJed Brown 
43084df9cb4SJed Brown   PetscFunctionBegin;
4314782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
43284df9cb4SJed Brown   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
4331566a47fSLisandro Dalcin    * function can only be called from inside TSSetFromOptions()  */
434e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");CHKERRQ(ierr);
435a264d7a6SBarry Smith   ierr = PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList,
4368caf3d72SBarry Smith                           ((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);CHKERRQ(ierr);
43784df9cb4SJed Brown   if (flg || !((PetscObject)adapt)->type_name) {
43884df9cb4SJed Brown     ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr);
43984df9cb4SJed Brown   }
4400298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",adapt->dt_min,&adapt->dt_min,NULL);CHKERRQ(ierr);
4410298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",adapt->dt_max,&adapt->dt_max,NULL);CHKERRQ(ierr);
4420298fd71SBarry 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);
4431c3436cfSJed Brown   ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
4447619abb3SShri   ierr = PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);CHKERRQ(ierr);
4457619abb3SShri   if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported");
4461c3436cfSJed Brown   if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);}
447e55864a3SBarry Smith   if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);}
44884df9cb4SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
44984df9cb4SJed Brown   PetscFunctionReturn(0);
45084df9cb4SJed Brown }
45184df9cb4SJed Brown 
45284df9cb4SJed Brown #undef __FUNCT__
45384df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidatesClear"
45484df9cb4SJed Brown /*@
45584df9cb4SJed Brown    TSAdaptCandidatesClear - clear any previously set candidate schemes
45684df9cb4SJed Brown 
45784df9cb4SJed Brown    Logically Collective
45884df9cb4SJed Brown 
45984df9cb4SJed Brown    Input Argument:
46084df9cb4SJed Brown .  adapt - adaptive controller
46184df9cb4SJed Brown 
46284df9cb4SJed Brown    Level: developer
46384df9cb4SJed Brown 
46484df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose()
46584df9cb4SJed Brown @*/
46684df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
46784df9cb4SJed Brown {
46884df9cb4SJed Brown   PetscErrorCode ierr;
46984df9cb4SJed Brown 
47084df9cb4SJed Brown   PetscFunctionBegin;
4714782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
47284df9cb4SJed Brown   ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr);
47384df9cb4SJed Brown   PetscFunctionReturn(0);
47484df9cb4SJed Brown }
47584df9cb4SJed Brown 
47684df9cb4SJed Brown #undef __FUNCT__
47784df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidateAdd"
47884df9cb4SJed Brown /*@C
47984df9cb4SJed Brown    TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from
48084df9cb4SJed Brown 
48184df9cb4SJed Brown    Logically Collective
48284df9cb4SJed Brown 
48384df9cb4SJed Brown    Input Arguments:
484552698daSJed Brown +  adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate()
48584df9cb4SJed Brown .  name - name of the candidate scheme to add
48684df9cb4SJed Brown .  order - order of the candidate scheme
48784df9cb4SJed Brown .  stageorder - stage order of the candidate scheme
4888d59e960SJed Brown .  ccfl - stability coefficient relative to explicit Euler, used for CFL constraints
48984df9cb4SJed Brown .  cost - relative measure of the amount of work required for the candidate scheme
49084df9cb4SJed Brown -  inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme
49184df9cb4SJed Brown 
49284df9cb4SJed Brown    Note:
49384df9cb4SJed Brown    This routine is not available in Fortran.
49484df9cb4SJed Brown 
49584df9cb4SJed Brown    Level: developer
49684df9cb4SJed Brown 
49784df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose()
49884df9cb4SJed Brown @*/
4998d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse)
50084df9cb4SJed Brown {
50184df9cb4SJed Brown   PetscInt c;
50284df9cb4SJed Brown 
50384df9cb4SJed Brown   PetscFunctionBegin;
50484df9cb4SJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
505ce94432eSBarry Smith   if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order);
50684df9cb4SJed Brown   if (inuse) {
507ce94432eSBarry 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()");
50884df9cb4SJed Brown     adapt->candidates.inuse_set = PETSC_TRUE;
50984df9cb4SJed Brown   }
5101c3436cfSJed Brown   /* first slot if this is the current scheme, otherwise the next available slot */
5111c3436cfSJed Brown   c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;
512bbd56ea5SKarl Rupp 
51384df9cb4SJed Brown   adapt->candidates.name[c]       = name;
51484df9cb4SJed Brown   adapt->candidates.order[c]      = order;
51584df9cb4SJed Brown   adapt->candidates.stageorder[c] = stageorder;
5168d59e960SJed Brown   adapt->candidates.ccfl[c]       = ccfl;
51784df9cb4SJed Brown   adapt->candidates.cost[c]       = cost;
51884df9cb4SJed Brown   adapt->candidates.n++;
51984df9cb4SJed Brown   PetscFunctionReturn(0);
52084df9cb4SJed Brown }
52184df9cb4SJed Brown 
52284df9cb4SJed Brown #undef __FUNCT__
5238d59e960SJed Brown #define __FUNCT__ "TSAdaptCandidatesGet"
5248d59e960SJed Brown /*@C
5258d59e960SJed Brown    TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost
5268d59e960SJed Brown 
5278d59e960SJed Brown    Not Collective
5288d59e960SJed Brown 
5298d59e960SJed Brown    Input Arguments:
5308d59e960SJed Brown .  adapt - time step adaptivity context
5318d59e960SJed Brown 
5328d59e960SJed Brown    Output Arguments:
5338d59e960SJed Brown +  n - number of candidate schemes, always at least 1
5348d59e960SJed Brown .  order - the order of each candidate scheme
5358d59e960SJed Brown .  stageorder - the stage order of each candidate scheme
5368d59e960SJed Brown .  ccfl - the CFL coefficient of each scheme
5378d59e960SJed Brown -  cost - the relative cost of each scheme
5388d59e960SJed Brown 
5398d59e960SJed Brown    Level: developer
5408d59e960SJed Brown 
5418d59e960SJed Brown    Note:
5428d59e960SJed Brown    The current scheme is always returned in the first slot
5438d59e960SJed Brown 
5448d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose()
5458d59e960SJed Brown @*/
5468d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost)
5478d59e960SJed Brown {
5488d59e960SJed Brown   PetscFunctionBegin;
5498d59e960SJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
5508d59e960SJed Brown   if (n) *n = adapt->candidates.n;
5518d59e960SJed Brown   if (order) *order = adapt->candidates.order;
5528d59e960SJed Brown   if (stageorder) *stageorder = adapt->candidates.stageorder;
5538d59e960SJed Brown   if (ccfl) *ccfl = adapt->candidates.ccfl;
5548d59e960SJed Brown   if (cost) *cost = adapt->candidates.cost;
5558d59e960SJed Brown   PetscFunctionReturn(0);
5568d59e960SJed Brown }
5578d59e960SJed Brown 
5588d59e960SJed Brown #undef __FUNCT__
55984df9cb4SJed Brown #define __FUNCT__ "TSAdaptChoose"
56084df9cb4SJed Brown /*@C
56184df9cb4SJed Brown    TSAdaptChoose - choose which method and step size to use for the next step
56284df9cb4SJed Brown 
56384df9cb4SJed Brown    Logically Collective
56484df9cb4SJed Brown 
56584df9cb4SJed Brown    Input Arguments:
56684df9cb4SJed Brown +  adapt - adaptive contoller
56784df9cb4SJed Brown -  h - current step size
56884df9cb4SJed Brown 
56984df9cb4SJed Brown    Output Arguments:
5701566a47fSLisandro Dalcin +  next_sc - optional, scheme to use for the next step
57184df9cb4SJed Brown .  next_h - step size to use for the next step
57284df9cb4SJed Brown -  accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size
57384df9cb4SJed Brown 
5741c3436cfSJed Brown    Note:
5751c3436cfSJed Brown    The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is
5761c3436cfSJed Brown    being retried after an initial rejection.
5771c3436cfSJed Brown 
57884df9cb4SJed Brown    Level: developer
57984df9cb4SJed Brown 
58084df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd()
58184df9cb4SJed Brown @*/
58284df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept)
58384df9cb4SJed Brown {
58484df9cb4SJed Brown   PetscErrorCode ierr;
5851566a47fSLisandro Dalcin   PetscInt       ncandidates = adapt->candidates.n;
5861566a47fSLisandro Dalcin   PetscInt       scheme = 0;
5870b99f514SJed Brown   PetscReal      wlte = -1.0;
58884df9cb4SJed Brown 
58984df9cb4SJed Brown   PetscFunctionBegin;
59084df9cb4SJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
59184df9cb4SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
5921566a47fSLisandro Dalcin   if (next_sc) PetscValidIntPointer(next_sc,4);
59384df9cb4SJed Brown   PetscValidPointer(next_h,5);
59484df9cb4SJed Brown   PetscValidIntPointer(accept,6);
5951566a47fSLisandro Dalcin   if (next_sc) *next_sc = 0;
5961566a47fSLisandro Dalcin 
5971566a47fSLisandro Dalcin   /* Do not mess with adaptivity while handling events*/
5981566a47fSLisandro Dalcin   if (ts->event && ts->event->status != TSEVENT_NONE) {
5991566a47fSLisandro Dalcin     *next_h = h;
6001566a47fSLisandro Dalcin     *accept = PETSC_TRUE;
6011566a47fSLisandro Dalcin     PetscFunctionReturn(0);
6021566a47fSLisandro Dalcin   }
6031566a47fSLisandro Dalcin 
6041566a47fSLisandro Dalcin   ierr = (*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte);CHKERRQ(ierr);
6051566a47fSLisandro 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);
6061566a47fSLisandro Dalcin   if (*next_h < 0) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h);
6071566a47fSLisandro Dalcin   if (next_sc) *next_sc = scheme;
6081566a47fSLisandro Dalcin 
60949354f04SShri Abhyankar   if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
61049354f04SShri Abhyankar     /* Reduce time step if it overshoots max time */
6111566a47fSLisandro Dalcin     if (ts->ptime + ts->time_step + *next_h >= ts->max_time) {
6121566a47fSLisandro Dalcin       PetscReal next_dt = ts->max_time - (ts->ptime + ts->time_step);
6138709b12bSShri Abhyankar       if (next_dt > PETSC_SMALL) *next_h = next_dt;
61449354f04SShri Abhyankar       else ts->reason = TS_CONVERGED_TIME;
61549354f04SShri Abhyankar     }
61649354f04SShri Abhyankar   }
6171c3436cfSJed Brown 
6181c3436cfSJed Brown   if (adapt->monitor) {
6191566a47fSLisandro Dalcin     const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : "";
6201c3436cfSJed Brown     ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
6210b99f514SJed Brown     if (wlte < 0) {
6221566a47fSLisandro 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);
6230b99f514SJed Brown     } else {
6241566a47fSLisandro 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);
6250b99f514SJed Brown     }
6261c3436cfSJed Brown     ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
6271c3436cfSJed Brown   }
62884df9cb4SJed Brown   PetscFunctionReturn(0);
62984df9cb4SJed Brown }
63084df9cb4SJed Brown 
63184df9cb4SJed Brown #undef __FUNCT__
63297335746SJed Brown #define __FUNCT__ "TSAdaptCheckStage"
63397335746SJed Brown /*@
63497335746SJed Brown    TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails)
63597335746SJed Brown 
63697335746SJed Brown    Collective
63797335746SJed Brown 
63897335746SJed Brown    Input Arguments:
63997335746SJed Brown +  adapt - adaptive controller context
640b295832fSPierre Barbier de Reuille .  ts - time stepper
641b295832fSPierre Barbier de Reuille .  t - Current simulation time
642b295832fSPierre Barbier de Reuille -  Y - Current solution vector
64397335746SJed Brown 
64497335746SJed Brown    Output Arguments:
64597335746SJed Brown .  accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject
64697335746SJed Brown 
64797335746SJed Brown    Level: developer
64897335746SJed Brown 
64997335746SJed Brown .seealso:
65097335746SJed Brown @*/
651b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept)
65297335746SJed Brown {
65397335746SJed Brown   PetscErrorCode      ierr;
6541566a47fSLisandro Dalcin   SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING;
65597335746SJed Brown 
65697335746SJed Brown   PetscFunctionBegin;
6574782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
6584782b174SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
6594782b174SLisandro Dalcin   PetscValidIntPointer(accept,3);
6601566a47fSLisandro Dalcin 
6611566a47fSLisandro Dalcin   if (ts->snes) {ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);}
66297335746SJed Brown   if (snesreason < 0) {
66397335746SJed Brown     *accept = PETSC_FALSE;
6646de24e2aSJed Brown     if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) {
66597335746SJed Brown       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
66648c19aefSLisandro 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);
66797335746SJed Brown       if (adapt->monitor) {
66897335746SJed Brown         ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
6691566a47fSLisandro 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);
67097335746SJed Brown         ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
67197335746SJed Brown       }
672cb9d8021SPierre Barbier de Reuille     }
673cb9d8021SPierre Barbier de Reuille   } else {
6741566a47fSLisandro Dalcin     *accept = PETSC_TRUE;
675b295832fSPierre Barbier de Reuille     ierr = TSFunctionDomainError(ts,t,Y,accept);CHKERRQ(ierr);
676cb9d8021SPierre Barbier de Reuille     if(*accept && adapt->checkstage) {
677b295832fSPierre Barbier de Reuille       ierr = (*adapt->checkstage)(adapt,ts,t,Y,accept);CHKERRQ(ierr);
678cb9d8021SPierre Barbier de Reuille     }
679cb9d8021SPierre Barbier de Reuille   }
680cb9d8021SPierre Barbier de Reuille 
6811566a47fSLisandro Dalcin   if(!(*accept) && !ts->reason) {
6821566a47fSLisandro Dalcin     PetscReal dt,new_dt;
6831566a47fSLisandro Dalcin     ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
684cb9d8021SPierre Barbier de Reuille     new_dt = dt * adapt->scale_solve_failed;
68597335746SJed Brown     ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr);
68697335746SJed Brown     if (adapt->monitor) {
68797335746SJed Brown       ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
68897335746SJed 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);
68997335746SJed Brown       ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
69097335746SJed Brown     }
69197335746SJed Brown   }
69297335746SJed Brown   PetscFunctionReturn(0);
69397335746SJed Brown }
69497335746SJed Brown 
69597335746SJed Brown #undef __FUNCT__
69684df9cb4SJed Brown #define __FUNCT__ "TSAdaptCreate"
69784df9cb4SJed Brown /*@
69884df9cb4SJed Brown   TSAdaptCreate - create an adaptive controller context for time stepping
69984df9cb4SJed Brown 
70084df9cb4SJed Brown   Collective on MPI_Comm
70184df9cb4SJed Brown 
70284df9cb4SJed Brown   Input Parameter:
70384df9cb4SJed Brown . comm - The communicator
70484df9cb4SJed Brown 
70584df9cb4SJed Brown   Output Parameter:
70684df9cb4SJed Brown . adapt - new TSAdapt object
70784df9cb4SJed Brown 
70884df9cb4SJed Brown   Level: developer
70984df9cb4SJed Brown 
71084df9cb4SJed Brown   Notes:
71184df9cb4SJed Brown   TSAdapt creation is handled by TS, so users should not need to call this function.
71284df9cb4SJed Brown 
71384df9cb4SJed Brown .keywords: TSAdapt, create
714552698daSJed Brown .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy()
71584df9cb4SJed Brown @*/
71684df9cb4SJed Brown PetscErrorCode  TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt)
71784df9cb4SJed Brown {
71884df9cb4SJed Brown   PetscErrorCode ierr;
71984df9cb4SJed Brown   TSAdapt        adapt;
72084df9cb4SJed Brown 
72184df9cb4SJed Brown   PetscFunctionBegin;
7223b3bcf4cSLisandro Dalcin   PetscValidPointer(inadapt,1);
7233b3bcf4cSLisandro Dalcin   *inadapt = NULL;
7243b3bcf4cSLisandro Dalcin   ierr = TSAdaptInitializePackage();CHKERRQ(ierr);
7253b3bcf4cSLisandro Dalcin 
72673107ff1SLisandro Dalcin   ierr = PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr);
7271c3436cfSJed Brown 
7281c3436cfSJed Brown   adapt->dt_min             = 1e-20;
7291c3436cfSJed Brown   adapt->dt_max             = 1e50;
73097335746SJed Brown   adapt->scale_solve_failed = 0.25;
7317619abb3SShri   adapt->wnormtype          = NORM_2;
732*7eef6055SBarry Smith   ierr = TSAdaptSetType(adapt,TSADAPTBASIC);CHKERRQ(ierr);
7331c3436cfSJed Brown 
73484df9cb4SJed Brown   *inadapt = adapt;
73584df9cb4SJed Brown   PetscFunctionReturn(0);
73684df9cb4SJed Brown }
737