xref: /petsc/src/ts/adapt/interface/tsadapt.c (revision 140e18c1ab5d6e82cb97419c1b69bc894c9508d2)
184df9cb4SJed Brown 
2b45d2f2cSJed Brown #include <petsc-private/tsimpl.h> /*I  "petscts.h" I*/
384df9cb4SJed Brown 
4*140e18c1SBarry Smith static PetscFunctionList TSAdaptList;
584df9cb4SJed Brown static PetscBool         TSAdaptPackageInitialized;
684df9cb4SJed Brown static PetscBool         TSAdaptRegisterAllCalled;
784df9cb4SJed Brown static PetscClassId      TSADAPT_CLASSID;
884df9cb4SJed Brown 
984df9cb4SJed Brown EXTERN_C_BEGIN
1084df9cb4SJed Brown PetscErrorCode  TSAdaptCreate_Basic(TSAdapt);
11b0f836d7SJed Brown PetscErrorCode  TSAdaptCreate_None(TSAdapt);
128d59e960SJed Brown PetscErrorCode  TSAdaptCreate_CFL(TSAdapt);
1384df9cb4SJed Brown EXTERN_C_END
1484df9cb4SJed Brown 
1584df9cb4SJed Brown #undef __FUNCT__
1684df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegister"
1784df9cb4SJed Brown /*@C
1884df9cb4SJed Brown    TSAdaptRegister - see TSAdaptRegisterDynamic()
1984df9cb4SJed Brown 
2084df9cb4SJed Brown    Level: advanced
2184df9cb4SJed Brown @*/
2284df9cb4SJed Brown PetscErrorCode  TSAdaptRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(TSAdapt))
2384df9cb4SJed Brown {
2484df9cb4SJed Brown   PetscErrorCode ierr;
2584df9cb4SJed Brown   char           fullname[PETSC_MAX_PATH_LEN];
2684df9cb4SJed Brown 
2784df9cb4SJed Brown   PetscFunctionBegin;
28*140e18c1SBarry Smith   ierr = PetscFunctionListConcat(path,name,fullname);CHKERRQ(ierr);
29*140e18c1SBarry Smith   ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&TSAdaptList,sname,fullname,(void(*)(void))function);CHKERRQ(ierr);
3084df9cb4SJed Brown   PetscFunctionReturn(0);
3184df9cb4SJed Brown }
3284df9cb4SJed Brown 
3384df9cb4SJed Brown #undef __FUNCT__
3484df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegisterAll"
3584df9cb4SJed Brown /*@C
3684df9cb4SJed Brown   TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt
3784df9cb4SJed Brown 
3884df9cb4SJed Brown   Not Collective
3984df9cb4SJed Brown 
4084df9cb4SJed Brown   Level: advanced
4184df9cb4SJed Brown 
4284df9cb4SJed Brown .keywords: TSAdapt, register, all
4384df9cb4SJed Brown 
4484df9cb4SJed Brown .seealso: TSAdaptRegisterDestroy()
4584df9cb4SJed Brown @*/
4684df9cb4SJed Brown PetscErrorCode  TSAdaptRegisterAll(const char path[])
4784df9cb4SJed Brown {
4884df9cb4SJed Brown   PetscErrorCode ierr;
4984df9cb4SJed Brown 
5084df9cb4SJed Brown   PetscFunctionBegin;
5184df9cb4SJed Brown   ierr = TSAdaptRegisterDynamic(TSADAPTBASIC,path,"TSAdaptCreate_Basic",TSAdaptCreate_Basic);CHKERRQ(ierr);
52b0f836d7SJed Brown   ierr = TSAdaptRegisterDynamic(TSADAPTNONE, path,"TSAdaptCreate_None", TSAdaptCreate_None);CHKERRQ(ierr);
538d59e960SJed Brown   ierr = TSAdaptRegisterDynamic(TSADAPTCFL,  path,"TSAdaptCreate_CFL",  TSAdaptCreate_CFL);CHKERRQ(ierr);
5484df9cb4SJed Brown   PetscFunctionReturn(0);
5584df9cb4SJed Brown }
5684df9cb4SJed Brown 
5784df9cb4SJed Brown #undef __FUNCT__
5884df9cb4SJed Brown #define __FUNCT__ "TSAdaptFinalizePackage"
5984df9cb4SJed Brown /*@C
6084df9cb4SJed Brown   TSFinalizePackage - This function destroys everything in the TS package. It is
6184df9cb4SJed Brown   called from PetscFinalize().
6284df9cb4SJed Brown 
6384df9cb4SJed Brown   Level: developer
6484df9cb4SJed Brown 
6584df9cb4SJed Brown .keywords: Petsc, destroy, package
6684df9cb4SJed Brown .seealso: PetscFinalize()
6784df9cb4SJed Brown @*/
6884df9cb4SJed Brown PetscErrorCode  TSAdaptFinalizePackage(void)
6984df9cb4SJed Brown {
7084df9cb4SJed Brown   PetscFunctionBegin;
7184df9cb4SJed Brown   TSAdaptPackageInitialized = PETSC_FALSE;
7284df9cb4SJed Brown   TSAdaptRegisterAllCalled  = PETSC_FALSE;
7384df9cb4SJed Brown   TSAdaptList               = PETSC_NULL;
7484df9cb4SJed Brown   PetscFunctionReturn(0);
7584df9cb4SJed Brown }
7684df9cb4SJed Brown 
7784df9cb4SJed Brown #undef __FUNCT__
7884df9cb4SJed Brown #define __FUNCT__ "TSAdaptInitializePackage"
7984df9cb4SJed Brown /*@C
8084df9cb4SJed Brown   TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is
8184df9cb4SJed Brown   called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to
8284df9cb4SJed Brown   TSCreate_GL() when using static libraries.
8384df9cb4SJed Brown 
8484df9cb4SJed Brown   Input Parameter:
8584df9cb4SJed Brown   path - The dynamic library path, or PETSC_NULL
8684df9cb4SJed Brown 
8784df9cb4SJed Brown   Level: developer
8884df9cb4SJed Brown 
8984df9cb4SJed Brown .keywords: TSAdapt, initialize, package
9084df9cb4SJed Brown .seealso: PetscInitialize()
9184df9cb4SJed Brown @*/
9284df9cb4SJed Brown PetscErrorCode  TSAdaptInitializePackage(const char path[])
9384df9cb4SJed Brown {
9484df9cb4SJed Brown   PetscErrorCode ierr;
9584df9cb4SJed Brown 
9684df9cb4SJed Brown   PetscFunctionBegin;
9784df9cb4SJed Brown   if (TSAdaptPackageInitialized) PetscFunctionReturn(0);
9884df9cb4SJed Brown   TSAdaptPackageInitialized = PETSC_TRUE;
9984df9cb4SJed Brown   ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr);
10084df9cb4SJed Brown   ierr = TSAdaptRegisterAll(path);CHKERRQ(ierr);
10184df9cb4SJed Brown   ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr);
10284df9cb4SJed Brown   PetscFunctionReturn(0);
10384df9cb4SJed Brown }
10484df9cb4SJed Brown 
10584df9cb4SJed Brown #undef __FUNCT__
10684df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegisterDestroy"
10784df9cb4SJed Brown /*@C
10884df9cb4SJed Brown    TSAdaptRegisterDestroy - Frees the list of adaptivity schemes that were registered by TSAdaptRegister()/TSAdaptRegisterDynamic().
10984df9cb4SJed Brown 
11084df9cb4SJed Brown    Not Collective
11184df9cb4SJed Brown 
11284df9cb4SJed Brown    Level: advanced
11384df9cb4SJed Brown 
11484df9cb4SJed Brown .keywords: TSAdapt, register, destroy
11584df9cb4SJed Brown .seealso: TSAdaptRegister(), TSAdaptRegisterAll(), TSAdaptRegisterDynamic()
11684df9cb4SJed Brown @*/
11784df9cb4SJed Brown PetscErrorCode  TSAdaptRegisterDestroy(void)
11884df9cb4SJed Brown {
11984df9cb4SJed Brown   PetscErrorCode ierr;
12084df9cb4SJed Brown 
12184df9cb4SJed Brown   PetscFunctionBegin;
122*140e18c1SBarry Smith   ierr = PetscFunctionListDestroy(&TSAdaptList);CHKERRQ(ierr);
12384df9cb4SJed Brown   TSAdaptRegisterAllCalled = PETSC_FALSE;
12484df9cb4SJed Brown   PetscFunctionReturn(0);
12584df9cb4SJed Brown }
12684df9cb4SJed Brown 
12784df9cb4SJed Brown 
12884df9cb4SJed Brown #undef __FUNCT__
12984df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetType"
13019fd82e9SBarry Smith PetscErrorCode  TSAdaptSetType(TSAdapt adapt,TSAdaptType type)
13184df9cb4SJed Brown {
13284df9cb4SJed Brown   PetscErrorCode ierr,(*r)(TSAdapt);
13384df9cb4SJed Brown 
13484df9cb4SJed Brown   PetscFunctionBegin;
135*140e18c1SBarry Smith   ierr = PetscFunctionListFind(((PetscObject)adapt)->comm,TSAdaptList,type,PETSC_TRUE,(void(**)(void))&r);CHKERRQ(ierr);
13684df9cb4SJed Brown   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type);
13784df9cb4SJed Brown   if (((PetscObject)adapt)->type_name) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);}
13884df9cb4SJed Brown   ierr = (*r)(adapt);CHKERRQ(ierr);
13984df9cb4SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr);
14084df9cb4SJed Brown   PetscFunctionReturn(0);
14184df9cb4SJed Brown }
14284df9cb4SJed Brown 
14384df9cb4SJed Brown #undef __FUNCT__
14484df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetOptionsPrefix"
14584df9cb4SJed Brown PetscErrorCode  TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[])
14684df9cb4SJed Brown {
14784df9cb4SJed Brown   PetscErrorCode ierr;
14884df9cb4SJed Brown 
14984df9cb4SJed Brown   PetscFunctionBegin;
15084df9cb4SJed Brown   ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr);
15184df9cb4SJed Brown   PetscFunctionReturn(0);
15284df9cb4SJed Brown }
15384df9cb4SJed Brown 
15484df9cb4SJed Brown #undef __FUNCT__
155ad6bc421SBarry Smith #define __FUNCT__ "TSAdaptLoad"
156ad6bc421SBarry Smith /*@C
157ad6bc421SBarry Smith   TSAdaptLoad - Loads a TSAdapt that has been stored in binary  with TSAdaptView().
158ad6bc421SBarry Smith 
159ad6bc421SBarry Smith   Collective on PetscViewer
160ad6bc421SBarry Smith 
161ad6bc421SBarry Smith   Input Parameters:
162ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or
163ad6bc421SBarry Smith            some related function before a call to TSAdaptLoad().
164ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
165ad6bc421SBarry Smith            HDF5 file viewer, obtained from PetscViewerHDF5Open()
166ad6bc421SBarry Smith 
167ad6bc421SBarry Smith    Level: intermediate
168ad6bc421SBarry Smith 
169ad6bc421SBarry Smith   Notes:
170ad6bc421SBarry Smith    The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored.
171ad6bc421SBarry Smith 
172ad6bc421SBarry Smith   Notes for advanced users:
173ad6bc421SBarry Smith   Most users should not need to know the details of the binary storage
174ad6bc421SBarry Smith   format, since TSAdaptLoad() and TSAdaptView() completely hide these details.
175ad6bc421SBarry Smith   But for anyone who's interested, the standard binary matrix storage
176ad6bc421SBarry Smith   format is
177ad6bc421SBarry Smith .vb
178ad6bc421SBarry Smith      has not yet been determined
179ad6bc421SBarry Smith .ve
180ad6bc421SBarry Smith 
181ad6bc421SBarry Smith .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad()
182ad6bc421SBarry Smith @*/
183ad6bc421SBarry Smith PetscErrorCode  TSAdaptLoad(TSAdapt tsadapt, PetscViewer viewer)
184ad6bc421SBarry Smith {
185ad6bc421SBarry Smith   PetscErrorCode ierr;
186ad6bc421SBarry Smith   PetscBool      isbinary;
187ad6bc421SBarry Smith   char           type[256];
188ad6bc421SBarry Smith 
189ad6bc421SBarry Smith   PetscFunctionBegin;
190ad6bc421SBarry Smith   PetscValidHeaderSpecific(tsadapt,TSADAPT_CLASSID,1);
191ad6bc421SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
192ad6bc421SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
193ad6bc421SBarry Smith   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
194ad6bc421SBarry Smith 
195ad6bc421SBarry Smith   ierr = PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);CHKERRQ(ierr);
196ad6bc421SBarry Smith   ierr = TSAdaptSetType(tsadapt, type);CHKERRQ(ierr);
197ad6bc421SBarry Smith   if (tsadapt->ops->load) {
198ad6bc421SBarry Smith     ierr = (*tsadapt->ops->load)(tsadapt,viewer);CHKERRQ(ierr);
199ad6bc421SBarry Smith   }
200ad6bc421SBarry Smith   PetscFunctionReturn(0);
201ad6bc421SBarry Smith }
202ad6bc421SBarry Smith 
203ad6bc421SBarry Smith #undef __FUNCT__
20484df9cb4SJed Brown #define __FUNCT__ "TSAdaptView"
20584df9cb4SJed Brown PetscErrorCode  TSAdaptView(TSAdapt adapt,PetscViewer viewer)
20684df9cb4SJed Brown {
20784df9cb4SJed Brown   PetscErrorCode ierr;
208ad6bc421SBarry Smith   PetscBool      iascii,isbinary;
20984df9cb4SJed Brown 
21084df9cb4SJed Brown   PetscFunctionBegin;
211251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
212ad6bc421SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
21384df9cb4SJed Brown   if (iascii) {
21484df9cb4SJed Brown     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer,"TSAdapt Object");CHKERRQ(ierr);
21584df9cb4SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  number of candidates %D\n",adapt->candidates.n);CHKERRQ(ierr);
21684df9cb4SJed Brown     if (adapt->ops->view) {
21784df9cb4SJed Brown       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
21884df9cb4SJed Brown       ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
21984df9cb4SJed Brown       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
22084df9cb4SJed Brown     }
221ad6bc421SBarry Smith   } else if (isbinary) {
222ad6bc421SBarry Smith     char type[256];
223ad6bc421SBarry Smith 
224ad6bc421SBarry Smith     /* need to save FILE_CLASS_ID for adapt class */
225ad6bc421SBarry Smith     ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr);
226ad6bc421SBarry Smith     ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr);
22784df9cb4SJed Brown   } else {
228f2c2a1b9SBarry Smith     if (adapt->ops->view) {
229f2c2a1b9SBarry Smith       ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
230f2c2a1b9SBarry Smith     }
23184df9cb4SJed Brown   }
23284df9cb4SJed Brown   PetscFunctionReturn(0);
23384df9cb4SJed Brown }
23484df9cb4SJed Brown 
23584df9cb4SJed Brown #undef __FUNCT__
23684df9cb4SJed Brown #define __FUNCT__ "TSAdaptDestroy"
23784df9cb4SJed Brown PetscErrorCode  TSAdaptDestroy(TSAdapt *adapt)
23884df9cb4SJed Brown {
23984df9cb4SJed Brown   PetscErrorCode ierr;
24084df9cb4SJed Brown 
24184df9cb4SJed Brown   PetscFunctionBegin;
24284df9cb4SJed Brown   if (!*adapt) PetscFunctionReturn(0);
24384df9cb4SJed Brown   PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1);
24484df9cb4SJed Brown   if (--((PetscObject)(*adapt))->refct > 0) {*adapt = 0; PetscFunctionReturn(0);}
24584df9cb4SJed Brown   if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);}
2461c3436cfSJed Brown   ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr);
24784df9cb4SJed Brown   ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr);
24884df9cb4SJed Brown   PetscFunctionReturn(0);
24984df9cb4SJed Brown }
25084df9cb4SJed Brown 
25184df9cb4SJed Brown #undef __FUNCT__
2521c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetMonitor"
2531c3436cfSJed Brown /*@
2541c3436cfSJed Brown    TSAdaptSetMonitor - Monitor the choices made by the adaptive controller
2551c3436cfSJed Brown 
2561c3436cfSJed Brown    Collective on TSAdapt
2571c3436cfSJed Brown 
2581c3436cfSJed Brown    Input Arguments:
2591c3436cfSJed Brown +  adapt - adaptive controller context
2601c3436cfSJed Brown -  flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable
2611c3436cfSJed Brown 
2621c3436cfSJed Brown    Level: intermediate
2631c3436cfSJed Brown 
2641c3436cfSJed Brown .seealso: TSAdaptChoose()
2651c3436cfSJed Brown @*/
2661c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg)
2671c3436cfSJed Brown {
2681c3436cfSJed Brown   PetscErrorCode ierr;
2691c3436cfSJed Brown 
2701c3436cfSJed Brown   PetscFunctionBegin;
2711c3436cfSJed Brown   if (flg) {
2721c3436cfSJed Brown     if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(((PetscObject)adapt)->comm,"stdout",&adapt->monitor);CHKERRQ(ierr);}
2731c3436cfSJed Brown   } else {
2741c3436cfSJed Brown     ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr);
2751c3436cfSJed Brown   }
2761c3436cfSJed Brown   PetscFunctionReturn(0);
2771c3436cfSJed Brown }
2781c3436cfSJed Brown 
2791c3436cfSJed Brown #undef __FUNCT__
2800873808bSJed Brown #define __FUNCT__ "TSAdaptSetCheckStage"
2810873808bSJed Brown /*@C
2820873808bSJed Brown    TSAdaptSetCheckStage - set a callback to check convergence for a stage
2830873808bSJed Brown 
2840873808bSJed Brown    Logically collective on TSAdapt
2850873808bSJed Brown 
2860873808bSJed Brown    Input Arguments:
2870873808bSJed Brown +  adapt - adaptive controller context
2880873808bSJed Brown -  func - stage check function
2890873808bSJed Brown 
2900873808bSJed Brown    Arguments of func:
2910873808bSJed Brown $  PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept)
2920873808bSJed Brown 
2930873808bSJed Brown +  adapt - adaptive controller context
2940873808bSJed Brown .  ts - time stepping context
2950873808bSJed Brown -  accept - pending choice of whether to accept, can be modified by this routine
2960873808bSJed Brown 
2970873808bSJed Brown    Level: advanced
2980873808bSJed Brown 
2990873808bSJed Brown .seealso: TSAdaptChoose()
3000873808bSJed Brown @*/
3010873808bSJed Brown PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscBool*))
3020873808bSJed Brown {
3030873808bSJed Brown 
3040873808bSJed Brown   PetscFunctionBegin;
3050873808bSJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
3060873808bSJed Brown   adapt->ops->checkstage = func;
3070873808bSJed Brown   PetscFunctionReturn(0);
3080873808bSJed Brown }
3090873808bSJed Brown 
3100873808bSJed Brown #undef __FUNCT__
3111c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetStepLimits"
3121c3436cfSJed Brown /*@
3131c3436cfSJed Brown    TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller
3141c3436cfSJed Brown 
3151c3436cfSJed Brown    Logically Collective
3161c3436cfSJed Brown 
3171c3436cfSJed Brown    Input Arguments:
318ad6bc421SBarry Smith +  adapt - time step adaptivity context, usually gotten with TSGetTSAdapt()
3191c3436cfSJed Brown .  hmin - minimum time step
3201c3436cfSJed Brown -  hmax - maximum time step
3211c3436cfSJed Brown 
3221c3436cfSJed Brown    Options Database Keys:
3231c3436cfSJed Brown +  -ts_adapt_dt_min - minimum time step
3241c3436cfSJed Brown -  -ts_adapt_dt_max - maximum time step
3251c3436cfSJed Brown 
3261c3436cfSJed Brown    Level: intermediate
3271c3436cfSJed Brown 
3281c3436cfSJed Brown .seealso: TSAdapt
3291c3436cfSJed Brown @*/
3301c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax)
3311c3436cfSJed Brown {
3321c3436cfSJed Brown 
3331c3436cfSJed Brown   PetscFunctionBegin;
3341c3436cfSJed Brown   if (hmin != PETSC_DECIDE) adapt->dt_min = hmin;
3351c3436cfSJed Brown   if (hmax != PETSC_DECIDE) adapt->dt_max = hmax;
3361c3436cfSJed Brown   PetscFunctionReturn(0);
3371c3436cfSJed Brown }
3381c3436cfSJed Brown 
3391c3436cfSJed Brown #undef __FUNCT__
34084df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetFromOptions"
34184df9cb4SJed Brown /*@
34284df9cb4SJed Brown    TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options.
34384df9cb4SJed Brown 
34484df9cb4SJed Brown    Collective on TSAdapt
34584df9cb4SJed Brown 
34684df9cb4SJed Brown    Input Parameter:
34784df9cb4SJed Brown .  adapt - the TSAdapt context
34884df9cb4SJed Brown 
34984df9cb4SJed Brown    Options Database Keys:
35084df9cb4SJed Brown .  -ts_adapt_type <type> - basic
35184df9cb4SJed Brown 
35284df9cb4SJed Brown    Level: advanced
35384df9cb4SJed Brown 
35484df9cb4SJed Brown    Notes:
35584df9cb4SJed Brown    This function is automatically called by TSSetFromOptions()
35684df9cb4SJed Brown 
357ad6bc421SBarry Smith .keywords: TS, TSGetTSAdapt(), TSAdaptSetType()
35884df9cb4SJed Brown 
35984df9cb4SJed Brown .seealso: TSGetType()
36084df9cb4SJed Brown @*/
36184df9cb4SJed Brown PetscErrorCode  TSAdaptSetFromOptions(TSAdapt adapt)
36284df9cb4SJed Brown {
36384df9cb4SJed Brown   PetscErrorCode ierr;
36484df9cb4SJed Brown   char           type[256] = TSADAPTBASIC;
3651c3436cfSJed Brown   PetscBool      set,flg;
36684df9cb4SJed Brown 
36784df9cb4SJed Brown   PetscFunctionBegin;
36884df9cb4SJed Brown   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
36984df9cb4SJed Brown   * function can only be called from inside TSSetFromOptions_GL()  */
37084df9cb4SJed Brown   ierr = PetscOptionsHead("TS Adaptivity options");CHKERRQ(ierr);
37184df9cb4SJed Brown   ierr = PetscOptionsList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList,
3728caf3d72SBarry Smith                           ((PetscObject)adapt)->type_name?((PetscObject)adapt)->type_name:type,type,sizeof(type),&flg);CHKERRQ(ierr);
37384df9cb4SJed Brown   if (flg || !((PetscObject)adapt)->type_name) {
37484df9cb4SJed Brown     ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr);
37584df9cb4SJed Brown   }
37684df9cb4SJed Brown   if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(adapt);CHKERRQ(ierr);}
3771c3436cfSJed Brown   ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",adapt->dt_min,&adapt->dt_min,PETSC_NULL);CHKERRQ(ierr);
3781c3436cfSJed Brown   ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",adapt->dt_max,&adapt->dt_max,PETSC_NULL);CHKERRQ(ierr);
37997335746SJed Brown   ierr = PetscOptionsReal("-ts_adapt_scale_solve_failed","Scale step by this factor if solve fails","",adapt->scale_solve_failed,&adapt->scale_solve_failed,PETSC_NULL);CHKERRQ(ierr);
3801c3436cfSJed Brown   ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
3811c3436cfSJed Brown   if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);}
38284df9cb4SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
38384df9cb4SJed Brown   PetscFunctionReturn(0);
38484df9cb4SJed Brown }
38584df9cb4SJed Brown 
38684df9cb4SJed Brown #undef __FUNCT__
38784df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidatesClear"
38884df9cb4SJed Brown /*@
38984df9cb4SJed Brown    TSAdaptCandidatesClear - clear any previously set candidate schemes
39084df9cb4SJed Brown 
39184df9cb4SJed Brown    Logically Collective
39284df9cb4SJed Brown 
39384df9cb4SJed Brown    Input Argument:
39484df9cb4SJed Brown .  adapt - adaptive controller
39584df9cb4SJed Brown 
39684df9cb4SJed Brown    Level: developer
39784df9cb4SJed Brown 
39884df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose()
39984df9cb4SJed Brown @*/
40084df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
40184df9cb4SJed Brown {
40284df9cb4SJed Brown   PetscErrorCode ierr;
40384df9cb4SJed Brown 
40484df9cb4SJed Brown   PetscFunctionBegin;
40584df9cb4SJed Brown   ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr);
40684df9cb4SJed Brown   PetscFunctionReturn(0);
40784df9cb4SJed Brown }
40884df9cb4SJed Brown 
40984df9cb4SJed Brown #undef __FUNCT__
41084df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidateAdd"
41184df9cb4SJed Brown /*@C
41284df9cb4SJed Brown    TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from
41384df9cb4SJed Brown 
41484df9cb4SJed Brown    Logically Collective
41584df9cb4SJed Brown 
41684df9cb4SJed Brown    Input Arguments:
417ad6bc421SBarry Smith +  adapt - time step adaptivity context, obtained with TSGetTSAdapt() or TSAdaptCreate()
41884df9cb4SJed Brown .  name - name of the candidate scheme to add
41984df9cb4SJed Brown .  order - order of the candidate scheme
42084df9cb4SJed Brown .  stageorder - stage order of the candidate scheme
4218d59e960SJed Brown .  ccfl - stability coefficient relative to explicit Euler, used for CFL constraints
42284df9cb4SJed Brown .  cost - relative measure of the amount of work required for the candidate scheme
42384df9cb4SJed Brown -  inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme
42484df9cb4SJed Brown 
42584df9cb4SJed Brown    Note:
42684df9cb4SJed Brown    This routine is not available in Fortran.
42784df9cb4SJed Brown 
42884df9cb4SJed Brown    Level: developer
42984df9cb4SJed Brown 
43084df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose()
43184df9cb4SJed Brown @*/
4328d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse)
43384df9cb4SJed Brown {
43484df9cb4SJed Brown   PetscInt c;
43584df9cb4SJed Brown 
43684df9cb4SJed Brown   PetscFunctionBegin;
43784df9cb4SJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
43884df9cb4SJed Brown   if (order < 1) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order);
43984df9cb4SJed Brown   if (inuse) {
44084df9cb4SJed Brown     if (adapt->candidates.inuse_set) SETERRQ(((PetscObject)adapt)->comm,PETSC_ERR_ARG_WRONGSTATE,"Cannot set the inuse method twice, maybe forgot to call TSAdaptCandidatesClear()");
44184df9cb4SJed Brown     adapt->candidates.inuse_set = PETSC_TRUE;
44284df9cb4SJed Brown   }
4431c3436cfSJed Brown   /* first slot if this is the current scheme, otherwise the next available slot */
4441c3436cfSJed Brown   c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;
44584df9cb4SJed Brown   adapt->candidates.name[c]         = name;
44684df9cb4SJed Brown   adapt->candidates.order[c]        = order;
44784df9cb4SJed Brown   adapt->candidates.stageorder[c]   = stageorder;
4488d59e960SJed Brown   adapt->candidates.ccfl[c]         = ccfl;
44984df9cb4SJed Brown   adapt->candidates.cost[c]         = cost;
45084df9cb4SJed Brown   adapt->candidates.n++;
45184df9cb4SJed Brown   PetscFunctionReturn(0);
45284df9cb4SJed Brown }
45384df9cb4SJed Brown 
45484df9cb4SJed Brown #undef __FUNCT__
4558d59e960SJed Brown #define __FUNCT__ "TSAdaptCandidatesGet"
4568d59e960SJed Brown /*@C
4578d59e960SJed Brown    TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost
4588d59e960SJed Brown 
4598d59e960SJed Brown    Not Collective
4608d59e960SJed Brown 
4618d59e960SJed Brown    Input Arguments:
4628d59e960SJed Brown .  adapt - time step adaptivity context
4638d59e960SJed Brown 
4648d59e960SJed Brown    Output Arguments:
4658d59e960SJed Brown +  n - number of candidate schemes, always at least 1
4668d59e960SJed Brown .  order - the order of each candidate scheme
4678d59e960SJed Brown .  stageorder - the stage order of each candidate scheme
4688d59e960SJed Brown .  ccfl - the CFL coefficient of each scheme
4698d59e960SJed Brown -  cost - the relative cost of each scheme
4708d59e960SJed Brown 
4718d59e960SJed Brown    Level: developer
4728d59e960SJed Brown 
4738d59e960SJed Brown    Note:
4748d59e960SJed Brown    The current scheme is always returned in the first slot
4758d59e960SJed Brown 
4768d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose()
4778d59e960SJed Brown @*/
4788d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost)
4798d59e960SJed Brown {
4808d59e960SJed Brown   PetscFunctionBegin;
4818d59e960SJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
4828d59e960SJed Brown   if (n) *n = adapt->candidates.n;
4838d59e960SJed Brown   if (order) *order = adapt->candidates.order;
4848d59e960SJed Brown   if (stageorder) *stageorder = adapt->candidates.stageorder;
4858d59e960SJed Brown   if (ccfl) *ccfl = adapt->candidates.ccfl;
4868d59e960SJed Brown   if (cost) *cost = adapt->candidates.cost;
4878d59e960SJed Brown   PetscFunctionReturn(0);
4888d59e960SJed Brown }
4898d59e960SJed Brown 
4908d59e960SJed Brown #undef __FUNCT__
49184df9cb4SJed Brown #define __FUNCT__ "TSAdaptChoose"
49284df9cb4SJed Brown /*@C
49384df9cb4SJed Brown    TSAdaptChoose - choose which method and step size to use for the next step
49484df9cb4SJed Brown 
49584df9cb4SJed Brown    Logically Collective
49684df9cb4SJed Brown 
49784df9cb4SJed Brown    Input Arguments:
49884df9cb4SJed Brown +  adapt - adaptive contoller
49984df9cb4SJed Brown -  h - current step size
50084df9cb4SJed Brown 
50184df9cb4SJed Brown    Output Arguments:
50284df9cb4SJed Brown +  next_sc - scheme to use for the next step
50384df9cb4SJed Brown .  next_h - step size to use for the next step
50484df9cb4SJed Brown -  accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size
50584df9cb4SJed Brown 
5061c3436cfSJed Brown    Note:
5071c3436cfSJed Brown    The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is
5081c3436cfSJed Brown    being retried after an initial rejection.
5091c3436cfSJed Brown 
51084df9cb4SJed Brown    Level: developer
51184df9cb4SJed Brown 
51284df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd()
51384df9cb4SJed Brown @*/
51484df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept)
51584df9cb4SJed Brown {
51684df9cb4SJed Brown   PetscErrorCode ierr;
5170b99f514SJed Brown   PetscReal wlte = -1.0;
51884df9cb4SJed Brown 
51984df9cb4SJed Brown   PetscFunctionBegin;
52084df9cb4SJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
52184df9cb4SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
52284df9cb4SJed Brown   PetscValidIntPointer(next_sc,4);
52384df9cb4SJed Brown   PetscValidPointer(next_h,5);
52484df9cb4SJed Brown   PetscValidIntPointer(accept,6);
52584df9cb4SJed Brown   if (adapt->candidates.n < 1) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_WRONGSTATE,"%D candidates have been registered",adapt->candidates.n);
52684df9cb4SJed Brown   if (!adapt->candidates.inuse_set) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_WRONGSTATE,"The current in-use scheme is not among the %D candidates",adapt->candidates.n);
5270b99f514SJed Brown   ierr = (*adapt->ops->choose)(adapt,ts,h,next_sc,next_h,accept,&wlte);CHKERRQ(ierr);
52849354f04SShri Abhyankar   if(*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
52949354f04SShri Abhyankar     /* Reduce time step if it overshoots max time */
53049354f04SShri Abhyankar     PetscReal max_time=ts->max_time;
53149354f04SShri Abhyankar     PetscReal next_dt=0.0;
53249354f04SShri Abhyankar     if(ts->ptime + ts->time_step + *next_h >= max_time) {
53349354f04SShri Abhyankar       next_dt = max_time - (ts->ptime + ts->time_step);
53449354f04SShri Abhyankar       if(next_dt != 0.0) *next_h = next_dt;
53549354f04SShri Abhyankar       else ts->reason = TS_CONVERGED_TIME;
53649354f04SShri Abhyankar     }
53749354f04SShri Abhyankar   }
5381c3436cfSJed Brown   if (*next_sc < 0 || adapt->candidates.n <= *next_sc) SETERRQ2(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Chosen scheme %D not in valid range 0..%D",*next_sc,adapt->candidates.n-1);
5391c3436cfSJed Brown   if (!(*next_h > 0.)) SETERRQ1(((PetscObject)adapt)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %G must be positive",*next_h);
5401c3436cfSJed Brown 
5411c3436cfSJed Brown   if (adapt->monitor) {
5421c3436cfSJed Brown     ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
5430b99f514SJed Brown     if (wlte < 0) {
5446de24e2aSJed Brown       ierr = PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt '%s': step %3D %s t=%-11g+%10.3e family='%s' scheme=%D:'%s' dt=%-10g\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);
5450b99f514SJed Brown     } else {
5466de24e2aSJed 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);
5470b99f514SJed Brown     }
5481c3436cfSJed Brown     ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
5491c3436cfSJed Brown   }
55084df9cb4SJed Brown   PetscFunctionReturn(0);
55184df9cb4SJed Brown }
55284df9cb4SJed Brown 
55384df9cb4SJed Brown #undef __FUNCT__
55497335746SJed Brown #define __FUNCT__ "TSAdaptCheckStage"
55597335746SJed Brown /*@
55697335746SJed Brown    TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails)
55797335746SJed Brown 
55897335746SJed Brown    Collective
55997335746SJed Brown 
56097335746SJed Brown    Input Arguments:
56197335746SJed Brown +  adapt - adaptive controller context
56297335746SJed Brown -  ts - time stepper
56397335746SJed Brown 
56497335746SJed Brown    Output Arguments:
56597335746SJed Brown .  accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject
56697335746SJed Brown 
56797335746SJed Brown    Level: developer
56897335746SJed Brown 
56997335746SJed Brown .seealso:
57097335746SJed Brown @*/
57197335746SJed Brown PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscBool *accept)
57297335746SJed Brown {
57397335746SJed Brown   PetscErrorCode      ierr;
57497335746SJed Brown   SNES                snes;
57597335746SJed Brown   SNESConvergedReason snesreason;
57697335746SJed Brown 
57797335746SJed Brown   PetscFunctionBegin;
57897335746SJed Brown   *accept = PETSC_TRUE;
57997335746SJed Brown   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
58097335746SJed Brown   ierr = SNESGetConvergedReason(snes,&snesreason);CHKERRQ(ierr);
58197335746SJed Brown   if (snesreason < 0) {
58297335746SJed Brown     PetscReal dt,new_dt;
58397335746SJed Brown     *accept = PETSC_FALSE;
58497335746SJed Brown     ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
5856de24e2aSJed Brown     if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) {
58697335746SJed Brown       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
58797335746SJed Brown       ierr = PetscInfo2(ts,"Step=%D, nonlinear solve solve failures %D greater than current TS allowed, stopping solve\n",ts->steps,ts->num_snes_failures);CHKERRQ(ierr);
58897335746SJed Brown       if (adapt->monitor) {
58997335746SJed Brown         ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
59097335746SJed Brown         ierr = PetscViewerASCIIPrintf(adapt->monitor,"    TSAdapt '%s': step %3D stage rejected t=%-11g+%10.3e, %D failures exceeds current TS allowed\n",((PetscObject)adapt)->type_name,ts->steps,(double)ts->ptime,dt,ts->num_snes_failures);CHKERRQ(ierr);
59197335746SJed Brown         ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
59297335746SJed Brown       }
59397335746SJed Brown     } else {
59497335746SJed Brown       new_dt = dt*adapt->scale_solve_failed;
59597335746SJed Brown       ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr);
59697335746SJed Brown       if (adapt->monitor) {
59797335746SJed Brown         ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
59897335746SJed 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);
59997335746SJed Brown         ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
60097335746SJed Brown       }
60197335746SJed Brown     }
60297335746SJed Brown   }
6030873808bSJed Brown   if (adapt->ops->checkstage) {ierr = (*adapt->ops->checkstage)(adapt,ts,accept);CHKERRQ(ierr);}
60497335746SJed Brown   PetscFunctionReturn(0);
60597335746SJed Brown }
60697335746SJed Brown 
6070873808bSJed Brown 
6080873808bSJed Brown 
60997335746SJed Brown #undef __FUNCT__
61084df9cb4SJed Brown #define __FUNCT__ "TSAdaptCreate"
61184df9cb4SJed Brown /*@
61284df9cb4SJed Brown   TSAdaptCreate - create an adaptive controller context for time stepping
61384df9cb4SJed Brown 
61484df9cb4SJed Brown   Collective on MPI_Comm
61584df9cb4SJed Brown 
61684df9cb4SJed Brown   Input Parameter:
61784df9cb4SJed Brown . comm - The communicator
61884df9cb4SJed Brown 
61984df9cb4SJed Brown   Output Parameter:
62084df9cb4SJed Brown . adapt - new TSAdapt object
62184df9cb4SJed Brown 
62284df9cb4SJed Brown   Level: developer
62384df9cb4SJed Brown 
62484df9cb4SJed Brown   Notes:
62584df9cb4SJed Brown   TSAdapt creation is handled by TS, so users should not need to call this function.
62684df9cb4SJed Brown 
62784df9cb4SJed Brown .keywords: TSAdapt, create
628ad6bc421SBarry Smith .seealso: TSGetTSAdapt(), TSAdaptSetType(), TSAdaptDestroy()
62984df9cb4SJed Brown @*/
63084df9cb4SJed Brown PetscErrorCode  TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt)
63184df9cb4SJed Brown {
63284df9cb4SJed Brown   PetscErrorCode ierr;
63384df9cb4SJed Brown   TSAdapt adapt;
63484df9cb4SJed Brown 
63584df9cb4SJed Brown   PetscFunctionBegin;
63684df9cb4SJed Brown   *inadapt = 0;
63784df9cb4SJed Brown   ierr = PetscHeaderCreate(adapt,_p_TSAdapt,struct _TSAdaptOps,TSADAPT_CLASSID,0,"TSAdapt","General Linear adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr);
6381c3436cfSJed Brown 
6391c3436cfSJed Brown   adapt->dt_min             = 1e-20;
6401c3436cfSJed Brown   adapt->dt_max             = 1e50;
64197335746SJed Brown   adapt->scale_solve_failed = 0.25;
6421c3436cfSJed Brown 
64384df9cb4SJed Brown   *inadapt = adapt;
64484df9cb4SJed Brown   PetscFunctionReturn(0);
64584df9cb4SJed Brown }
646