1*24989b8cSPeter Brune #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 2*24989b8cSPeter Brune #include <petsc-private/dmimpl.h> /*I "petscdm.h" I*/ 3*24989b8cSPeter Brune 4*24989b8cSPeter Brune 5*24989b8cSPeter Brune #undef __FUNCT__ 6*24989b8cSPeter Brune #define __FUNCT__ "DMCoarsenHook_TSDM" 7*24989b8cSPeter Brune /* Attaches the SNESDM to the coarse level. 8*24989b8cSPeter Brune * Under what conditions should we copy versus duplicate? 9*24989b8cSPeter Brune */ 10*24989b8cSPeter Brune static PetscErrorCode DMCoarsenHook_TSDM(DM dm,DM dmc,void *ctx) 11*24989b8cSPeter Brune { 12*24989b8cSPeter Brune PetscErrorCode ierr; 13*24989b8cSPeter Brune 14*24989b8cSPeter Brune PetscFunctionBegin; 15*24989b8cSPeter Brune ierr = DMTSCopyContext(dm,dmc);CHKERRQ(ierr); 16*24989b8cSPeter Brune PetscFunctionReturn(0); 17*24989b8cSPeter Brune } 18*24989b8cSPeter Brune 19*24989b8cSPeter Brune #undef __FUNCT__ 20*24989b8cSPeter Brune #define __FUNCT__ "DMRestrictHook_SNESDM" 21*24989b8cSPeter Brune /* This could restrict auxiliary information to the coarse level. 22*24989b8cSPeter Brune */ 23*24989b8cSPeter Brune static PetscErrorCode DMRestrictHook_TSDM(DM dm,Mat Restrict,Vec rscale,Mat Inject,DM dmc,void *ctx) 24*24989b8cSPeter Brune { 25*24989b8cSPeter Brune 26*24989b8cSPeter Brune PetscFunctionBegin; 27*24989b8cSPeter Brune PetscFunctionReturn(0); 28*24989b8cSPeter Brune } 29*24989b8cSPeter Brune 30*24989b8cSPeter Brune 31*24989b8cSPeter Brune #undef __FUNCT__ 32*24989b8cSPeter Brune #define __FUNCT__ "PetscContainerDestroy_TSDM" 33*24989b8cSPeter Brune static PetscErrorCode PetscContainerDestroy_TSDM(void *ctx) 34*24989b8cSPeter Brune { 35*24989b8cSPeter Brune PetscErrorCode ierr; 36*24989b8cSPeter Brune TSDM tsdm = (TSDM)ctx; 37*24989b8cSPeter Brune 38*24989b8cSPeter Brune PetscFunctionBegin; 39*24989b8cSPeter Brune if (tsdm->destroy) {ierr = (*tsdm->destroy)(tsdm);CHKERRQ(ierr);} 40*24989b8cSPeter Brune ierr = PetscFree(tsdm);CHKERRQ(ierr); 41*24989b8cSPeter Brune PetscFunctionReturn(0); 42*24989b8cSPeter Brune } 43*24989b8cSPeter Brune 44*24989b8cSPeter Brune #undef __FUNCT__ 45*24989b8cSPeter Brune #define __FUNCT__ "DMTSGetContext" 46*24989b8cSPeter Brune /*@C 47*24989b8cSPeter Brune DMTSGetContext - get read-only private TSDM context from a DM 48*24989b8cSPeter Brune 49*24989b8cSPeter Brune Not Collective 50*24989b8cSPeter Brune 51*24989b8cSPeter Brune Input Argument: 52*24989b8cSPeter Brune . dm - DM to be used with TS 53*24989b8cSPeter Brune 54*24989b8cSPeter Brune Output Argument: 55*24989b8cSPeter Brune . tsdm - private TSDM context 56*24989b8cSPeter Brune 57*24989b8cSPeter Brune Level: developer 58*24989b8cSPeter Brune 59*24989b8cSPeter Brune Notes: 60*24989b8cSPeter Brune Use DMTSGetContextWrite() if write access is needed. The DMTSSetXXX API should be used wherever possible. 61*24989b8cSPeter Brune 62*24989b8cSPeter Brune .seealso: DMTSGetContextWrite() 63*24989b8cSPeter Brune @*/ 64*24989b8cSPeter Brune PetscErrorCode DMTSGetContext(DM dm,TSDM *tsdm) 65*24989b8cSPeter Brune { 66*24989b8cSPeter Brune PetscErrorCode ierr; 67*24989b8cSPeter Brune PetscContainer container; 68*24989b8cSPeter Brune TSDM tsdmnew; 69*24989b8cSPeter Brune 70*24989b8cSPeter Brune 71*24989b8cSPeter Brune PetscFunctionBegin; 72*24989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 73*24989b8cSPeter Brune ierr = PetscObjectQuery((PetscObject)dm,"TSDM",(PetscObject*)&container);CHKERRQ(ierr); 74*24989b8cSPeter Brune if (container) { 75*24989b8cSPeter Brune ierr = PetscContainerGetPointer(container,(void**)tsdm);CHKERRQ(ierr); 76*24989b8cSPeter Brune } else { 77*24989b8cSPeter Brune ierr = PetscInfo(dm,"Creating new TSDM\n");CHKERRQ(ierr); 78*24989b8cSPeter Brune ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr); 79*24989b8cSPeter Brune ierr = PetscNewLog(dm,struct _n_TSDM,&tsdmnew);CHKERRQ(ierr); 80*24989b8cSPeter Brune ierr = PetscContainerSetPointer(container,tsdmnew);CHKERRQ(ierr); 81*24989b8cSPeter Brune ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_TSDM);CHKERRQ(ierr); 82*24989b8cSPeter Brune ierr = PetscObjectCompose((PetscObject)dm,"TSDM",(PetscObject)container);CHKERRQ(ierr); 83*24989b8cSPeter Brune ierr = PetscContainerGetPointer(container,(void**)tsdm);CHKERRQ(ierr); 84*24989b8cSPeter Brune ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 85*24989b8cSPeter Brune } 86*24989b8cSPeter Brune PetscFunctionReturn(0); 87*24989b8cSPeter Brune } 88*24989b8cSPeter Brune 89*24989b8cSPeter Brune #undef __FUNCT__ 90*24989b8cSPeter Brune #define __FUNCT__ "DMTSGetContextWrite" 91*24989b8cSPeter Brune /*@C 92*24989b8cSPeter Brune DMTSGetContextWrite - get write access to private TSDM context from a DM 93*24989b8cSPeter Brune 94*24989b8cSPeter Brune Not Collective 95*24989b8cSPeter Brune 96*24989b8cSPeter Brune Input Argument: 97*24989b8cSPeter Brune . dm - DM to be used with TS 98*24989b8cSPeter Brune 99*24989b8cSPeter Brune Output Argument: 100*24989b8cSPeter Brune . tsdm - private TSDM context 101*24989b8cSPeter Brune 102*24989b8cSPeter Brune Level: developer 103*24989b8cSPeter Brune 104*24989b8cSPeter Brune .seealso: DMTSGetContext() 105*24989b8cSPeter Brune @*/ 106*24989b8cSPeter Brune PetscErrorCode DMTSGetContextWrite(DM dm,TSDM *tsdm) 107*24989b8cSPeter Brune { 108*24989b8cSPeter Brune PetscErrorCode ierr; 109*24989b8cSPeter Brune TSDM sdm; 110*24989b8cSPeter Brune 111*24989b8cSPeter Brune PetscFunctionBegin; 112*24989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 113*24989b8cSPeter Brune ierr = DMTSGetContext(dm,&sdm);CHKERRQ(ierr); 114*24989b8cSPeter Brune if (!sdm->originaldm) sdm->originaldm = dm; 115*24989b8cSPeter Brune if (sdm->originaldm != dm) { /* Copy on write */ 116*24989b8cSPeter Brune PetscContainer container; 117*24989b8cSPeter Brune TSDM oldsdm = sdm; 118*24989b8cSPeter Brune ierr = PetscInfo(dm,"Copying TSDM due to write\n");CHKERRQ(ierr); 119*24989b8cSPeter Brune ierr = PetscContainerCreate(((PetscObject)dm)->comm,&container);CHKERRQ(ierr); 120*24989b8cSPeter Brune ierr = PetscNewLog(dm,struct _n_TSDM,&sdm);CHKERRQ(ierr); 121*24989b8cSPeter Brune ierr = PetscMemcpy(sdm,oldsdm,sizeof *sdm);CHKERRQ(ierr); 122*24989b8cSPeter Brune ierr = PetscContainerSetPointer(container,sdm);CHKERRQ(ierr); 123*24989b8cSPeter Brune ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_TSDM);CHKERRQ(ierr); 124*24989b8cSPeter Brune ierr = PetscObjectCompose((PetscObject)dm,"TSDM",(PetscObject)container);CHKERRQ(ierr); 125*24989b8cSPeter Brune ierr = PetscContainerDestroy(&container);CHKERRQ(ierr); 126*24989b8cSPeter Brune } 127*24989b8cSPeter Brune *tsdm = sdm; 128*24989b8cSPeter Brune PetscFunctionReturn(0); 129*24989b8cSPeter Brune } 130*24989b8cSPeter Brune 131*24989b8cSPeter Brune #undef __FUNCT__ 132*24989b8cSPeter Brune #define __FUNCT__ "DMTSCopyContext" 133*24989b8cSPeter Brune /*@C 134*24989b8cSPeter Brune DMTSCopyContext - copies a DM context to a new DM 135*24989b8cSPeter Brune 136*24989b8cSPeter Brune Logically Collective 137*24989b8cSPeter Brune 138*24989b8cSPeter Brune Input Arguments: 139*24989b8cSPeter Brune + dmsrc - DM to obtain context from 140*24989b8cSPeter Brune - dmdest - DM to add context to 141*24989b8cSPeter Brune 142*24989b8cSPeter Brune Level: developer 143*24989b8cSPeter Brune 144*24989b8cSPeter Brune Note: 145*24989b8cSPeter Brune The context is copied by reference. This function does not ensure that a context exists. 146*24989b8cSPeter Brune 147*24989b8cSPeter Brune .seealso: DMTSGetContext(), TSSetDM() 148*24989b8cSPeter Brune @*/ 149*24989b8cSPeter Brune PetscErrorCode DMTSCopyContext(DM dmsrc,DM dmdest) 150*24989b8cSPeter Brune { 151*24989b8cSPeter Brune PetscErrorCode ierr; 152*24989b8cSPeter Brune PetscContainer container; 153*24989b8cSPeter Brune 154*24989b8cSPeter Brune PetscFunctionBegin; 155*24989b8cSPeter Brune PetscValidHeaderSpecific(dmsrc,DM_CLASSID,1); 156*24989b8cSPeter Brune PetscValidHeaderSpecific(dmdest,DM_CLASSID,2); 157*24989b8cSPeter Brune ierr = PetscObjectQuery((PetscObject)dmsrc,"TSDM",(PetscObject*)&container);CHKERRQ(ierr); 158*24989b8cSPeter Brune if (container) { 159*24989b8cSPeter Brune ierr = PetscObjectCompose((PetscObject)dmdest,"TSDM",(PetscObject)container);CHKERRQ(ierr); 160*24989b8cSPeter Brune ierr = DMCoarsenHookAdd(dmdest,DMCoarsenHook_TSDM,DMRestrictHook_TSDM,PETSC_NULL);CHKERRQ(ierr); 161*24989b8cSPeter Brune } 162*24989b8cSPeter Brune PetscFunctionReturn(0); 163*24989b8cSPeter Brune } 164*24989b8cSPeter Brune 165*24989b8cSPeter Brune #undef __FUNCT__ 166*24989b8cSPeter Brune #define __FUNCT__ "DMTSSetIFunction" 167*24989b8cSPeter Brune /*@C 168*24989b8cSPeter Brune DMTSSetIFunction - set TS implicit function evaluation function 169*24989b8cSPeter Brune 170*24989b8cSPeter Brune Not Collective 171*24989b8cSPeter Brune 172*24989b8cSPeter Brune Input Arguments: 173*24989b8cSPeter Brune + dm - DM to be used with TS 174*24989b8cSPeter Brune . func - function evaluation function, see TSSetIFunction() for calling sequence 175*24989b8cSPeter Brune - ctx - context for residual evaluation 176*24989b8cSPeter Brune 177*24989b8cSPeter Brune Level: advanced 178*24989b8cSPeter Brune 179*24989b8cSPeter Brune Note: 180*24989b8cSPeter Brune TSSetFunction() is normally used, but it calls this function internally because the user context is actually 181*24989b8cSPeter Brune associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 182*24989b8cSPeter Brune not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 183*24989b8cSPeter Brune 184*24989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 185*24989b8cSPeter Brune @*/ 186*24989b8cSPeter Brune PetscErrorCode DMTSSetIFunction(DM dm,TSIFunction func,void *ctx) 187*24989b8cSPeter Brune { 188*24989b8cSPeter Brune PetscErrorCode ierr; 189*24989b8cSPeter Brune TSDM tsdm; 190*24989b8cSPeter Brune 191*24989b8cSPeter Brune PetscFunctionBegin; 192*24989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 193*24989b8cSPeter Brune ierr = DMTSGetContextWrite(dm,&tsdm);CHKERRQ(ierr); 194*24989b8cSPeter Brune if (func) tsdm->ifunction = func; 195*24989b8cSPeter Brune if (ctx) tsdm->ifunctionctx = ctx; 196*24989b8cSPeter Brune PetscFunctionReturn(0); 197*24989b8cSPeter Brune } 198*24989b8cSPeter Brune 199*24989b8cSPeter Brune #undef __FUNCT__ 200*24989b8cSPeter Brune #define __FUNCT__ "DMTSGetIFunction" 201*24989b8cSPeter Brune /*@C 202*24989b8cSPeter Brune DMTSGetIFunction - get TS implicit residual evaluation function 203*24989b8cSPeter Brune 204*24989b8cSPeter Brune Not Collective 205*24989b8cSPeter Brune 206*24989b8cSPeter Brune Input Argument: 207*24989b8cSPeter Brune . dm - DM to be used with TS 208*24989b8cSPeter Brune 209*24989b8cSPeter Brune Output Arguments: 210*24989b8cSPeter Brune + func - function evaluation function, see TSSetIFunction() for calling sequence 211*24989b8cSPeter Brune - ctx - context for residual evaluation 212*24989b8cSPeter Brune 213*24989b8cSPeter Brune Level: advanced 214*24989b8cSPeter Brune 215*24989b8cSPeter Brune Note: 216*24989b8cSPeter Brune TSGetFunction() is normally used, but it calls this function internally because the user context is actually 217*24989b8cSPeter Brune associated with the DM. 218*24989b8cSPeter Brune 219*24989b8cSPeter Brune .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction() 220*24989b8cSPeter Brune @*/ 221*24989b8cSPeter Brune PetscErrorCode DMTSGetIFunction(DM dm,TSIFunction *func,void **ctx) 222*24989b8cSPeter Brune { 223*24989b8cSPeter Brune PetscErrorCode ierr; 224*24989b8cSPeter Brune TSDM tsdm; 225*24989b8cSPeter Brune 226*24989b8cSPeter Brune PetscFunctionBegin; 227*24989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 228*24989b8cSPeter Brune ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr); 229*24989b8cSPeter Brune if (func) *func = tsdm->ifunction; 230*24989b8cSPeter Brune if (ctx) *ctx = tsdm->ifunctionctx; 231*24989b8cSPeter Brune PetscFunctionReturn(0); 232*24989b8cSPeter Brune } 233*24989b8cSPeter Brune 234*24989b8cSPeter Brune 235*24989b8cSPeter Brune #undef __FUNCT__ 236*24989b8cSPeter Brune #define __FUNCT__ "DMTSSetRHSFunction" 237*24989b8cSPeter Brune /*@C 238*24989b8cSPeter Brune DMTSSetRHSFunction - set TS explicit residual evaluation function 239*24989b8cSPeter Brune 240*24989b8cSPeter Brune Not Collective 241*24989b8cSPeter Brune 242*24989b8cSPeter Brune Input Arguments: 243*24989b8cSPeter Brune + dm - DM to be used with TS 244*24989b8cSPeter Brune . func - RHS function evaluation function, see TSSetRHSFunction() for calling sequence 245*24989b8cSPeter Brune - ctx - context for residual evaluation 246*24989b8cSPeter Brune 247*24989b8cSPeter Brune Level: advanced 248*24989b8cSPeter Brune 249*24989b8cSPeter Brune Note: 250*24989b8cSPeter Brune TSSetFunction() is normally used, but it calls this function internally because the user context is actually 251*24989b8cSPeter Brune associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 252*24989b8cSPeter Brune not. If DM took a more central role at some later date, this could become the primary method of setting the residual. 253*24989b8cSPeter Brune 254*24989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 255*24989b8cSPeter Brune @*/ 256*24989b8cSPeter Brune PetscErrorCode DMTSSetRHSFunction(DM dm,TSRHSFunction func,void *ctx) 257*24989b8cSPeter Brune { 258*24989b8cSPeter Brune PetscErrorCode ierr; 259*24989b8cSPeter Brune TSDM tsdm; 260*24989b8cSPeter Brune 261*24989b8cSPeter Brune PetscFunctionBegin; 262*24989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 263*24989b8cSPeter Brune ierr = DMTSGetContextWrite(dm,&tsdm);CHKERRQ(ierr); 264*24989b8cSPeter Brune if (func) tsdm->rhsfunction = func; 265*24989b8cSPeter Brune if (ctx) tsdm->rhsfunctionctx = ctx; 266*24989b8cSPeter Brune PetscFunctionReturn(0); 267*24989b8cSPeter Brune } 268*24989b8cSPeter Brune 269*24989b8cSPeter Brune #undef __FUNCT__ 270*24989b8cSPeter Brune #define __FUNCT__ "DMTSGetRHSFunction" 271*24989b8cSPeter Brune /*@C 272*24989b8cSPeter Brune DMTSGetRHSFunction - get TS explicit residual evaluation function 273*24989b8cSPeter Brune 274*24989b8cSPeter Brune Not Collective 275*24989b8cSPeter Brune 276*24989b8cSPeter Brune Input Argument: 277*24989b8cSPeter Brune . dm - DM to be used with TS 278*24989b8cSPeter Brune 279*24989b8cSPeter Brune Output Arguments: 280*24989b8cSPeter Brune + func - residual evaluation function, see TSSetRHSFunction() for calling sequence 281*24989b8cSPeter Brune - ctx - context for residual evaluation 282*24989b8cSPeter Brune 283*24989b8cSPeter Brune Level: advanced 284*24989b8cSPeter Brune 285*24989b8cSPeter Brune Note: 286*24989b8cSPeter Brune TSGetFunction() is normally used, but it calls this function internally because the user context is actually 287*24989b8cSPeter Brune associated with the DM. 288*24989b8cSPeter Brune 289*24989b8cSPeter Brune .seealso: DMTSSetContext(), DMTSSetFunction(), TSSetFunction() 290*24989b8cSPeter Brune @*/ 291*24989b8cSPeter Brune PetscErrorCode DMTSGetRHSFunction(DM dm,TSRHSFunction *func,void **ctx) 292*24989b8cSPeter Brune { 293*24989b8cSPeter Brune PetscErrorCode ierr; 294*24989b8cSPeter Brune TSDM tsdm; 295*24989b8cSPeter Brune 296*24989b8cSPeter Brune PetscFunctionBegin; 297*24989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 298*24989b8cSPeter Brune ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr); 299*24989b8cSPeter Brune if (func) *func = tsdm->rhsfunction; 300*24989b8cSPeter Brune if (ctx) *ctx = tsdm->rhsfunctionctx; 301*24989b8cSPeter Brune PetscFunctionReturn(0); 302*24989b8cSPeter Brune } 303*24989b8cSPeter Brune 304*24989b8cSPeter Brune #undef __FUNCT__ 305*24989b8cSPeter Brune #define __FUNCT__ "DMTSSetIJacobian" 306*24989b8cSPeter Brune /*@C 307*24989b8cSPeter Brune DMTSSetIJacobian - set TS Jacobian evaluation function 308*24989b8cSPeter Brune 309*24989b8cSPeter Brune Not Collective 310*24989b8cSPeter Brune 311*24989b8cSPeter Brune Input Argument: 312*24989b8cSPeter Brune + dm - DM to be used with TS 313*24989b8cSPeter Brune . func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence 314*24989b8cSPeter Brune - ctx - context for residual evaluation 315*24989b8cSPeter Brune 316*24989b8cSPeter Brune Level: advanced 317*24989b8cSPeter Brune 318*24989b8cSPeter Brune Note: 319*24989b8cSPeter Brune TSSetJacobian() is normally used, but it calls this function internally because the user context is actually 320*24989b8cSPeter Brune associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 321*24989b8cSPeter Brune not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 322*24989b8cSPeter Brune 323*24989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian() 324*24989b8cSPeter Brune @*/ 325*24989b8cSPeter Brune PetscErrorCode DMTSSetIJacobian(DM dm,TSIJacobian func,void *ctx) 326*24989b8cSPeter Brune { 327*24989b8cSPeter Brune PetscErrorCode ierr; 328*24989b8cSPeter Brune TSDM sdm; 329*24989b8cSPeter Brune 330*24989b8cSPeter Brune PetscFunctionBegin; 331*24989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 332*24989b8cSPeter Brune ierr = DMTSGetContextWrite(dm,&sdm);CHKERRQ(ierr); 333*24989b8cSPeter Brune if (func) sdm->ijacobian = func; 334*24989b8cSPeter Brune if (ctx) sdm->ijacobianctx = ctx; 335*24989b8cSPeter Brune PetscFunctionReturn(0); 336*24989b8cSPeter Brune } 337*24989b8cSPeter Brune 338*24989b8cSPeter Brune #undef __FUNCT__ 339*24989b8cSPeter Brune #define __FUNCT__ "DMTSGetIJacobian" 340*24989b8cSPeter Brune /*@C 341*24989b8cSPeter Brune DMTSGetIJacobian - get TS Jacobian evaluation function 342*24989b8cSPeter Brune 343*24989b8cSPeter Brune Not Collective 344*24989b8cSPeter Brune 345*24989b8cSPeter Brune Input Argument: 346*24989b8cSPeter Brune . dm - DM to be used with TS 347*24989b8cSPeter Brune 348*24989b8cSPeter Brune Output Arguments: 349*24989b8cSPeter Brune + func - Jacobian evaluation function, see TSSetIJacobian() for calling sequence 350*24989b8cSPeter Brune - ctx - context for residual evaluation 351*24989b8cSPeter Brune 352*24989b8cSPeter Brune Level: advanced 353*24989b8cSPeter Brune 354*24989b8cSPeter Brune Note: 355*24989b8cSPeter Brune TSGetJacobian() is normally used, but it calls this function internally because the user context is actually 356*24989b8cSPeter Brune associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 357*24989b8cSPeter Brune not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 358*24989b8cSPeter Brune 359*24989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 360*24989b8cSPeter Brune @*/ 361*24989b8cSPeter Brune PetscErrorCode DMTSGetIJacobian(DM dm,TSIJacobian *func,void **ctx) 362*24989b8cSPeter Brune { 363*24989b8cSPeter Brune PetscErrorCode ierr; 364*24989b8cSPeter Brune TSDM tsdm; 365*24989b8cSPeter Brune 366*24989b8cSPeter Brune PetscFunctionBegin; 367*24989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 368*24989b8cSPeter Brune ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr); 369*24989b8cSPeter Brune if (func) *func = tsdm->ijacobian; 370*24989b8cSPeter Brune if (ctx) *ctx = tsdm->ijacobianctx; 371*24989b8cSPeter Brune PetscFunctionReturn(0); 372*24989b8cSPeter Brune } 373*24989b8cSPeter Brune 374*24989b8cSPeter Brune 375*24989b8cSPeter Brune #undef __FUNCT__ 376*24989b8cSPeter Brune #define __FUNCT__ "DMTSSetRHSJacobian" 377*24989b8cSPeter Brune /*@C 378*24989b8cSPeter Brune DMTSSetRHSJacobian - set TS Jacobian evaluation function 379*24989b8cSPeter Brune 380*24989b8cSPeter Brune Not Collective 381*24989b8cSPeter Brune 382*24989b8cSPeter Brune Input Argument: 383*24989b8cSPeter Brune + dm - DM to be used with TS 384*24989b8cSPeter Brune . func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence 385*24989b8cSPeter Brune - ctx - context for residual evaluation 386*24989b8cSPeter Brune 387*24989b8cSPeter Brune Level: advanced 388*24989b8cSPeter Brune 389*24989b8cSPeter Brune Note: 390*24989b8cSPeter Brune TSSetJacobian() is normally used, but it calls this function internally because the user context is actually 391*24989b8cSPeter Brune associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 392*24989b8cSPeter Brune not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 393*24989b8cSPeter Brune 394*24989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSGetJacobian(), TSSetJacobian() 395*24989b8cSPeter Brune @*/ 396*24989b8cSPeter Brune PetscErrorCode DMTSSetRHSJacobian(DM dm,TSRHSJacobian func,void *ctx) 397*24989b8cSPeter Brune { 398*24989b8cSPeter Brune PetscErrorCode ierr; 399*24989b8cSPeter Brune TSDM tsdm; 400*24989b8cSPeter Brune 401*24989b8cSPeter Brune PetscFunctionBegin; 402*24989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 403*24989b8cSPeter Brune ierr = DMTSGetContextWrite(dm,&tsdm);CHKERRQ(ierr); 404*24989b8cSPeter Brune if (func) tsdm->rhsjacobian = func; 405*24989b8cSPeter Brune if (ctx) tsdm->rhsjacobianctx = ctx; 406*24989b8cSPeter Brune PetscFunctionReturn(0); 407*24989b8cSPeter Brune } 408*24989b8cSPeter Brune 409*24989b8cSPeter Brune #undef __FUNCT__ 410*24989b8cSPeter Brune #define __FUNCT__ "DMTSGetRHSJacobian" 411*24989b8cSPeter Brune /*@C 412*24989b8cSPeter Brune DMTSGetRHSJacobian - get TS Jacobian evaluation function 413*24989b8cSPeter Brune 414*24989b8cSPeter Brune Not Collective 415*24989b8cSPeter Brune 416*24989b8cSPeter Brune Input Argument: 417*24989b8cSPeter Brune . dm - DM to be used with TS 418*24989b8cSPeter Brune 419*24989b8cSPeter Brune Output Arguments: 420*24989b8cSPeter Brune + func - Jacobian evaluation function, see TSSetRHSJacobian() for calling sequence 421*24989b8cSPeter Brune - ctx - context for residual evaluation 422*24989b8cSPeter Brune 423*24989b8cSPeter Brune Level: advanced 424*24989b8cSPeter Brune 425*24989b8cSPeter Brune Note: 426*24989b8cSPeter Brune TSGetJacobian() is normally used, but it calls this function internally because the user context is actually 427*24989b8cSPeter Brune associated with the DM. This makes the interface consistent regardless of whether the user interacts with a DM or 428*24989b8cSPeter Brune not. If DM took a more central role at some later date, this could become the primary method of setting the Jacobian. 429*24989b8cSPeter Brune 430*24989b8cSPeter Brune .seealso: DMTSSetContext(), TSSetFunction(), DMTSSetJacobian() 431*24989b8cSPeter Brune @*/ 432*24989b8cSPeter Brune PetscErrorCode DMTSGetRHSJacobian(DM dm,TSRHSJacobian *func,void **ctx) 433*24989b8cSPeter Brune { 434*24989b8cSPeter Brune PetscErrorCode ierr; 435*24989b8cSPeter Brune TSDM tsdm; 436*24989b8cSPeter Brune 437*24989b8cSPeter Brune PetscFunctionBegin; 438*24989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 439*24989b8cSPeter Brune ierr = DMTSGetContext(dm,&tsdm);CHKERRQ(ierr); 440*24989b8cSPeter Brune if (func) *func = tsdm->rhsjacobian; 441*24989b8cSPeter Brune if (ctx) *ctx = tsdm->rhsjacobianctx; 442*24989b8cSPeter Brune PetscFunctionReturn(0); 443*24989b8cSPeter Brune } 444*24989b8cSPeter Brune 445*24989b8cSPeter Brune #if 0 446*24989b8cSPeter Brune 447*24989b8cSPeter Brune #undef __FUNCT__ 448*24989b8cSPeter Brune #define __FUNCT__ "TSDefaultComputeFunction_DMLegacy" 449*24989b8cSPeter Brune static PetscErrorCode TSDefaultComputeFunction_DMLegacy(TS ts,Vec X,Vec F,void *ctx) 450*24989b8cSPeter Brune { 451*24989b8cSPeter Brune PetscErrorCode ierr; 452*24989b8cSPeter Brune DM dm; 453*24989b8cSPeter Brune 454*24989b8cSPeter Brune PetscFunctionBegin; 455*24989b8cSPeter Brune ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 456*24989b8cSPeter Brune ierr = DMComputeFunction(dm,X,F);CHKERRQ(ierr); 457*24989b8cSPeter Brune PetscFunctionReturn(0); 458*24989b8cSPeter Brune } 459*24989b8cSPeter Brune 460*24989b8cSPeter Brune #undef __FUNCT__ 461*24989b8cSPeter Brune #define __FUNCT__ "TSDefaultComputeJacobian_DMLegacy" 462*24989b8cSPeter Brune static PetscErrorCode TSDefaultComputeJacobian_DMLegacy(TS ts,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx) 463*24989b8cSPeter Brune { 464*24989b8cSPeter Brune PetscErrorCode ierr; 465*24989b8cSPeter Brune DM dm; 466*24989b8cSPeter Brune 467*24989b8cSPeter Brune PetscFunctionBegin; 468*24989b8cSPeter Brune ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 469*24989b8cSPeter Brune ierr = DMComputeJacobian(dm,X,*A,*B,mstr);CHKERRQ(ierr); 470*24989b8cSPeter Brune PetscFunctionReturn(0); 471*24989b8cSPeter Brune } 472*24989b8cSPeter Brune 473*24989b8cSPeter Brune #undef __FUNCT__ 474*24989b8cSPeter Brune #define __FUNCT__ "DMTSSetUpLegacy" 475*24989b8cSPeter Brune /* Sets up calling of legacy DM routines */ 476*24989b8cSPeter Brune PetscErrorCode DMTSSetUpLegacy(DM dm) 477*24989b8cSPeter Brune { 478*24989b8cSPeter Brune PetscErrorCode ierr; 479*24989b8cSPeter Brune TSDM sdm; 480*24989b8cSPeter Brune 481*24989b8cSPeter Brune PetscFunctionBegin; 482*24989b8cSPeter Brune ierr = DMTSGetContext(dm,&sdm);CHKERRQ(ierr); 483*24989b8cSPeter Brune if (!sdm->computefunction) {ierr = DMTSSetFunction(dm,TSDefaultComputeFunction_DMLegacy,PETSC_NULL);CHKERRQ(ierr);} 484*24989b8cSPeter Brune ierr = DMTSGetContext(dm,&sdm);CHKERRQ(ierr); 485*24989b8cSPeter Brune if (!sdm->computejacobian) { 486*24989b8cSPeter Brune if (dm->ops->functionj) { 487*24989b8cSPeter Brune ierr = DMTSSetJacobian(dm,TSDefaultComputeJacobian_DMLegacy,PETSC_NULL);CHKERRQ(ierr); 488*24989b8cSPeter Brune } else { 489*24989b8cSPeter Brune ierr = DMTSSetJacobian(dm,TSDefaultComputeJacobianColor,PETSC_NULL);CHKERRQ(ierr); 490*24989b8cSPeter Brune } 491*24989b8cSPeter Brune } 492*24989b8cSPeter Brune PetscFunctionReturn(0); 493*24989b8cSPeter Brune } 494*24989b8cSPeter Brune 495*24989b8cSPeter Brune #endif 496