184df9cb4SJed Brown 2af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 384df9cb4SJed Brown 41b9b13dfSLisandro Dalcin PetscClassId TSADAPT_CLASSID; 51b9b13dfSLisandro Dalcin 6140e18c1SBarry Smith static PetscFunctionList TSAdaptList; 784df9cb4SJed Brown static PetscBool TSAdaptPackageInitialized; 884df9cb4SJed Brown static PetscBool TSAdaptRegisterAllCalled; 984df9cb4SJed Brown 10726095e4SEmil Constantinescu PETSC_EXTERN PetscErrorCode TSAdaptCreate_GLEE(TSAdapt); 118cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_None(TSAdapt); 121566a47fSLisandro Dalcin PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt); 138cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt); 1484df9cb4SJed Brown 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 /*@C 5384df9cb4SJed Brown TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt 5484df9cb4SJed Brown 5584df9cb4SJed Brown Not Collective 5684df9cb4SJed Brown 5784df9cb4SJed Brown Level: advanced 5884df9cb4SJed Brown 5984df9cb4SJed Brown .keywords: TSAdapt, register, all 6084df9cb4SJed Brown 6184df9cb4SJed Brown .seealso: TSAdaptRegisterDestroy() 6284df9cb4SJed Brown @*/ 63607a6623SBarry Smith PetscErrorCode TSAdaptRegisterAll(void) 6484df9cb4SJed Brown { 6584df9cb4SJed Brown PetscErrorCode ierr; 6684df9cb4SJed Brown 6784df9cb4SJed Brown PetscFunctionBegin; 680f51fdf8SToby Isaac if (TSAdaptRegisterAllCalled) PetscFunctionReturn(0); 690f51fdf8SToby Isaac TSAdaptRegisterAllCalled = PETSC_TRUE; 70726095e4SEmil Constantinescu ierr = TSAdaptRegister(TSADAPTGLEE ,TSAdaptCreate_GLEE);CHKERRQ(ierr); 71bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None);CHKERRQ(ierr); 721566a47fSLisandro Dalcin ierr = TSAdaptRegister(TSADAPTBASIC,TSAdaptCreate_Basic);CHKERRQ(ierr); 73bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTCFL, TSAdaptCreate_CFL);CHKERRQ(ierr); 7484df9cb4SJed Brown PetscFunctionReturn(0); 7584df9cb4SJed Brown } 7684df9cb4SJed Brown 7784df9cb4SJed Brown /*@C 78560360afSLisandro Dalcin TSAdaptFinalizePackage - This function destroys everything in the TS package. It is 7984df9cb4SJed Brown called from PetscFinalize(). 8084df9cb4SJed Brown 8184df9cb4SJed Brown Level: developer 8284df9cb4SJed Brown 8384df9cb4SJed Brown .keywords: Petsc, destroy, package 8484df9cb4SJed Brown .seealso: PetscFinalize() 8584df9cb4SJed Brown @*/ 8684df9cb4SJed Brown PetscErrorCode TSAdaptFinalizePackage(void) 8784df9cb4SJed Brown { 8837e93019SBarry Smith PetscErrorCode ierr; 8937e93019SBarry Smith 9084df9cb4SJed Brown PetscFunctionBegin; 9137e93019SBarry Smith ierr = PetscFunctionListDestroy(&TSAdaptList);CHKERRQ(ierr); 9284df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_FALSE; 9384df9cb4SJed Brown TSAdaptRegisterAllCalled = PETSC_FALSE; 9484df9cb4SJed Brown PetscFunctionReturn(0); 9584df9cb4SJed Brown } 9684df9cb4SJed Brown 9784df9cb4SJed Brown /*@C 9884df9cb4SJed Brown TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is 9984df9cb4SJed Brown called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to 10026d28e4eSEmil Constantinescu TSCreate_GLLE() when using static libraries. 10184df9cb4SJed Brown 10284df9cb4SJed Brown Level: developer 10384df9cb4SJed Brown 10484df9cb4SJed Brown .keywords: TSAdapt, initialize, package 10584df9cb4SJed Brown .seealso: PetscInitialize() 10684df9cb4SJed Brown @*/ 107607a6623SBarry Smith PetscErrorCode TSAdaptInitializePackage(void) 10884df9cb4SJed Brown { 10984df9cb4SJed Brown PetscErrorCode ierr; 11084df9cb4SJed Brown 11184df9cb4SJed Brown PetscFunctionBegin; 11284df9cb4SJed Brown if (TSAdaptPackageInitialized) PetscFunctionReturn(0); 11384df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_TRUE; 11484df9cb4SJed Brown ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr); 115607a6623SBarry Smith ierr = TSAdaptRegisterAll();CHKERRQ(ierr); 11684df9cb4SJed Brown ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr); 11784df9cb4SJed Brown PetscFunctionReturn(0); 11884df9cb4SJed Brown } 11984df9cb4SJed Brown 1207eef6055SBarry Smith /*@C 1217eef6055SBarry Smith TSAdaptSetType - sets the approach used for the error adapter, currently there is only TSADAPTBASIC and TSADAPTNONE 1227eef6055SBarry Smith 1237eef6055SBarry Smith Logicially Collective on TSAdapt 1247eef6055SBarry Smith 1257eef6055SBarry Smith Input Parameter: 1267eef6055SBarry Smith + adapt - the TS error adapter, most likely obtained with TSGetAdapt() 1277eef6055SBarry Smith - type - either TSADAPTBASIC or TSADAPTNONE 1287eef6055SBarry Smith 1297eef6055SBarry Smith Options Database: 1307eef6055SBarry Smith . -ts_adapt_type basic or none - to setting the adapter type 1317eef6055SBarry Smith 1327eef6055SBarry Smith Level: intermediate 1337eef6055SBarry Smith 1347eef6055SBarry Smith .keywords: TSAdapt, create 1357eef6055SBarry Smith 1367eef6055SBarry Smith .seealso: TSGetAdapt(), TSAdaptDestroy(), TSAdaptType, TSAdaptGetType() 1377eef6055SBarry Smith @*/ 13819fd82e9SBarry Smith PetscErrorCode TSAdaptSetType(TSAdapt adapt,TSAdaptType type) 13984df9cb4SJed Brown { 140ef749922SLisandro Dalcin PetscBool match; 14184df9cb4SJed Brown PetscErrorCode ierr,(*r)(TSAdapt); 14284df9cb4SJed Brown 14384df9cb4SJed Brown PetscFunctionBegin; 1444782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 145ef749922SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)adapt,type,&match);CHKERRQ(ierr); 146ef749922SLisandro Dalcin if (match) PetscFunctionReturn(0); 1471c9cd337SJed Brown ierr = PetscFunctionListFind(TSAdaptList,type,&r);CHKERRQ(ierr); 14884df9cb4SJed Brown if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type); 149ef749922SLisandro Dalcin if (adapt->ops->destroy) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);} 1504cefc2ffSBarry Smith ierr = PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));CHKERRQ(ierr); 15184df9cb4SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr); 15268ae941aSLisandro Dalcin ierr = (*r)(adapt);CHKERRQ(ierr); 15384df9cb4SJed Brown PetscFunctionReturn(0); 15484df9cb4SJed Brown } 15584df9cb4SJed Brown 15684df9cb4SJed Brown PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[]) 15784df9cb4SJed Brown { 15884df9cb4SJed Brown PetscErrorCode ierr; 15984df9cb4SJed Brown 16084df9cb4SJed Brown PetscFunctionBegin; 1614782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 16284df9cb4SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr); 16384df9cb4SJed Brown PetscFunctionReturn(0); 16484df9cb4SJed Brown } 16584df9cb4SJed Brown 166ad6bc421SBarry Smith /*@C 167ad6bc421SBarry Smith TSAdaptLoad - Loads a TSAdapt that has been stored in binary with TSAdaptView(). 168ad6bc421SBarry Smith 169ad6bc421SBarry Smith Collective on PetscViewer 170ad6bc421SBarry Smith 171ad6bc421SBarry Smith Input Parameters: 172ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or 173ad6bc421SBarry Smith some related function before a call to TSAdaptLoad(). 174ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 175ad6bc421SBarry Smith HDF5 file viewer, obtained from PetscViewerHDF5Open() 176ad6bc421SBarry Smith 177ad6bc421SBarry Smith Level: intermediate 178ad6bc421SBarry Smith 179ad6bc421SBarry Smith Notes: 180ad6bc421SBarry Smith The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored. 181ad6bc421SBarry Smith 182ad6bc421SBarry Smith Notes for advanced users: 183ad6bc421SBarry Smith Most users should not need to know the details of the binary storage 184ad6bc421SBarry Smith format, since TSAdaptLoad() and TSAdaptView() completely hide these details. 185ad6bc421SBarry Smith But for anyone who's interested, the standard binary matrix storage 186ad6bc421SBarry Smith format is 187ad6bc421SBarry Smith .vb 188ad6bc421SBarry Smith has not yet been determined 189ad6bc421SBarry Smith .ve 190ad6bc421SBarry Smith 191ad6bc421SBarry Smith .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad() 192ad6bc421SBarry Smith @*/ 1934782b174SLisandro Dalcin PetscErrorCode TSAdaptLoad(TSAdapt adapt,PetscViewer viewer) 194ad6bc421SBarry Smith { 195ad6bc421SBarry Smith PetscErrorCode ierr; 196ad6bc421SBarry Smith PetscBool isbinary; 197ad6bc421SBarry Smith char type[256]; 198ad6bc421SBarry Smith 199ad6bc421SBarry Smith PetscFunctionBegin; 2004782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 201ad6bc421SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 202ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 203ad6bc421SBarry Smith if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 204ad6bc421SBarry Smith 205060da220SMatthew G. Knepley ierr = PetscViewerBinaryRead(viewer,type,256,NULL,PETSC_CHAR);CHKERRQ(ierr); 2064782b174SLisandro Dalcin ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 2074782b174SLisandro Dalcin if (adapt->ops->load) { 2084782b174SLisandro Dalcin ierr = (*adapt->ops->load)(adapt,viewer);CHKERRQ(ierr); 209ad6bc421SBarry Smith } 210ad6bc421SBarry Smith PetscFunctionReturn(0); 211ad6bc421SBarry Smith } 212ad6bc421SBarry Smith 21384df9cb4SJed Brown PetscErrorCode TSAdaptView(TSAdapt adapt,PetscViewer viewer) 21484df9cb4SJed Brown { 21584df9cb4SJed Brown PetscErrorCode ierr; 216ad6bc421SBarry Smith PetscBool iascii,isbinary; 21784df9cb4SJed Brown 21884df9cb4SJed Brown PetscFunctionBegin; 2194782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 2204782b174SLisandro Dalcin if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer);CHKERRQ(ierr);} 2214782b174SLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2224782b174SLisandro Dalcin PetscCheckSameComm(adapt,1,viewer,2); 223251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 224ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 22584df9cb4SJed Brown if (iascii) { 226dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr); 22784df9cb4SJed Brown ierr = PetscViewerASCIIPrintf(viewer," number of candidates %D\n",adapt->candidates.n);CHKERRQ(ierr); 228bf997491SLisandro Dalcin if (adapt->always_accept) {ierr = PetscViewerASCIIPrintf(viewer," always accepting steps\n");CHKERRQ(ierr);} 229*bc0b8257SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," maximum allowed timestep %g\n",(double)adapt->dt_max);CHKERRQ(ierr); 230*bc0b8257SSatish Balay ierr = PetscViewerASCIIPrintf(viewer," minimum allowed timestep %g\n",(double)adapt->dt_min);CHKERRQ(ierr); 23184df9cb4SJed Brown if (adapt->ops->view) { 23284df9cb4SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 23384df9cb4SJed Brown ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 23484df9cb4SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 23584df9cb4SJed Brown } 236ad6bc421SBarry Smith } else if (isbinary) { 237ad6bc421SBarry Smith char type[256]; 238ad6bc421SBarry Smith 239ad6bc421SBarry Smith /* need to save FILE_CLASS_ID for adapt class */ 240ad6bc421SBarry Smith ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr); 241ad6bc421SBarry Smith ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr); 242bbd56ea5SKarl Rupp } else if (adapt->ops->view) { 243f2c2a1b9SBarry Smith ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 244f2c2a1b9SBarry Smith } 24584df9cb4SJed Brown PetscFunctionReturn(0); 24684df9cb4SJed Brown } 24784df9cb4SJed Brown 248099af0a3SLisandro Dalcin /*@ 249099af0a3SLisandro Dalcin TSAdaptReset - Resets a TSAdapt context. 250099af0a3SLisandro Dalcin 251099af0a3SLisandro Dalcin Collective on TS 252099af0a3SLisandro Dalcin 253099af0a3SLisandro Dalcin Input Parameter: 254099af0a3SLisandro Dalcin . adapt - the TSAdapt context obtained from TSAdaptCreate() 255099af0a3SLisandro Dalcin 256099af0a3SLisandro Dalcin Level: developer 257099af0a3SLisandro Dalcin 258099af0a3SLisandro Dalcin .seealso: TSAdaptCreate(), TSAdaptDestroy() 259099af0a3SLisandro Dalcin @*/ 260099af0a3SLisandro Dalcin PetscErrorCode TSAdaptReset(TSAdapt adapt) 261099af0a3SLisandro Dalcin { 262099af0a3SLisandro Dalcin PetscErrorCode ierr; 263099af0a3SLisandro Dalcin 264099af0a3SLisandro Dalcin PetscFunctionBegin; 265099af0a3SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 266099af0a3SLisandro Dalcin if (adapt->ops->reset) {ierr = (*adapt->ops->reset)(adapt);CHKERRQ(ierr);} 267099af0a3SLisandro Dalcin PetscFunctionReturn(0); 268099af0a3SLisandro Dalcin } 269099af0a3SLisandro Dalcin 27084df9cb4SJed Brown PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 27184df9cb4SJed Brown { 27284df9cb4SJed Brown PetscErrorCode ierr; 27384df9cb4SJed Brown 27484df9cb4SJed Brown PetscFunctionBegin; 27584df9cb4SJed Brown if (!*adapt) PetscFunctionReturn(0); 27684df9cb4SJed Brown PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 2774782b174SLisandro Dalcin if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);} 278099af0a3SLisandro Dalcin 279099af0a3SLisandro Dalcin ierr = TSAdaptReset(*adapt);CHKERRQ(ierr); 280099af0a3SLisandro Dalcin 28184df9cb4SJed Brown if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);} 2821c3436cfSJed Brown ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr); 28384df9cb4SJed Brown ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr); 28484df9cb4SJed Brown PetscFunctionReturn(0); 28584df9cb4SJed Brown } 28684df9cb4SJed Brown 2871c3436cfSJed Brown /*@ 2881c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 2891c3436cfSJed Brown 2901c3436cfSJed Brown Collective on TSAdapt 2911c3436cfSJed Brown 2921c3436cfSJed Brown Input Arguments: 2931c3436cfSJed Brown + adapt - adaptive controller context 2941c3436cfSJed Brown - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 2951c3436cfSJed Brown 296bf997491SLisandro Dalcin Options Database Keys: 297bf997491SLisandro Dalcin . -ts_adapt_monitor 298bf997491SLisandro Dalcin 2991c3436cfSJed Brown Level: intermediate 3001c3436cfSJed Brown 3011c3436cfSJed Brown .seealso: TSAdaptChoose() 3021c3436cfSJed Brown @*/ 3031c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 3041c3436cfSJed Brown { 3051c3436cfSJed Brown PetscErrorCode ierr; 3061c3436cfSJed Brown 3071c3436cfSJed Brown PetscFunctionBegin; 3084782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 3094782b174SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flg,2); 3101c3436cfSJed Brown if (flg) { 311ce94432eSBarry Smith if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);} 3121c3436cfSJed Brown } else { 3131c3436cfSJed Brown ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr); 3141c3436cfSJed Brown } 3151c3436cfSJed Brown PetscFunctionReturn(0); 3161c3436cfSJed Brown } 3171c3436cfSJed Brown 3180873808bSJed Brown /*@C 319bf997491SLisandro Dalcin TSAdaptSetCheckStage - Set a callback to check convergence for a stage 3200873808bSJed Brown 3210873808bSJed Brown Logically collective on TSAdapt 3220873808bSJed Brown 3230873808bSJed Brown Input Arguments: 3240873808bSJed Brown + adapt - adaptive controller context 3250873808bSJed Brown - func - stage check function 3260873808bSJed Brown 3270873808bSJed Brown Arguments of func: 3280873808bSJed Brown $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 3290873808bSJed Brown 3300873808bSJed Brown + adapt - adaptive controller context 3310873808bSJed Brown . ts - time stepping context 3320873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine 3330873808bSJed Brown 3340873808bSJed Brown Level: advanced 3350873808bSJed Brown 3360873808bSJed Brown .seealso: TSAdaptChoose() 3370873808bSJed Brown @*/ 338b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscReal,Vec,PetscBool*)) 3390873808bSJed Brown { 3400873808bSJed Brown 3410873808bSJed Brown PetscFunctionBegin; 3420873808bSJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 34368ae941aSLisandro Dalcin adapt->checkstage = func; 3440873808bSJed Brown PetscFunctionReturn(0); 3450873808bSJed Brown } 3460873808bSJed Brown 3471c3436cfSJed Brown /*@ 348bf997491SLisandro Dalcin TSAdaptSetAlwaysAccept - Set whether to always accept steps regardless of 349bf997491SLisandro Dalcin any error or stability condition not meeting the prescribed goal. 350bf997491SLisandro Dalcin 351bf997491SLisandro Dalcin 352bf997491SLisandro Dalcin Input Arguments: 353bf997491SLisandro Dalcin + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 354bf997491SLisandro Dalcin - flag - whether to always accept steps 355bf997491SLisandro Dalcin 356bf997491SLisandro Dalcin Options Database Keys: 357bf997491SLisandro Dalcin . -ts_adapt_always_accept 358bf997491SLisandro Dalcin 359bf997491SLisandro Dalcin Level: intermediate 360bf997491SLisandro Dalcin 361bf997491SLisandro Dalcin .seealso: TSAdapt, TSAdaptChoose() 362bf997491SLisandro Dalcin @*/ 363bf997491SLisandro Dalcin PetscErrorCode TSAdaptSetAlwaysAccept(TSAdapt adapt,PetscBool flag) 364bf997491SLisandro Dalcin { 365bf997491SLisandro Dalcin PetscFunctionBegin; 366bf997491SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 367bf997491SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flag,2); 368bf997491SLisandro Dalcin adapt->always_accept = flag; 369bf997491SLisandro Dalcin PetscFunctionReturn(0); 370bf997491SLisandro Dalcin } 371bf997491SLisandro Dalcin 372bf997491SLisandro Dalcin /*@ 3731c3436cfSJed Brown TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller 3741c3436cfSJed Brown 3751c3436cfSJed Brown Logically Collective 3761c3436cfSJed Brown 3771c3436cfSJed Brown Input Arguments: 378552698daSJed Brown + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 3791c3436cfSJed Brown . hmin - minimum time step 3801c3436cfSJed Brown - hmax - maximum time step 3811c3436cfSJed Brown 3821c3436cfSJed Brown Options Database Keys: 3831c3436cfSJed Brown + -ts_adapt_dt_min - minimum time step 3841c3436cfSJed Brown - -ts_adapt_dt_max - maximum time step 3851c3436cfSJed Brown 3861c3436cfSJed Brown Level: intermediate 3871c3436cfSJed Brown 388bf997491SLisandro Dalcin .seealso: TSAdapt, TSAdaptChoose() 3891c3436cfSJed Brown @*/ 3901c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 3911c3436cfSJed Brown { 3921c3436cfSJed Brown 3931c3436cfSJed Brown PetscFunctionBegin; 3944782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 395b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,hmin,2); 396b1f5bccaSLisandro Dalcin PetscValidLogicalCollectiveReal(adapt,hmax,3); 397b1f5bccaSLisandro Dalcin if (hmin != PETSC_DEFAULT && hmin < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmin); 398b1f5bccaSLisandro Dalcin if (hmax != PETSC_DEFAULT && hmax < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Minimum time step %g must be non negative",(double)hmax); 399b1f5bccaSLisandro Dalcin if (hmin != PETSC_DEFAULT) adapt->dt_min = hmin; 400b1f5bccaSLisandro Dalcin if (hmax != PETSC_DEFAULT) adapt->dt_max = hmax; 401b1f5bccaSLisandro Dalcin #if defined(PETSC_USE_DEBUG) 402b1f5bccaSLisandro Dalcin hmin = adapt->dt_min; 403b1f5bccaSLisandro Dalcin hmax = adapt->dt_max; 404b1f5bccaSLisandro Dalcin if (hmax <= hmin) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Maximum time step %g must geather than minimum time step %g",(double)hmax,(double)hmin); 405b1f5bccaSLisandro Dalcin #endif 4061c3436cfSJed Brown PetscFunctionReturn(0); 4071c3436cfSJed Brown } 4081c3436cfSJed Brown 409e55864a3SBarry Smith /* 41084df9cb4SJed Brown TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 41184df9cb4SJed Brown 41284df9cb4SJed Brown Collective on TSAdapt 41384df9cb4SJed Brown 41484df9cb4SJed Brown Input Parameter: 41584df9cb4SJed Brown . adapt - the TSAdapt context 41684df9cb4SJed Brown 41784df9cb4SJed Brown Options Database Keys: 41823de1d84SBarry Smith + -ts_adapt_type <type> - basic 419bf997491SLisandro Dalcin . -ts_adapt_always_accept - always accept steps regardless of error/stability goals 42023de1d84SBarry Smith . -ts_adapt_dt_min <min> - minimum timestep to use 42123de1d84SBarry Smith . -ts_adapt_dt_max <max> - maximum timestep to use 42223de1d84SBarry Smith . -ts_adapt_scale_solve_failed <scale> - scale timestep by this factor if a solve fails 42323de1d84SBarry Smith - -ts_adapt_wnormtype <2 or infinity> - type of norm for computing error estimates 42484df9cb4SJed Brown 42584df9cb4SJed Brown Level: advanced 42684df9cb4SJed Brown 42784df9cb4SJed Brown Notes: 42884df9cb4SJed Brown This function is automatically called by TSSetFromOptions() 42984df9cb4SJed Brown 43023de1d84SBarry Smith .keywords: TS, TSGetAdapt(), TSAdaptSetType(), TSAdaptSetStepLimits() 43184df9cb4SJed Brown 43284df9cb4SJed Brown .seealso: TSGetType() 433e55864a3SBarry Smith */ 4344416b707SBarry Smith PetscErrorCode TSAdaptSetFromOptions(PetscOptionItems *PetscOptionsObject,TSAdapt adapt) 43584df9cb4SJed Brown { 43684df9cb4SJed Brown PetscErrorCode ierr; 43784df9cb4SJed Brown char type[256] = TSADAPTBASIC; 438bf997491SLisandro Dalcin PetscReal hmin,hmax; 4391c3436cfSJed Brown PetscBool set,flg; 44084df9cb4SJed Brown 44184df9cb4SJed Brown PetscFunctionBegin; 4424782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 44384df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 4441566a47fSLisandro Dalcin * function can only be called from inside TSSetFromOptions() */ 445e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");CHKERRQ(ierr); 44623de1d84SBarry Smith ierr = PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList,((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);CHKERRQ(ierr); 44784df9cb4SJed Brown if (flg || !((PetscObject)adapt)->type_name) { 44884df9cb4SJed Brown ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 44984df9cb4SJed Brown } 450bf997491SLisandro Dalcin ierr = PetscOptionsBool("-ts_adapt_always_accept","Always accept the step","TSAdaptSetAlwaysAccept",adapt->always_accept,&flg,&set);CHKERRQ(ierr); 451bf997491SLisandro Dalcin if (set) {ierr = TSAdaptSetAlwaysAccept(adapt,flg);CHKERRQ(ierr);} 452bf997491SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",hmin=adapt->dt_min,&hmin,&set);CHKERRQ(ierr); 453bf997491SLisandro Dalcin ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",hmax=adapt->dt_max,&hmax,&flg);CHKERRQ(ierr); 454bf997491SLisandro Dalcin if (set || flg) {ierr = TSAdaptSetStepLimits(adapt,hmin,hmax);CHKERRQ(ierr);} 4550298fd71SBarry 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); 4561c3436cfSJed Brown ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 457bf997491SLisandro Dalcin if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);} 4587619abb3SShri ierr = PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);CHKERRQ(ierr); 4597619abb3SShri if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported"); 460e55864a3SBarry Smith if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);} 46184df9cb4SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 46284df9cb4SJed Brown PetscFunctionReturn(0); 46384df9cb4SJed Brown } 46484df9cb4SJed Brown 46584df9cb4SJed Brown /*@ 46684df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes 46784df9cb4SJed Brown 46884df9cb4SJed Brown Logically Collective 46984df9cb4SJed Brown 47084df9cb4SJed Brown Input Argument: 47184df9cb4SJed Brown . adapt - adaptive controller 47284df9cb4SJed Brown 47384df9cb4SJed Brown Level: developer 47484df9cb4SJed Brown 47584df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose() 47684df9cb4SJed Brown @*/ 47784df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 47884df9cb4SJed Brown { 47984df9cb4SJed Brown PetscErrorCode ierr; 48084df9cb4SJed Brown 48184df9cb4SJed Brown PetscFunctionBegin; 4824782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 48384df9cb4SJed Brown ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr); 48484df9cb4SJed Brown PetscFunctionReturn(0); 48584df9cb4SJed Brown } 48684df9cb4SJed Brown 48784df9cb4SJed Brown /*@C 48884df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 48984df9cb4SJed Brown 49084df9cb4SJed Brown Logically Collective 49184df9cb4SJed Brown 49284df9cb4SJed Brown Input Arguments: 493552698daSJed Brown + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 49484df9cb4SJed Brown . name - name of the candidate scheme to add 49584df9cb4SJed Brown . order - order of the candidate scheme 49684df9cb4SJed Brown . stageorder - stage order of the candidate scheme 4978d59e960SJed Brown . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 49884df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme 49984df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 50084df9cb4SJed Brown 50184df9cb4SJed Brown Note: 50284df9cb4SJed Brown This routine is not available in Fortran. 50384df9cb4SJed Brown 50484df9cb4SJed Brown Level: developer 50584df9cb4SJed Brown 50684df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose() 50784df9cb4SJed Brown @*/ 5088d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse) 50984df9cb4SJed Brown { 51084df9cb4SJed Brown PetscInt c; 51184df9cb4SJed Brown 51284df9cb4SJed Brown PetscFunctionBegin; 51384df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 514ce94432eSBarry Smith if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order); 51584df9cb4SJed Brown if (inuse) { 516ce94432eSBarry 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()"); 51784df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE; 51884df9cb4SJed Brown } 5191c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */ 5201c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 521bbd56ea5SKarl Rupp 52284df9cb4SJed Brown adapt->candidates.name[c] = name; 52384df9cb4SJed Brown adapt->candidates.order[c] = order; 52484df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder; 5258d59e960SJed Brown adapt->candidates.ccfl[c] = ccfl; 52684df9cb4SJed Brown adapt->candidates.cost[c] = cost; 52784df9cb4SJed Brown adapt->candidates.n++; 52884df9cb4SJed Brown PetscFunctionReturn(0); 52984df9cb4SJed Brown } 53084df9cb4SJed Brown 5318d59e960SJed Brown /*@C 5328d59e960SJed Brown TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 5338d59e960SJed Brown 5348d59e960SJed Brown Not Collective 5358d59e960SJed Brown 5368d59e960SJed Brown Input Arguments: 5378d59e960SJed Brown . adapt - time step adaptivity context 5388d59e960SJed Brown 5398d59e960SJed Brown Output Arguments: 5408d59e960SJed Brown + n - number of candidate schemes, always at least 1 5418d59e960SJed Brown . order - the order of each candidate scheme 5428d59e960SJed Brown . stageorder - the stage order of each candidate scheme 5438d59e960SJed Brown . ccfl - the CFL coefficient of each scheme 5448d59e960SJed Brown - cost - the relative cost of each scheme 5458d59e960SJed Brown 5468d59e960SJed Brown Level: developer 5478d59e960SJed Brown 5488d59e960SJed Brown Note: 5498d59e960SJed Brown The current scheme is always returned in the first slot 5508d59e960SJed Brown 5518d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose() 5528d59e960SJed Brown @*/ 5538d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost) 5548d59e960SJed Brown { 5558d59e960SJed Brown PetscFunctionBegin; 5568d59e960SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 5578d59e960SJed Brown if (n) *n = adapt->candidates.n; 5588d59e960SJed Brown if (order) *order = adapt->candidates.order; 5598d59e960SJed Brown if (stageorder) *stageorder = adapt->candidates.stageorder; 5608d59e960SJed Brown if (ccfl) *ccfl = adapt->candidates.ccfl; 5618d59e960SJed Brown if (cost) *cost = adapt->candidates.cost; 5628d59e960SJed Brown PetscFunctionReturn(0); 5638d59e960SJed Brown } 5648d59e960SJed Brown 56584df9cb4SJed Brown /*@C 56684df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step 56784df9cb4SJed Brown 56884df9cb4SJed Brown Logically Collective 56984df9cb4SJed Brown 57084df9cb4SJed Brown Input Arguments: 57184df9cb4SJed Brown + adapt - adaptive contoller 57284df9cb4SJed Brown - h - current step size 57384df9cb4SJed Brown 57484df9cb4SJed Brown Output Arguments: 5751566a47fSLisandro Dalcin + next_sc - optional, scheme to use for the next step 57684df9cb4SJed Brown . next_h - step size to use for the next step 57784df9cb4SJed Brown - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 57884df9cb4SJed Brown 5791c3436cfSJed Brown Note: 5801c3436cfSJed Brown The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 5811c3436cfSJed Brown being retried after an initial rejection. 5821c3436cfSJed Brown 58384df9cb4SJed Brown Level: developer 58484df9cb4SJed Brown 58584df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd() 58684df9cb4SJed Brown @*/ 58784df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 58884df9cb4SJed Brown { 58984df9cb4SJed Brown PetscErrorCode ierr; 5901566a47fSLisandro Dalcin PetscInt ncandidates = adapt->candidates.n; 5911566a47fSLisandro Dalcin PetscInt scheme = 0; 5920b99f514SJed Brown PetscReal wlte = -1.0; 5935e4ed32fSEmil Constantinescu PetscReal wltea = -1.0; 5945e4ed32fSEmil Constantinescu PetscReal wlter = -1.0; 59584df9cb4SJed Brown 59684df9cb4SJed Brown PetscFunctionBegin; 59784df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 59884df9cb4SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,2); 5991566a47fSLisandro Dalcin if (next_sc) PetscValidIntPointer(next_sc,4); 60084df9cb4SJed Brown PetscValidPointer(next_h,5); 60184df9cb4SJed Brown PetscValidIntPointer(accept,6); 6021566a47fSLisandro Dalcin if (next_sc) *next_sc = 0; 6031566a47fSLisandro Dalcin 6041566a47fSLisandro Dalcin /* Do not mess with adaptivity while handling events*/ 6051566a47fSLisandro Dalcin if (ts->event && ts->event->status != TSEVENT_NONE) { 6061566a47fSLisandro Dalcin *next_h = h; 6071566a47fSLisandro Dalcin *accept = PETSC_TRUE; 6081566a47fSLisandro Dalcin PetscFunctionReturn(0); 6091566a47fSLisandro Dalcin } 6101566a47fSLisandro Dalcin 6115e4ed32fSEmil Constantinescu ierr = (*adapt->ops->choose)(adapt,ts,h,&scheme,next_h,accept,&wlte,&wltea,&wlter);CHKERRQ(ierr); 6121566a47fSLisandro Dalcin if (scheme < 0 || (ncandidates > 0 && scheme >= ncandidates)) SETERRQ2(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Chosen scheme %D not in valid range 0..%D",scheme,ncandidates-1); 6131566a47fSLisandro Dalcin if (*next_h < 0) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h); 6141566a47fSLisandro Dalcin if (next_sc) *next_sc = scheme; 6151566a47fSLisandro Dalcin 61649354f04SShri Abhyankar if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 61749354f04SShri Abhyankar /* Reduce time step if it overshoots max time */ 6181566a47fSLisandro Dalcin if (ts->ptime + ts->time_step + *next_h >= ts->max_time) { 6191566a47fSLisandro Dalcin PetscReal next_dt = ts->max_time - (ts->ptime + ts->time_step); 6208709b12bSShri Abhyankar if (next_dt > PETSC_SMALL) *next_h = next_dt; 62149354f04SShri Abhyankar else ts->reason = TS_CONVERGED_TIME; 62249354f04SShri Abhyankar } 62349354f04SShri Abhyankar } 6241c3436cfSJed Brown 6251c3436cfSJed Brown if (adapt->monitor) { 6261566a47fSLisandro Dalcin const char *sc_name = (scheme < ncandidates) ? adapt->candidates.name[scheme] : ""; 6271c3436cfSJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 6280b99f514SJed Brown if (wlte < 0) { 6290dea241dSBarry Smith ierr = PetscViewerASCIIPrintf(adapt->monitor," TSAdapt %s %s %D:%s step %3D %s t=%-11g+%10.3e dt=%-10.3e\n",((PetscObject)adapt)->type_name,((PetscObject)ts)->type_name,scheme,sc_name,ts->steps,*accept ? "accepted" : "rejected",(double)ts->ptime,(double)h,(double)*next_h);CHKERRQ(ierr); 6300b99f514SJed Brown } else { 6310dea241dSBarry Smith ierr = PetscViewerASCIIPrintf(adapt->monitor," TSAdapt %s %s %D:%s step %3D %s t=%-11g+%10.3e dt=%-10.3e wlte=%5.3g wltea=%5.3g wlter=%5.3g\n",((PetscObject)adapt)->type_name,((PetscObject)ts)->type_name,scheme,sc_name,ts->steps,*accept ? "accepted" : "rejected",(double)ts->ptime,(double)h,(double)*next_h,(double)wlte,(double)wltea,(double)wlter);CHKERRQ(ierr); 6320b99f514SJed Brown } 6331c3436cfSJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 6341c3436cfSJed Brown } 63584df9cb4SJed Brown PetscFunctionReturn(0); 63684df9cb4SJed Brown } 63784df9cb4SJed Brown 63897335746SJed Brown /*@ 63997335746SJed Brown TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails) 64097335746SJed Brown 64197335746SJed Brown Collective 64297335746SJed Brown 64397335746SJed Brown Input Arguments: 64497335746SJed Brown + adapt - adaptive controller context 645b295832fSPierre Barbier de Reuille . ts - time stepper 646b295832fSPierre Barbier de Reuille . t - Current simulation time 647b295832fSPierre Barbier de Reuille - Y - Current solution vector 64897335746SJed Brown 64997335746SJed Brown Output Arguments: 65097335746SJed Brown . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 65197335746SJed Brown 65297335746SJed Brown Level: developer 65397335746SJed Brown 65497335746SJed Brown .seealso: 65597335746SJed Brown @*/ 656b295832fSPierre Barbier de Reuille PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,Vec Y,PetscBool *accept) 65797335746SJed Brown { 65897335746SJed Brown PetscErrorCode ierr; 6591566a47fSLisandro Dalcin SNESConvergedReason snesreason = SNES_CONVERGED_ITERATING; 66097335746SJed Brown 66197335746SJed Brown PetscFunctionBegin; 6624782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 6634782b174SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,2); 6644782b174SLisandro Dalcin PetscValidIntPointer(accept,3); 6651566a47fSLisandro Dalcin 6661566a47fSLisandro Dalcin if (ts->snes) {ierr = SNESGetConvergedReason(ts->snes,&snesreason);CHKERRQ(ierr);} 66797335746SJed Brown if (snesreason < 0) { 66897335746SJed Brown *accept = PETSC_FALSE; 6696de24e2aSJed Brown if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 67097335746SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 67148c19aefSLisandro 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); 67297335746SJed Brown if (adapt->monitor) { 67397335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 6740dea241dSBarry Smith ierr = PetscViewerASCIIPrintf(adapt->monitor," TSAdapt %s step %3D stage rejected t=%-11g+%10.3e, nonlinear solve failures %D greater than current TS allowed\n",((PetscObject)adapt)->type_name,ts->steps,(double)ts->ptime,(double)ts->time_step,ts->num_snes_failures);CHKERRQ(ierr); 67597335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 67697335746SJed Brown } 677cb9d8021SPierre Barbier de Reuille } 678cb9d8021SPierre Barbier de Reuille } else { 6791566a47fSLisandro Dalcin *accept = PETSC_TRUE; 680b295832fSPierre Barbier de Reuille ierr = TSFunctionDomainError(ts,t,Y,accept);CHKERRQ(ierr); 681cb9d8021SPierre Barbier de Reuille if(*accept && adapt->checkstage) { 682b295832fSPierre Barbier de Reuille ierr = (*adapt->checkstage)(adapt,ts,t,Y,accept);CHKERRQ(ierr); 683cb9d8021SPierre Barbier de Reuille } 684cb9d8021SPierre Barbier de Reuille } 685cb9d8021SPierre Barbier de Reuille 6861566a47fSLisandro Dalcin if(!(*accept) && !ts->reason) { 6871566a47fSLisandro Dalcin PetscReal dt,new_dt; 6881566a47fSLisandro Dalcin ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 689cb9d8021SPierre Barbier de Reuille new_dt = dt * adapt->scale_solve_failed; 69097335746SJed Brown ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr); 69197335746SJed Brown if (adapt->monitor) { 69297335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 6930dea241dSBarry Smith 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); 69497335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 69597335746SJed Brown } 69697335746SJed Brown } 69797335746SJed Brown PetscFunctionReturn(0); 69897335746SJed Brown } 69997335746SJed Brown 70084df9cb4SJed Brown /*@ 70184df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping 70284df9cb4SJed Brown 70384df9cb4SJed Brown Collective on MPI_Comm 70484df9cb4SJed Brown 70584df9cb4SJed Brown Input Parameter: 70684df9cb4SJed Brown . comm - The communicator 70784df9cb4SJed Brown 70884df9cb4SJed Brown Output Parameter: 70984df9cb4SJed Brown . adapt - new TSAdapt object 71084df9cb4SJed Brown 71184df9cb4SJed Brown Level: developer 71284df9cb4SJed Brown 71384df9cb4SJed Brown Notes: 71484df9cb4SJed Brown TSAdapt creation is handled by TS, so users should not need to call this function. 71584df9cb4SJed Brown 71684df9cb4SJed Brown .keywords: TSAdapt, create 717552698daSJed Brown .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy() 71884df9cb4SJed Brown @*/ 71984df9cb4SJed Brown PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 72084df9cb4SJed Brown { 72184df9cb4SJed Brown PetscErrorCode ierr; 72284df9cb4SJed Brown TSAdapt adapt; 72384df9cb4SJed Brown 72484df9cb4SJed Brown PetscFunctionBegin; 7253b3bcf4cSLisandro Dalcin PetscValidPointer(inadapt,1); 7263b3bcf4cSLisandro Dalcin *inadapt = NULL; 7273b3bcf4cSLisandro Dalcin ierr = TSAdaptInitializePackage();CHKERRQ(ierr); 7283b3bcf4cSLisandro Dalcin 72973107ff1SLisandro Dalcin ierr = PetscHeaderCreate(adapt,TSADAPT_CLASSID,"TSAdapt","Time stepping adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr); 7301c3436cfSJed Brown 731bf997491SLisandro Dalcin adapt->always_accept = PETSC_FALSE; 7321c3436cfSJed Brown adapt->dt_min = 1e-20; 733*bc0b8257SSatish Balay #if defined(PETSC_USE_REAL_SINGLE) 734*bc0b8257SSatish Balay adapt->dt_max = 1e35; 735*bc0b8257SSatish Balay #else 7361c3436cfSJed Brown adapt->dt_max = 1e50; 737*bc0b8257SSatish Balay #endif 73897335746SJed Brown adapt->scale_solve_failed = 0.25; 7397619abb3SShri adapt->wnormtype = NORM_2; 7407eef6055SBarry Smith ierr = TSAdaptSetType(adapt,TSADAPTBASIC);CHKERRQ(ierr); 7411c3436cfSJed Brown 74284df9cb4SJed Brown *inadapt = adapt; 74384df9cb4SJed Brown PetscFunctionReturn(0); 74484df9cb4SJed Brown } 745