124989b8cSPeter Brune #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 224989b8cSPeter Brune #include <petsc-private/dmimpl.h> /*I "petscdm.h" I*/ 324989b8cSPeter Brune 4d74926cbSBarry Smith #undef __FUNCT__ 5d74926cbSBarry Smith #define __FUNCT__ "DMTSDestroy" 6d74926cbSBarry Smith static PetscErrorCode DMTSDestroy(DMTS *kdm) 7d74926cbSBarry Smith { 8d74926cbSBarry Smith PetscErrorCode ierr; 9d74926cbSBarry Smith 10d74926cbSBarry Smith PetscFunctionBegin; 11d74926cbSBarry Smith if (!*kdm) PetscFunctionReturn(0); 12d74926cbSBarry Smith PetscValidHeaderSpecific((*kdm),DMTS_CLASSID,1); 13d74926cbSBarry Smith if (--((PetscObject)(*kdm))->refct > 0) {*kdm = 0; PetscFunctionReturn(0);} 14d74926cbSBarry Smith if ((*kdm)->ops->destroy) {ierr = ((*kdm)->ops->destroy)(*kdm);CHKERRQ(ierr);} 15d74926cbSBarry Smith ierr = PetscHeaderDestroy(kdm);CHKERRQ(ierr); 16d74926cbSBarry Smith PetscFunctionReturn(0); 17d74926cbSBarry Smith } 18d74926cbSBarry Smith 19d74926cbSBarry Smith #undef __FUNCT__ 20d74926cbSBarry Smith #define __FUNCT__ "DMTSCreate" 21d74926cbSBarry Smith static PetscErrorCode DMTSCreate(MPI_Comm comm,DMTS *kdm) 22d74926cbSBarry Smith { 23d74926cbSBarry Smith PetscErrorCode ierr; 24d74926cbSBarry Smith 25d74926cbSBarry Smith PetscFunctionBegin; 26d74926cbSBarry Smith #ifndef PETSC_USE_DYNAMIC_LIBRARIES 27d74926cbSBarry Smith ierr = TSInitializePackage(PETSC_NULL);CHKERRQ(ierr); 28d74926cbSBarry Smith #endif 29b4615a05SBarry Smith ierr = PetscHeaderCreate(*kdm, _p_DMTS, struct _DMTSOps, DMTS_CLASSID, -1, "DMTS", "DMTS", "DMTS", comm, DMTSDestroy, PETSC_NULL);CHKERRQ(ierr); 30d74926cbSBarry Smith ierr = PetscMemzero((*kdm)->ops, sizeof(struct _DMTSOps));CHKERRQ(ierr); 31d74926cbSBarry Smith PetscFunctionReturn(0); 32d74926cbSBarry Smith } 3324989b8cSPeter Brune 3424989b8cSPeter Brune #undef __FUNCT__ 35942e3340SBarry Smith #define __FUNCT__ "DMCoarsenHook_DMTS" 362a34c10cSBarry Smith /* Attaches the DMTS to the coarse level. 3724989b8cSPeter Brune * Under what conditions should we copy versus duplicate? 3824989b8cSPeter Brune */ 39942e3340SBarry Smith static PetscErrorCode DMCoarsenHook_DMTS(DM dm,DM dmc,void *ctx) 4024989b8cSPeter Brune { 4124989b8cSPeter Brune PetscErrorCode ierr; 4224989b8cSPeter Brune 4324989b8cSPeter Brune PetscFunctionBegin; 44942e3340SBarry Smith ierr = DMCopyDMTS(dm,dmc);CHKERRQ(ierr); 4524989b8cSPeter Brune PetscFunctionReturn(0); 4624989b8cSPeter Brune } 4724989b8cSPeter Brune 4824989b8cSPeter Brune #undef __FUNCT__ 49942e3340SBarry Smith #define __FUNCT__ "DMRestrictHook_DMTS" 5024989b8cSPeter Brune /* This could restrict auxiliary information to the coarse level. 5124989b8cSPeter Brune */ 52942e3340SBarry Smith static PetscErrorCode DMRestrictHook_DMTS(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx) 5324989b8cSPeter Brune { 5424989b8cSPeter Brune 5524989b8cSPeter Brune PetscFunctionBegin; 5624989b8cSPeter Brune PetscFunctionReturn(0); 5724989b8cSPeter Brune } 5824989b8cSPeter Brune 5924989b8cSPeter Brune #undef __FUNCT__ 60*258e1594SPeter Brune #define __FUNCT__ "DMSubDomainHook_DMTS" 61*258e1594SPeter Brune static PetscErrorCode DMSubDomainHook_DMTS(DM dm,DM subdm,void *ctx) 62*258e1594SPeter Brune { 63*258e1594SPeter Brune PetscErrorCode ierr; 64*258e1594SPeter Brune 65*258e1594SPeter Brune PetscFunctionBegin; 66*258e1594SPeter Brune ierr = DMCopyDMTS(dm,subdm);CHKERRQ(ierr); 67*258e1594SPeter Brune PetscFunctionReturn(0); 68*258e1594SPeter Brune } 69*258e1594SPeter Brune 70*258e1594SPeter Brune #undef __FUNCT__ 71*258e1594SPeter Brune #define __FUNCT__ "DMSubDomainRestrictHook_DMTS" 72*258e1594SPeter Brune /* This could restrict auxiliary information to the coarse level. 73*258e1594SPeter Brune */ 74*258e1594SPeter Brune static PetscErrorCode DMSubDomainRestrictHook_DMTS(DM dm,VecScatter gscat,VecScatter lscat,DM subdm,void *ctx) 75*258e1594SPeter Brune { 76*258e1594SPeter Brune PetscFunctionBegin; 77*258e1594SPeter Brune PetscFunctionReturn(0); 78*258e1594SPeter Brune } 79*258e1594SPeter Brune 80*258e1594SPeter Brune #undef __FUNCT__ 81d74926cbSBarry Smith #define __FUNCT__ "DMTSCopy" 82d74926cbSBarry Smith /*@C 83d74926cbSBarry Smith DMTSCopy - copies the information in a DMTS to another DMTS 84d74926cbSBarry Smith 85d74926cbSBarry Smith Not Collective 86d74926cbSBarry Smith 87d74926cbSBarry Smith Input Argument: 88d74926cbSBarry Smith + kdm - Original DMTS 89d74926cbSBarry Smith - nkdm - DMTS to receive the data, should have been created with DMTSCreate() 90d74926cbSBarry Smith 91d74926cbSBarry Smith Level: developer 92d74926cbSBarry Smith 93d74926cbSBarry Smith .seealso: DMTSCreate(), DMTSDestroy() 94d74926cbSBarry Smith @*/ 95d74926cbSBarry Smith PetscErrorCode DMTSCopy(DMTS kdm,DMTS nkdm) 9624989b8cSPeter Brune { 9724989b8cSPeter Brune PetscErrorCode ierr; 9824989b8cSPeter Brune 9924989b8cSPeter Brune PetscFunctionBegin; 100d74926cbSBarry Smith PetscValidHeaderSpecific(kdm,DMTS_CLASSID,1); 101d74926cbSBarry Smith PetscValidHeaderSpecific(nkdm,DMTS_CLASSID,2); 102d74926cbSBarry Smith nkdm->ops->rhsfunction = kdm->ops->rhsfunction; 103d74926cbSBarry Smith nkdm->ops->rhsjacobian = kdm->ops->rhsjacobian; 104d74926cbSBarry Smith nkdm->ops->ifunction = kdm->ops->ifunction; 105d74926cbSBarry Smith nkdm->ops->ijacobian = kdm->ops->ijacobian; 106d74926cbSBarry Smith nkdm->ops->solution = kdm->ops->solution; 107d74926cbSBarry Smith nkdm->ops->destroy = kdm->ops->destroy; 108d74926cbSBarry Smith nkdm->ops->duplicate = kdm->ops->duplicate; 109d74926cbSBarry Smith 110d74926cbSBarry Smith nkdm->rhsfunctionctx = kdm->rhsfunctionctx; 111d74926cbSBarry Smith nkdm->rhsjacobianctx = kdm->rhsjacobianctx; 112d74926cbSBarry Smith nkdm->ifunctionctx = kdm->ifunctionctx; 113d74926cbSBarry Smith nkdm->ijacobianctx = kdm->ijacobianctx; 114d74926cbSBarry Smith nkdm->solutionctx = kdm->solutionctx; 115d74926cbSBarry Smith 116d74926cbSBarry Smith nkdm->data = kdm->data; 117d74926cbSBarry Smith 118d74926cbSBarry Smith /* 119d74926cbSBarry Smith nkdm->fortran_func_pointers[0] = kdm->fortran_func_pointers[0]; 120d74926cbSBarry Smith nkdm->fortran_func_pointers[1] = kdm->fortran_func_pointers[1]; 121d74926cbSBarry Smith nkdm->fortran_func_pointers[2] = kdm->fortran_func_pointers[2]; 122d74926cbSBarry Smith */ 123d74926cbSBarry Smith 124d74926cbSBarry Smith /* implementation specific copy hooks */ 125d74926cbSBarry Smith if (kdm->ops->duplicate) {ierr = (*kdm->ops->duplicate)(kdm,nkdm);CHKERRQ(ierr);} 12624989b8cSPeter Brune PetscFunctionReturn(0); 12724989b8cSPeter Brune } 12824989b8cSPeter Brune 12924989b8cSPeter Brune #undef __FUNCT__ 130942e3340SBarry Smith #define __FUNCT__ "DMGetDMTS" 13124989b8cSPeter Brune /*@C 132942e3340SBarry Smith DMGetDMTS - get read-only private DMTS context from a DM 13324989b8cSPeter Brune 13424989b8cSPeter Brune Not Collective 13524989b8cSPeter Brune 13624989b8cSPeter Brune Input Argument: 13724989b8cSPeter Brune . dm - DM to be used with TS 13824989b8cSPeter Brune 13924989b8cSPeter Brune Output Argument: 140942e3340SBarry Smith . tsdm - private DMTS context 14124989b8cSPeter Brune 14224989b8cSPeter Brune Level: developer 14324989b8cSPeter Brune 14424989b8cSPeter Brune Notes: 145942e3340SBarry Smith Use DMGetDMTSWrite() if write access is needed. The DMTSSetXXX API should be used wherever possible. 14624989b8cSPeter Brune 147942e3340SBarry Smith .seealso: DMGetDMTSWrite() 14824989b8cSPeter Brune @*/ 149942e3340SBarry Smith PetscErrorCode DMGetDMTS(DM dm,DMTS *tsdm) 15024989b8cSPeter Brune { 15124989b8cSPeter Brune PetscErrorCode ierr; 15224989b8cSPeter Brune 15324989b8cSPeter Brune PetscFunctionBegin; 15424989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 1552a34c10cSBarry Smith *tsdm = (DMTS) dm->dmts; 156d74926cbSBarry Smith if (!*tsdm) { 157942e3340SBarry Smith ierr = PetscInfo(dm,"Creating new DMTS\n");CHKERRQ(ierr); 158d74926cbSBarry Smith ierr = DMTSCreate(((PetscObject)dm)->comm,tsdm);CHKERRQ(ierr); 1592a34c10cSBarry Smith dm->dmts = (PetscObject) *tsdm; 160942e3340SBarry Smith ierr = DMCoarsenHookAdd(dm,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,PETSC_NULL);CHKERRQ(ierr); 161*258e1594SPeter Brune ierr = DMSubDomainHookAdd(dm,DMSubDomainHook_DMTS,DMSubDomainRestrictHook_DMTS,PETSC_NULL);CHKERRQ(ierr); 16224989b8cSPeter Brune } 16324989b8cSPeter Brune PetscFunctionReturn(0); 16424989b8cSPeter Brune } 16524989b8cSPeter Brune 16624989b8cSPeter Brune #undef __FUNCT__ 167942e3340SBarry Smith #define __FUNCT__ "DMGetDMTSWrite" 16824989b8cSPeter Brune /*@C 169942e3340SBarry Smith DMGetDMTSWrite - get write access to private DMTS context from a DM 17024989b8cSPeter Brune 17124989b8cSPeter Brune Not Collective 17224989b8cSPeter Brune 17324989b8cSPeter Brune Input Argument: 17424989b8cSPeter Brune . dm - DM to be used with TS 17524989b8cSPeter Brune 17624989b8cSPeter Brune Output Argument: 177942e3340SBarry Smith . tsdm - private DMTS context 17824989b8cSPeter Brune 17924989b8cSPeter Brune Level: developer 18024989b8cSPeter Brune 181942e3340SBarry Smith .seealso: DMGetDMTS() 18224989b8cSPeter Brune @*/ 183942e3340SBarry Smith PetscErrorCode DMGetDMTSWrite(DM dm,DMTS *tsdm) 18424989b8cSPeter Brune { 18524989b8cSPeter Brune PetscErrorCode ierr; 186942e3340SBarry Smith DMTS sdm; 18724989b8cSPeter Brune 18824989b8cSPeter Brune PetscFunctionBegin; 18924989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 190942e3340SBarry Smith ierr = DMGetDMTS(dm,&sdm);CHKERRQ(ierr); 19124989b8cSPeter Brune if (!sdm->originaldm) sdm->originaldm = dm; 19224989b8cSPeter Brune if (sdm->originaldm != dm) { /* Copy on write */ 1932a34c10cSBarry Smith DMTS oldsdm = sdm; 194942e3340SBarry Smith ierr = PetscInfo(dm,"Copying DMTS due to write\n");CHKERRQ(ierr); 195d74926cbSBarry Smith ierr = DMTSCreate(((PetscObject)dm)->comm,&sdm);CHKERRQ(ierr); 196d74926cbSBarry Smith ierr = DMTSCopy(oldsdm,sdm);CHKERRQ(ierr); 1972a34c10cSBarry Smith ierr = DMTSDestroy((DMTS*)&dm->dmts);CHKERRQ(ierr); 1982a34c10cSBarry Smith dm->dmts = (PetscObject) sdm; 19924989b8cSPeter Brune } 20024989b8cSPeter Brune *tsdm = sdm; 20124989b8cSPeter Brune PetscFunctionReturn(0); 20224989b8cSPeter Brune } 20324989b8cSPeter Brune 20424989b8cSPeter Brune #undef __FUNCT__ 205942e3340SBarry Smith #define __FUNCT__ "DMCopyDMTS" 20624989b8cSPeter Brune /*@C 207942e3340SBarry Smith DMCopyDMTS - copies a DM context to a new DM 20824989b8cSPeter Brune 20924989b8cSPeter Brune Logically Collective 21024989b8cSPeter Brune 21124989b8cSPeter Brune Input Arguments: 21224989b8cSPeter Brune + dmsrc - DM to obtain context from 21324989b8cSPeter Brune - dmdest - DM to add context to 21424989b8cSPeter Brune 21524989b8cSPeter Brune Level: developer 21624989b8cSPeter Brune 21724989b8cSPeter Brune Note: 21824989b8cSPeter Brune The context is copied by reference. This function does not ensure that a context exists. 21924989b8cSPeter Brune 220942e3340SBarry Smith .seealso: DMGetDMTS(), TSSetDM() 22124989b8cSPeter Brune @*/ 222942e3340SBarry Smith PetscErrorCode DMCopyDMTS(DM dmsrc,DM dmdest) 22324989b8cSPeter Brune { 22424989b8cSPeter Brune PetscErrorCode ierr; 22524989b8cSPeter Brune 22624989b8cSPeter Brune PetscFunctionBegin; 22724989b8cSPeter Brune PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1); 22824989b8cSPeter Brune PetscValidHeaderSpecific(dmdest,DM_CLASSID,2); 2292a34c10cSBarry Smith ierr = DMTSDestroy((DMTS*)&dmdest->dmts);CHKERRQ(ierr); 2302a34c10cSBarry Smith dmdest->dmts = dmsrc->dmts; 2312a34c10cSBarry Smith ierr = PetscObjectReference(dmdest->dmts);CHKERRQ(ierr); 232942e3340SBarry Smith ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_DMTS,DMRestrictHook_DMTS,PETSC_NULL);CHKERRQ(ierr); 23324989b8cSPeter Brune PetscFunctionReturn(0); 23424989b8cSPeter Brune } 23524989b8cSPeter Brune 23624989b8cSPeter Brune #undef __FUNCT__ 23724989b8cSPeter Brune #define __FUNCT__ "DMTSSetIFunction" 23824989b8cSPeter Brune /*@C 23924989b8cSPeter Brune DMTSSetIFunction - set TS implicit function evaluation function 24024989b8cSPeter Brune 24124989b8cSPeter Brune Not Collective 24224989b8cSPeter Brune 24324989b8cSPeter Brune Input Arguments: 24424989b8cSPeter Brune + dm - DM to be used with TS 24524989b8cSPeter Brune . func - function evaluation function, see TSSetIFunction() for calling sequence 24624989b8cSPeter Brune - ctx - context for residual evaluation 24724989b8cSPeter Brune 24824989b8cSPeter Brune Level: advanced 24924989b8cSPeter Brune 25024989b8cSPeter Brune Note: 25124989b8cSPeter Brune TSSetFunction() is normally used, but it calls this function internally because the user context is actually 25224989b8cSPeter Brune associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 25324989b8cSPeter Brune not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 25424989b8cSPeter Brune 25524989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 25624989b8cSPeter Brune @*/ 25724989b8cSPeter Brune PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx) 25824989b8cSPeter Brune { 25924989b8cSPeter Brune PetscErrorCode ierr; 260942e3340SBarry Smith DMTS tsdm; 26124989b8cSPeter Brune 26224989b8cSPeter Brune PetscFunctionBegin; 26324989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 264942e3340SBarry Smith ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 265d74926cbSBarry Smith if (func) tsdm->ops->ifunction = func; 26624989b8cSPeter Brune if (ctx) tsdm->ifunctionctx = ctx; 26724989b8cSPeter Brune PetscFunctionReturn(0); 26824989b8cSPeter Brune } 26924989b8cSPeter Brune 27024989b8cSPeter Brune #undef __FUNCT__ 27124989b8cSPeter Brune #define __FUNCT__ "DMTSGetIFunction" 27224989b8cSPeter Brune /*@C 27324989b8cSPeter Brune DMTSGetIFunction - get TS implicit residual evaluation function 27424989b8cSPeter Brune 27524989b8cSPeter Brune Not Collective 27624989b8cSPeter Brune 27724989b8cSPeter Brune Input Argument: 27824989b8cSPeter Brune . dm - DM to be used with TS 27924989b8cSPeter Brune 28024989b8cSPeter Brune Output Arguments: 28124989b8cSPeter Brune + func - function evaluation function, see TSSetIFunction() for calling sequence 28224989b8cSPeter Brune - ctx - context for residual evaluation 28324989b8cSPeter Brune 28424989b8cSPeter Brune Level: advanced 28524989b8cSPeter Brune 28624989b8cSPeter Brune Note: 28724989b8cSPeter Brune TSGetFunction() is normally used, but it calls this function internally because the user context is actually 28824989b8cSPeter Brune associated with the DM. 28924989b8cSPeter Brune 29024989b8cSPeter Brune .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction() 29124989b8cSPeter Brune @*/ 29224989b8cSPeter Brune PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx) 29324989b8cSPeter Brune { 29424989b8cSPeter Brune PetscErrorCode ierr; 295942e3340SBarry Smith DMTS tsdm; 29624989b8cSPeter Brune 29724989b8cSPeter Brune PetscFunctionBegin; 29824989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 299942e3340SBarry Smith ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 300d74926cbSBarry Smith if (func) *func = tsdm->ops->ifunction; 30124989b8cSPeter Brune if (ctx) *ctx = tsdm->ifunctionctx; 30224989b8cSPeter Brune PetscFunctionReturn(0); 30324989b8cSPeter Brune } 30424989b8cSPeter Brune 30524989b8cSPeter Brune 30624989b8cSPeter Brune #undef __FUNCT__ 30724989b8cSPeter Brune #define __FUNCT__ "DMTSSetRHSFunction" 30824989b8cSPeter Brune /*@C 30924989b8cSPeter Brune DMTSSetRHSFunction - set TS explicit residual evaluation function 31024989b8cSPeter Brune 31124989b8cSPeter Brune Not Collective 31224989b8cSPeter Brune 31324989b8cSPeter Brune Input Arguments: 31424989b8cSPeter Brune + dm - DM to be used with TS 31524989b8cSPeter Brune . func - RHS function evaluation function, see TSSetRHSFunction() for calling sequence 31624989b8cSPeter Brune - ctx - context for residual evaluation 31724989b8cSPeter Brune 31824989b8cSPeter Brune Level: advanced 31924989b8cSPeter Brune 32024989b8cSPeter Brune Note: 321ef20d060SBarry Smith TSSetRSHFunction() is normally used, but it calls this function internally because the user context is actually 32224989b8cSPeter Brune associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 32324989b8cSPeter Brune not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 32424989b8cSPeter Brune 32524989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 32624989b8cSPeter Brune @*/ 32724989b8cSPeter Brune PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx) 32824989b8cSPeter Brune { 32924989b8cSPeter Brune PetscErrorCode ierr; 330942e3340SBarry Smith DMTS tsdm; 33124989b8cSPeter Brune 33224989b8cSPeter Brune PetscFunctionBegin; 33324989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 334942e3340SBarry Smith ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 335d74926cbSBarry Smith if (func) tsdm->ops->rhsfunction = func; 33624989b8cSPeter Brune if (ctx) tsdm->rhsfunctionctx = ctx; 33724989b8cSPeter Brune PetscFunctionReturn(0); 33824989b8cSPeter Brune } 33924989b8cSPeter Brune 34024989b8cSPeter Brune #undef __FUNCT__ 341ef20d060SBarry Smith #define __FUNCT__ "DMTSGetSolutionFunction" 342ef20d060SBarry Smith /*@C 343ef20d060SBarry Smith DMTSGetSolutionFunction - gets the TS solution evaluation function 344ef20d060SBarry Smith 345ef20d060SBarry Smith Not Collective 346ef20d060SBarry Smith 347ef20d060SBarry Smith Input Arguments: 348ef20d060SBarry Smith . dm - DM to be used with TS 349ef20d060SBarry Smith 350ef20d060SBarry Smith Output Parameters: 351ef20d060SBarry Smith + func - solution function evaluation function, see TSSetSolution() for calling sequence 352ef20d060SBarry Smith - ctx - context for solution evaluation 353ef20d060SBarry Smith 354ef20d060SBarry Smith Level: advanced 355ef20d060SBarry Smith 356ef20d060SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 357ef20d060SBarry Smith @*/ 358ef20d060SBarry Smith PetscErrorCode DMTSGetSolutionFunction(DM dm,TSSolutionFunction *func,void **ctx) 359ef20d060SBarry Smith { 360ef20d060SBarry Smith PetscErrorCode ierr; 361942e3340SBarry Smith DMTS tsdm; 362ef20d060SBarry Smith 363ef20d060SBarry Smith PetscFunctionBegin; 364ef20d060SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 365942e3340SBarry Smith ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 366d74926cbSBarry Smith if (func) *func = tsdm->ops->solution; 367ef20d060SBarry Smith if (ctx) *ctx = tsdm->solutionctx; 368ef20d060SBarry Smith PetscFunctionReturn(0); 369ef20d060SBarry Smith } 370ef20d060SBarry Smith 371ef20d060SBarry Smith #undef __FUNCT__ 372ef20d060SBarry Smith #define __FUNCT__ "DMTSSetSolutionFunction" 373ef20d060SBarry Smith /*@C 374ef20d060SBarry Smith DMTSSetSolutionFunction - set TS solution evaluation function 375ef20d060SBarry Smith 376ef20d060SBarry Smith Not Collective 377ef20d060SBarry Smith 378ef20d060SBarry Smith Input Arguments: 379ef20d060SBarry Smith + dm - DM to be used with TS 380ef20d060SBarry Smith . func - solution function evaluation function, see TSSetSolution() for calling sequence 381ef20d060SBarry Smith - ctx - context for solution evaluation 382ef20d060SBarry Smith 383ef20d060SBarry Smith Level: advanced 384ef20d060SBarry Smith 385ef20d060SBarry Smith Note: 386ef20d060SBarry Smith TSSetSolutionFunction() is normally used, but it calls this function internally because the user context is actually 387ef20d060SBarry Smith associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 388ef20d060SBarry Smith not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 389ef20d060SBarry Smith 390ef20d060SBarry Smith .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 391ef20d060SBarry Smith @*/ 392ef20d060SBarry Smith PetscErrorCode DMTSSetSolutionFunction(DM dm,TSSolutionFunction func,void *ctx) 393ef20d060SBarry Smith { 394ef20d060SBarry Smith PetscErrorCode ierr; 395942e3340SBarry Smith DMTS tsdm; 396ef20d060SBarry Smith 397ef20d060SBarry Smith PetscFunctionBegin; 398ef20d060SBarry Smith PetscValidHeaderSpecific(dm,DM_CLASSID,1); 399942e3340SBarry Smith ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 400d74926cbSBarry Smith if (func) tsdm->ops->solution = func; 401ef20d060SBarry Smith if (ctx) tsdm->solutionctx = ctx; 402ef20d060SBarry Smith PetscFunctionReturn(0); 403ef20d060SBarry Smith } 404ef20d060SBarry Smith 405ef20d060SBarry Smith #undef __FUNCT__ 40624989b8cSPeter Brune #define __FUNCT__ "DMTSGetRHSFunction" 40724989b8cSPeter Brune /*@C 40824989b8cSPeter Brune DMTSGetRHSFunction - get TS explicit residual evaluation function 40924989b8cSPeter Brune 41024989b8cSPeter Brune Not Collective 41124989b8cSPeter Brune 41224989b8cSPeter Brune Input Argument: 41324989b8cSPeter Brune . dm - DM to be used with TS 41424989b8cSPeter Brune 41524989b8cSPeter Brune Output Arguments: 41624989b8cSPeter Brune + func - residual evaluation function, see TSSetRHSFunction() for calling sequence 41724989b8cSPeter Brune - ctx - context for residual evaluation 41824989b8cSPeter Brune 41924989b8cSPeter Brune Level: advanced 42024989b8cSPeter Brune 42124989b8cSPeter Brune Note: 42224989b8cSPeter Brune TSGetFunction() is normally used, but it calls this function internally because the user context is actually 42324989b8cSPeter Brune associated with the DM. 42424989b8cSPeter Brune 42524989b8cSPeter Brune .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction() 42624989b8cSPeter Brune @*/ 42724989b8cSPeter Brune PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx) 42824989b8cSPeter Brune { 42924989b8cSPeter Brune PetscErrorCode ierr; 430942e3340SBarry Smith DMTS tsdm; 43124989b8cSPeter Brune 43224989b8cSPeter Brune PetscFunctionBegin; 43324989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 434942e3340SBarry Smith ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 435d74926cbSBarry Smith if (func) *func = tsdm->ops->rhsfunction; 43624989b8cSPeter Brune if (ctx) *ctx = tsdm->rhsfunctionctx; 43724989b8cSPeter Brune PetscFunctionReturn(0); 43824989b8cSPeter Brune } 43924989b8cSPeter Brune 44024989b8cSPeter Brune #undef __FUNCT__ 44124989b8cSPeter Brune #define __FUNCT__ "DMTSSetIJacobian" 44224989b8cSPeter Brune /*@C 44324989b8cSPeter Brune DMTSSetIJacobian - set TS Jacobian evaluation function 44424989b8cSPeter Brune 44524989b8cSPeter Brune Not Collective 44624989b8cSPeter Brune 44724989b8cSPeter Brune Input Argument: 44824989b8cSPeter Brune + dm - DM to be used with TS 44924989b8cSPeter Brune . func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence 45024989b8cSPeter Brune - ctx - context for residual evaluation 45124989b8cSPeter Brune 45224989b8cSPeter Brune Level: advanced 45324989b8cSPeter Brune 45424989b8cSPeter Brune Note: 45524989b8cSPeter Brune TSSetJacobian() is normally used, but it calls this function internally because the user context is actually 45624989b8cSPeter Brune associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 45724989b8cSPeter Brune not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 45824989b8cSPeter Brune 45924989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian() 46024989b8cSPeter Brune @*/ 46124989b8cSPeter Brune PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx) 46224989b8cSPeter Brune { 46324989b8cSPeter Brune PetscErrorCode ierr; 464942e3340SBarry Smith DMTS sdm; 46524989b8cSPeter Brune 46624989b8cSPeter Brune PetscFunctionBegin; 46724989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 468942e3340SBarry Smith ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr); 469d74926cbSBarry Smith if (func) sdm->ops->ijacobian = func; 47024989b8cSPeter Brune if (ctx) sdm->ijacobianctx = ctx; 47124989b8cSPeter Brune PetscFunctionReturn(0); 47224989b8cSPeter Brune } 47324989b8cSPeter Brune 47424989b8cSPeter Brune #undef __FUNCT__ 47524989b8cSPeter Brune #define __FUNCT__ "DMTSGetIJacobian" 47624989b8cSPeter Brune /*@C 47724989b8cSPeter Brune DMTSGetIJacobian - get TS Jacobian evaluation function 47824989b8cSPeter Brune 47924989b8cSPeter Brune Not Collective 48024989b8cSPeter Brune 48124989b8cSPeter Brune Input Argument: 48224989b8cSPeter Brune . dm - DM to be used with TS 48324989b8cSPeter Brune 48424989b8cSPeter Brune Output Arguments: 48524989b8cSPeter Brune + func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence 48624989b8cSPeter Brune - ctx - context for residual evaluation 48724989b8cSPeter Brune 48824989b8cSPeter Brune Level: advanced 48924989b8cSPeter Brune 49024989b8cSPeter Brune Note: 49124989b8cSPeter Brune TSGetJacobian() is normally used, but it calls this function internally because the user context is actually 49224989b8cSPeter Brune associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 49324989b8cSPeter Brune not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 49424989b8cSPeter Brune 49524989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 49624989b8cSPeter Brune @*/ 49724989b8cSPeter Brune PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx) 49824989b8cSPeter Brune { 49924989b8cSPeter Brune PetscErrorCode ierr; 500942e3340SBarry Smith DMTS tsdm; 50124989b8cSPeter Brune 50224989b8cSPeter Brune PetscFunctionBegin; 50324989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 504942e3340SBarry Smith ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 505d74926cbSBarry Smith if (func) *func = tsdm->ops->ijacobian; 50624989b8cSPeter Brune if (ctx) *ctx = tsdm->ijacobianctx; 50724989b8cSPeter Brune PetscFunctionReturn(0); 50824989b8cSPeter Brune } 50924989b8cSPeter Brune 51024989b8cSPeter Brune 51124989b8cSPeter Brune #undef __FUNCT__ 51224989b8cSPeter Brune #define __FUNCT__ "DMTSSetRHSJacobian" 51324989b8cSPeter Brune /*@C 51424989b8cSPeter Brune DMTSSetRHSJacobian - set TS Jacobian evaluation function 51524989b8cSPeter Brune 51624989b8cSPeter Brune Not Collective 51724989b8cSPeter Brune 51824989b8cSPeter Brune Input Argument: 51924989b8cSPeter Brune + dm - DM to be used with TS 52024989b8cSPeter Brune . func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence 52124989b8cSPeter Brune - ctx - context for residual evaluation 52224989b8cSPeter Brune 52324989b8cSPeter Brune Level: advanced 52424989b8cSPeter Brune 52524989b8cSPeter Brune Note: 52624989b8cSPeter Brune TSSetJacobian() is normally used, but it calls this function internally because the user context is actually 52724989b8cSPeter Brune associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 52824989b8cSPeter Brune not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 52924989b8cSPeter Brune 53024989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian() 53124989b8cSPeter Brune @*/ 53224989b8cSPeter Brune PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx) 53324989b8cSPeter Brune { 53424989b8cSPeter Brune PetscErrorCode ierr; 535942e3340SBarry Smith DMTS tsdm; 53624989b8cSPeter Brune 53724989b8cSPeter Brune PetscFunctionBegin; 53824989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 539942e3340SBarry Smith ierr = DMGetDMTSWrite(dm,&tsdm);CHKERRQ(ierr); 540d74926cbSBarry Smith if (func) tsdm->ops->rhsjacobian = func; 54124989b8cSPeter Brune if (ctx) tsdm->rhsjacobianctx = ctx; 54224989b8cSPeter Brune PetscFunctionReturn(0); 54324989b8cSPeter Brune } 54424989b8cSPeter Brune 54524989b8cSPeter Brune #undef __FUNCT__ 54624989b8cSPeter Brune #define __FUNCT__ "DMTSGetRHSJacobian" 54724989b8cSPeter Brune /*@C 54824989b8cSPeter Brune DMTSGetRHSJacobian - get TS Jacobian evaluation function 54924989b8cSPeter Brune 55024989b8cSPeter Brune Not Collective 55124989b8cSPeter Brune 55224989b8cSPeter Brune Input Argument: 55324989b8cSPeter Brune . dm - DM to be used with TS 55424989b8cSPeter Brune 55524989b8cSPeter Brune Output Arguments: 55624989b8cSPeter Brune + func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence 55724989b8cSPeter Brune - ctx - context for residual evaluation 55824989b8cSPeter Brune 55924989b8cSPeter Brune Level: advanced 56024989b8cSPeter Brune 56124989b8cSPeter Brune Note: 56224989b8cSPeter Brune TSGetJacobian() is normally used, but it calls this function internally because the user context is actually 56324989b8cSPeter Brune associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 56424989b8cSPeter Brune not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 56524989b8cSPeter Brune 56624989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 56724989b8cSPeter Brune @*/ 56824989b8cSPeter Brune PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx) 56924989b8cSPeter Brune { 57024989b8cSPeter Brune PetscErrorCode ierr; 571942e3340SBarry Smith DMTS tsdm; 57224989b8cSPeter Brune 57324989b8cSPeter Brune PetscFunctionBegin; 57424989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 575942e3340SBarry Smith ierr = DMGetDMTS(dm,&tsdm);CHKERRQ(ierr); 576d74926cbSBarry Smith if (func) *func = tsdm->ops->rhsjacobian; 57724989b8cSPeter Brune if (ctx) *ctx = tsdm->rhsjacobianctx; 57824989b8cSPeter Brune PetscFunctionReturn(0); 57924989b8cSPeter Brune } 580