1af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2af0996ceSBarry Smith #include <petsc/private/dmimpl.h> 324989b8cSPeter Brune 4d74926cbSBarry Smith static PetscErrorCode DMTSDestroy(DMTS *kdm) 5d74926cbSBarry Smith { 6d74926cbSBarry Smith PetscErrorCode ierr; 7d74926cbSBarry Smith 8d74926cbSBarry Smith PetscFunctionBegin; 9d74926cbSBarry Smith if (!*kdm) PetscFunctionReturn(0); 10d74926cbSBarry Smith PetscValidHeaderSpecific((*kdm),DMTS_CLASSID,1); 11c793f718SLisandro Dalcin if (--((PetscObject)(*kdm))->refct > 0) {*kdm = NULL; PetscFunctionReturn(0);} 12d74926cbSBarry Smith if ((*kdm)->ops->destroy) {ierr = ((*kdm)->ops->destroy)(*kdm);CHKERRQ(ierr);} 13d74926cbSBarry Smith ierr = PetscHeaderDestroy(kdm);CHKERRQ(ierr); 14d74926cbSBarry Smith PetscFunctionReturn(0); 15d74926cbSBarry Smith } 16d74926cbSBarry Smith 172d53ad75SBarry Smith PetscErrorCode DMTSLoad(DMTS kdm,PetscViewer viewer) 182d53ad75SBarry Smith { 192d53ad75SBarry Smith PetscErrorCode ierr; 202d53ad75SBarry Smith 212d53ad75SBarry Smith PetscFunctionBegin; 22060da220SMatthew G. Knepley ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunction,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr); 23060da220SMatthew G. Knepley ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionview,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr); 24060da220SMatthew G. Knepley ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ifunctionload,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr); 25ad6bc421SBarry Smith if (kdm->ops->ifunctionload) { 26ad6bc421SBarry Smith ierr = (*kdm->ops->ifunctionload)(&kdm->ifunctionctx,viewer);CHKERRQ(ierr); 27ad6bc421SBarry Smith } 28060da220SMatthew G. Knepley ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobian,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr); 29060da220SMatthew G. Knepley ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianview,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr); 30060da220SMatthew G. Knepley ierr = PetscViewerBinaryRead(viewer,&kdm->ops->ijacobianload,1,NULL,PETSC_FUNCTION);CHKERRQ(ierr); 31ad6bc421SBarry Smith if (kdm->ops->ijacobianload) { 32ad6bc421SBarry Smith ierr = (*kdm->ops->ijacobianload)(&kdm->ijacobianctx,viewer);CHKERRQ(ierr); 33ad6bc421SBarry Smith } 342d53ad75SBarry Smith PetscFunctionReturn(0); 352d53ad75SBarry Smith } 362d53ad75SBarry Smith 372d53ad75SBarry Smith PetscErrorCode DMTSView(DMTS kdm,PetscViewer viewer) 382d53ad75SBarry Smith { 392d53ad75SBarry Smith PetscErrorCode ierr; 402d53ad75SBarry Smith PetscBool isascii,isbinary; 412d53ad75SBarry Smith 422d53ad75SBarry Smith PetscFunctionBegin; 432d53ad75SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr); 442d53ad75SBarry Smith ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERBINARY,&isbinary);CHKERRQ(ierr); 452d53ad75SBarry Smith if (isascii) { 46c7a10e08SBarry Smith #if defined(PETSC_SERIALIZE_FUNCTIONS) 472d53ad75SBarry Smith const char *fname; 482d53ad75SBarry Smith 492d53ad75SBarry Smith ierr = PetscFPTFind(kdm->ops->ifunction,&fname);CHKERRQ(ierr); 502d53ad75SBarry Smith if (fname) { 512d53ad75SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," IFunction used by TS: %s\n",fname);CHKERRQ(ierr); 522d53ad75SBarry Smith } 532d53ad75SBarry Smith ierr = PetscFPTFind(kdm->ops->ijacobian,&fname);CHKERRQ(ierr); 542d53ad75SBarry Smith if (fname) { 552d53ad75SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," IJacobian function used by TS: %s\n",fname);CHKERRQ(ierr); 562d53ad75SBarry Smith } 57c7a10e08SBarry Smith #endif 582d53ad75SBarry Smith } else if (isbinary) { 593964eb88SJed Brown struct { 603964eb88SJed Brown TSIFunction ifunction; 619200755eSBarry Smith } funcstruct; 629200755eSBarry Smith struct { 633964eb88SJed Brown PetscErrorCode (*ifunctionview)(void*,PetscViewer); 649200755eSBarry Smith } funcviewstruct; 659200755eSBarry Smith struct { 663964eb88SJed Brown PetscErrorCode (*ifunctionload)(void**,PetscViewer); 679200755eSBarry Smith } funcloadstruct; 683964eb88SJed Brown struct { 693964eb88SJed Brown TSIJacobian ijacobian; 709200755eSBarry Smith } jacstruct; 719200755eSBarry Smith struct { 723964eb88SJed Brown PetscErrorCode (*ijacobianview)(void*,PetscViewer); 739200755eSBarry Smith } jacviewstruct; 749200755eSBarry Smith struct { 753964eb88SJed Brown PetscErrorCode (*ijacobianload)(void**,PetscViewer); 769200755eSBarry Smith } jacloadstruct; 773964eb88SJed Brown 789200755eSBarry Smith funcstruct.ifunction = kdm->ops->ifunction; 799200755eSBarry Smith funcviewstruct.ifunctionview = kdm->ops->ifunctionview; 809200755eSBarry Smith funcloadstruct.ifunctionload = kdm->ops->ifunctionload; 81f253e43cSLisandro Dalcin ierr = PetscViewerBinaryWrite(viewer,&funcstruct,1,PETSC_FUNCTION);CHKERRQ(ierr); 82f253e43cSLisandro Dalcin ierr = PetscViewerBinaryWrite(viewer,&funcviewstruct,1,PETSC_FUNCTION);CHKERRQ(ierr); 83f253e43cSLisandro Dalcin ierr = PetscViewerBinaryWrite(viewer,&funcloadstruct,1,PETSC_FUNCTION);CHKERRQ(ierr); 84ad6bc421SBarry Smith if (kdm->ops->ifunctionview) { 85ad6bc421SBarry Smith ierr = (*kdm->ops->ifunctionview)(kdm->ifunctionctx,viewer);CHKERRQ(ierr); 86ad6bc421SBarry Smith } 879200755eSBarry Smith jacstruct.ijacobian = kdm->ops->ijacobian; 889200755eSBarry Smith jacviewstruct.ijacobianview = kdm->ops->ijacobianview; 899200755eSBarry Smith jacloadstruct.ijacobianload = kdm->ops->ijacobianload; 90f253e43cSLisandro Dalcin ierr = PetscViewerBinaryWrite(viewer,&jacstruct,1,PETSC_FUNCTION);CHKERRQ(ierr); 91f253e43cSLisandro Dalcin ierr = PetscViewerBinaryWrite(viewer,&jacviewstruct,1,PETSC_FUNCTION);CHKERRQ(ierr); 92f253e43cSLisandro Dalcin ierr = PetscViewerBinaryWrite(viewer,&jacloadstruct,1,PETSC_FUNCTION);CHKERRQ(ierr); 93ad6bc421SBarry Smith if (kdm->ops->ijacobianview) { 94ad6bc421SBarry Smith ierr = (*kdm->ops->ijacobianview)(kdm->ijacobianctx,viewer);CHKERRQ(ierr); 95ad6bc421SBarry Smith } 962d53ad75SBarry Smith } 972d53ad75SBarry Smith PetscFunctionReturn(0); 982d53ad75SBarry Smith } 992d53ad75SBarry Smith 100d74926cbSBarry Smith static PetscErrorCode DMTSCreate(MPI_Comm comm,DMTS *kdm) 101d74926cbSBarry Smith { 102d74926cbSBarry Smith PetscErrorCode ierr; 103d74926cbSBarry Smith 104d74926cbSBarry Smith PetscFunctionBegin; 105607a6623SBarry Smith ierr = TSInitializePackage();CHKERRQ(ierr); 10673107ff1SLisandro Dalcin ierr = PetscHeaderCreate(*kdm, DMTS_CLASSID, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, DMTSView);CHKERRQ(ierr); 107d74926cbSBarry Smith PetscFunctionReturn(0); 108d74926cbSBarry Smith } 10924989b8cSPeter Brune 1102a34c10cSBarry Smith /* Attaches the DMTS to the coarse level. 11124989b8cSPeter Brune * Under what conditions should we copy versus duplicate? 11224989b8cSPeter Brune */ 113942e3340SBarry Smith static PetscErrorCode DMCoarsenHook_DMTS(DM dm,DM dmc,void *ctx) 11424989b8cSPeter Brune { 11524989b8cSPeter Brune PetscErrorCode ierr; 11624989b8cSPeter Brune 11724989b8cSPeter Brune PetscFunctionBegin; 118942e3340SBarry Smith ierr = DMCopyDMTS(dm,dmc);CHKERRQ(ierr); 11924989b8cSPeter Brune PetscFunctionReturn(0); 12024989b8cSPeter Brune } 12124989b8cSPeter Brune 12224989b8cSPeter Brune /* This could restrict auxiliary information to the coarse level. 12324989b8cSPeter Brune */ 124942e3340SBarry Smith static PetscErrorCode DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx) 12524989b8cSPeter Brune { 12624989b8cSPeter Brune 12724989b8cSPeter Brune PetscFunctionBegin; 12824989b8cSPeter Brune PetscFunctionReturn(0); 12924989b8cSPeter Brune } 13024989b8cSPeter Brune 131258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_DMTS(DM dm,DM subdm,void *ctx) 132258e1594SPeter Brune { 133258e1594SPeter Brune PetscErrorCode ierr; 134258e1594SPeter Brune 135258e1594SPeter Brune PetscFunctionBegin; 136258e1594SPeter Brune ierr = DMCopyDMTS(dm,subdm);CHKERRQ(ierr); 137258e1594SPeter Brune PetscFunctionReturn(0); 138258e1594SPeter Brune } 139258e1594SPeter Brune 140258e1594SPeter Brune /* This could restrict auxiliary information to the coarse level. 141258e1594SPeter Brune */ 142258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 143258e1594SPeter Brune { 144258e1594SPeter Brune PetscFunctionBegin; 145258e1594SPeter Brune PetscFunctionReturn(0); 146258e1594SPeter Brune } 147258e1594SPeter Brune 148d74926cbSBarry Smith /*@C 149d74926cbSBarry Smith DMTSCopy - copies the information in a DMTS to another DMTS 150d74926cbSBarry Smith 151d74926cbSBarry Smith Not Collective 152d74926cbSBarry Smith 153*4165533cSJose E. Roman Input Parameters: 154d74926cbSBarry Smith + kdm - Original DMTS 155d74926cbSBarry Smith - nkdm - DMTS to receive the data, should have been created with DMTSCreate() 156d74926cbSBarry Smith 157d74926cbSBarry Smith Level: developer 158d74926cbSBarry Smith 159d74926cbSBarry Smith .seealso: DMTSCreate(), DMTSDestroy() 160d74926cbSBarry Smith @*/ 161d74926cbSBarry Smith PetscErrorCode DMTSCopy(DMTS kdm,DMTS nkdm) 16224989b8cSPeter Brune { 16324989b8cSPeter Brune PetscErrorCode ierr; 16424989b8cSPeter Brune 16524989b8cSPeter Brune PetscFunctionBegin; 166d74926cbSBarry Smith PetscValidHeaderSpecific(kdm,DMTS_CLASSID,1); 167d74926cbSBarry Smith PetscValidHeaderSpecific(nkdm,DMTS_CLASSID,2); 168d74926cbSBarry Smith nkdm->ops->rhsfunction = kdm->ops->rhsfunction; 169d74926cbSBarry Smith nkdm->ops->rhsjacobian = kdm->ops->rhsjacobian; 170d74926cbSBarry Smith nkdm->ops->ifunction = kdm->ops->ifunction; 171d74926cbSBarry Smith nkdm->ops->ijacobian = kdm->ops->ijacobian; 172efe9872eSLisandro Dalcin nkdm->ops->i2function = kdm->ops->i2function; 173efe9872eSLisandro Dalcin nkdm->ops->i2jacobian = kdm->ops->i2jacobian; 174d74926cbSBarry Smith nkdm->ops->solution = kdm->ops->solution; 175d74926cbSBarry Smith nkdm->ops->destroy = kdm->ops->destroy; 176d74926cbSBarry Smith nkdm->ops->duplicate = kdm->ops->duplicate; 177d74926cbSBarry Smith 178d74926cbSBarry Smith nkdm->rhsfunctionctx = kdm->rhsfunctionctx; 179d74926cbSBarry Smith nkdm->rhsjacobianctx = kdm->rhsjacobianctx; 180d74926cbSBarry Smith nkdm->ifunctionctx = kdm->ifunctionctx; 181d74926cbSBarry Smith nkdm->ijacobianctx = kdm->ijacobianctx; 182efe9872eSLisandro Dalcin nkdm->i2functionctx = kdm->i2functionctx; 183efe9872eSLisandro Dalcin nkdm->i2jacobianctx = kdm->i2jacobianctx; 184d74926cbSBarry Smith nkdm->solutionctx = kdm->solutionctx; 185d74926cbSBarry Smith 186d74926cbSBarry Smith nkdm->data = kdm->data; 187d74926cbSBarry Smith 188d74926cbSBarry Smith /* 189d74926cbSBarry Smith nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0]; 190d74926cbSBarry Smith nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1]; 191d74926cbSBarry Smith nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2]; 192d74926cbSBarry Smith */ 193d74926cbSBarry Smith 194d74926cbSBarry Smith /* implementation specific copy hooks */ 195d74926cbSBarry Smith if (kdm->ops->duplicate) {ierr = (*kdm->ops->duplicate)(kdm,nkdm);CHKERRQ(ierr);} 19624989b8cSPeter Brune PetscFunctionReturn(0); 19724989b8cSPeter Brune } 19824989b8cSPeter Brune 19924989b8cSPeter Brune /*@C 200942e3340SBarry Smith DMGetDMTS - get read-only private DMTS context from a DM 20124989b8cSPeter Brune 20224989b8cSPeter Brune Not Collective 20324989b8cSPeter Brune 204*4165533cSJose E. Roman Input Parameter: 20524989b8cSPeter Brune . dm - DM to be used with TS 20624989b8cSPeter Brune 207*4165533cSJose E. Roman Output Parameter: 208942e3340SBarry Smith . tsdm - private DMTS context 20924989b8cSPeter Brune 21024989b8cSPeter Brune Level: developer 21124989b8cSPeter Brune 21224989b8cSPeter Brune Notes: 213942e3340SBarry Smith Use DMGetDMTSWrite() if write access is needed. The DMTSSetXXX API should be used wherever possible. 21424989b8cSPeter Brune 215942e3340SBarry Smith .seealso: DMGetDMTSWrite() 21624989b8cSPeter Brune @*/ 217942e3340SBarry Smith PetscErrorCode DMGetDMTS(DM dm,DMTS *tsdm) 21824989b8cSPeter Brune { 21924989b8cSPeter Brune PetscErrorCode ierr; 22024989b8cSPeter Brune 22124989b8cSPeter Brune PetscFunctionBegin; 22224989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2232a34c10cSBarry Smith *tsdm = (DMTS) dm->dmts; 224d74926cbSBarry Smith if (!*tsdm) { 225942e3340SBarry Smith ierr = PetscInfo(dm,"Creating new DMTS\n");CHKERRQ(ierr); 226ce94432eSBarry Smith ierr = DMTSCreate(PetscObjectComm((PetscObject)dm),tsdm);CHKERRQ(ierr); 2272a34c10cSBarry Smith dm->dmts = (PetscObject) *tsdm; 2285c87d4f4SJunchao Zhang (*tsdm)->originaldm = dm; 2290298fd71SBarry Smith ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);CHKERRQ(ierr); 2300298fd71SBarry Smith ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);CHKERRQ(ierr); 23124989b8cSPeter Brune } 23224989b8cSPeter Brune PetscFunctionReturn(0); 23324989b8cSPeter Brune } 23424989b8cSPeter Brune 23524989b8cSPeter Brune /*@C 236942e3340SBarry Smith DMGetDMTSWrite - get write access to private DMTS context from a DM 23724989b8cSPeter Brune 23824989b8cSPeter Brune Not Collective 23924989b8cSPeter Brune 240*4165533cSJose E. Roman Input Parameter: 24124989b8cSPeter Brune . dm - DM to be used with TS 24224989b8cSPeter Brune 243*4165533cSJose E. Roman Output Parameter: 244942e3340SBarry Smith . tsdm - private DMTS context 24524989b8cSPeter Brune 24624989b8cSPeter Brune Level: developer 24724989b8cSPeter Brune 248942e3340SBarry Smith .seealso: DMGetDMTS() 24924989b8cSPeter Brune @*/ 250942e3340SBarry Smith PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm) 25124989b8cSPeter Brune { 25224989b8cSPeter Brune PetscErrorCode ierr; 253942e3340SBarry Smith DMTS sdm; 25424989b8cSPeter Brune 25524989b8cSPeter Brune PetscFunctionBegin; 25624989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 257942e3340SBarry Smith ierr = DMGetDMTS(dm,&sdm);CHKERRQ(ierr); 2585c87d4f4SJunchao Zhang if (!sdm->originaldm) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"DMTS has a NULL originaldm"); 25924989b8cSPeter Brune if (sdm->originaldm != dm) { /* Copy on write */ 2602a34c10cSBarry Smith DMTS oldsdm = sdm; 261942e3340SBarry Smith ierr = PetscInfo(dm,"Copying DMTS due to write\n");CHKERRQ(ierr); 262ce94432eSBarry Smith ierr = DMTSCreate(PetscObjectComm((PetscObject)dm),&sdm);CHKERRQ(ierr); 263d74926cbSBarry Smith ierr = DMTSCopy(oldsdm,sdm);CHKERRQ(ierr); 2642a34c10cSBarry Smith ierr = DMTSDestroy((DMTS*)&dm->dmts);CHKERRQ(ierr); 2652a34c10cSBarry Smith dm->dmts = (PetscObject) sdm; 2665c87d4f4SJunchao Zhang sdm->originaldm = dm; 26724989b8cSPeter Brune } 26824989b8cSPeter Brune *tsdm = sdm; 26924989b8cSPeter Brune PetscFunctionReturn(0); 27024989b8cSPeter Brune } 27124989b8cSPeter Brune 27224989b8cSPeter Brune /*@C 273942e3340SBarry Smith DMCopyDMTS - copies a DM context to a new DM 27424989b8cSPeter Brune 27524989b8cSPeter Brune Logically Collective 27624989b8cSPeter Brune 277*4165533cSJose E. Roman Input Parameters: 27824989b8cSPeter Brune + dmsrc - DM to obtain context from 27924989b8cSPeter Brune - dmdest - DM to add context to 28024989b8cSPeter Brune 28124989b8cSPeter Brune Level: developer 28224989b8cSPeter Brune 28324989b8cSPeter Brune Note: 28424989b8cSPeter Brune The context is copied by reference. This function does not ensure that a context exists. 28524989b8cSPeter Brune 286942e3340SBarry Smith .seealso: DMGetDMTS(), TSSetDM() 28724989b8cSPeter Brune @*/ 288942e3340SBarry Smith PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest) 28924989b8cSPeter Brune { 29024989b8cSPeter Brune PetscErrorCode ierr; 29124989b8cSPeter Brune 29224989b8cSPeter Brune PetscFunctionBegin; 29324989b8cSPeter Brune PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1); 29424989b8cSPeter Brune PetscValidHeaderSpecific(dmdest,DM_CLASSID,2); 2952a34c10cSBarry Smith ierr = DMTSDestroy((DMTS*)&dmdest->dmts);CHKERRQ(ierr); 2962a34c10cSBarry Smith dmdest->dmts = dmsrc->dmts; 2972a34c10cSBarry Smith ierr = PetscObjectReference(dmdest->dmts);CHKERRQ(ierr); 2980298fd71SBarry Smith ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,NULL);CHKERRQ(ierr); 2990298fd71SBarry Smith ierr = DMSubDomainHookAdd(dmdest,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,NULL);CHKERRQ(ierr); 30024989b8cSPeter Brune PetscFunctionReturn(0); 30124989b8cSPeter Brune } 30224989b8cSPeter Brune 30324989b8cSPeter Brune /*@C 30424989b8cSPeter Brune DMTSSetIFunction - set TS implicit function evaluation function 30524989b8cSPeter Brune 30624989b8cSPeter Brune Not Collective 30724989b8cSPeter Brune 308*4165533cSJose E. Roman Input Parameters: 30924989b8cSPeter Brune + dm - DM to be used with TS 310a96d6ef6SBarry Smith . func - function evaluating f(t,u,u_t) 31124989b8cSPeter Brune - ctx - context for residual evaluation 31224989b8cSPeter Brune 313a96d6ef6SBarry Smith Calling sequence of func: 314a96d6ef6SBarry Smith $ PetscErrorCode func(TS ts,PetscReal t,Vec u,Vec u_t,Vec F,ctx); 315a96d6ef6SBarry Smith 316a96d6ef6SBarry Smith + t - time at step/stage being solved 317a96d6ef6SBarry Smith . u - state vector 318a96d6ef6SBarry Smith . u_t - time derivative of state vector 319a96d6ef6SBarry Smith . F - function vector 320a96d6ef6SBarry Smith - ctx - [optional] user-defined context for matrix evaluation routine 321a96d6ef6SBarry Smith 32224989b8cSPeter Brune Level: advanced 32324989b8cSPeter Brune 32424989b8cSPeter Brune Note: 32524989b8cSPeter Brune TSSetFunction() is normally used, but it calls this function internally because the user context is actually 32624989b8cSPeter Brune associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 32724989b8cSPeter Brune not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 32824989b8cSPeter Brune 329a96d6ef6SBarry Smith .seealso: DMTSSetContext(), TSSetIFunction(), DMTSSetJacobian() 33024989b8cSPeter Brune @*/ 33124989b8cSPeter Brune PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx) 33224989b8cSPeter Brune { 33324989b8cSPeter Brune PetscErrorCode ierr; 334942e3340SBarry Smith DMTS tsdm; 33524989b8cSPeter Brune 33624989b8cSPeter Brune PetscFunctionBegin; 33724989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 338942e3340SBarry Smith ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 339d74926cbSBarry Smith if (func) tsdm->ops->ifunction = func; 34024989b8cSPeter Brune if (ctx) tsdm->ifunctionctx = ctx; 34124989b8cSPeter Brune PetscFunctionReturn(0); 34224989b8cSPeter Brune } 34324989b8cSPeter Brune 34424989b8cSPeter Brune /*@C 34524989b8cSPeter Brune DMTSGetIFunction - get TS implicit residual evaluation function 34624989b8cSPeter Brune 34724989b8cSPeter Brune Not Collective 34824989b8cSPeter Brune 349*4165533cSJose E. Roman Input Parameter: 35024989b8cSPeter Brune . dm - DM to be used with TS 35124989b8cSPeter Brune 352*4165533cSJose E. Roman Output Parameters: 35324989b8cSPeter Brune + func - function evaluation function, see TSSetIFunction() for calling sequence 35424989b8cSPeter Brune - ctx - context for residual evaluation 35524989b8cSPeter Brune 35624989b8cSPeter Brune Level: advanced 35724989b8cSPeter Brune 35824989b8cSPeter Brune Note: 35924989b8cSPeter Brune TSGetFunction() is normally used, but it calls this function internally because the user context is actually 36024989b8cSPeter Brune associated with the DM. 36124989b8cSPeter Brune 36224989b8cSPeter Brune .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction() 36324989b8cSPeter Brune @*/ 36424989b8cSPeter Brune PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx) 36524989b8cSPeter Brune { 36624989b8cSPeter Brune PetscErrorCode ierr; 367942e3340SBarry Smith DMTS tsdm; 36824989b8cSPeter Brune 36924989b8cSPeter Brune PetscFunctionBegin; 37024989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 371942e3340SBarry Smith ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 372d74926cbSBarry Smith if (func) *func = tsdm->ops->ifunction; 37324989b8cSPeter Brune if (ctx) *ctx = tsdm->ifunctionctx; 37424989b8cSPeter Brune PetscFunctionReturn(0); 37524989b8cSPeter Brune } 37624989b8cSPeter Brune 377efe9872eSLisandro Dalcin /*@C 378efe9872eSLisandro Dalcin DMTSSetI2Function - set TS implicit function evaluation function for 2nd order systems 379efe9872eSLisandro Dalcin 380efe9872eSLisandro Dalcin Not Collective 381efe9872eSLisandro Dalcin 382*4165533cSJose E. Roman Input Parameters: 383efe9872eSLisandro Dalcin + dm - DM to be used with TS 384a96d6ef6SBarry Smith . fun - function evaluation routine 385efe9872eSLisandro Dalcin - ctx - context for residual evaluation 386efe9872eSLisandro Dalcin 387a96d6ef6SBarry Smith Calling sequence of fun: 388a96d6ef6SBarry Smith $ PetscErrorCode fun(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,Vec F,ctx); 389a96d6ef6SBarry Smith 390a96d6ef6SBarry Smith + t - time at step/stage being solved 391a96d6ef6SBarry Smith . U - state vector 392a96d6ef6SBarry Smith . U_t - time derivative of state vector 393a96d6ef6SBarry Smith . U_tt - second time derivative of state vector 394a96d6ef6SBarry Smith . F - function vector 395a96d6ef6SBarry Smith - ctx - [optional] user-defined context for matrix evaluation routine (may be NULL) 396a96d6ef6SBarry Smith 397efe9872eSLisandro Dalcin Level: advanced 398efe9872eSLisandro Dalcin 399efe9872eSLisandro Dalcin Note: 400efe9872eSLisandro Dalcin TSSetI2Function() is normally used, but it calls this function internally because the user context is actually 401efe9872eSLisandro Dalcin associated with the DM. 402efe9872eSLisandro Dalcin 403efe9872eSLisandro Dalcin .seealso: TSSetI2Function() 404efe9872eSLisandro Dalcin @*/ 405efe9872eSLisandro Dalcin PetscErrorCode DMTSSetI2Function(DM dm,TSI2Function fun,void *ctx) 406efe9872eSLisandro Dalcin { 407efe9872eSLisandro Dalcin DMTS tsdm; 408efe9872eSLisandro Dalcin PetscErrorCode ierr; 409efe9872eSLisandro Dalcin 410efe9872eSLisandro Dalcin PetscFunctionBegin; 411efe9872eSLisandro Dalcin PetscValidHeaderSpecific(dm,DM_CLASSID,1); 412efe9872eSLisandro Dalcin ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 413efe9872eSLisandro Dalcin if (fun) tsdm->ops->i2function = fun; 414efe9872eSLisandro Dalcin if (ctx) tsdm->i2functionctx = ctx; 415efe9872eSLisandro Dalcin PetscFunctionReturn(0); 416efe9872eSLisandro Dalcin } 417efe9872eSLisandro Dalcin 418efe9872eSLisandro Dalcin /*@C 419efe9872eSLisandro Dalcin DMTSGetI2Function - get TS implicit residual evaluation function for 2nd order systems 420efe9872eSLisandro Dalcin 421efe9872eSLisandro Dalcin Not Collective 422efe9872eSLisandro Dalcin 423*4165533cSJose E. Roman Input Parameter: 424efe9872eSLisandro Dalcin . dm - DM to be used with TS 425efe9872eSLisandro Dalcin 426*4165533cSJose E. Roman Output Parameters: 427efe9872eSLisandro Dalcin + fun - function evaluation function, see TSSetI2Function() for calling sequence 428efe9872eSLisandro Dalcin - ctx - context for residual evaluation 429efe9872eSLisandro Dalcin 430efe9872eSLisandro Dalcin Level: advanced 431efe9872eSLisandro Dalcin 432efe9872eSLisandro Dalcin Note: 433efe9872eSLisandro Dalcin TSGetI2Function() is normally used, but it calls this function internally because the user context is actually 434efe9872eSLisandro Dalcin associated with the DM. 435efe9872eSLisandro Dalcin 436efe9872eSLisandro Dalcin .seealso: DMTSSetI2Function(),TSGetI2Function() 437efe9872eSLisandro Dalcin @*/ 438efe9872eSLisandro Dalcin PetscErrorCode DMTSGetI2Function(DM dm,TSI2Function *fun,void **ctx) 439efe9872eSLisandro Dalcin { 440efe9872eSLisandro Dalcin DMTS tsdm; 441efe9872eSLisandro Dalcin PetscErrorCode ierr; 442efe9872eSLisandro Dalcin 443efe9872eSLisandro Dalcin PetscFunctionBegin; 444efe9872eSLisandro Dalcin PetscValidHeaderSpecific(dm,DM_CLASSID,1); 445efe9872eSLisandro Dalcin ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 446efe9872eSLisandro Dalcin if (fun) *fun = tsdm->ops->i2function; 447efe9872eSLisandro Dalcin if (ctx) *ctx = tsdm->i2functionctx; 448efe9872eSLisandro Dalcin PetscFunctionReturn(0); 449efe9872eSLisandro Dalcin } 450efe9872eSLisandro Dalcin 451efe9872eSLisandro Dalcin /*@C 452efe9872eSLisandro Dalcin DMTSSetI2Jacobian - set TS implicit Jacobian evaluation function for 2nd order systems 453efe9872eSLisandro Dalcin 454efe9872eSLisandro Dalcin Not Collective 455efe9872eSLisandro Dalcin 456*4165533cSJose E. Roman Input Parameters: 457efe9872eSLisandro Dalcin + dm - DM to be used with TS 458a96d6ef6SBarry Smith . fun - Jacobian evaluation routine 459efe9872eSLisandro Dalcin - ctx - context for Jacobian evaluation 460efe9872eSLisandro Dalcin 461a96d6ef6SBarry Smith Calling sequence of jac: 462a96d6ef6SBarry Smith $ PetscErrorCode jac(TS ts,PetscReal t,Vec U,Vec U_t,Vec U_tt,PetscReal v,PetscReal a,Mat J,Mat P,void *ctx); 463a96d6ef6SBarry Smith 464a96d6ef6SBarry Smith + t - time at step/stage being solved 465a96d6ef6SBarry Smith . U - state vector 466a96d6ef6SBarry Smith . U_t - time derivative of state vector 467a96d6ef6SBarry Smith . U_tt - second time derivative of state vector 468a96d6ef6SBarry Smith . v - shift for U_t 469a96d6ef6SBarry Smith . a - shift for U_tt 470a96d6ef6SBarry Smith . J - Jacobian of G(U) = F(t,U,W+v*U,W'+a*U), equivalent to dF/dU + v*dF/dU_t + a*dF/dU_tt 471a96d6ef6SBarry Smith . P - preconditioning matrix for J, may be same as J 472a96d6ef6SBarry Smith - ctx - [optional] user-defined context for matrix evaluation routine 473a96d6ef6SBarry Smith 474efe9872eSLisandro Dalcin Level: advanced 475efe9872eSLisandro Dalcin 476efe9872eSLisandro Dalcin Note: 477efe9872eSLisandro Dalcin TSSetI2Jacobian() is normally used, but it calls this function internally because the user context is actually 478efe9872eSLisandro Dalcin associated with the DM. 479efe9872eSLisandro Dalcin 480efe9872eSLisandro Dalcin .seealso: TSSetI2Jacobian() 481efe9872eSLisandro Dalcin @*/ 482efe9872eSLisandro Dalcin PetscErrorCode DMTSSetI2Jacobian(DM dm,TSI2Jacobian jac,void *ctx) 483efe9872eSLisandro Dalcin { 484efe9872eSLisandro Dalcin DMTS tsdm; 485efe9872eSLisandro Dalcin PetscErrorCode ierr; 486efe9872eSLisandro Dalcin 487efe9872eSLisandro Dalcin PetscFunctionBegin; 488efe9872eSLisandro Dalcin PetscValidHeaderSpecific(dm,DM_CLASSID,1); 489efe9872eSLisandro Dalcin ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 490efe9872eSLisandro Dalcin if (jac) tsdm->ops->i2jacobian = jac; 491efe9872eSLisandro Dalcin if (ctx) tsdm->i2jacobianctx = ctx; 492efe9872eSLisandro Dalcin PetscFunctionReturn(0); 493efe9872eSLisandro Dalcin } 494efe9872eSLisandro Dalcin 495efe9872eSLisandro Dalcin /*@C 496efe9872eSLisandro Dalcin DMTSGetI2Jacobian - get TS implicit Jacobian evaluation function for 2nd order systems 497efe9872eSLisandro Dalcin 498efe9872eSLisandro Dalcin Not Collective 499efe9872eSLisandro Dalcin 500*4165533cSJose E. Roman Input Parameter: 501efe9872eSLisandro Dalcin . dm - DM to be used with TS 502efe9872eSLisandro Dalcin 503*4165533cSJose E. Roman Output Parameters: 504efe9872eSLisandro Dalcin + jac - Jacobian evaluation function, see TSSetI2Jacobian() for calling sequence 505efe9872eSLisandro Dalcin - ctx - context for Jacobian evaluation 506efe9872eSLisandro Dalcin 507efe9872eSLisandro Dalcin Level: advanced 508efe9872eSLisandro Dalcin 509efe9872eSLisandro Dalcin Note: 510efe9872eSLisandro Dalcin TSGetI2Jacobian() is normally used, but it calls this function internally because the user context is actually 511efe9872eSLisandro Dalcin associated with the DM. 512efe9872eSLisandro Dalcin 513efe9872eSLisandro Dalcin .seealso: DMTSSetI2Jacobian(),TSGetI2Jacobian() 514efe9872eSLisandro Dalcin @*/ 515efe9872eSLisandro Dalcin PetscErrorCode DMTSGetI2Jacobian(DM dm,TSI2Jacobian *jac,void **ctx) 516efe9872eSLisandro Dalcin { 517efe9872eSLisandro Dalcin DMTS tsdm; 518efe9872eSLisandro Dalcin PetscErrorCode ierr; 519efe9872eSLisandro Dalcin 520efe9872eSLisandro Dalcin PetscFunctionBegin; 521efe9872eSLisandro Dalcin PetscValidHeaderSpecific(dm,DM_CLASSID,1); 522efe9872eSLisandro Dalcin ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 523efe9872eSLisandro Dalcin if (jac) *jac = tsdm->ops->i2jacobian; 524efe9872eSLisandro Dalcin if (ctx) *ctx = tsdm->i2jacobianctx; 525efe9872eSLisandro Dalcin PetscFunctionReturn(0); 526efe9872eSLisandro Dalcin } 52724989b8cSPeter Brune 52824989b8cSPeter Brune /*@C 52924989b8cSPeter Brune DMTSSetRHSFunction - set TS explicit residual evaluation function 53024989b8cSPeter Brune 53124989b8cSPeter Brune Not Collective 53224989b8cSPeter Brune 533*4165533cSJose E. Roman Input Parameters: 53424989b8cSPeter Brune + dm - DM to be used with TS 535a96d6ef6SBarry Smith . func - RHS function evaluation routine 53624989b8cSPeter Brune - ctx - context for residual evaluation 53724989b8cSPeter Brune 538a96d6ef6SBarry Smith Calling sequence of func: 539a96d6ef6SBarry Smith $ PetscErrorCode func(TS ts,PetscReal t,Vec u,Vec F,void *ctx); 540a96d6ef6SBarry Smith 541a96d6ef6SBarry Smith + ts - timestep context 542a96d6ef6SBarry Smith . t - current timestep 543a96d6ef6SBarry Smith . u - input vector 544a96d6ef6SBarry Smith . F - function vector 545a96d6ef6SBarry Smith - ctx - [optional] user-defined function context 546a96d6ef6SBarry Smith 54724989b8cSPeter Brune Level: advanced 54824989b8cSPeter Brune 54924989b8cSPeter Brune Note: 550ef20d060SBarry Smith TSSetRSHFunction() is normally used, but it calls this function internally because the user context is actually 55124989b8cSPeter Brune associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 55224989b8cSPeter Brune not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 55324989b8cSPeter Brune 554a96d6ef6SBarry Smith .seealso: DMTSSetContext(), TSSetRHSFunction(), DMTSSetJacobian() 55524989b8cSPeter Brune @*/ 55624989b8cSPeter Brune PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx) 55724989b8cSPeter Brune { 55824989b8cSPeter Brune PetscErrorCode ierr; 559942e3340SBarry Smith DMTS tsdm; 56024989b8cSPeter Brune 56124989b8cSPeter Brune PetscFunctionBegin; 56224989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 563942e3340SBarry Smith ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 564d74926cbSBarry Smith if (func) tsdm->ops->rhsfunction = func; 56524989b8cSPeter Brune if (ctx) tsdm->rhsfunctionctx = ctx; 56624989b8cSPeter Brune PetscFunctionReturn(0); 56724989b8cSPeter Brune } 56824989b8cSPeter Brune 569ef20d060SBarry Smith /*@C 570e3c11fc1SJed Brown DMTSSetTransientVariable - sets function to transform from state to transient variables 571e3c11fc1SJed Brown 572e3c11fc1SJed Brown Logically Collective 573e3c11fc1SJed Brown 574*4165533cSJose E. Roman Input Parameters: 575e3c11fc1SJed Brown + dm - DM to be used with TS 576a96d6ef6SBarry Smith . tvar - a function that transforms to transient variables 577e3c11fc1SJed Brown - ctx - a context for tvar 578e3c11fc1SJed Brown 579a96d6ef6SBarry Smith Calling sequence of tvar: 580a96d6ef6SBarry Smith $ PetscErrorCode tvar(TS ts,Vec p,Vec c,void *ctx); 581a96d6ef6SBarry Smith 582a96d6ef6SBarry Smith + ts - timestep context 583a96d6ef6SBarry Smith . p - input vector (primative form) 584a96d6ef6SBarry Smith . c - output vector, transient variables (conservative form) 585a96d6ef6SBarry Smith - ctx - [optional] user-defined function context 586a96d6ef6SBarry Smith 587e3c11fc1SJed Brown Level: advanced 588e3c11fc1SJed Brown 589e3c11fc1SJed Brown Notes: 590e3c11fc1SJed Brown This is typically used to transform from primitive to conservative variables so that a time integrator (e.g., TSBDF) 591e3c11fc1SJed Brown can be conservative. In this context, primitive variables P are used to model the state (e.g., because they lead to 592e3c11fc1SJed Brown well-conditioned formulations even in limiting cases such as low-Mach or zero porosity). The transient variable is 593e3c11fc1SJed Brown C(P), specified by calling this function. An IFunction thus receives arguments (P, Cdot) and the IJacobian must be 594e3c11fc1SJed Brown evaluated via the chain rule, as in 595e3c11fc1SJed Brown 596e3c11fc1SJed Brown dF/dP + shift * dF/dCdot dC/dP. 597e3c11fc1SJed Brown 598438f35afSJed Brown .seealso: TSSetTransientVariable(), DMTSGetTransientVariable(), DMTSSetIFunction(), DMTSSetIJacobian() 599e3c11fc1SJed Brown @*/ 600e3c11fc1SJed Brown PetscErrorCode DMTSSetTransientVariable(DM dm,TSTransientVariable tvar,void *ctx) 601e3c11fc1SJed Brown { 602e3c11fc1SJed Brown PetscErrorCode ierr; 603e3c11fc1SJed Brown DMTS dmts; 604e3c11fc1SJed Brown 605e3c11fc1SJed Brown PetscFunctionBegin; 606e3c11fc1SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 607e3c11fc1SJed Brown ierr = DMGetDMTSWrite(dm,&dmts);CHKERRQ(ierr); 608e3c11fc1SJed Brown dmts->ops->transientvar = tvar; 609e3c11fc1SJed Brown dmts->transientvarctx = ctx; 610e3c11fc1SJed Brown PetscFunctionReturn(0); 611e3c11fc1SJed Brown } 612e3c11fc1SJed Brown 613e3c11fc1SJed Brown /*@C 614e3c11fc1SJed Brown DMTSGetTransientVariable - gets function to transform from state to transient variables 615e3c11fc1SJed Brown 616e3c11fc1SJed Brown Logically Collective 617e3c11fc1SJed Brown 618*4165533cSJose E. Roman Input Parameter: 619e3c11fc1SJed Brown . dm - DM to be used with TS 620e3c11fc1SJed Brown 621*4165533cSJose E. Roman Output Parameters: 622a96d6ef6SBarry Smith + tvar - a function that transforms to transient variables 623e3c11fc1SJed Brown - ctx - a context for tvar 624e3c11fc1SJed Brown 625e3c11fc1SJed Brown Level: advanced 626e3c11fc1SJed Brown 627e3c11fc1SJed Brown .seealso: DMTSSetTransientVariable(), DMTSGetIFunction(), DMTSGetIJacobian() 628e3c11fc1SJed Brown @*/ 629e3c11fc1SJed Brown PetscErrorCode DMTSGetTransientVariable(DM dm,TSTransientVariable *tvar,void *ctx) 630e3c11fc1SJed Brown { 631e3c11fc1SJed Brown PetscErrorCode ierr; 632e3c11fc1SJed Brown DMTS dmts; 633e3c11fc1SJed Brown 634e3c11fc1SJed Brown PetscFunctionBegin; 635e3c11fc1SJed Brown PetscValidHeaderSpecific(dm,DM_CLASSID,1); 636e3c11fc1SJed Brown ierr = DMGetDMTS(dm,&dmts);CHKERRQ(ierr); 637e3c11fc1SJed Brown if (tvar) *tvar = dmts->ops->transientvar; 638e3c11fc1SJed Brown if (ctx) *(void**)ctx = dmts->transientvarctx; 639e3c11fc1SJed Brown PetscFunctionReturn(0); 640e3c11fc1SJed Brown } 641e3c11fc1SJed Brown 642e3c11fc1SJed Brown /*@C 643ef20d060SBarry Smith DMTSGetSolutionFunction - gets the TS solution evaluation function 644ef20d060SBarry Smith 645ef20d060SBarry Smith Not Collective 646ef20d060SBarry Smith 647*4165533cSJose E. Roman Input Parameter: 648ef20d060SBarry Smith . dm - DM to be used with TS 649ef20d060SBarry Smith 650ef20d060SBarry Smith Output Parameters: 651ef20d060SBarry Smith + func - solution function evaluation function, see TSSetSolution() for calling sequence 652ef20d060SBarry Smith - ctx - context for solution evaluation 653ef20d060SBarry Smith 654ef20d060SBarry Smith Level: advanced 655ef20d060SBarry Smith 6560ed3bfb6SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), DMTSSetSolutionFunction() 657ef20d060SBarry Smith @*/ 658ef20d060SBarry Smith PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx) 659ef20d060SBarry Smith { 660ef20d060SBarry Smith PetscErrorCode ierr; 661942e3340SBarry Smith DMTS tsdm; 662ef20d060SBarry Smith 663ef20d060SBarry Smith PetscFunctionBegin; 664ef20d060SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 665942e3340SBarry Smith ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 666d74926cbSBarry Smith if (func) *func = tsdm->ops->solution; 667ef20d060SBarry Smith if (ctx) *ctx = tsdm->solutionctx; 668ef20d060SBarry Smith PetscFunctionReturn(0); 669ef20d060SBarry Smith } 670ef20d060SBarry Smith 671ef20d060SBarry Smith /*@C 672ef20d060SBarry Smith DMTSSetSolutionFunction - set TS solution evaluation function 673ef20d060SBarry Smith 674ef20d060SBarry Smith Not Collective 675ef20d060SBarry Smith 676*4165533cSJose E. Roman Input Parameters: 677ef20d060SBarry Smith + dm - DM to be used with TS 678a96d6ef6SBarry Smith . func - solution function evaluation routine 679ef20d060SBarry Smith - ctx - context for solution evaluation 680ef20d060SBarry Smith 681a96d6ef6SBarry Smith Calling sequence of f: 682a96d6ef6SBarry Smith $ PetscErrorCode f(TS ts,PetscReal t,Vec u,void *ctx); 683a96d6ef6SBarry Smith 684a96d6ef6SBarry Smith + ts - timestep context 685a96d6ef6SBarry Smith . t - current timestep 686a96d6ef6SBarry Smith . u - output vector 687a96d6ef6SBarry Smith - ctx - [optional] user-defined function context 688a96d6ef6SBarry Smith 689ef20d060SBarry Smith Level: advanced 690ef20d060SBarry Smith 691ef20d060SBarry Smith Note: 692ef20d060SBarry Smith TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually 693ef20d060SBarry Smith associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 694ef20d060SBarry Smith not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 695ef20d060SBarry Smith 6960ed3bfb6SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), DMTSGetSolutionFunction() 697ef20d060SBarry Smith @*/ 698ef20d060SBarry Smith PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx) 699ef20d060SBarry Smith { 700ef20d060SBarry Smith PetscErrorCode ierr; 701942e3340SBarry Smith DMTS tsdm; 702ef20d060SBarry Smith 703ef20d060SBarry Smith PetscFunctionBegin; 704ef20d060SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 705942e3340SBarry Smith ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 706d74926cbSBarry Smith if (func) tsdm->ops->solution = func; 707ef20d060SBarry Smith if (ctx) tsdm->solutionctx = ctx; 708ef20d060SBarry Smith PetscFunctionReturn(0); 709ef20d060SBarry Smith } 710ef20d060SBarry Smith 7119b7cd975SBarry Smith /*@C 7129b7cd975SBarry Smith DMTSSetForcingFunction - set TS forcing function evaluation function 7139b7cd975SBarry Smith 7149b7cd975SBarry Smith Not Collective 7159b7cd975SBarry Smith 716*4165533cSJose E. Roman Input Parameters: 7179b7cd975SBarry Smith + dm - DM to be used with TS 718a96d6ef6SBarry Smith . f - forcing function evaluation routine 7199b7cd975SBarry Smith - ctx - context for solution evaluation 7209b7cd975SBarry Smith 721a96d6ef6SBarry Smith Calling sequence of func: 722a96d6ef6SBarry Smith $ PetscErrorCode func (TS ts,PetscReal t,Vec f,void *ctx); 723a96d6ef6SBarry Smith 724a96d6ef6SBarry Smith + ts - timestep context 725a96d6ef6SBarry Smith . t - current timestep 726a96d6ef6SBarry Smith . f - output vector 727a96d6ef6SBarry Smith - ctx - [optional] user-defined function context 728a96d6ef6SBarry Smith 7299b7cd975SBarry Smith Level: advanced 7309b7cd975SBarry Smith 7319b7cd975SBarry Smith Note: 7329b7cd975SBarry Smith TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually 7339b7cd975SBarry Smith associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 7349b7cd975SBarry Smith not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 7359b7cd975SBarry Smith 7369b7cd975SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction() 7379b7cd975SBarry Smith @*/ 738d56366bfSLisandro Dalcin PetscErrorCode DMTSSetForcingFunction(DM dm,TSForcingFunction f,void *ctx) 7399b7cd975SBarry Smith { 7409b7cd975SBarry Smith PetscErrorCode ierr; 7419b7cd975SBarry Smith DMTS tsdm; 7429b7cd975SBarry Smith 7439b7cd975SBarry Smith PetscFunctionBegin; 7449b7cd975SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 7459b7cd975SBarry Smith ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 746f8b49ee9SBarry Smith if (f) tsdm->ops->forcing = f; 7479b7cd975SBarry Smith if (ctx) tsdm->forcingctx = ctx; 7489b7cd975SBarry Smith PetscFunctionReturn(0); 7499b7cd975SBarry Smith } 7509b7cd975SBarry Smith 7519b7cd975SBarry Smith /*@C 7529b7cd975SBarry Smith DMTSGetForcingFunction - get TS forcing function evaluation function 7539b7cd975SBarry Smith 7549b7cd975SBarry Smith Not Collective 7559b7cd975SBarry Smith 756*4165533cSJose E. Roman Input Parameter: 7579b7cd975SBarry Smith . dm - DM to be used with TS 7589b7cd975SBarry Smith 759*4165533cSJose E. Roman Output Parameters: 760f8b49ee9SBarry Smith + f - forcing function evaluation function; see TSForcingFunction for details 7619b7cd975SBarry Smith - ctx - context for solution evaluation 7629b7cd975SBarry Smith 7639b7cd975SBarry Smith Level: advanced 7649b7cd975SBarry Smith 7659b7cd975SBarry Smith Note: 7669b7cd975SBarry Smith TSSetForcingFunction() is normally used, but it calls this function internally because the user context is actually 7679b7cd975SBarry Smith associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 7689b7cd975SBarry Smith not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 7699b7cd975SBarry Smith 7709b7cd975SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian(), TSSetForcingFunction(), DMTSGetForcingFunction() 7719b7cd975SBarry Smith @*/ 772d56366bfSLisandro Dalcin PetscErrorCode DMTSGetForcingFunction(DM dm,TSForcingFunction *f,void **ctx) 7739b7cd975SBarry Smith { 7749b7cd975SBarry Smith PetscErrorCode ierr; 7759b7cd975SBarry Smith DMTS tsdm; 7769b7cd975SBarry Smith 7779b7cd975SBarry Smith PetscFunctionBegin; 7789b7cd975SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 7799b7cd975SBarry Smith ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 780f8b49ee9SBarry Smith if (f) *f = tsdm->ops->forcing; 7819b7cd975SBarry Smith if (ctx) *ctx = tsdm->forcingctx; 7829b7cd975SBarry Smith PetscFunctionReturn(0); 7839b7cd975SBarry Smith } 7849b7cd975SBarry Smith 78524989b8cSPeter Brune /*@C 78624989b8cSPeter Brune DMTSGetRHSFunction - get TS explicit residual evaluation function 78724989b8cSPeter Brune 78824989b8cSPeter Brune Not Collective 78924989b8cSPeter Brune 790*4165533cSJose E. Roman Input Parameter: 79124989b8cSPeter Brune . dm - DM to be used with TS 79224989b8cSPeter Brune 793*4165533cSJose E. Roman Output Parameters: 79424989b8cSPeter Brune + func - residual evaluation function, see TSSetRHSFunction() for calling sequence 79524989b8cSPeter Brune - ctx - context for residual evaluation 79624989b8cSPeter Brune 79724989b8cSPeter Brune Level: advanced 79824989b8cSPeter Brune 79924989b8cSPeter Brune Note: 80024989b8cSPeter Brune TSGetFunction() is normally used, but it calls this function internally because the user context is actually 80124989b8cSPeter Brune associated with the DM. 80224989b8cSPeter Brune 803a96d6ef6SBarry Smith .seealso: DMTSSetContext(), DMTSSetRHSFunction(), TSSetRHSFunction() 80424989b8cSPeter Brune @*/ 80524989b8cSPeter Brune PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx) 80624989b8cSPeter Brune { 80724989b8cSPeter Brune PetscErrorCode ierr; 808942e3340SBarry Smith DMTS tsdm; 80924989b8cSPeter Brune 81024989b8cSPeter Brune PetscFunctionBegin; 81124989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 812942e3340SBarry Smith ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 813d74926cbSBarry Smith if (func) *func = tsdm->ops->rhsfunction; 81424989b8cSPeter Brune if (ctx) *ctx = tsdm->rhsfunctionctx; 81524989b8cSPeter Brune PetscFunctionReturn(0); 81624989b8cSPeter Brune } 81724989b8cSPeter Brune 81824989b8cSPeter Brune /*@C 81924989b8cSPeter Brune DMTSSetIJacobian - set TS Jacobian evaluation function 82024989b8cSPeter Brune 82124989b8cSPeter Brune Not Collective 82224989b8cSPeter Brune 823*4165533cSJose E. Roman Input Parameters: 82424989b8cSPeter Brune + dm - DM to be used with TS 825a96d6ef6SBarry Smith . func - Jacobian evaluation routine 82624989b8cSPeter Brune - ctx - context for residual evaluation 82724989b8cSPeter Brune 828a96d6ef6SBarry Smith Calling sequence of f: 829a96d6ef6SBarry Smith $ PetscErrorCode f(TS ts,PetscReal t,Vec U,Vec U_t,PetscReal a,Mat Amat,Mat Pmat,void *ctx); 830a96d6ef6SBarry Smith 831a96d6ef6SBarry Smith + t - time at step/stage being solved 832a96d6ef6SBarry Smith . U - state vector 833a96d6ef6SBarry Smith . U_t - time derivative of state vector 834a96d6ef6SBarry Smith . a - shift 835a96d6ef6SBarry Smith . Amat - (approximate) Jacobian of F(t,U,W+a*U), equivalent to dF/dU + a*dF/dU_t 836a96d6ef6SBarry Smith . Pmat - matrix used for constructing preconditioner, usually the same as Amat 837a96d6ef6SBarry Smith - ctx - [optional] user-defined context for matrix evaluation routine 838a96d6ef6SBarry Smith 83924989b8cSPeter Brune Level: advanced 84024989b8cSPeter Brune 84124989b8cSPeter Brune Note: 84224989b8cSPeter Brune TSSetJacobian() is normally used, but it calls this function internally because the user context is actually 84324989b8cSPeter Brune associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 84424989b8cSPeter Brune not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 84524989b8cSPeter Brune 846a96d6ef6SBarry Smith .seealso: DMTSSetContext(), TSSetRHSFunction(), DMTSGetJacobian(), TSSetIJacobian(), TSSetIFunction() 84724989b8cSPeter Brune @*/ 84824989b8cSPeter Brune PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx) 84924989b8cSPeter Brune { 85024989b8cSPeter Brune PetscErrorCode ierr; 851942e3340SBarry Smith DMTS sdm; 85224989b8cSPeter Brune 85324989b8cSPeter Brune PetscFunctionBegin; 85424989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 855942e3340SBarry Smith ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr); 856d74926cbSBarry Smith if (func) sdm->ops->ijacobian = func; 85724989b8cSPeter Brune if (ctx) sdm->ijacobianctx = ctx; 85824989b8cSPeter Brune PetscFunctionReturn(0); 85924989b8cSPeter Brune } 86024989b8cSPeter Brune 86124989b8cSPeter Brune /*@C 86224989b8cSPeter Brune DMTSGetIJacobian - get TS Jacobian evaluation function 86324989b8cSPeter Brune 86424989b8cSPeter Brune Not Collective 86524989b8cSPeter Brune 866*4165533cSJose E. Roman Input Parameter: 86724989b8cSPeter Brune . dm - DM to be used with TS 86824989b8cSPeter Brune 869*4165533cSJose E. Roman Output Parameters: 87024989b8cSPeter Brune + func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence 87124989b8cSPeter Brune - ctx - context for residual evaluation 87224989b8cSPeter Brune 87324989b8cSPeter Brune Level: advanced 87424989b8cSPeter Brune 87524989b8cSPeter Brune Note: 87624989b8cSPeter Brune TSGetJacobian() is normally used, but it calls this function internally because the user context is actually 87724989b8cSPeter Brune associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 87824989b8cSPeter Brune not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 87924989b8cSPeter Brune 88024989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 88124989b8cSPeter Brune @*/ 88224989b8cSPeter Brune PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx) 88324989b8cSPeter Brune { 88424989b8cSPeter Brune PetscErrorCode ierr; 885942e3340SBarry Smith DMTS tsdm; 88624989b8cSPeter Brune 88724989b8cSPeter Brune PetscFunctionBegin; 88824989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 889942e3340SBarry Smith ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 890d74926cbSBarry Smith if (func) *func = tsdm->ops->ijacobian; 89124989b8cSPeter Brune if (ctx) *ctx = tsdm->ijacobianctx; 89224989b8cSPeter Brune PetscFunctionReturn(0); 89324989b8cSPeter Brune } 89424989b8cSPeter Brune 89524989b8cSPeter Brune /*@C 89624989b8cSPeter Brune DMTSSetRHSJacobian - set TS Jacobian evaluation function 89724989b8cSPeter Brune 89824989b8cSPeter Brune Not Collective 89924989b8cSPeter Brune 900*4165533cSJose E. Roman Input Parameters: 90124989b8cSPeter Brune + dm - DM to be used with TS 902a96d6ef6SBarry Smith . func - Jacobian evaluation routine 90324989b8cSPeter Brune - ctx - context for residual evaluation 90424989b8cSPeter Brune 905a96d6ef6SBarry Smith Calling sequence of func: 906a96d6ef6SBarry Smith $ PetscErrorCode func(TS ts,PetscReal t,Vec u,Mat A,Mat B,void *ctx); 907a96d6ef6SBarry Smith 908a96d6ef6SBarry Smith + t - current timestep 909a96d6ef6SBarry Smith . u - input vector 910a96d6ef6SBarry Smith . Amat - (approximate) Jacobian matrix 911a96d6ef6SBarry Smith . Pmat - matrix from which preconditioner is to be constructed (usually the same as Amat) 912a96d6ef6SBarry Smith - ctx - [optional] user-defined context for matrix evaluation routine 913a96d6ef6SBarry Smith 91424989b8cSPeter Brune Level: advanced 91524989b8cSPeter Brune 91624989b8cSPeter Brune Note: 91724989b8cSPeter Brune TSSetJacobian() is normally used, but it calls this function internally because the user context is actually 91824989b8cSPeter Brune associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 91924989b8cSPeter Brune not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 92024989b8cSPeter Brune 921a96d6ef6SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetRHSJacobian() 92224989b8cSPeter Brune @*/ 92324989b8cSPeter Brune PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx) 92424989b8cSPeter Brune { 92524989b8cSPeter Brune PetscErrorCode ierr; 926942e3340SBarry Smith DMTS tsdm; 92724989b8cSPeter Brune 92824989b8cSPeter Brune PetscFunctionBegin; 92924989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 930942e3340SBarry Smith ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 931d74926cbSBarry Smith if (func) tsdm->ops->rhsjacobian = func; 93224989b8cSPeter Brune if (ctx) tsdm->rhsjacobianctx = ctx; 93324989b8cSPeter Brune PetscFunctionReturn(0); 93424989b8cSPeter Brune } 93524989b8cSPeter Brune 93624989b8cSPeter Brune /*@C 93724989b8cSPeter Brune DMTSGetRHSJacobian - get TS Jacobian evaluation function 93824989b8cSPeter Brune 93924989b8cSPeter Brune Not Collective 94024989b8cSPeter Brune 941*4165533cSJose E. Roman Input Parameter: 94224989b8cSPeter Brune . dm - DM to be used with TS 94324989b8cSPeter Brune 944*4165533cSJose E. Roman Output Parameters: 94524989b8cSPeter Brune + func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence 94624989b8cSPeter Brune - ctx - context for residual evaluation 94724989b8cSPeter Brune 94824989b8cSPeter Brune Level: advanced 94924989b8cSPeter Brune 95024989b8cSPeter Brune Note: 95124989b8cSPeter Brune TSGetJacobian() is normally used, but it calls this function internally because the user context is actually 95224989b8cSPeter Brune associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 95324989b8cSPeter Brune not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 95424989b8cSPeter Brune 955a96d6ef6SBarry Smith .seealso: DMTSSetContext(), TSSetRHSFunction(), DMTSSetRHSJacobian(), TSSetRHSJacobian() 95624989b8cSPeter Brune @*/ 95724989b8cSPeter Brune PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx) 95824989b8cSPeter Brune { 95924989b8cSPeter Brune PetscErrorCode ierr; 960942e3340SBarry Smith DMTS tsdm; 96124989b8cSPeter Brune 96224989b8cSPeter Brune PetscFunctionBegin; 96324989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 964942e3340SBarry Smith ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 965d74926cbSBarry Smith if (func) *func = tsdm->ops->rhsjacobian; 96624989b8cSPeter Brune if (ctx) *ctx = tsdm->rhsjacobianctx; 96724989b8cSPeter Brune PetscFunctionReturn(0); 96824989b8cSPeter Brune } 969ad6bc421SBarry Smith 970ad6bc421SBarry Smith /*@C 971ad6bc421SBarry Smith DMTSSetIFunctionSerialize - sets functions used to view and load a IFunction context 972ad6bc421SBarry Smith 973ad6bc421SBarry Smith Not Collective 974ad6bc421SBarry Smith 975*4165533cSJose E. Roman Input Parameters: 976ad6bc421SBarry Smith + dm - DM to be used with TS 977ad6bc421SBarry Smith . view - viewer function 978ad6bc421SBarry Smith - load - loading function 979ad6bc421SBarry Smith 980ad6bc421SBarry Smith Level: advanced 981ad6bc421SBarry Smith 982ad6bc421SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 983ad6bc421SBarry Smith @*/ 984ad6bc421SBarry Smith PetscErrorCode DMTSSetIFunctionSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer)) 985ad6bc421SBarry Smith { 986ad6bc421SBarry Smith PetscErrorCode ierr; 987ad6bc421SBarry Smith DMTS tsdm; 988ad6bc421SBarry Smith 989ad6bc421SBarry Smith PetscFunctionBegin; 990ad6bc421SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 991ad6bc421SBarry Smith ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 992ad6bc421SBarry Smith tsdm->ops->ifunctionview = view; 993ad6bc421SBarry Smith tsdm->ops->ifunctionload = load; 994ad6bc421SBarry Smith PetscFunctionReturn(0); 995ad6bc421SBarry Smith } 996ad6bc421SBarry Smith 997ad6bc421SBarry Smith /*@C 998ad6bc421SBarry Smith DMTSSetIJacobianSerialize - sets functions used to view and load a IJacobian context 999ad6bc421SBarry Smith 1000ad6bc421SBarry Smith Not Collective 1001ad6bc421SBarry Smith 1002*4165533cSJose E. Roman Input Parameters: 1003ad6bc421SBarry Smith + dm - DM to be used with TS 1004ad6bc421SBarry Smith . view - viewer function 1005ad6bc421SBarry Smith - load - loading function 1006ad6bc421SBarry Smith 1007ad6bc421SBarry Smith Level: advanced 1008ad6bc421SBarry Smith 1009ad6bc421SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 1010ad6bc421SBarry Smith @*/ 1011ad6bc421SBarry Smith PetscErrorCode DMTSSetIJacobianSerialize(DM dm,PetscErrorCode (*view)(void*,PetscViewer),PetscErrorCode (*load)(void**,PetscViewer)) 1012ad6bc421SBarry Smith { 1013ad6bc421SBarry Smith PetscErrorCode ierr; 1014ad6bc421SBarry Smith DMTS tsdm; 1015ad6bc421SBarry Smith 1016ad6bc421SBarry Smith PetscFunctionBegin; 1017ad6bc421SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1018ad6bc421SBarry Smith ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 1019ad6bc421SBarry Smith tsdm->ops->ijacobianview = view; 1020ad6bc421SBarry Smith tsdm->ops->ijacobianload = load; 1021ad6bc421SBarry Smith PetscFunctionReturn(0); 1022ad6bc421SBarry Smith } 1023