184df9cb4SJed Brown 2b45d2f2cSJed Brown #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 384df9cb4SJed Brown 4140e18c1SBarry Smith static PetscFunctionList TSAdaptList; 584df9cb4SJed Brown static PetscBool TSAdaptPackageInitialized; 684df9cb4SJed Brown static PetscBool TSAdaptRegisterAllCalled; 784df9cb4SJed Brown static PetscClassId TSADAPT_CLASSID; 884df9cb4SJed Brown 98cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_Basic(TSAdapt); 108cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_None(TSAdapt); 118cc058d9SJed Brown PETSC_EXTERN PetscErrorCode TSAdaptCreate_CFL(TSAdapt); 1284df9cb4SJed Brown 1384df9cb4SJed Brown #undef __FUNCT__ 1484df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegister" 1584df9cb4SJed Brown /*@C 161c84c290SBarry Smith TSAdaptRegister - adds a TSAdapt implementation 171c84c290SBarry Smith 181c84c290SBarry Smith Not Collective 191c84c290SBarry Smith 201c84c290SBarry Smith Input Parameters: 211c84c290SBarry Smith + name_scheme - name of user-defined adaptivity scheme 221c84c290SBarry Smith - routine_create - routine to create method context 231c84c290SBarry Smith 241c84c290SBarry Smith Notes: 251c84c290SBarry Smith TSAdaptRegister() may be called multiple times to add several user-defined families. 261c84c290SBarry Smith 271c84c290SBarry Smith Sample usage: 281c84c290SBarry Smith .vb 29bdf89e91SBarry Smith TSAdaptRegister("my_scheme",MySchemeCreate); 301c84c290SBarry Smith .ve 311c84c290SBarry Smith 321c84c290SBarry Smith Then, your scheme can be chosen with the procedural interface via 331c84c290SBarry Smith $ TSAdaptSetType(ts,"my_scheme") 341c84c290SBarry Smith or at runtime via the option 351c84c290SBarry Smith $ -ts_adapt_type my_scheme 3684df9cb4SJed Brown 3784df9cb4SJed Brown Level: advanced 381c84c290SBarry Smith 391c84c290SBarry Smith .keywords: TSAdapt, register 401c84c290SBarry Smith 411c84c290SBarry Smith .seealso: TSAdaptRegisterAll() 4284df9cb4SJed Brown @*/ 43bdf89e91SBarry Smith PetscErrorCode TSAdaptRegister(const char sname[],PetscErrorCode (*function)(TSAdapt)) 4484df9cb4SJed Brown { 4584df9cb4SJed Brown PetscErrorCode ierr; 4684df9cb4SJed Brown 4784df9cb4SJed Brown PetscFunctionBegin; 48a240a19fSJed Brown ierr = PetscFunctionListAdd(&TSAdaptList,sname,function);CHKERRQ(ierr); 4984df9cb4SJed Brown PetscFunctionReturn(0); 5084df9cb4SJed Brown } 5184df9cb4SJed Brown 5284df9cb4SJed Brown #undef __FUNCT__ 5384df9cb4SJed Brown #define __FUNCT__ "TSAdaptRegisterAll" 5484df9cb4SJed Brown /*@C 5584df9cb4SJed Brown TSAdaptRegisterAll - Registers all of the adaptivity schemes in TSAdapt 5684df9cb4SJed Brown 5784df9cb4SJed Brown Not Collective 5884df9cb4SJed Brown 5984df9cb4SJed Brown Level: advanced 6084df9cb4SJed Brown 6184df9cb4SJed Brown .keywords: TSAdapt, register, all 6284df9cb4SJed Brown 6384df9cb4SJed Brown .seealso: TSAdaptRegisterDestroy() 6484df9cb4SJed Brown @*/ 65607a6623SBarry Smith PetscErrorCode TSAdaptRegisterAll(void) 6684df9cb4SJed Brown { 6784df9cb4SJed Brown PetscErrorCode ierr; 6884df9cb4SJed Brown 6984df9cb4SJed Brown PetscFunctionBegin; 700f51fdf8SToby Isaac if (TSAdaptRegisterAllCalled) PetscFunctionReturn(0); 710f51fdf8SToby Isaac TSAdaptRegisterAllCalled = PETSC_TRUE; 72bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTBASIC,TSAdaptCreate_Basic);CHKERRQ(ierr); 73bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTNONE, TSAdaptCreate_None);CHKERRQ(ierr); 74bdf89e91SBarry Smith ierr = TSAdaptRegister(TSADAPTCFL, TSAdaptCreate_CFL);CHKERRQ(ierr); 7584df9cb4SJed Brown PetscFunctionReturn(0); 7684df9cb4SJed Brown } 7784df9cb4SJed Brown 7884df9cb4SJed Brown #undef __FUNCT__ 7984df9cb4SJed Brown #define __FUNCT__ "TSAdaptFinalizePackage" 8084df9cb4SJed Brown /*@C 8184df9cb4SJed Brown TSFinalizePackage - This function destroys everything in the TS package. It is 8284df9cb4SJed Brown called from PetscFinalize(). 8384df9cb4SJed Brown 8484df9cb4SJed Brown Level: developer 8584df9cb4SJed Brown 8684df9cb4SJed Brown .keywords: Petsc, destroy, package 8784df9cb4SJed Brown .seealso: PetscFinalize() 8884df9cb4SJed Brown @*/ 8984df9cb4SJed Brown PetscErrorCode TSAdaptFinalizePackage(void) 9084df9cb4SJed Brown { 9137e93019SBarry Smith PetscErrorCode ierr; 9237e93019SBarry Smith 9384df9cb4SJed Brown PetscFunctionBegin; 9437e93019SBarry Smith ierr = PetscFunctionListDestroy(&TSAdaptList);CHKERRQ(ierr); 9584df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_FALSE; 9684df9cb4SJed Brown TSAdaptRegisterAllCalled = PETSC_FALSE; 9784df9cb4SJed Brown PetscFunctionReturn(0); 9884df9cb4SJed Brown } 9984df9cb4SJed Brown 10084df9cb4SJed Brown #undef __FUNCT__ 10184df9cb4SJed Brown #define __FUNCT__ "TSAdaptInitializePackage" 10284df9cb4SJed Brown /*@C 10384df9cb4SJed Brown TSAdaptInitializePackage - This function initializes everything in the TSAdapt package. It is 10484df9cb4SJed Brown called from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to 10584df9cb4SJed Brown TSCreate_GL() when using static libraries. 10684df9cb4SJed Brown 10784df9cb4SJed Brown Level: developer 10884df9cb4SJed Brown 10984df9cb4SJed Brown .keywords: TSAdapt, initialize, package 11084df9cb4SJed Brown .seealso: PetscInitialize() 11184df9cb4SJed Brown @*/ 112607a6623SBarry Smith PetscErrorCode TSAdaptInitializePackage(void) 11384df9cb4SJed Brown { 11484df9cb4SJed Brown PetscErrorCode ierr; 11584df9cb4SJed Brown 11684df9cb4SJed Brown PetscFunctionBegin; 11784df9cb4SJed Brown if (TSAdaptPackageInitialized) PetscFunctionReturn(0); 11884df9cb4SJed Brown TSAdaptPackageInitialized = PETSC_TRUE; 11984df9cb4SJed Brown ierr = PetscClassIdRegister("TSAdapt",&TSADAPT_CLASSID);CHKERRQ(ierr); 120607a6623SBarry Smith ierr = TSAdaptRegisterAll();CHKERRQ(ierr); 12184df9cb4SJed Brown ierr = PetscRegisterFinalize(TSAdaptFinalizePackage);CHKERRQ(ierr); 12284df9cb4SJed Brown PetscFunctionReturn(0); 12384df9cb4SJed Brown } 12484df9cb4SJed Brown 12584df9cb4SJed Brown #undef __FUNCT__ 12684df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetType" 12719fd82e9SBarry Smith PetscErrorCode TSAdaptSetType(TSAdapt adapt,TSAdaptType type) 12884df9cb4SJed Brown { 129ef749922SLisandro Dalcin PetscBool match; 13084df9cb4SJed Brown PetscErrorCode ierr,(*r)(TSAdapt); 13184df9cb4SJed Brown 13284df9cb4SJed Brown PetscFunctionBegin; 1334782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 134ef749922SLisandro Dalcin ierr = PetscObjectTypeCompare((PetscObject)adapt,type,&match);CHKERRQ(ierr); 135ef749922SLisandro Dalcin if (match) PetscFunctionReturn(0); 1361c9cd337SJed Brown ierr = PetscFunctionListFind(TSAdaptList,type,&r);CHKERRQ(ierr); 13784df9cb4SJed Brown if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown TSAdapt type \"%s\" given",type); 138ef749922SLisandro Dalcin if (adapt->ops->destroy) {ierr = (*adapt->ops->destroy)(adapt);CHKERRQ(ierr);} 1394cefc2ffSBarry Smith ierr = PetscMemzero(adapt->ops,sizeof(struct _TSAdaptOps));CHKERRQ(ierr); 14084df9cb4SJed Brown ierr = PetscObjectChangeTypeName((PetscObject)adapt,type);CHKERRQ(ierr); 14168ae941aSLisandro Dalcin ierr = (*r)(adapt);CHKERRQ(ierr); 14284df9cb4SJed Brown PetscFunctionReturn(0); 14384df9cb4SJed Brown } 14484df9cb4SJed Brown 14584df9cb4SJed Brown #undef __FUNCT__ 14684df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetOptionsPrefix" 14784df9cb4SJed Brown PetscErrorCode TSAdaptSetOptionsPrefix(TSAdapt adapt,const char prefix[]) 14884df9cb4SJed Brown { 14984df9cb4SJed Brown PetscErrorCode ierr; 15084df9cb4SJed Brown 15184df9cb4SJed Brown PetscFunctionBegin; 1524782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 15384df9cb4SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)adapt,prefix);CHKERRQ(ierr); 15484df9cb4SJed Brown PetscFunctionReturn(0); 15584df9cb4SJed Brown } 15684df9cb4SJed Brown 15784df9cb4SJed Brown #undef __FUNCT__ 158ad6bc421SBarry Smith #define __FUNCT__ "TSAdaptLoad" 159ad6bc421SBarry Smith /*@C 160ad6bc421SBarry Smith TSAdaptLoad - Loads a TSAdapt that has been stored in binary with TSAdaptView(). 161ad6bc421SBarry Smith 162ad6bc421SBarry Smith Collective on PetscViewer 163ad6bc421SBarry Smith 164ad6bc421SBarry Smith Input Parameters: 165ad6bc421SBarry Smith + newdm - the newly loaded TSAdapt, this needs to have been created with TSAdaptCreate() or 166ad6bc421SBarry Smith some related function before a call to TSAdaptLoad(). 167ad6bc421SBarry Smith - viewer - binary file viewer, obtained from PetscViewerBinaryOpen() or 168ad6bc421SBarry Smith HDF5 file viewer, obtained from PetscViewerHDF5Open() 169ad6bc421SBarry Smith 170ad6bc421SBarry Smith Level: intermediate 171ad6bc421SBarry Smith 172ad6bc421SBarry Smith Notes: 173ad6bc421SBarry Smith The type is determined by the data in the file, any type set into the TSAdapt before this call is ignored. 174ad6bc421SBarry Smith 175ad6bc421SBarry Smith Notes for advanced users: 176ad6bc421SBarry Smith Most users should not need to know the details of the binary storage 177ad6bc421SBarry Smith format, since TSAdaptLoad() and TSAdaptView() completely hide these details. 178ad6bc421SBarry Smith But for anyone who's interested, the standard binary matrix storage 179ad6bc421SBarry Smith format is 180ad6bc421SBarry Smith .vb 181ad6bc421SBarry Smith has not yet been determined 182ad6bc421SBarry Smith .ve 183ad6bc421SBarry Smith 184ad6bc421SBarry Smith .seealso: PetscViewerBinaryOpen(), TSAdaptView(), MatLoad(), VecLoad() 185ad6bc421SBarry Smith @*/ 1864782b174SLisandro Dalcin PetscErrorCode TSAdaptLoad(TSAdapt adapt,PetscViewer viewer) 187ad6bc421SBarry Smith { 188ad6bc421SBarry Smith PetscErrorCode ierr; 189ad6bc421SBarry Smith PetscBool isbinary; 190ad6bc421SBarry Smith char type[256]; 191ad6bc421SBarry Smith 192ad6bc421SBarry Smith PetscFunctionBegin; 1934782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 194ad6bc421SBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 195ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 196ad6bc421SBarry Smith if (!isbinary) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Invalid viewer; open viewer with PetscViewerBinaryOpen()"); 197ad6bc421SBarry Smith 198ad6bc421SBarry Smith ierr = PetscViewerBinaryRead(viewer,type,256,PETSC_CHAR);CHKERRQ(ierr); 1994782b174SLisandro Dalcin ierr = TSAdaptSetType(adapt, type);CHKERRQ(ierr); 2004782b174SLisandro Dalcin if (adapt->ops->load) { 2014782b174SLisandro Dalcin ierr = (*adapt->ops->load)(adapt,viewer);CHKERRQ(ierr); 202ad6bc421SBarry Smith } 203ad6bc421SBarry Smith PetscFunctionReturn(0); 204ad6bc421SBarry Smith } 205ad6bc421SBarry Smith 206ad6bc421SBarry Smith #undef __FUNCT__ 20784df9cb4SJed Brown #define __FUNCT__ "TSAdaptView" 20884df9cb4SJed Brown PetscErrorCode TSAdaptView(TSAdapt adapt,PetscViewer viewer) 20984df9cb4SJed Brown { 21084df9cb4SJed Brown PetscErrorCode ierr; 211ad6bc421SBarry Smith PetscBool iascii,isbinary; 21284df9cb4SJed Brown 21384df9cb4SJed Brown PetscFunctionBegin; 2144782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 2154782b174SLisandro Dalcin if (!viewer) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)adapt),&viewer);CHKERRQ(ierr);} 2164782b174SLisandro Dalcin PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2174782b174SLisandro Dalcin PetscCheckSameComm(adapt,1,viewer,2); 218251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 219ad6bc421SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 22084df9cb4SJed Brown if (iascii) { 221dae58748SBarry Smith ierr = PetscObjectPrintClassNamePrefixType((PetscObject)adapt,viewer);CHKERRQ(ierr); 22284df9cb4SJed Brown ierr = PetscViewerASCIIPrintf(viewer," number of candidates %D\n",adapt->candidates.n);CHKERRQ(ierr); 22384df9cb4SJed Brown if (adapt->ops->view) { 22484df9cb4SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 22584df9cb4SJed Brown ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 22684df9cb4SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 22784df9cb4SJed Brown } 228ad6bc421SBarry Smith } else if (isbinary) { 229ad6bc421SBarry Smith char type[256]; 230ad6bc421SBarry Smith 231ad6bc421SBarry Smith /* need to save FILE_CLASS_ID for adapt class */ 232ad6bc421SBarry Smith ierr = PetscStrncpy(type,((PetscObject)adapt)->type_name,256);CHKERRQ(ierr); 233ad6bc421SBarry Smith ierr = PetscViewerBinaryWrite(viewer,type,256,PETSC_CHAR,PETSC_FALSE);CHKERRQ(ierr); 234bbd56ea5SKarl Rupp } else if (adapt->ops->view) { 235f2c2a1b9SBarry Smith ierr = (*adapt->ops->view)(adapt,viewer);CHKERRQ(ierr); 236f2c2a1b9SBarry Smith } 23784df9cb4SJed Brown PetscFunctionReturn(0); 23884df9cb4SJed Brown } 23984df9cb4SJed Brown 24084df9cb4SJed Brown #undef __FUNCT__ 241099af0a3SLisandro Dalcin #define __FUNCT__ "TSAdaptReset" 242099af0a3SLisandro Dalcin /*@ 243099af0a3SLisandro Dalcin TSAdaptReset - Resets a TSAdapt context. 244099af0a3SLisandro Dalcin 245099af0a3SLisandro Dalcin Collective on TS 246099af0a3SLisandro Dalcin 247099af0a3SLisandro Dalcin Input Parameter: 248099af0a3SLisandro Dalcin . adapt - the TSAdapt context obtained from TSAdaptCreate() 249099af0a3SLisandro Dalcin 250099af0a3SLisandro Dalcin Level: developer 251099af0a3SLisandro Dalcin 252099af0a3SLisandro Dalcin .seealso: TSAdaptCreate(), TSAdaptDestroy() 253099af0a3SLisandro Dalcin @*/ 254099af0a3SLisandro Dalcin PetscErrorCode TSAdaptReset(TSAdapt adapt) 255099af0a3SLisandro Dalcin { 256099af0a3SLisandro Dalcin PetscErrorCode ierr; 257099af0a3SLisandro Dalcin 258099af0a3SLisandro Dalcin PetscFunctionBegin; 259099af0a3SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 260099af0a3SLisandro Dalcin if (adapt->ops->reset) {ierr = (*adapt->ops->reset)(adapt);CHKERRQ(ierr);} 261099af0a3SLisandro Dalcin PetscFunctionReturn(0); 262099af0a3SLisandro Dalcin } 263099af0a3SLisandro Dalcin 264099af0a3SLisandro Dalcin #undef __FUNCT__ 26584df9cb4SJed Brown #define __FUNCT__ "TSAdaptDestroy" 26684df9cb4SJed Brown PetscErrorCode TSAdaptDestroy(TSAdapt *adapt) 26784df9cb4SJed Brown { 26884df9cb4SJed Brown PetscErrorCode ierr; 26984df9cb4SJed Brown 27084df9cb4SJed Brown PetscFunctionBegin; 27184df9cb4SJed Brown if (!*adapt) PetscFunctionReturn(0); 27284df9cb4SJed Brown PetscValidHeaderSpecific(*adapt,TSADAPT_CLASSID,1); 2734782b174SLisandro Dalcin if (--((PetscObject)(*adapt))->refct > 0) {*adapt = NULL; PetscFunctionReturn(0);} 274099af0a3SLisandro Dalcin 275099af0a3SLisandro Dalcin ierr = TSAdaptReset(*adapt);CHKERRQ(ierr); 276099af0a3SLisandro Dalcin 27784df9cb4SJed Brown if ((*adapt)->ops->destroy) {ierr = (*(*adapt)->ops->destroy)(*adapt);CHKERRQ(ierr);} 2781c3436cfSJed Brown ierr = PetscViewerDestroy(&(*adapt)->monitor);CHKERRQ(ierr); 27984df9cb4SJed Brown ierr = PetscHeaderDestroy(adapt);CHKERRQ(ierr); 28084df9cb4SJed Brown PetscFunctionReturn(0); 28184df9cb4SJed Brown } 28284df9cb4SJed Brown 28384df9cb4SJed Brown #undef __FUNCT__ 2841c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetMonitor" 2851c3436cfSJed Brown /*@ 2861c3436cfSJed Brown TSAdaptSetMonitor - Monitor the choices made by the adaptive controller 2871c3436cfSJed Brown 2881c3436cfSJed Brown Collective on TSAdapt 2891c3436cfSJed Brown 2901c3436cfSJed Brown Input Arguments: 2911c3436cfSJed Brown + adapt - adaptive controller context 2921c3436cfSJed Brown - flg - PETSC_TRUE to active a monitor, PETSC_FALSE to disable 2931c3436cfSJed Brown 2941c3436cfSJed Brown Level: intermediate 2951c3436cfSJed Brown 2961c3436cfSJed Brown .seealso: TSAdaptChoose() 2971c3436cfSJed Brown @*/ 2981c3436cfSJed Brown PetscErrorCode TSAdaptSetMonitor(TSAdapt adapt,PetscBool flg) 2991c3436cfSJed Brown { 3001c3436cfSJed Brown PetscErrorCode ierr; 3011c3436cfSJed Brown 3021c3436cfSJed Brown PetscFunctionBegin; 3034782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 3044782b174SLisandro Dalcin PetscValidLogicalCollectiveBool(adapt,flg,2); 3051c3436cfSJed Brown if (flg) { 306ce94432eSBarry Smith if (!adapt->monitor) {ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)adapt),"stdout",&adapt->monitor);CHKERRQ(ierr);} 3071c3436cfSJed Brown } else { 3081c3436cfSJed Brown ierr = PetscViewerDestroy(&adapt->monitor);CHKERRQ(ierr); 3091c3436cfSJed Brown } 3101c3436cfSJed Brown PetscFunctionReturn(0); 3111c3436cfSJed Brown } 3121c3436cfSJed Brown 3131c3436cfSJed Brown #undef __FUNCT__ 3140873808bSJed Brown #define __FUNCT__ "TSAdaptSetCheckStage" 3150873808bSJed Brown /*@C 3160873808bSJed Brown TSAdaptSetCheckStage - set a callback to check convergence for a stage 3170873808bSJed Brown 3180873808bSJed Brown Logically collective on TSAdapt 3190873808bSJed Brown 3200873808bSJed Brown Input Arguments: 3210873808bSJed Brown + adapt - adaptive controller context 3220873808bSJed Brown - func - stage check function 3230873808bSJed Brown 3240873808bSJed Brown Arguments of func: 3250873808bSJed Brown $ PetscErrorCode func(TSAdapt adapt,TS ts,PetscBool *accept) 3260873808bSJed Brown 3270873808bSJed Brown + adapt - adaptive controller context 3280873808bSJed Brown . ts - time stepping context 3290873808bSJed Brown - accept - pending choice of whether to accept, can be modified by this routine 3300873808bSJed Brown 3310873808bSJed Brown Level: advanced 3320873808bSJed Brown 3330873808bSJed Brown .seealso: TSAdaptChoose() 3340873808bSJed Brown @*/ 3350873808bSJed Brown PetscErrorCode TSAdaptSetCheckStage(TSAdapt adapt,PetscErrorCode (*func)(TSAdapt,TS,PetscBool*)) 3360873808bSJed Brown { 3370873808bSJed Brown 3380873808bSJed Brown PetscFunctionBegin; 3390873808bSJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 34068ae941aSLisandro Dalcin adapt->checkstage = func; 3410873808bSJed Brown PetscFunctionReturn(0); 3420873808bSJed Brown } 3430873808bSJed Brown 3440873808bSJed Brown #undef __FUNCT__ 3451c3436cfSJed Brown #define __FUNCT__ "TSAdaptSetStepLimits" 3461c3436cfSJed Brown /*@ 3471c3436cfSJed Brown TSAdaptSetStepLimits - Set minimum and maximum step sizes to be considered by the controller 3481c3436cfSJed Brown 3491c3436cfSJed Brown Logically Collective 3501c3436cfSJed Brown 3511c3436cfSJed Brown Input Arguments: 352552698daSJed Brown + adapt - time step adaptivity context, usually gotten with TSGetAdapt() 3531c3436cfSJed Brown . hmin - minimum time step 3541c3436cfSJed Brown - hmax - maximum time step 3551c3436cfSJed Brown 3561c3436cfSJed Brown Options Database Keys: 3571c3436cfSJed Brown + -ts_adapt_dt_min - minimum time step 3581c3436cfSJed Brown - -ts_adapt_dt_max - maximum time step 3591c3436cfSJed Brown 3601c3436cfSJed Brown Level: intermediate 3611c3436cfSJed Brown 3621c3436cfSJed Brown .seealso: TSAdapt 3631c3436cfSJed Brown @*/ 3641c3436cfSJed Brown PetscErrorCode TSAdaptSetStepLimits(TSAdapt adapt,PetscReal hmin,PetscReal hmax) 3651c3436cfSJed Brown { 3661c3436cfSJed Brown 3671c3436cfSJed Brown PetscFunctionBegin; 3684782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 3691c3436cfSJed Brown if (hmin != PETSC_DECIDE) adapt->dt_min = hmin; 3701c3436cfSJed Brown if (hmax != PETSC_DECIDE) adapt->dt_max = hmax; 3711c3436cfSJed Brown PetscFunctionReturn(0); 3721c3436cfSJed Brown } 3731c3436cfSJed Brown 3741c3436cfSJed Brown #undef __FUNCT__ 37584df9cb4SJed Brown #define __FUNCT__ "TSAdaptSetFromOptions" 376e55864a3SBarry Smith /* 37784df9cb4SJed Brown TSAdaptSetFromOptions - Sets various TSAdapt parameters from user options. 37884df9cb4SJed Brown 37984df9cb4SJed Brown Collective on TSAdapt 38084df9cb4SJed Brown 38184df9cb4SJed Brown Input Parameter: 38284df9cb4SJed Brown . adapt - the TSAdapt context 38384df9cb4SJed Brown 38484df9cb4SJed Brown Options Database Keys: 38584df9cb4SJed Brown . -ts_adapt_type <type> - basic 38684df9cb4SJed Brown 38784df9cb4SJed Brown Level: advanced 38884df9cb4SJed Brown 38984df9cb4SJed Brown Notes: 39084df9cb4SJed Brown This function is automatically called by TSSetFromOptions() 39184df9cb4SJed Brown 392552698daSJed Brown .keywords: TS, TSGetAdapt(), TSAdaptSetType() 39384df9cb4SJed Brown 39484df9cb4SJed Brown .seealso: TSGetType() 395e55864a3SBarry Smith */ 3968c34d3f5SBarry Smith PetscErrorCode TSAdaptSetFromOptions(PetscOptions *PetscOptionsObject,TSAdapt adapt) 39784df9cb4SJed Brown { 39884df9cb4SJed Brown PetscErrorCode ierr; 39984df9cb4SJed Brown char type[256] = TSADAPTBASIC; 4001c3436cfSJed Brown PetscBool set,flg; 40184df9cb4SJed Brown 40284df9cb4SJed Brown PetscFunctionBegin; 4034782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 40484df9cb4SJed Brown /* This should use PetscOptionsBegin() if/when this becomes an object used outside of TS, but currently this 40584df9cb4SJed Brown * function can only be called from inside TSSetFromOptions_GL() */ 406e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"TS Adaptivity options");CHKERRQ(ierr); 407a264d7a6SBarry Smith ierr = PetscOptionsFList("-ts_adapt_type","Algorithm to use for adaptivity","TSAdaptSetType",TSAdaptList, 4088caf3d72SBarry Smith ((PetscObject)adapt)->type_name ? ((PetscObject)adapt)->type_name : type,type,sizeof(type),&flg);CHKERRQ(ierr); 40984df9cb4SJed Brown if (flg || !((PetscObject)adapt)->type_name) { 41084df9cb4SJed Brown ierr = TSAdaptSetType(adapt,type);CHKERRQ(ierr); 41184df9cb4SJed Brown } 4120298fd71SBarry Smith ierr = PetscOptionsReal("-ts_adapt_dt_min","Minimum time step considered","TSAdaptSetStepLimits",adapt->dt_min,&adapt->dt_min,NULL);CHKERRQ(ierr); 4130298fd71SBarry Smith ierr = PetscOptionsReal("-ts_adapt_dt_max","Maximum time step considered","TSAdaptSetStepLimits",adapt->dt_max,&adapt->dt_max,NULL);CHKERRQ(ierr); 4140298fd71SBarry 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); 4151c3436cfSJed Brown ierr = PetscOptionsBool("-ts_adapt_monitor","Print choices made by adaptive controller","TSAdaptSetMonitor",adapt->monitor ? PETSC_TRUE : PETSC_FALSE,&flg,&set);CHKERRQ(ierr); 4167619abb3SShri ierr = PetscOptionsEnum("-ts_adapt_wnormtype","Type of norm computed for error estimation","",NormTypes,(PetscEnum)adapt->wnormtype,(PetscEnum*)&adapt->wnormtype,NULL);CHKERRQ(ierr); 4177619abb3SShri if (adapt->wnormtype != NORM_2 && adapt->wnormtype != NORM_INFINITY) SETERRQ(PetscObjectComm((PetscObject)adapt),PETSC_ERR_SUP,"Only 2-norm and infinite norm supported"); 4181c3436cfSJed Brown if (set) {ierr = TSAdaptSetMonitor(adapt,flg);CHKERRQ(ierr);} 419e55864a3SBarry Smith if (adapt->ops->setfromoptions) {ierr = (*adapt->ops->setfromoptions)(PetscOptionsObject,adapt);CHKERRQ(ierr);} 42084df9cb4SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 42184df9cb4SJed Brown PetscFunctionReturn(0); 42284df9cb4SJed Brown } 42384df9cb4SJed Brown 42484df9cb4SJed Brown #undef __FUNCT__ 42584df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidatesClear" 42684df9cb4SJed Brown /*@ 42784df9cb4SJed Brown TSAdaptCandidatesClear - clear any previously set candidate schemes 42884df9cb4SJed Brown 42984df9cb4SJed Brown Logically Collective 43084df9cb4SJed Brown 43184df9cb4SJed Brown Input Argument: 43284df9cb4SJed Brown . adapt - adaptive controller 43384df9cb4SJed Brown 43484df9cb4SJed Brown Level: developer 43584df9cb4SJed Brown 43684df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCreate(), TSAdaptCandidateAdd(), TSAdaptChoose() 43784df9cb4SJed Brown @*/ 43884df9cb4SJed Brown PetscErrorCode TSAdaptCandidatesClear(TSAdapt adapt) 43984df9cb4SJed Brown { 44084df9cb4SJed Brown PetscErrorCode ierr; 44184df9cb4SJed Brown 44284df9cb4SJed Brown PetscFunctionBegin; 4434782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 44484df9cb4SJed Brown ierr = PetscMemzero(&adapt->candidates,sizeof(adapt->candidates));CHKERRQ(ierr); 44584df9cb4SJed Brown PetscFunctionReturn(0); 44684df9cb4SJed Brown } 44784df9cb4SJed Brown 44884df9cb4SJed Brown #undef __FUNCT__ 44984df9cb4SJed Brown #define __FUNCT__ "TSAdaptCandidateAdd" 45084df9cb4SJed Brown /*@C 45184df9cb4SJed Brown TSAdaptCandidateAdd - add a candidate scheme for the adaptive controller to select from 45284df9cb4SJed Brown 45384df9cb4SJed Brown Logically Collective 45484df9cb4SJed Brown 45584df9cb4SJed Brown Input Arguments: 456552698daSJed Brown + adapt - time step adaptivity context, obtained with TSGetAdapt() or TSAdaptCreate() 45784df9cb4SJed Brown . name - name of the candidate scheme to add 45884df9cb4SJed Brown . order - order of the candidate scheme 45984df9cb4SJed Brown . stageorder - stage order of the candidate scheme 4608d59e960SJed Brown . ccfl - stability coefficient relative to explicit Euler, used for CFL constraints 46184df9cb4SJed Brown . cost - relative measure of the amount of work required for the candidate scheme 46284df9cb4SJed Brown - inuse - indicates that this scheme is the one currently in use, this flag can only be set for one scheme 46384df9cb4SJed Brown 46484df9cb4SJed Brown Note: 46584df9cb4SJed Brown This routine is not available in Fortran. 46684df9cb4SJed Brown 46784df9cb4SJed Brown Level: developer 46884df9cb4SJed Brown 46984df9cb4SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptChoose() 47084df9cb4SJed Brown @*/ 4718d59e960SJed Brown PetscErrorCode TSAdaptCandidateAdd(TSAdapt adapt,const char name[],PetscInt order,PetscInt stageorder,PetscReal ccfl,PetscReal cost,PetscBool inuse) 47284df9cb4SJed Brown { 47384df9cb4SJed Brown PetscInt c; 47484df9cb4SJed Brown 47584df9cb4SJed Brown PetscFunctionBegin; 47684df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 477ce94432eSBarry Smith if (order < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Classical order %D must be a positive integer",order); 47884df9cb4SJed Brown if (inuse) { 479ce94432eSBarry 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()"); 48084df9cb4SJed Brown adapt->candidates.inuse_set = PETSC_TRUE; 48184df9cb4SJed Brown } 4821c3436cfSJed Brown /* first slot if this is the current scheme, otherwise the next available slot */ 4831c3436cfSJed Brown c = inuse ? 0 : !adapt->candidates.inuse_set + adapt->candidates.n; 484bbd56ea5SKarl Rupp 48584df9cb4SJed Brown adapt->candidates.name[c] = name; 48684df9cb4SJed Brown adapt->candidates.order[c] = order; 48784df9cb4SJed Brown adapt->candidates.stageorder[c] = stageorder; 4888d59e960SJed Brown adapt->candidates.ccfl[c] = ccfl; 48984df9cb4SJed Brown adapt->candidates.cost[c] = cost; 49084df9cb4SJed Brown adapt->candidates.n++; 49184df9cb4SJed Brown PetscFunctionReturn(0); 49284df9cb4SJed Brown } 49384df9cb4SJed Brown 49484df9cb4SJed Brown #undef __FUNCT__ 4958d59e960SJed Brown #define __FUNCT__ "TSAdaptCandidatesGet" 4968d59e960SJed Brown /*@C 4978d59e960SJed Brown TSAdaptCandidatesGet - Get the list of candidate orders of accuracy and cost 4988d59e960SJed Brown 4998d59e960SJed Brown Not Collective 5008d59e960SJed Brown 5018d59e960SJed Brown Input Arguments: 5028d59e960SJed Brown . adapt - time step adaptivity context 5038d59e960SJed Brown 5048d59e960SJed Brown Output Arguments: 5058d59e960SJed Brown + n - number of candidate schemes, always at least 1 5068d59e960SJed Brown . order - the order of each candidate scheme 5078d59e960SJed Brown . stageorder - the stage order of each candidate scheme 5088d59e960SJed Brown . ccfl - the CFL coefficient of each scheme 5098d59e960SJed Brown - cost - the relative cost of each scheme 5108d59e960SJed Brown 5118d59e960SJed Brown Level: developer 5128d59e960SJed Brown 5138d59e960SJed Brown Note: 5148d59e960SJed Brown The current scheme is always returned in the first slot 5158d59e960SJed Brown 5168d59e960SJed Brown .seealso: TSAdaptCandidatesClear(), TSAdaptCandidateAdd(), TSAdaptChoose() 5178d59e960SJed Brown @*/ 5188d59e960SJed Brown PetscErrorCode TSAdaptCandidatesGet(TSAdapt adapt,PetscInt *n,const PetscInt **order,const PetscInt **stageorder,const PetscReal **ccfl,const PetscReal **cost) 5198d59e960SJed Brown { 5208d59e960SJed Brown PetscFunctionBegin; 5218d59e960SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 5228d59e960SJed Brown if (n) *n = adapt->candidates.n; 5238d59e960SJed Brown if (order) *order = adapt->candidates.order; 5248d59e960SJed Brown if (stageorder) *stageorder = adapt->candidates.stageorder; 5258d59e960SJed Brown if (ccfl) *ccfl = adapt->candidates.ccfl; 5268d59e960SJed Brown if (cost) *cost = adapt->candidates.cost; 5278d59e960SJed Brown PetscFunctionReturn(0); 5288d59e960SJed Brown } 5298d59e960SJed Brown 5308d59e960SJed Brown #undef __FUNCT__ 53184df9cb4SJed Brown #define __FUNCT__ "TSAdaptChoose" 53284df9cb4SJed Brown /*@C 53384df9cb4SJed Brown TSAdaptChoose - choose which method and step size to use for the next step 53484df9cb4SJed Brown 53584df9cb4SJed Brown Logically Collective 53684df9cb4SJed Brown 53784df9cb4SJed Brown Input Arguments: 53884df9cb4SJed Brown + adapt - adaptive contoller 53984df9cb4SJed Brown - h - current step size 54084df9cb4SJed Brown 54184df9cb4SJed Brown Output Arguments: 54284df9cb4SJed Brown + next_sc - scheme to use for the next step 54384df9cb4SJed Brown . next_h - step size to use for the next step 54484df9cb4SJed Brown - accept - PETSC_TRUE to accept the current step, PETSC_FALSE to repeat the current step with the new step size 54584df9cb4SJed Brown 5461c3436cfSJed Brown Note: 5471c3436cfSJed Brown The input value of parameter accept is retained from the last time step, so it will be PETSC_FALSE if the step is 5481c3436cfSJed Brown being retried after an initial rejection. 5491c3436cfSJed Brown 55084df9cb4SJed Brown Level: developer 55184df9cb4SJed Brown 55284df9cb4SJed Brown .seealso: TSAdapt, TSAdaptCandidatesClear(), TSAdaptCandidateAdd() 55384df9cb4SJed Brown @*/ 55484df9cb4SJed Brown PetscErrorCode TSAdaptChoose(TSAdapt adapt,TS ts,PetscReal h,PetscInt *next_sc,PetscReal *next_h,PetscBool *accept) 55584df9cb4SJed Brown { 55684df9cb4SJed Brown PetscErrorCode ierr; 5570b99f514SJed Brown PetscReal wlte = -1.0; 55884df9cb4SJed Brown 55984df9cb4SJed Brown PetscFunctionBegin; 56084df9cb4SJed Brown PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 56184df9cb4SJed Brown PetscValidHeaderSpecific(ts,TS_CLASSID,2); 56284df9cb4SJed Brown PetscValidIntPointer(next_sc,4); 56384df9cb4SJed Brown PetscValidPointer(next_h,5); 56484df9cb4SJed Brown PetscValidIntPointer(accept,6); 565ce94432eSBarry Smith if (adapt->candidates.n < 1) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_WRONGSTATE,"%D candidates have been registered",adapt->candidates.n); 566ce94432eSBarry 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); 5670b99f514SJed Brown ierr = (*adapt->ops->choose)(adapt,ts,h,next_sc,next_h,accept,&wlte);CHKERRQ(ierr); 56849354f04SShri Abhyankar if (*accept && ts->exact_final_time == TS_EXACTFINALTIME_MATCHSTEP) { 56949354f04SShri Abhyankar /* Reduce time step if it overshoots max time */ 57049354f04SShri Abhyankar PetscReal max_time = ts->max_time; 57149354f04SShri Abhyankar PetscReal next_dt = 0.0; 57249354f04SShri Abhyankar if (ts->ptime + ts->time_step + *next_h >= max_time) { 57349354f04SShri Abhyankar next_dt = max_time - (ts->ptime + ts->time_step); 5748709b12bSShri Abhyankar if (next_dt > PETSC_SMALL) *next_h = next_dt; 57549354f04SShri Abhyankar else ts->reason = TS_CONVERGED_TIME; 57649354f04SShri Abhyankar } 57749354f04SShri Abhyankar } 578ce94432eSBarry 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); 57957622a8eSBarry Smith if (!(*next_h > 0.)) SETERRQ1(PetscObjectComm((PetscObject)adapt),PETSC_ERR_ARG_OUTOFRANGE,"Computed step size %g must be positive",(double)*next_h); 5801c3436cfSJed Brown 5811c3436cfSJed Brown if (adapt->monitor) { 5821c3436cfSJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 5830b99f514SJed Brown if (wlte < 0) { 58448c19aefSLisandro Dalcin ierr = PetscViewerASCIIPrintf(adapt->monitor," TSAdapt '%s': step %3D %s t=%-11g+%10.3e family='%s' scheme=%D:'%s' dt=%-10.3e\n",((PetscObject)adapt)->type_name,ts->steps,*accept ? "accepted" : "rejected",(double)ts->ptime,(double)h,((PetscObject)ts)->type_name,*next_sc,adapt->candidates.name[*next_sc],(double)*next_h);CHKERRQ(ierr); 5850b99f514SJed Brown } else { 5866de24e2aSJed 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); 5870b99f514SJed Brown } 5881c3436cfSJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 5891c3436cfSJed Brown } 59084df9cb4SJed Brown PetscFunctionReturn(0); 59184df9cb4SJed Brown } 59284df9cb4SJed Brown 59384df9cb4SJed Brown #undef __FUNCT__ 59497335746SJed Brown #define __FUNCT__ "TSAdaptCheckStage" 59597335746SJed Brown /*@ 59697335746SJed Brown TSAdaptCheckStage - checks whether to accept a stage, (e.g. reject and change time step size if nonlinear solve fails) 59797335746SJed Brown 59897335746SJed Brown Collective 59997335746SJed Brown 60097335746SJed Brown Input Arguments: 60197335746SJed Brown + adapt - adaptive controller context 60297335746SJed Brown - ts - time stepper 60397335746SJed Brown 60497335746SJed Brown Output Arguments: 60597335746SJed Brown . accept - PETSC_TRUE to accept the stage, PETSC_FALSE to reject 60697335746SJed Brown 60797335746SJed Brown Level: developer 60897335746SJed Brown 60997335746SJed Brown .seealso: 61097335746SJed Brown @*/ 611*cb9d8021SPierre Barbier de Reuille PetscErrorCode TSAdaptCheckStage(TSAdapt adapt,TS ts,PetscReal t,PetscInt stage,Vec* Y,PetscBool *accept) 61297335746SJed Brown { 61397335746SJed Brown PetscErrorCode ierr; 61497335746SJed Brown SNES snes; 61597335746SJed Brown SNESConvergedReason snesreason; 616*cb9d8021SPierre Barbier de Reuille PetscReal dt,new_dt; 61797335746SJed Brown 61897335746SJed Brown PetscFunctionBegin; 6194782b174SLisandro Dalcin PetscValidHeaderSpecific(adapt,TSADAPT_CLASSID,1); 6204782b174SLisandro Dalcin PetscValidHeaderSpecific(ts,TS_CLASSID,2); 6214782b174SLisandro Dalcin PetscValidIntPointer(accept,3); 62297335746SJed Brown *accept = PETSC_TRUE; 62397335746SJed Brown ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr); 62497335746SJed Brown ierr = SNESGetConvergedReason(snes,&snesreason);CHKERRQ(ierr); 62597335746SJed Brown if (snesreason < 0) { 62697335746SJed Brown *accept = PETSC_FALSE; 62797335746SJed Brown ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 6286de24e2aSJed Brown if (++ts->num_snes_failures >= ts->max_snes_failures && ts->max_snes_failures > 0) { 62997335746SJed Brown ts->reason = TS_DIVERGED_NONLINEAR_SOLVE; 63048c19aefSLisandro 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); 63197335746SJed Brown if (adapt->monitor) { 63297335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 63348c19aefSLisandro Dalcin ierr = PetscViewerASCIIPrintf(adapt->monitor," TSAdapt '%s': step %3D stage rejected t=%-11g+%10.3e, nonlinear solve failures %D greater than current TS allowed\n",((PetscObject)adapt)->type_name,ts->steps,(double)ts->ptime,dt,ts->num_snes_failures);CHKERRQ(ierr); 63497335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 63597335746SJed Brown } 63697335746SJed Brown } else { 63797335746SJed Brown new_dt = dt*adapt->scale_solve_failed; 638*cb9d8021SPierre Barbier de Reuille } 639*cb9d8021SPierre Barbier de Reuille } else { 640*cb9d8021SPierre Barbier de Reuille ierr = TSFunctionDomainError(ts,t,stage,Y,accept);CHKERRQ(ierr); 641*cb9d8021SPierre Barbier de Reuille if(*accept && adapt->checkstage) { 642*cb9d8021SPierre Barbier de Reuille ierr = (*adapt->checkstage)(adapt,ts,accept);CHKERRQ(ierr); 643*cb9d8021SPierre Barbier de Reuille } 644*cb9d8021SPierre Barbier de Reuille } 645*cb9d8021SPierre Barbier de Reuille 646*cb9d8021SPierre Barbier de Reuille if(!(*accept)) { 647*cb9d8021SPierre Barbier de Reuille new_dt = dt*adapt->scale_solve_failed; 64897335746SJed Brown ierr = TSSetTimeStep(ts,new_dt);CHKERRQ(ierr); 64997335746SJed Brown if (adapt->monitor) { 65097335746SJed Brown ierr = PetscViewerASCIIAddTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 65197335746SJed 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); 65297335746SJed Brown ierr = PetscViewerASCIISubtractTab(adapt->monitor,((PetscObject)adapt)->tablevel);CHKERRQ(ierr); 65397335746SJed Brown } 65497335746SJed Brown } 655*cb9d8021SPierre Barbier de Reuille 65697335746SJed Brown PetscFunctionReturn(0); 65797335746SJed Brown } 65897335746SJed Brown 6590873808bSJed Brown 6600873808bSJed Brown 66197335746SJed Brown #undef __FUNCT__ 66284df9cb4SJed Brown #define __FUNCT__ "TSAdaptCreate" 66384df9cb4SJed Brown /*@ 66484df9cb4SJed Brown TSAdaptCreate - create an adaptive controller context for time stepping 66584df9cb4SJed Brown 66684df9cb4SJed Brown Collective on MPI_Comm 66784df9cb4SJed Brown 66884df9cb4SJed Brown Input Parameter: 66984df9cb4SJed Brown . comm - The communicator 67084df9cb4SJed Brown 67184df9cb4SJed Brown Output Parameter: 67284df9cb4SJed Brown . adapt - new TSAdapt object 67384df9cb4SJed Brown 67484df9cb4SJed Brown Level: developer 67584df9cb4SJed Brown 67684df9cb4SJed Brown Notes: 67784df9cb4SJed Brown TSAdapt creation is handled by TS, so users should not need to call this function. 67884df9cb4SJed Brown 67984df9cb4SJed Brown .keywords: TSAdapt, create 680552698daSJed Brown .seealso: TSGetAdapt(), TSAdaptSetType(), TSAdaptDestroy() 68184df9cb4SJed Brown @*/ 68284df9cb4SJed Brown PetscErrorCode TSAdaptCreate(MPI_Comm comm,TSAdapt *inadapt) 68384df9cb4SJed Brown { 68484df9cb4SJed Brown PetscErrorCode ierr; 68584df9cb4SJed Brown TSAdapt adapt; 68684df9cb4SJed Brown 68784df9cb4SJed Brown PetscFunctionBegin; 6883b3bcf4cSLisandro Dalcin PetscValidPointer(inadapt,1); 6893b3bcf4cSLisandro Dalcin *inadapt = NULL; 6903b3bcf4cSLisandro Dalcin ierr = TSAdaptInitializePackage();CHKERRQ(ierr); 6913b3bcf4cSLisandro Dalcin 69267c2884eSBarry Smith ierr = PetscHeaderCreate(adapt,_p_TSAdapt,struct _TSAdaptOps,TSADAPT_CLASSID,"TSAdapt","General Linear adaptivity","TS",comm,TSAdaptDestroy,TSAdaptView);CHKERRQ(ierr); 6931c3436cfSJed Brown 6941c3436cfSJed Brown adapt->dt_min = 1e-20; 6951c3436cfSJed Brown adapt->dt_max = 1e50; 69697335746SJed Brown adapt->scale_solve_failed = 0.25; 6977619abb3SShri adapt->wnormtype = NORM_2; 6981c3436cfSJed Brown 69984df9cb4SJed Brown *inadapt = adapt; 70084df9cb4SJed Brown PetscFunctionReturn(0); 70184df9cb4SJed Brown } 702