xref: /petsc/src/ts/adapt/interface/tsadapt.c (revision 8cc058d9cd56c1ccb3be12a47760ddfc446aaffc)
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 
9*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt);
10*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_None(TSAdapt);
11*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt);
1284df9cb4SJed Brown 
1384df9cb4SJed Brown #undef __FUNCT__
1484df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegister"
1584df9cb4SJed Brown /*@C
1684df9cb4SJed Brown    TSAdaptRegister - see TSAdaptRegisterDynamic()
1784df9cb4SJed Brown 
1884df9cb4SJed Brown    Level: advanced
1984df9cb4SJed Brown @*/
2084df9cb4SJed Brown PetscErrorCode  TSAdaptRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(TSAdapt))
2184df9cb4SJed Brown {
2284df9cb4SJed Brown   PetscErrorCode ierr;
2384df9cb4SJed Brown   char           fullname[PETSC_MAX_PATH_LEN];
2484df9cb4SJed Brown 
2584df9cb4SJed Brown   PetscFunctionBegin;
26140e18c1SBarry Smith   ierr = PetscFunctionListConcat(path,name,fullname);CHKERRQ(ierr);
27140e18c1SBarry Smith   ierr = PetscFunctionListAdd(PETSC_COMM_WORLD,&TSAdaptList,sname,fullname,(void(*)(void))function);CHKERRQ(ierr);
2884df9cb4SJed Brown   PetscFunctionReturn(0);
2984df9cb4SJed Brown }
3084df9cb4SJed Brown 
3184df9cb4SJed Brown #undef __FUNCT__
3284df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegisterAll"
3384df9cb4SJed Brown /*@C
3484df9cb4SJed Brown   TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt
3584df9cb4SJed Brown 
3684df9cb4SJed Brown   Not Collective
3784df9cb4SJed Brown 
3884df9cb4SJed Brown   Level: advanced
3984df9cb4SJed Brown 
4084df9cb4SJed Brown .keywords: TSAdapt, register, all
4184df9cb4SJed Brown 
4284df9cb4SJed Brown .seealso: TSAdaptRegisterDestroy()
4384df9cb4SJed Brown @*/
4484df9cb4SJed Brown PetscErrorCode  TSAdaptRegisterAll(const char path[])
4584df9cb4SJed Brown {
4684df9cb4SJed Brown   PetscErrorCode ierr;
4784df9cb4SJed Brown 
4884df9cb4SJed Brown   PetscFunctionBegin;
4984df9cb4SJed Brown   ierr = TSAdaptRegisterDynamic(TSADAPTBASIC,path,"TSAdaptCreate_Basic",TSAdaptCreate_Basic);CHKERRQ(ierr);
50b0f836d7SJed Brown   ierr = TSAdaptRegisterDynamic(TSADAPTNONE, path,"TSAdaptCreate_None", TSAdaptCreate_None);CHKERRQ(ierr);
518d59e960SJed Brown   ierr = TSAdaptRegisterDynamic(TSADAPTCFL,  path,"TSAdaptCreate_CFL",  TSAdaptCreate_CFL);CHKERRQ(ierr);
5284df9cb4SJed Brown   PetscFunctionReturn(0);
5384df9cb4SJed Brown }
5484df9cb4SJed Brown 
5584df9cb4SJed Brown #undef __FUNCT__
5684df9cb4SJed Brown #define __FUNCT__ "TSAdaptFinalizePackage"
5784df9cb4SJed Brown /*@C
5884df9cb4SJed Brown   TSFinalizePackage - This function destroys everything in the TS package. It is
5984df9cb4SJed Brown   called from PetscFinalize().
6084df9cb4SJed Brown 
6184df9cb4SJed Brown   Level: developer
6284df9cb4SJed Brown 
6384df9cb4SJed Brown .keywords: Petsc, destroy, package
6484df9cb4SJed Brown .seealso: PetscFinalize()
6584df9cb4SJed Brown @*/
6684df9cb4SJed Brown PetscErrorCode  TSAdaptFinalizePackage(void)
6784df9cb4SJed Brown {
6884df9cb4SJed Brown   PetscFunctionBegin;
6984df9cb4SJed Brown   TSAdaptPackageInitialized = PETSC_FALSE;
7084df9cb4SJed Brown   TSAdaptRegisterAllCalled  = PETSC_FALSE;
710298fd71SBarry Smith   TSAdaptList               = NULL;
7284df9cb4SJed Brown   PetscFunctionReturn(0);
7384df9cb4SJed Brown }
7484df9cb4SJed Brown 
7584df9cb4SJed Brown #undef __FUNCT__
7684df9cb4SJed Brown #define __FUNCT__ "TSAdaptInitializePackage"
7784df9cb4SJed Brown /*@C
7884df9cb4SJed Brown   TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is
7984df9cb4SJed Brown   called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to
8084df9cb4SJed Brown   TSCreate_GL() when using static libraries.
8184df9cb4SJed Brown 
8284df9cb4SJed Brown   Input Parameter:
830298fd71SBarry Smith   path - The dynamic library path, or NULL
8484df9cb4SJed Brown 
8584df9cb4SJed Brown   Level: developer
8684df9cb4SJed Brown 
8784df9cb4SJed Brown .keywords: TSAdapt, initialize, package
8884df9cb4SJed Brown .seealso: PetscInitialize()
8984df9cb4SJed Brown @*/
9084df9cb4SJed Brown PetscErrorCode  TSAdaptInitializePackage(const char path[])
9184df9cb4SJed Brown {
9284df9cb4SJed Brown   PetscErrorCode ierr;
9384df9cb4SJed Brown 
9484df9cb4SJed Brown   PetscFunctionBegin;
9584df9cb4SJed Brown   if (TSAdaptPackageInitialized) PetscFunctionReturn(0);
9684df9cb4SJed Brown   TSAdaptPackageInitialized = PETSC_TRUE;
9784df9cb4SJed Brown   ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr);
9884df9cb4SJed Brown   ierr = TSAdaptRegisterAll(path);CHKERRQ(ierr);
9984df9cb4SJed Brown   ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr);
10084df9cb4SJed Brown   PetscFunctionReturn(0);
10184df9cb4SJed Brown }
10284df9cb4SJed Brown 
10384df9cb4SJed Brown #undef __FUNCT__
10484df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegisterDestroy"
10584df9cb4SJed Brown /*@C
10684df9cb4SJed Brown    TSAdaptRegisterDestroy - Frees the list of adaptivity schemes that were registered by TSAdaptRegister()/TSAdaptRegisterDynamic().
10784df9cb4SJed Brown 
10884df9cb4SJed Brown    Not Collective
10984df9cb4SJed Brown 
11084df9cb4SJed Brown    Level: advanced
11184df9cb4SJed Brown 
11284df9cb4SJed Brown .keywords: TSAdapt, register, destroy
11384df9cb4SJed Brown .seealso: TSAdaptRegister(), TSAdaptRegisterAll(), TSAdaptRegisterDynamic()
11484df9cb4SJed Brown @*/
11584df9cb4SJed Brown PetscErrorCode  TSAdaptRegisterDestroy(void)
11684df9cb4SJed Brown {
11784df9cb4SJed Brown   PetscErrorCode ierr;
11884df9cb4SJed Brown 
11984df9cb4SJed Brown   PetscFunctionBegin;
120140e18c1SBarry Smith   ierr = PetscFunctionListDestroy(&TSAdaptList);CHKERRQ(ierr);
12184df9cb4SJed Brown   TSAdaptRegisterAllCalled = PETSC_FALSE;
12284df9cb4SJed Brown   PetscFunctionReturn(0);
12384df9cb4SJed Brown }
12484df9cb4SJed Brown 
12584df9cb4SJed Brown 
12684df9cb4SJed Brown #undef __FUNCT__
12784df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetType"
12819fd82e9SBarry Smith PetscErrorCode  TSAdaptSetType(TSAdapt adapt,TSAdaptType type)
12984df9cb4SJed Brown {
13084df9cb4SJed Brown   PetscErrorCode ierr,(*r)(TSAdapt);
13184df9cb4SJed Brown 
13284df9cb4SJed Brown   PetscFunctionBegin;
133ce94432eSBarry Smith   ierr = PetscFunctionListFind(PetscObjectComm((PetscObject)adapt),TSAdaptList,type,PETSC_TRUE,(void(**)(void))&r);CHKERRQ(ierr);
13484df9cb4SJed Brown   if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type);
13584df9cb4SJed Brown   if (((PetscObject)adapt)->type_name) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);}
13684df9cb4SJed Brown   ierr = (*r)(adapt);CHKERRQ(ierr);
13784df9cb4SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr);
13884df9cb4SJed Brown   PetscFunctionReturn(0);
13984df9cb4SJed Brown }
14084df9cb4SJed Brown 
14184df9cb4SJed Brown #undef __FUNCT__
14284df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetOptionsPrefix"
14384df9cb4SJed Brown PetscErrorCode  TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[])
14484df9cb4SJed Brown {
14584df9cb4SJed Brown   PetscErrorCode ierr;
14684df9cb4SJed Brown 
14784df9cb4SJed Brown   PetscFunctionBegin;
14884df9cb4SJed Brown   ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr);
14984df9cb4SJed Brown   PetscFunctionReturn(0);
15084df9cb4SJed Brown }
15184df9cb4SJed Brown 
15284df9cb4SJed Brown #undef __FUNCT__
153ad6bc421SBarry Smith #define __FUNCT__ "TSAdaptLoad"
154ad6bc421SBarry Smith /*@C
155ad6bc421SBarry Smith   TSAdaptLoad - Loads a TSAdapt that has been stored in binary  with TSAdaptView().
156ad6bc421SBarry Smith 
157ad6bc421SBarry Smith   Collective on PetscViewer
158ad6bc421SBarry Smith 
159ad6bc421SBarry Smith   Input Parameters:
160ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or
161ad6bc421SBarry Smith            some related function before a call to TSAdaptLoad().
162ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or
163ad6bc421SBarry Smith            HDF5 file viewer, obtained from PetscViewerHDF5Open()
164ad6bc421SBarry Smith 
165ad6bc421SBarry Smith    Level: intermediate
166ad6bc421SBarry Smith 
167ad6bc421SBarry Smith   Notes:
168ad6bc421SBarry Smith    The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored.
169ad6bc421SBarry Smith 
170ad6bc421SBarry Smith   Notes for advanced users:
171ad6bc421SBarry Smith   Most users should not need to know the details of the binary storage
172ad6bc421SBarry Smith   format, since TSAdaptLoad() and TSAdaptView() completely hide these details.
173ad6bc421SBarry Smith   But for anyone who's interested, the standard binary matrix storage
174ad6bc421SBarry Smith   format is
175ad6bc421SBarry Smith .vb
176ad6bc421SBarry Smith      has not yet been determined
177ad6bc421SBarry Smith .ve
178ad6bc421SBarry Smith 
179ad6bc421SBarry Smith .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad()
180ad6bc421SBarry Smith @*/
181ad6bc421SBarry Smith PetscErrorCode  TSAdaptLoad(TSAdapt tsadapt, PetscViewer viewer)
182ad6bc421SBarry Smith {
183ad6bc421SBarry Smith   PetscErrorCode ierr;
184ad6bc421SBarry Smith   PetscBool      isbinary;
185ad6bc421SBarry Smith   char           type[256];
186ad6bc421SBarry Smith 
187ad6bc421SBarry Smith   PetscFunctionBegin;
188ad6bc421SBarry Smith   PetscValidHeaderSpecific(tsadapt,TSADAPT_CLASSID,1);
189ad6bc421SBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2);
190ad6bc421SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
191ad6bc421SBarry Smith   if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()");
192ad6bc421SBarry Smith 
193ad6bc421SBarry Smith   ierr = PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);CHKERRQ(ierr);
194ad6bc421SBarry Smith   ierr = TSAdaptSetType(tsadapt, type);CHKERRQ(ierr);
195ad6bc421SBarry Smith   if (tsadapt->ops->load) {
196ad6bc421SBarry Smith     ierr = (*tsadapt->ops->load)(tsadapt,viewer);CHKERRQ(ierr);
197ad6bc421SBarry Smith   }
198ad6bc421SBarry Smith   PetscFunctionReturn(0);
199ad6bc421SBarry Smith }
200ad6bc421SBarry Smith 
201ad6bc421SBarry Smith #undef __FUNCT__
20284df9cb4SJed Brown #define __FUNCT__ "TSAdaptView"
20384df9cb4SJed Brown PetscErrorCode  TSAdaptView(TSAdapt adapt,PetscViewer viewer)
20484df9cb4SJed Brown {
20584df9cb4SJed Brown   PetscErrorCode ierr;
206ad6bc421SBarry Smith   PetscBool      iascii,isbinary;
20784df9cb4SJed Brown 
20884df9cb4SJed Brown   PetscFunctionBegin;
209251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
210ad6bc421SBarry Smith   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr);
21184df9cb4SJed Brown   if (iascii) {
21284df9cb4SJed Brown     ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer,"TSAdapt Object");CHKERRQ(ierr);
21384df9cb4SJed Brown     ierr = PetscViewerASCIIPrintf(viewer,"  number of candidates %D\n",adapt->candidates.n);CHKERRQ(ierr);
21484df9cb4SJed Brown     if (adapt->ops->view) {
21584df9cb4SJed Brown       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
21684df9cb4SJed Brown       ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
21784df9cb4SJed Brown       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
21884df9cb4SJed Brown     }
219ad6bc421SBarry Smith   } else if (isbinary) {
220ad6bc421SBarry Smith     char type[256];
221ad6bc421SBarry Smith 
222ad6bc421SBarry Smith     /* need to save FILE_CLASS_ID for adapt class */
223ad6bc421SBarry Smith     ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr);
224ad6bc421SBarry Smith     ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr);
225bbd56ea5SKarl Rupp   } else if (adapt->ops->view) {
226f2c2a1b9SBarry Smith     ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr);
227f2c2a1b9SBarry Smith   }
22884df9cb4SJed Brown   PetscFunctionReturn(0);
22984df9cb4SJed Brown }
23084df9cb4SJed Brown 
23184df9cb4SJed Brown #undef __FUNCT__
23284df9cb4SJed Brown #define __FUNCT__ "TSAdaptDestroy"
23384df9cb4SJed Brown PetscErrorCode  TSAdaptDestroy(TSAdapt *adapt)
23484df9cb4SJed Brown {
23584df9cb4SJed Brown   PetscErrorCode ierr;
23684df9cb4SJed Brown 
23784df9cb4SJed Brown   PetscFunctionBegin;
23884df9cb4SJed Brown   if (!*adapt) PetscFunctionReturn(0);
23984df9cb4SJed Brown   PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1);
24084df9cb4SJed Brown   if (--((PetscObject)(*adapt))->refct > 0) {*adapt = 0; PetscFunctionReturn(0);}
24184df9cb4SJed Brown   if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);}
2421c3436cfSJed Brown   ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr);
24384df9cb4SJed Brown   ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr);
24484df9cb4SJed Brown   PetscFunctionReturn(0);
24584df9cb4SJed Brown }
24684df9cb4SJed Brown 
24784df9cb4SJed Brown #undef __FUNCT__
2481c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetMonitor"
2491c3436cfSJed Brown /*@
2501c3436cfSJed Brown    TSAdaptSetMonitor - Monitor the choices made by the adaptive controller
2511c3436cfSJed Brown 
2521c3436cfSJed Brown    Collective on TSAdapt
2531c3436cfSJed Brown 
2541c3436cfSJed Brown    Input Arguments:
2551c3436cfSJed Brown +  adapt - adaptive controller context
2561c3436cfSJed Brown -  flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable
2571c3436cfSJed Brown 
2581c3436cfSJed Brown    Level: intermediate
2591c3436cfSJed Brown 
2601c3436cfSJed Brown .seealso: TSAdaptChoose()
2611c3436cfSJed Brown @*/
2621c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg)
2631c3436cfSJed Brown {
2641c3436cfSJed Brown   PetscErrorCode ierr;
2651c3436cfSJed Brown 
2661c3436cfSJed Brown   PetscFunctionBegin;
2671c3436cfSJed Brown   if (flg) {
268ce94432eSBarry Smith     if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);}
2691c3436cfSJed Brown   } else {
2701c3436cfSJed Brown     ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr);
2711c3436cfSJed Brown   }
2721c3436cfSJed Brown   PetscFunctionReturn(0);
2731c3436cfSJed Brown }
2741c3436cfSJed Brown 
2751c3436cfSJed Brown #undef __FUNCT__
2760873808bSJed Brown #define __FUNCT__ "TSAdaptSetCheckStage"
2770873808bSJed Brown /*@C
2780873808bSJed Brown    TSAdaptSetCheckStage - set a callback to check convergence for a stage
2790873808bSJed Brown 
2800873808bSJed Brown    Logically collective on TSAdapt
2810873808bSJed Brown 
2820873808bSJed Brown    Input Arguments:
2830873808bSJed Brown +  adapt - adaptive controller context
2840873808bSJed Brown -  func - stage check function
2850873808bSJed Brown 
2860873808bSJed Brown    Arguments of func:
2870873808bSJed Brown $  PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept)
2880873808bSJed Brown 
2890873808bSJed Brown +  adapt - adaptive controller context
2900873808bSJed Brown .  ts - time stepping context
2910873808bSJed Brown -  accept - pending choice of whether to accept, can be modified by this routine
2920873808bSJed Brown 
2930873808bSJed Brown    Level: advanced
2940873808bSJed Brown 
2950873808bSJed Brown .seealso: TSAdaptChoose()
2960873808bSJed Brown @*/
2970873808bSJed Brown PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscBool*))
2980873808bSJed Brown {
2990873808bSJed Brown 
3000873808bSJed Brown   PetscFunctionBegin;
3010873808bSJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
3020873808bSJed Brown   adapt->ops->checkstage = func;
3030873808bSJed Brown   PetscFunctionReturn(0);
3040873808bSJed Brown }
3050873808bSJed Brown 
3060873808bSJed Brown #undef __FUNCT__
3071c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetStepLimits"
3081c3436cfSJed Brown /*@
3091c3436cfSJed Brown    TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller
3101c3436cfSJed Brown 
3111c3436cfSJed Brown    Logically Collective
3121c3436cfSJed Brown 
3131c3436cfSJed Brown    Input Arguments:
314ad6bc421SBarry Smith +  adapt - time step adaptivity context, usually gotten with TSGetTSAdapt()
3151c3436cfSJed Brown .  hmin - minimum time step
3161c3436cfSJed Brown -  hmax - maximum time step
3171c3436cfSJed Brown 
3181c3436cfSJed Brown    Options Database Keys:
3191c3436cfSJed Brown +  -ts_adapt_dt_min - minimum time step
3201c3436cfSJed Brown -  -ts_adapt_dt_max - maximum time step
3211c3436cfSJed Brown 
3221c3436cfSJed Brown    Level: intermediate
3231c3436cfSJed Brown 
3241c3436cfSJed Brown .seealso: TSAdapt
3251c3436cfSJed Brown @*/
3261c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax)
3271c3436cfSJed Brown {
3281c3436cfSJed Brown 
3291c3436cfSJed Brown   PetscFunctionBegin;
3301c3436cfSJed Brown   if (hmin != PETSC_DECIDE) adapt->dt_min = hmin;
3311c3436cfSJed Brown   if (hmax != PETSC_DECIDE) adapt->dt_max = hmax;
3321c3436cfSJed Brown   PetscFunctionReturn(0);
3331c3436cfSJed Brown }
3341c3436cfSJed Brown 
3351c3436cfSJed Brown #undef __FUNCT__
33684df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetFromOptions"
33784df9cb4SJed Brown /*@
33884df9cb4SJed Brown    TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options.
33984df9cb4SJed Brown 
34084df9cb4SJed Brown    Collective on TSAdapt
34184df9cb4SJed Brown 
34284df9cb4SJed Brown    Input Parameter:
34384df9cb4SJed Brown .  adapt - the TSAdapt context
34484df9cb4SJed Brown 
34584df9cb4SJed Brown    Options Database Keys:
34684df9cb4SJed Brown .  -ts_adapt_type <type> - basic
34784df9cb4SJed Brown 
34884df9cb4SJed Brown    Level: advanced
34984df9cb4SJed Brown 
35084df9cb4SJed Brown    Notes:
35184df9cb4SJed Brown    This function is automatically called by TSSetFromOptions()
35284df9cb4SJed Brown 
353ad6bc421SBarry Smith .keywords: TS, TSGetTSAdapt(), TSAdaptSetType()
35484df9cb4SJed Brown 
35584df9cb4SJed Brown .seealso: TSGetType()
35684df9cb4SJed Brown @*/
35784df9cb4SJed Brown PetscErrorCode  TSAdaptSetFromOptions(TSAdapt adapt)
35884df9cb4SJed Brown {
35984df9cb4SJed Brown   PetscErrorCode ierr;
36084df9cb4SJed Brown   char           type[256] = TSADAPTBASIC;
3611c3436cfSJed Brown   PetscBool      set,flg;
36284df9cb4SJed Brown 
36384df9cb4SJed Brown   PetscFunctionBegin;
36484df9cb4SJed Brown   /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this
36584df9cb4SJed Brown   * function can only be called from inside TSSetFromOptions_GL()  */
36684df9cb4SJed Brown   ierr = PetscOptionsHead("TS Adaptivity options");CHKERRQ(ierr);
36784df9cb4SJed Brown   ierr = PetscOptionsList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList,
3688caf3d72SBarry Smith                           ((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);CHKERRQ(ierr);
36984df9cb4SJed Brown   if (flg || !((PetscObject)adapt)->type_name) {
37084df9cb4SJed Brown     ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr);
37184df9cb4SJed Brown   }
37284df9cb4SJed Brown   if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(adapt);CHKERRQ(ierr);}
3730298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",adapt->dt_min,&adapt->dt_min,NULL);CHKERRQ(ierr);
3740298fd71SBarry Smith   ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",adapt->dt_max,&adapt->dt_max,NULL);CHKERRQ(ierr);
3750298fd71SBarry 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);
3761c3436cfSJed Brown   ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr);
3771c3436cfSJed Brown   if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);}
37884df9cb4SJed Brown   ierr = PetscOptionsTail();CHKERRQ(ierr);
37984df9cb4SJed Brown   PetscFunctionReturn(0);
38084df9cb4SJed Brown }
38184df9cb4SJed Brown 
38284df9cb4SJed Brown #undef __FUNCT__
38384df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidatesClear"
38484df9cb4SJed Brown /*@
38584df9cb4SJed Brown    TSAdaptCandidatesClear - clear any previously set candidate schemes
38684df9cb4SJed Brown 
38784df9cb4SJed Brown    Logically Collective
38884df9cb4SJed Brown 
38984df9cb4SJed Brown    Input Argument:
39084df9cb4SJed Brown .  adapt - adaptive controller
39184df9cb4SJed Brown 
39284df9cb4SJed Brown    Level: developer
39384df9cb4SJed Brown 
39484df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose()
39584df9cb4SJed Brown @*/
39684df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt)
39784df9cb4SJed Brown {
39884df9cb4SJed Brown   PetscErrorCode ierr;
39984df9cb4SJed Brown 
40084df9cb4SJed Brown   PetscFunctionBegin;
40184df9cb4SJed Brown   ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr);
40284df9cb4SJed Brown   PetscFunctionReturn(0);
40384df9cb4SJed Brown }
40484df9cb4SJed Brown 
40584df9cb4SJed Brown #undef __FUNCT__
40684df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidateAdd"
40784df9cb4SJed Brown /*@C
40884df9cb4SJed Brown    TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from
40984df9cb4SJed Brown 
41084df9cb4SJed Brown    Logically Collective
41184df9cb4SJed Brown 
41284df9cb4SJed Brown    Input Arguments:
413ad6bc421SBarry Smith +  adapt - time step adaptivity context, obtained with TSGetTSAdapt() or TSAdaptCreate()
41484df9cb4SJed Brown .  name - name of the candidate scheme to add
41584df9cb4SJed Brown .  order - order of the candidate scheme
41684df9cb4SJed Brown .  stageorder - stage order of the candidate scheme
4178d59e960SJed Brown .  ccfl - stability coefficient relative to explicit Euler, used for CFL constraints
41884df9cb4SJed Brown .  cost - relative measure of the amount of work required for the candidate scheme
41984df9cb4SJed Brown -  inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme
42084df9cb4SJed Brown 
42184df9cb4SJed Brown    Note:
42284df9cb4SJed Brown    This routine is not available in Fortran.
42384df9cb4SJed Brown 
42484df9cb4SJed Brown    Level: developer
42584df9cb4SJed Brown 
42684df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose()
42784df9cb4SJed Brown @*/
4288d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse)
42984df9cb4SJed Brown {
43084df9cb4SJed Brown   PetscInt c;
43184df9cb4SJed Brown 
43284df9cb4SJed Brown   PetscFunctionBegin;
43384df9cb4SJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
434ce94432eSBarry Smith   if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order);
43584df9cb4SJed Brown   if (inuse) {
436ce94432eSBarry 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()");
43784df9cb4SJed Brown     adapt->candidates.inuse_set = PETSC_TRUE;
43884df9cb4SJed Brown   }
4391c3436cfSJed Brown   /* first slot if this is the current scheme, otherwise the next available slot */
4401c3436cfSJed Brown   c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n;
441bbd56ea5SKarl Rupp 
44284df9cb4SJed Brown   adapt->candidates.name[c]       = name;
44384df9cb4SJed Brown   adapt->candidates.order[c]      = order;
44484df9cb4SJed Brown   adapt->candidates.stageorder[c] = stageorder;
4458d59e960SJed Brown   adapt->candidates.ccfl[c]       = ccfl;
44684df9cb4SJed Brown   adapt->candidates.cost[c]       = cost;
44784df9cb4SJed Brown   adapt->candidates.n++;
44884df9cb4SJed Brown   PetscFunctionReturn(0);
44984df9cb4SJed Brown }
45084df9cb4SJed Brown 
45184df9cb4SJed Brown #undef __FUNCT__
4528d59e960SJed Brown #define __FUNCT__ "TSAdaptCandidatesGet"
4538d59e960SJed Brown /*@C
4548d59e960SJed Brown    TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost
4558d59e960SJed Brown 
4568d59e960SJed Brown    Not Collective
4578d59e960SJed Brown 
4588d59e960SJed Brown    Input Arguments:
4598d59e960SJed Brown .  adapt - time step adaptivity context
4608d59e960SJed Brown 
4618d59e960SJed Brown    Output Arguments:
4628d59e960SJed Brown +  n - number of candidate schemes, always at least 1
4638d59e960SJed Brown .  order - the order of each candidate scheme
4648d59e960SJed Brown .  stageorder - the stage order of each candidate scheme
4658d59e960SJed Brown .  ccfl - the CFL coefficient of each scheme
4668d59e960SJed Brown -  cost - the relative cost of each scheme
4678d59e960SJed Brown 
4688d59e960SJed Brown    Level: developer
4698d59e960SJed Brown 
4708d59e960SJed Brown    Note:
4718d59e960SJed Brown    The current scheme is always returned in the first slot
4728d59e960SJed Brown 
4738d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose()
4748d59e960SJed Brown @*/
4758d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost)
4768d59e960SJed Brown {
4778d59e960SJed Brown   PetscFunctionBegin;
4788d59e960SJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
4798d59e960SJed Brown   if (n) *n = adapt->candidates.n;
4808d59e960SJed Brown   if (order) *order = adapt->candidates.order;
4818d59e960SJed Brown   if (stageorder) *stageorder = adapt->candidates.stageorder;
4828d59e960SJed Brown   if (ccfl) *ccfl = adapt->candidates.ccfl;
4838d59e960SJed Brown   if (cost) *cost = adapt->candidates.cost;
4848d59e960SJed Brown   PetscFunctionReturn(0);
4858d59e960SJed Brown }
4868d59e960SJed Brown 
4878d59e960SJed Brown #undef __FUNCT__
48884df9cb4SJed Brown #define __FUNCT__ "TSAdaptChoose"
48984df9cb4SJed Brown /*@C
49084df9cb4SJed Brown    TSAdaptChoose - choose which method and step size to use for the next step
49184df9cb4SJed Brown 
49284df9cb4SJed Brown    Logically Collective
49384df9cb4SJed Brown 
49484df9cb4SJed Brown    Input Arguments:
49584df9cb4SJed Brown +  adapt - adaptive contoller
49684df9cb4SJed Brown -  h - current step size
49784df9cb4SJed Brown 
49884df9cb4SJed Brown    Output Arguments:
49984df9cb4SJed Brown +  next_sc - scheme to use for the next step
50084df9cb4SJed Brown .  next_h - step size to use for the next step
50184df9cb4SJed Brown -  accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size
50284df9cb4SJed Brown 
5031c3436cfSJed Brown    Note:
5041c3436cfSJed Brown    The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is
5051c3436cfSJed Brown    being retried after an initial rejection.
5061c3436cfSJed Brown 
50784df9cb4SJed Brown    Level: developer
50884df9cb4SJed Brown 
50984df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd()
51084df9cb4SJed Brown @*/
51184df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept)
51284df9cb4SJed Brown {
51384df9cb4SJed Brown   PetscErrorCode ierr;
5140b99f514SJed Brown   PetscReal      wlte = -1.0;
51584df9cb4SJed Brown 
51684df9cb4SJed Brown   PetscFunctionBegin;
51784df9cb4SJed Brown   PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1);
51884df9cb4SJed Brown   PetscValidHeaderSpecific(ts,TS_CLASSID,2);
51984df9cb4SJed Brown   PetscValidIntPointer(next_sc,4);
52084df9cb4SJed Brown   PetscValidPointer(next_h,5);
52184df9cb4SJed Brown   PetscValidIntPointer(accept,6);
522ce94432eSBarry Smith   if (adapt->candidates.n < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"%D candidates have been registered",adapt->candidates.n);
523ce94432eSBarry 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);
5240b99f514SJed Brown   ierr = (*adapt->ops->choose)(adapt,ts,h,next_sc,next_h,accept,&wlte);CHKERRQ(ierr);
52549354f04SShri Abhyankar   if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) {
52649354f04SShri Abhyankar     /* Reduce time step if it overshoots max time */
52749354f04SShri Abhyankar     PetscReal max_time = ts->max_time;
52849354f04SShri Abhyankar     PetscReal next_dt  = 0.0;
52949354f04SShri Abhyankar     if (ts->ptime + ts->time_step + *next_h >= max_time) {
53049354f04SShri Abhyankar       next_dt = max_time - (ts->ptime + ts->time_step);
5318709b12bSShri Abhyankar       if (next_dt > PETSC_SMALL) *next_h = next_dt;
53249354f04SShri Abhyankar       else ts->reason = TS_CONVERGED_TIME;
53349354f04SShri Abhyankar     }
53449354f04SShri Abhyankar   }
535ce94432eSBarry 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);
536ce94432eSBarry Smith   if (!(*next_h > 0.)) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %G must be positive",*next_h);
5371c3436cfSJed Brown 
5381c3436cfSJed Brown   if (adapt->monitor) {
5391c3436cfSJed Brown     ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
5400b99f514SJed Brown     if (wlte < 0) {
5416de24e2aSJed 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);
5420b99f514SJed Brown     } else {
5436de24e2aSJed 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);
5440b99f514SJed Brown     }
5451c3436cfSJed Brown     ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
5461c3436cfSJed Brown   }
54784df9cb4SJed Brown   PetscFunctionReturn(0);
54884df9cb4SJed Brown }
54984df9cb4SJed Brown 
55084df9cb4SJed Brown #undef __FUNCT__
55197335746SJed Brown #define __FUNCT__ "TSAdaptCheckStage"
55297335746SJed Brown /*@
55397335746SJed Brown    TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails)
55497335746SJed Brown 
55597335746SJed Brown    Collective
55697335746SJed Brown 
55797335746SJed Brown    Input Arguments:
55897335746SJed Brown +  adapt - adaptive controller context
55997335746SJed Brown -  ts - time stepper
56097335746SJed Brown 
56197335746SJed Brown    Output Arguments:
56297335746SJed Brown .  accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject
56397335746SJed Brown 
56497335746SJed Brown    Level: developer
56597335746SJed Brown 
56697335746SJed Brown .seealso:
56797335746SJed Brown @*/
56897335746SJed Brown PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscBool *accept)
56997335746SJed Brown {
57097335746SJed Brown   PetscErrorCode      ierr;
57197335746SJed Brown   SNES                snes;
57297335746SJed Brown   SNESConvergedReason snesreason;
57397335746SJed Brown 
57497335746SJed Brown   PetscFunctionBegin;
57597335746SJed Brown   *accept = PETSC_TRUE;
57697335746SJed Brown   ierr    = TSGetSNES(ts,&snes);CHKERRQ(ierr);
57797335746SJed Brown   ierr    = SNESGetConvergedReason(snes,&snesreason);CHKERRQ(ierr);
57897335746SJed Brown   if (snesreason < 0) {
57997335746SJed Brown     PetscReal dt,new_dt;
58097335746SJed Brown     *accept = PETSC_FALSE;
58197335746SJed Brown     ierr    = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
5826de24e2aSJed Brown     if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) {
58397335746SJed Brown       ts->reason = TS_DIVERGED_NONLINEAR_SOLVE;
58497335746SJed 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);
58597335746SJed Brown       if (adapt->monitor) {
58697335746SJed Brown         ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
58797335746SJed 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);
58897335746SJed Brown         ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
58997335746SJed Brown       }
59097335746SJed Brown     } else {
59197335746SJed Brown       new_dt = dt*adapt->scale_solve_failed;
59297335746SJed Brown       ierr   = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr);
59397335746SJed Brown       if (adapt->monitor) {
59497335746SJed Brown         ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
59597335746SJed 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);
59697335746SJed Brown         ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr);
59797335746SJed Brown       }
59897335746SJed Brown     }
59997335746SJed Brown   }
6000873808bSJed Brown   if (adapt->ops->checkstage) {ierr = (*adapt->ops->checkstage)(adapt,ts,accept);CHKERRQ(ierr);}
60197335746SJed Brown   PetscFunctionReturn(0);
60297335746SJed Brown }
60397335746SJed Brown 
6040873808bSJed Brown 
6050873808bSJed Brown 
60697335746SJed Brown #undef __FUNCT__
60784df9cb4SJed Brown #define __FUNCT__ "TSAdaptCreate"
60884df9cb4SJed Brown /*@
60984df9cb4SJed Brown   TSAdaptCreate - create an adaptive controller context for time stepping
61084df9cb4SJed Brown 
61184df9cb4SJed Brown   Collective on MPI_Comm
61284df9cb4SJed Brown 
61384df9cb4SJed Brown   Input Parameter:
61484df9cb4SJed Brown . comm - The communicator
61584df9cb4SJed Brown 
61684df9cb4SJed Brown   Output Parameter:
61784df9cb4SJed Brown . adapt - new TSAdapt object
61884df9cb4SJed Brown 
61984df9cb4SJed Brown   Level: developer
62084df9cb4SJed Brown 
62184df9cb4SJed Brown   Notes:
62284df9cb4SJed Brown   TSAdapt creation is handled by TS, so users should not need to call this function.
62384df9cb4SJed Brown 
62484df9cb4SJed Brown .keywords: TSAdapt, create
625ad6bc421SBarry Smith .seealso: TSGetTSAdapt(), TSAdaptSetType(), TSAdaptDestroy()
62684df9cb4SJed Brown @*/
62784df9cb4SJed Brown PetscErrorCode  TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt)
62884df9cb4SJed Brown {
62984df9cb4SJed Brown   PetscErrorCode ierr;
63084df9cb4SJed Brown   TSAdapt        adapt;
63184df9cb4SJed Brown 
63284df9cb4SJed Brown   PetscFunctionBegin;
63384df9cb4SJed Brown   *inadapt = 0;
63467c2884eSBarry Smith   ierr     = PetscHeaderCreate(adapt,_p_TSAdapt,struct _TSAdaptOps,TSADAPT_CLASSID,"TSAdapt","General Linear adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr);
6351c3436cfSJed Brown 
6361c3436cfSJed Brown   adapt->dt_min             = 1e-20;
6371c3436cfSJed Brown   adapt->dt_max             = 1e50;
63897335746SJed Brown   adapt->scale_solve_failed = 0.25;
6391c3436cfSJed Brown 
64084df9cb4SJed Brown   *inadapt = adapt;
64184df9cb4SJed Brown   PetscFunctionReturn(0);
64284df9cb4SJed Brown }
643