xref: /petsc/src/ts/adapt/interface/tsadapt.c (revision 060da220a69df799affd088f0100c1f4dac789cb)
184df9cb4SJed Brown 
2b45d2f2cSJed Brown #include <petsc-private/tsimpl.h> /*I  "petscts.h" I*/
384df9cb4SJed Brown 
4140e18c1SBarry Smith static PetscFunctionList TSAdaptList;
584df9cb4SJed Brown static PetscBool         TSAdaptPackageInitialized;
684df9cb4SJed Brown static PetscBool         TSAdaptRegisterAllCalled;
784df9cb4SJed Brown static PetscClassId      TSADAPT_CLASSID;
884df9cb4SJed Brown 
98cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt);
108cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_None(TSAdapt);
118cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt);
1284df9cb4SJed Brown 
1384df9cb4SJed Brown #undef __FUNCT__
1484df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegister"
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 #undef __FUNCT__
5384df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegisterAll"
5484df9cb4SJed Brown /*@C
5584df9cb4SJed Brown   TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt
5684df9cb4SJed Brown 
5784df9cb4SJed Brown   Not Collective
5884df9cb4SJed Brown 
5984df9cb4SJed Brown   Level: advanced
6084df9cb4SJed Brown 
6184df9cb4SJed Brown .keywords: TSAdapt, register, all
6284df9cb4SJed Brown 
6384df9cb4SJed Brown .seealso: TSAdaptRegisterDestroy()
6484df9cb4SJed Brown @*/
65607a6623SBarry Smith PetscErrorCode  TSAdaptRegisterAll(void)
6684df9cb4SJed Brown {
6784df9cb4SJed Brown   PetscErrorCode ierr;
6884df9cb4SJed Brown 
6984df9cb4SJed Brown   PetscFunctionBegin;
700f51fdf8SToby Isaac   if (TSAdaptRegisterAllCalled) PetscFunctionReturn(0);
710f51fdf8SToby Isaac   TSAdaptRegisterAllCalled = PETSC_TRUE;
72bdf89e91SBarry Smith   ierr = TSAdaptRegister(TSADAPTBASIC,TSAdaptCreate_Basic);CHKERRQ(ierr);
73bdf89e91SBarry Smith   ierr = TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None);CHKERRQ(ierr);
74bdf89e91SBarry Smith   ierr = TSAdaptRegister(TSADAPTCFL,  TSAdaptCreate_CFL);CHKERRQ(ierr);
7584df9cb4SJed Brown   PetscFunctionReturn(0);
7684df9cb4SJed Brown }
7784df9cb4SJed Brown 
7884df9cb4SJed Brown #undef __FUNCT__
7984df9cb4SJed Brown #define __FUNCT__ "TSAdaptFinalizePackage"
8084df9cb4SJed Brown /*@C
8184df9cb4SJed Brown   TSFinalizePackage - This function destroys everything in the TS package. It is
8284df9cb4SJed Brown   called from PetscFinalize().
8384df9cb4SJed Brown 
8484df9cb4SJed Brown   Level: developer
8584df9cb4SJed Brown 
8684df9cb4SJed Brown .keywords: Petsc, destroy, package
8784df9cb4SJed Brown .seealso: PetscFinalize()
8884df9cb4SJed Brown @*/
8984df9cb4SJed Brown PetscErrorCode  TSAdaptFinalizePackage(void)
9084df9cb4SJed Brown {
9137e93019SBarry Smith   PetscErrorCode ierr;
9237e93019SBarry Smith 
9384df9cb4SJed Brown   PetscFunctionBegin;
9437e93019SBarry Smith   ierr = PetscFunctionListDestroy(&TSAdaptList);CHKERRQ(ierr);
9584df9cb4SJed Brown   TSAdaptPackageInitialized = PETSC_FALSE;
9684df9cb4SJed Brown   TSAdaptRegisterAllCalled  = PETSC_FALSE;
9784df9cb4SJed Brown   PetscFunctionReturn(0);
9884df9cb4SJed Brown }
9984df9cb4SJed Brown 
10084df9cb4SJed Brown #undef __FUNCT__
10184df9cb4SJed Brown #define __FUNCT__ "TSAdaptInitializePackage"
10284df9cb4SJed Brown /*@C
10384df9cb4SJed Brown   TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is
10484df9cb4SJed Brown   called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to
10584df9cb4SJed Brown   TSCreate_GL() when using static libraries.
10684df9cb4SJed Brown 
10784df9cb4SJed Brown   Level: developer
10884df9cb4SJed Brown 
10984df9cb4SJed Brown .keywords: TSAdapt, initialize, package
11084df9cb4SJed Brown .seealso: PetscInitialize()
11184df9cb4SJed Brown @*/
112607a6623SBarry Smith PetscErrorCode  TSAdaptInitializePackage(void)
11384df9cb4SJed Brown {
11484df9cb4SJed Brown   PetscErrorCode ierr;
11584df9cb4SJed Brown 
11684df9cb4SJed Brown   PetscFunctionBegin;
11784df9cb4SJed Brown   if (TSAdaptPackageInitialized) PetscFunctionReturn(0);
11884df9cb4SJed Brown   TSAdaptPackageInitialized = PETSC_TRUE;
11984df9cb4SJed Brown   ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr);
120607a6623SBarry Smith   ierr = TSAdaptRegisterAll();CHKERRQ(ierr);
12184df9cb4SJed Brown   ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr);
12284df9cb4SJed Brown   PetscFunctionReturn(0);
12384df9cb4SJed Brown }
12484df9cb4SJed Brown 
12584df9cb4SJed Brown #undef __FUNCT__
12684df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetType"
12719fd82e9SBarry Smith PetscErrorCode  TSAdaptSetType(TSAdapt adapt,TSAdaptType type)
12884df9cb4SJed Brown {
129ef749922SLisandro Dalcin   PetscBool      match;
1303403ddbcSLisandro Dalcin   PetscErrorCode (*checkstage)(TSAdapt,TS,PetscBool*);
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   /* Reinitialize function pointers in TSAdaptOps structure */
1413403ddbcSLisandro Dalcin   checkstage = adapt->ops->checkstage;
1424cefc2ffSBarry Smith   ierr = PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));CHKERRQ(ierr);
1433403ddbcSLisandro Dalcin   adapt->ops->checkstage = checkstage;
14484df9cb4SJed Brown   ierr = (*r)(adapt);CHKERRQ(ierr);
14584df9cb4SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr);
14684df9cb4SJed Brown   PetscFunctionReturn(0);
14784df9cb4SJed Brown }
14884df9cb4SJed Brown 
14984df9cb4SJed Brown #undef __FUNCT__
15084df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetOptionsPrefix"
15184df9cb4SJed Brown PetscErrorCode  TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[])
15284df9cb4SJed Brown {
15384df9cb4SJed Brown   PetscErrorCode ierr;
15484df9cb4SJed Brown 
15584df9cb4SJed Brown   PetscFunctionBegin;
1564782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
15784df9cb4SJed Brown   ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr);
15884df9cb4SJed Brown   PetscFunctionReturn(0);
15984df9cb4SJed Brown }
16084df9cb4SJed Brown 
16184df9cb4SJed Brown #undef __FUNCT__
162ad6bc421SBarry Smith #define __FUNCT__ "TSAdaptLoad"
163ad6bc421SBarry Smith /*@C
164ad6bc421SBarry Smith   TSAdaptLoad - Loads a TSAdapt that has been stored in binary  with TSAdaptView().
165ad6bc421SBarry Smith 
166ad6bc421SBarry Smith   Collective on PetscViewer
167ad6bc421SBarry Smith 
168ad6bc421SBarry Smith   Input Parameters:
169ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or
170ad6bc421SBarry Smith            some related function before a call to TSAdaptLoad().
171ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
172ad6bc421SBarry Smith            HDF5 file viewer, obtained from PetscViewerHDF5Open()
173ad6bc421SBarry Smith 
174ad6bc421SBarry Smith    Level: intermediate
175ad6bc421SBarry Smith 
176ad6bc421SBarry Smith   Notes:
177ad6bc421SBarry Smith    The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored.
178ad6bc421SBarry Smith 
179ad6bc421SBarry Smith   Notes for advanced users:
180ad6bc421SBarry Smith   Most users should not need to know the details of the binary storage
181ad6bc421SBarry Smith   format, since TSAdaptLoad() and TSAdaptView() completely hide these details.
182ad6bc421SBarry Smith   But for anyone who's interested, the standard binary matrix storage
183ad6bc421SBarry Smith   format is
184ad6bc421SBarry Smith .vb
185ad6bc421SBarry Smith      has not yet been determined
186ad6bc421SBarry Smith .ve
187ad6bc421SBarry Smith 
188ad6bc421SBarry Smith .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad()
189ad6bc421SBarry Smith @*/
1904782b174SLisandro Dalcin PetscErrorCode  TSAdaptLoad(TSAdapt adapt,PetscViewer viewer)
191ad6bc421SBarry Smith {
192ad6bc421SBarry Smith   PetscErrorCode ierr;
193ad6bc421SBarry Smith   PetscBool      isbinary;
194ad6bc421SBarry Smith   char           type[256];
195ad6bc421SBarry Smith 
196ad6bc421SBarry Smith   PetscFunctionBegin;
1974782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
198ad6bc421SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
199ad6bc421SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
200ad6bc421SBarry Smith   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
201ad6bc421SBarry Smith 
202*060da220SMatthew G. Knepley   ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr);
2034782b174SLisandro Dalcin   ierr = TSAdaptSetType(adapt, type);CHKERRQ(ierr);
2044782b174SLisandro Dalcin   if (adapt->ops->load) {
2054782b174SLisandro Dalcin     ierr = (*adapt->ops->load)(adapt,viewer);CHKERRQ(ierr);
206ad6bc421SBarry Smith   }
207ad6bc421SBarry Smith   PetscFunctionReturn(0);
208ad6bc421SBarry Smith }
209ad6bc421SBarry Smith 
210ad6bc421SBarry Smith #undef __FUNCT__
21184df9cb4SJed Brown #define __FUNCT__ "TSAdaptView"
21284df9cb4SJed Brown PetscErrorCode  TSAdaptView(TSAdapt adapt,PetscViewer viewer)
21384df9cb4SJed Brown {
21484df9cb4SJed Brown   PetscErrorCode ierr;
215ad6bc421SBarry Smith   PetscBool      iascii,isbinary;
21684df9cb4SJed Brown 
21784df9cb4SJed Brown   PetscFunctionBegin;
2184782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
2194782b174SLisandro Dalcin   if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer);CHKERRQ(ierr);}
2204782b174SLisandro Dalcin   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
2214782b174SLisandro Dalcin   PetscCheckSameComm(adapt,1,viewer,2);
222251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
223ad6bc421SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
22484df9cb4SJed Brown   if (iascii) {
225dae58748SBarry Smith     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr);
22684df9cb4SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  number of candidates %D\n",adapt->candidates.n);CHKERRQ(ierr);
22784df9cb4SJed Brown     if (adapt->ops->view) {
22884df9cb4SJed Brown       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
22984df9cb4SJed Brown       ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
23084df9cb4SJed Brown       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
23184df9cb4SJed Brown     }
232ad6bc421SBarry Smith   } else if (isbinary) {
233ad6bc421SBarry Smith     char type[256];
234ad6bc421SBarry Smith 
235ad6bc421SBarry Smith     /* need to save FILE_CLASS_ID for adapt class */
236ad6bc421SBarry Smith     ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr);
237ad6bc421SBarry Smith     ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr);
238bbd56ea5SKarl Rupp   } else if (adapt->ops->view) {
239f2c2a1b9SBarry Smith     ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
240f2c2a1b9SBarry Smith   }
24184df9cb4SJed Brown   PetscFunctionReturn(0);
24284df9cb4SJed Brown }
24384df9cb4SJed Brown 
24484df9cb4SJed Brown #undef __FUNCT__
24584df9cb4SJed Brown #define __FUNCT__ "TSAdaptDestroy"
24684df9cb4SJed Brown PetscErrorCode  TSAdaptDestroy(TSAdapt *adapt)
24784df9cb4SJed Brown {
24884df9cb4SJed Brown   PetscErrorCode ierr;
24984df9cb4SJed Brown 
25084df9cb4SJed Brown   PetscFunctionBegin;
25184df9cb4SJed Brown   if (!*adapt) PetscFunctionReturn(0);
25284df9cb4SJed Brown   PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1);
2534782b174SLisandro Dalcin   if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);}
25484df9cb4SJed Brown   if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);}
2551c3436cfSJed Brown   ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr);
25684df9cb4SJed Brown   ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr);
25784df9cb4SJed Brown   PetscFunctionReturn(0);
25884df9cb4SJed Brown }
25984df9cb4SJed Brown 
26084df9cb4SJed Brown #undef __FUNCT__
2611c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetMonitor"
2621c3436cfSJed Brown /*@
2631c3436cfSJed Brown    TSAdaptSetMonitor - Monitor the choices made by the adaptive controller
2641c3436cfSJed Brown 
2651c3436cfSJed Brown    Collective on TSAdapt
2661c3436cfSJed Brown 
2671c3436cfSJed Brown    Input Arguments:
2681c3436cfSJed Brown +  adapt - adaptive controller context
2691c3436cfSJed Brown -  flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable
2701c3436cfSJed Brown 
2711c3436cfSJed Brown    Level: intermediate
2721c3436cfSJed Brown 
2731c3436cfSJed Brown .seealso: TSAdaptChoose()
2741c3436cfSJed Brown @*/
2751c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg)
2761c3436cfSJed Brown {
2771c3436cfSJed Brown   PetscErrorCode ierr;
2781c3436cfSJed Brown 
2791c3436cfSJed Brown   PetscFunctionBegin;
2804782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
2814782b174SLisandro Dalcin   PetscValidLogicalCollectiveBool(adapt,flg,2);
2821c3436cfSJed Brown   if (flg) {
283ce94432eSBarry Smith     if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);}
2841c3436cfSJed Brown   } else {
2851c3436cfSJed Brown     ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr);
2861c3436cfSJed Brown   }
2871c3436cfSJed Brown   PetscFunctionReturn(0);
2881c3436cfSJed Brown }
2891c3436cfSJed Brown 
2901c3436cfSJed Brown #undef __FUNCT__
2910873808bSJed Brown #define __FUNCT__ "TSAdaptSetCheckStage"
2920873808bSJed Brown /*@C
2930873808bSJed Brown    TSAdaptSetCheckStage - set a callback to check convergence for a stage
2940873808bSJed Brown 
2950873808bSJed Brown    Logically collective on TSAdapt
2960873808bSJed Brown 
2970873808bSJed Brown    Input Arguments:
2980873808bSJed Brown +  adapt - adaptive controller context
2990873808bSJed Brown -  func - stage check function
3000873808bSJed Brown 
3010873808bSJed Brown    Arguments of func:
3020873808bSJed Brown $  PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept)
3030873808bSJed Brown 
3040873808bSJed Brown +  adapt - adaptive controller context
3050873808bSJed Brown .  ts - time stepping context
3060873808bSJed Brown -  accept - pending choice of whether to accept, can be modified by this routine
3070873808bSJed Brown 
3080873808bSJed Brown    Level: advanced
3090873808bSJed Brown 
3100873808bSJed Brown .seealso: TSAdaptChoose()
3110873808bSJed Brown @*/
3120873808bSJed Brown PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscBool*))
3130873808bSJed Brown {
3140873808bSJed Brown 
3150873808bSJed Brown   PetscFunctionBegin;
3160873808bSJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
3170873808bSJed Brown   adapt->ops->checkstage = func;
3180873808bSJed Brown   PetscFunctionReturn(0);
3190873808bSJed Brown }
3200873808bSJed Brown 
3210873808bSJed Brown #undef __FUNCT__
3221c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetStepLimits"
3231c3436cfSJed Brown /*@
3241c3436cfSJed Brown    TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller
3251c3436cfSJed Brown 
3261c3436cfSJed Brown    Logically Collective
3271c3436cfSJed Brown 
3281c3436cfSJed Brown    Input Arguments:
329552698daSJed Brown +  adapt - time step adaptivity context, usually gotten with TSGetAdapt()
3301c3436cfSJed Brown .  hmin - minimum time step
3311c3436cfSJed Brown -  hmax - maximum time step
3321c3436cfSJed Brown 
3331c3436cfSJed Brown    Options Database Keys:
3341c3436cfSJed Brown +  -ts_adapt_dt_min - minimum time step
3351c3436cfSJed Brown -  -ts_adapt_dt_max - maximum time step
3361c3436cfSJed Brown 
3371c3436cfSJed Brown    Level: intermediate
3381c3436cfSJed Brown 
3391c3436cfSJed Brown .seealso: TSAdapt
3401c3436cfSJed Brown @*/
3411c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax)
3421c3436cfSJed Brown {
3431c3436cfSJed Brown 
3441c3436cfSJed Brown   PetscFunctionBegin;
3454782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
3461c3436cfSJed Brown   if (hmin != PETSC_DECIDE) adapt->dt_min = hmin;
3471c3436cfSJed Brown   if (hmax != PETSC_DECIDE) adapt->dt_max = hmax;
3481c3436cfSJed Brown   PetscFunctionReturn(0);
3491c3436cfSJed Brown }
3501c3436cfSJed Brown 
3511c3436cfSJed Brown #undef __FUNCT__
35284df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetFromOptions"
353e55864a3SBarry Smith /*
35484df9cb4SJed Brown    TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options.
35584df9cb4SJed Brown 
35684df9cb4SJed Brown    Collective on TSAdapt
35784df9cb4SJed Brown 
35884df9cb4SJed Brown    Input Parameter:
35984df9cb4SJed Brown .  adapt - the TSAdapt context
36084df9cb4SJed Brown 
36184df9cb4SJed Brown    Options Database Keys:
36284df9cb4SJed Brown .  -ts_adapt_type <type> - basic
36384df9cb4SJed Brown 
36484df9cb4SJed Brown    Level: advanced
36584df9cb4SJed Brown 
36684df9cb4SJed Brown    Notes:
36784df9cb4SJed Brown    This function is automatically called by TSSetFromOptions()
36884df9cb4SJed Brown 
369552698daSJed Brown .keywords: TS, TSGetAdapt(), TSAdaptSetType()
37084df9cb4SJed Brown 
37184df9cb4SJed Brown .seealso: TSGetType()
372e55864a3SBarry Smith */
3738c34d3f5SBarry Smith PetscErrorCode  TSAdaptSetFromOptions(PetscOptions *PetscOptionsObject,TSAdapt adapt)
37484df9cb4SJed Brown {
37584df9cb4SJed Brown   PetscErrorCode ierr;
37684df9cb4SJed Brown   char           type[256] = TSADAPTBASIC;
3771c3436cfSJed Brown   PetscBool      set,flg;
37884df9cb4SJed Brown 
37984df9cb4SJed Brown   PetscFunctionBegin;
3804782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
38184df9cb4SJed Brown   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
38284df9cb4SJed Brown   * function can only be called from inside TSSetFromOptions_GL()  */
383e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");CHKERRQ(ierr);
384a264d7a6SBarry Smith   ierr = PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList,
3858caf3d72SBarry Smith                           ((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);CHKERRQ(ierr);
38684df9cb4SJed Brown   if (flg || !((PetscObject)adapt)->type_name) {
38784df9cb4SJed Brown     ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr);
38884df9cb4SJed Brown   }
3890298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",adapt->dt_min,&adapt->dt_min,NULL);CHKERRQ(ierr);
3900298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",adapt->dt_max,&adapt->dt_max,NULL);CHKERRQ(ierr);
3910298fd71SBarry 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);
3921c3436cfSJed Brown   ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
3937619abb3SShri   ierr = PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);CHKERRQ(ierr);
3947619abb3SShri   if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported");
3951c3436cfSJed Brown   if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);}
396e55864a3SBarry Smith   if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);}
39784df9cb4SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
39884df9cb4SJed Brown   PetscFunctionReturn(0);
39984df9cb4SJed Brown }
40084df9cb4SJed Brown 
40184df9cb4SJed Brown #undef __FUNCT__
40284df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidatesClear"
40384df9cb4SJed Brown /*@
40484df9cb4SJed Brown    TSAdaptCandidatesClear - clear any previously set candidate schemes
40584df9cb4SJed Brown 
40684df9cb4SJed Brown    Logically Collective
40784df9cb4SJed Brown 
40884df9cb4SJed Brown    Input Argument:
40984df9cb4SJed Brown .  adapt - adaptive controller
41084df9cb4SJed Brown 
41184df9cb4SJed Brown    Level: developer
41284df9cb4SJed Brown 
41384df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose()
41484df9cb4SJed Brown @*/
41584df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
41684df9cb4SJed Brown {
41784df9cb4SJed Brown   PetscErrorCode ierr;
41884df9cb4SJed Brown 
41984df9cb4SJed Brown   PetscFunctionBegin;
4204782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
42184df9cb4SJed Brown   ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr);
42284df9cb4SJed Brown   PetscFunctionReturn(0);
42384df9cb4SJed Brown }
42484df9cb4SJed Brown 
42584df9cb4SJed Brown #undef __FUNCT__
42684df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidateAdd"
42784df9cb4SJed Brown /*@C
42884df9cb4SJed Brown    TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from
42984df9cb4SJed Brown 
43084df9cb4SJed Brown    Logically Collective
43184df9cb4SJed Brown 
43284df9cb4SJed Brown    Input Arguments:
433552698daSJed Brown +  adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate()
43484df9cb4SJed Brown .  name - name of the candidate scheme to add
43584df9cb4SJed Brown .  order - order of the candidate scheme
43684df9cb4SJed Brown .  stageorder - stage order of the candidate scheme
4378d59e960SJed Brown .  ccfl - stability coefficient relative to explicit Euler, used for CFL constraints
43884df9cb4SJed Brown .  cost - relative measure of the amount of work required for the candidate scheme
43984df9cb4SJed Brown -  inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme
44084df9cb4SJed Brown 
44184df9cb4SJed Brown    Note:
44284df9cb4SJed Brown    This routine is not available in Fortran.
44384df9cb4SJed Brown 
44484df9cb4SJed Brown    Level: developer
44584df9cb4SJed Brown 
44684df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose()
44784df9cb4SJed Brown @*/
4488d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse)
44984df9cb4SJed Brown {
45084df9cb4SJed Brown   PetscInt c;
45184df9cb4SJed Brown 
45284df9cb4SJed Brown   PetscFunctionBegin;
45384df9cb4SJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
454ce94432eSBarry Smith   if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order);
45584df9cb4SJed Brown   if (inuse) {
456ce94432eSBarry 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()");
45784df9cb4SJed Brown     adapt->candidates.inuse_set = PETSC_TRUE;
45884df9cb4SJed Brown   }
4591c3436cfSJed Brown   /* first slot if this is the current scheme, otherwise the next available slot */
4601c3436cfSJed Brown   c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;
461bbd56ea5SKarl Rupp 
46284df9cb4SJed Brown   adapt->candidates.name[c]       = name;
46384df9cb4SJed Brown   adapt->candidates.order[c]      = order;
46484df9cb4SJed Brown   adapt->candidates.stageorder[c] = stageorder;
4658d59e960SJed Brown   adapt->candidates.ccfl[c]       = ccfl;
46684df9cb4SJed Brown   adapt->candidates.cost[c]       = cost;
46784df9cb4SJed Brown   adapt->candidates.n++;
46884df9cb4SJed Brown   PetscFunctionReturn(0);
46984df9cb4SJed Brown }
47084df9cb4SJed Brown 
47184df9cb4SJed Brown #undef __FUNCT__
4728d59e960SJed Brown #define __FUNCT__ "TSAdaptCandidatesGet"
4738d59e960SJed Brown /*@C
4748d59e960SJed Brown    TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost
4758d59e960SJed Brown 
4768d59e960SJed Brown    Not Collective
4778d59e960SJed Brown 
4788d59e960SJed Brown    Input Arguments:
4798d59e960SJed Brown .  adapt - time step adaptivity context
4808d59e960SJed Brown 
4818d59e960SJed Brown    Output Arguments:
4828d59e960SJed Brown +  n - number of candidate schemes, always at least 1
4838d59e960SJed Brown .  order - the order of each candidate scheme
4848d59e960SJed Brown .  stageorder - the stage order of each candidate scheme
4858d59e960SJed Brown .  ccfl - the CFL coefficient of each scheme
4868d59e960SJed Brown -  cost - the relative cost of each scheme
4878d59e960SJed Brown 
4888d59e960SJed Brown    Level: developer
4898d59e960SJed Brown 
4908d59e960SJed Brown    Note:
4918d59e960SJed Brown    The current scheme is always returned in the first slot
4928d59e960SJed Brown 
4938d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose()
4948d59e960SJed Brown @*/
4958d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost)
4968d59e960SJed Brown {
4978d59e960SJed Brown   PetscFunctionBegin;
4988d59e960SJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
4998d59e960SJed Brown   if (n) *n = adapt->candidates.n;
5008d59e960SJed Brown   if (order) *order = adapt->candidates.order;
5018d59e960SJed Brown   if (stageorder) *stageorder = adapt->candidates.stageorder;
5028d59e960SJed Brown   if (ccfl) *ccfl = adapt->candidates.ccfl;
5038d59e960SJed Brown   if (cost) *cost = adapt->candidates.cost;
5048d59e960SJed Brown   PetscFunctionReturn(0);
5058d59e960SJed Brown }
5068d59e960SJed Brown 
5078d59e960SJed Brown #undef __FUNCT__
50884df9cb4SJed Brown #define __FUNCT__ "TSAdaptChoose"
50984df9cb4SJed Brown /*@C
51084df9cb4SJed Brown    TSAdaptChoose - choose which method and step size to use for the next step
51184df9cb4SJed Brown 
51284df9cb4SJed Brown    Logically Collective
51384df9cb4SJed Brown 
51484df9cb4SJed Brown    Input Arguments:
51584df9cb4SJed Brown +  adapt - adaptive contoller
51684df9cb4SJed Brown -  h - current step size
51784df9cb4SJed Brown 
51884df9cb4SJed Brown    Output Arguments:
51984df9cb4SJed Brown +  next_sc - scheme to use for the next step
52084df9cb4SJed Brown .  next_h - step size to use for the next step
52184df9cb4SJed Brown -  accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size
52284df9cb4SJed Brown 
5231c3436cfSJed Brown    Note:
5241c3436cfSJed Brown    The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is
5251c3436cfSJed Brown    being retried after an initial rejection.
5261c3436cfSJed Brown 
52784df9cb4SJed Brown    Level: developer
52884df9cb4SJed Brown 
52984df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd()
53084df9cb4SJed Brown @*/
53184df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept)
53284df9cb4SJed Brown {
53384df9cb4SJed Brown   PetscErrorCode ierr;
5340b99f514SJed Brown   PetscReal      wlte = -1.0;
53584df9cb4SJed Brown 
53684df9cb4SJed Brown   PetscFunctionBegin;
53784df9cb4SJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
53884df9cb4SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
53984df9cb4SJed Brown   PetscValidIntPointer(next_sc,4);
54084df9cb4SJed Brown   PetscValidPointer(next_h,5);
54184df9cb4SJed Brown   PetscValidIntPointer(accept,6);
542ce94432eSBarry Smith   if (adapt->candidates.n < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"%D candidates have been registered",adapt->candidates.n);
543ce94432eSBarry Smith   if (!adapt->candidates.inuse_set) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"The current in-use scheme is not among the %D candidates",adapt->candidates.n);
5440b99f514SJed Brown   ierr = (*adapt->ops->choose)(adapt,ts,h,next_sc,next_h,accept,&wlte);CHKERRQ(ierr);
54549354f04SShri Abhyankar   if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
54649354f04SShri Abhyankar     /* Reduce time step if it overshoots max time */
54749354f04SShri Abhyankar     PetscReal max_time = ts->max_time;
54849354f04SShri Abhyankar     PetscReal next_dt  = 0.0;
54949354f04SShri Abhyankar     if (ts->ptime + ts->time_step + *next_h >= max_time) {
55049354f04SShri Abhyankar       next_dt = max_time - (ts->ptime + ts->time_step);
5518709b12bSShri Abhyankar       if (next_dt > PETSC_SMALL) *next_h = next_dt;
55249354f04SShri Abhyankar       else ts->reason = TS_CONVERGED_TIME;
55349354f04SShri Abhyankar     }
55449354f04SShri Abhyankar   }
555ce94432eSBarry Smith   if (*next_sc < 0 || adapt->candidates.n <= *next_sc) SETERRQ2(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Chosen scheme %D not in valid range 0..%D",*next_sc,adapt->candidates.n-1);
55657622a8eSBarry Smith   if (!(*next_h > 0.)) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h);
5571c3436cfSJed Brown 
5581c3436cfSJed Brown   if (adapt->monitor) {
5591c3436cfSJed Brown     ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
5600b99f514SJed Brown     if (wlte < 0) {
56148c19aefSLisandro 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,*next_sc,adapt->candidates.name[*next_sc],(double)*next_h);CHKERRQ(ierr);
5620b99f514SJed Brown     } else {
5636de24e2aSJed Brown       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,*next_sc,adapt->candidates.name[*next_sc],(double)*next_h);CHKERRQ(ierr);
5640b99f514SJed Brown     }
5651c3436cfSJed Brown     ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
5661c3436cfSJed Brown   }
56784df9cb4SJed Brown   PetscFunctionReturn(0);
56884df9cb4SJed Brown }
56984df9cb4SJed Brown 
57084df9cb4SJed Brown #undef __FUNCT__
57197335746SJed Brown #define __FUNCT__ "TSAdaptCheckStage"
57297335746SJed Brown /*@
57397335746SJed Brown    TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails)
57497335746SJed Brown 
57597335746SJed Brown    Collective
57697335746SJed Brown 
57797335746SJed Brown    Input Arguments:
57897335746SJed Brown +  adapt - adaptive controller context
57997335746SJed Brown -  ts - time stepper
58097335746SJed Brown 
58197335746SJed Brown    Output Arguments:
58297335746SJed Brown .  accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject
58397335746SJed Brown 
58497335746SJed Brown    Level: developer
58597335746SJed Brown 
58697335746SJed Brown .seealso:
58797335746SJed Brown @*/
58897335746SJed Brown PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscBool *accept)
58997335746SJed Brown {
59097335746SJed Brown   PetscErrorCode      ierr;
59197335746SJed Brown   SNES                snes;
59297335746SJed Brown   SNESConvergedReason snesreason;
59397335746SJed Brown 
59497335746SJed Brown   PetscFunctionBegin;
5954782b174SLisandro Dalcin   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
5964782b174SLisandro Dalcin   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
5974782b174SLisandro Dalcin   PetscValidIntPointer(accept,3);
59897335746SJed Brown   *accept = PETSC_TRUE;
59997335746SJed Brown   ierr    = TSGetSNES(ts,&snes);CHKERRQ(ierr);
60097335746SJed Brown   ierr    = SNESGetConvergedReason(snes,&snesreason);CHKERRQ(ierr);
60197335746SJed Brown   if (snesreason < 0) {
60297335746SJed Brown     PetscReal dt,new_dt;
60397335746SJed Brown     *accept = PETSC_FALSE;
60497335746SJed Brown     ierr    = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
6056de24e2aSJed Brown     if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) {
60697335746SJed Brown       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
60748c19aefSLisandro 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);
60897335746SJed Brown       if (adapt->monitor) {
60997335746SJed Brown         ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
61048c19aefSLisandro 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,dt,ts->num_snes_failures);CHKERRQ(ierr);
61197335746SJed Brown         ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
61297335746SJed Brown       }
61397335746SJed Brown     } else {
61497335746SJed Brown       new_dt = dt*adapt->scale_solve_failed;
61597335746SJed Brown       ierr   = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr);
61697335746SJed Brown       if (adapt->monitor) {
61797335746SJed Brown         ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
61897335746SJed 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);
61997335746SJed Brown         ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
62097335746SJed Brown       }
62197335746SJed Brown     }
62297335746SJed Brown   }
6230873808bSJed Brown   if (adapt->ops->checkstage) {ierr = (*adapt->ops->checkstage)(adapt,ts,accept);CHKERRQ(ierr);}
62497335746SJed Brown   PetscFunctionReturn(0);
62597335746SJed Brown }
62697335746SJed Brown 
6270873808bSJed Brown 
6280873808bSJed Brown 
62997335746SJed Brown #undef __FUNCT__
63084df9cb4SJed Brown #define __FUNCT__ "TSAdaptCreate"
63184df9cb4SJed Brown /*@
63284df9cb4SJed Brown   TSAdaptCreate - create an adaptive controller context for time stepping
63384df9cb4SJed Brown 
63484df9cb4SJed Brown   Collective on MPI_Comm
63584df9cb4SJed Brown 
63684df9cb4SJed Brown   Input Parameter:
63784df9cb4SJed Brown . comm - The communicator
63884df9cb4SJed Brown 
63984df9cb4SJed Brown   Output Parameter:
64084df9cb4SJed Brown . adapt - new TSAdapt object
64184df9cb4SJed Brown 
64284df9cb4SJed Brown   Level: developer
64384df9cb4SJed Brown 
64484df9cb4SJed Brown   Notes:
64584df9cb4SJed Brown   TSAdapt creation is handled by TS, so users should not need to call this function.
64684df9cb4SJed Brown 
64784df9cb4SJed Brown .keywords: TSAdapt, create
648552698daSJed Brown .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy()
64984df9cb4SJed Brown @*/
65084df9cb4SJed Brown PetscErrorCode  TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt)
65184df9cb4SJed Brown {
65284df9cb4SJed Brown   PetscErrorCode ierr;
65384df9cb4SJed Brown   TSAdapt        adapt;
65484df9cb4SJed Brown 
65584df9cb4SJed Brown   PetscFunctionBegin;
6563b3bcf4cSLisandro Dalcin   PetscValidPointer(inadapt,1);
6573b3bcf4cSLisandro Dalcin   *inadapt = NULL;
6583b3bcf4cSLisandro Dalcin   ierr = TSAdaptInitializePackage();CHKERRQ(ierr);
6593b3bcf4cSLisandro Dalcin 
66067c2884eSBarry Smith   ierr = PetscHeaderCreate(adapt,_p_TSAdapt,struct _TSAdaptOps,TSADAPT_CLASSID,"TSAdapt","General Linear adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr);
6611c3436cfSJed Brown 
6621c3436cfSJed Brown   adapt->dt_min             = 1e-20;
6631c3436cfSJed Brown   adapt->dt_max             = 1e50;
66497335746SJed Brown   adapt->scale_solve_failed = 0.25;
6657619abb3SShri   adapt->wnormtype          = NORM_2;
6661c3436cfSJed Brown 
66784df9cb4SJed Brown   *inadapt = adapt;
66884df9cb4SJed Brown   PetscFunctionReturn(0);
66984df9cb4SJed Brown }
670