124989b8cSPeter Brune #include <petscdmda.h> /*I "petscdmda.h" I*/ 224989b8cSPeter Brune #include <petsc-private/dmimpl.h> 324989b8cSPeter Brune #include <petsc-private/tsimpl.h> /*I "petscts.h" I*/ 424989b8cSPeter Brune 524989b8cSPeter Brune /* This structure holds the user-provided DMDA callbacks */ 624989b8cSPeter Brune typedef struct { 724989b8cSPeter Brune PetscErrorCode (*ifunctionlocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*); 824989b8cSPeter Brune PetscErrorCode (*rhsfunctionlocal)(DMDALocalInfo*,PetscReal,void*,void*,void*); 924989b8cSPeter Brune PetscErrorCode (*ijacobianlocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,MatStructure*,void*); 1024989b8cSPeter Brune PetscErrorCode (*rhsjacobianlocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,MatStructure*,void*); 1124989b8cSPeter Brune void *ifunctionlocalctx; 1224989b8cSPeter Brune void *ijacobianlocalctx; 1324989b8cSPeter Brune void *rhsfunctionlocalctx; 1424989b8cSPeter Brune void *rhsjacobianlocalctx; 1524989b8cSPeter Brune InsertMode ifunctionlocalimode; 1624989b8cSPeter Brune InsertMode rhsfunctionlocalimode; 17942e3340SBarry Smith } DMTS_DA; 1824989b8cSPeter Brune 1924989b8cSPeter Brune #undef __FUNCT__ 20942e3340SBarry Smith #define __FUNCT__ "DMTSDestroy_DMDA" 21942e3340SBarry Smith static PetscErrorCode DMTSDestroy_DMDA(DMTS sdm) 2224989b8cSPeter Brune { 2324989b8cSPeter Brune PetscErrorCode ierr; 2424989b8cSPeter Brune 2524989b8cSPeter Brune PetscFunctionBegin; 2624989b8cSPeter Brune ierr = PetscFree(sdm->data);CHKERRQ(ierr); 2724989b8cSPeter Brune PetscFunctionReturn(0); 2824989b8cSPeter Brune } 2924989b8cSPeter Brune 3024989b8cSPeter Brune #undef __FUNCT__ 31942e3340SBarry Smith #define __FUNCT__ "DMTSDuplicate_DMDA" 32*d74926cbSBarry Smith static PetscErrorCode DMTSDuplicate_DMDA(DMTS oldsdm,DMTS sdm) 330313d0e0SJed Brown { 340313d0e0SJed Brown PetscErrorCode ierr; 350313d0e0SJed Brown 360313d0e0SJed Brown PetscFunctionBegin; 37*d74926cbSBarry Smith ierr = PetscNewLog(sdm,DMTS_DA,&sdm->data);CHKERRQ(ierr); 38942e3340SBarry Smith if (oldsdm->data) {ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMTS_DA));CHKERRQ(ierr);} 390313d0e0SJed Brown PetscFunctionReturn(0); 400313d0e0SJed Brown } 410313d0e0SJed Brown 420313d0e0SJed Brown #undef __FUNCT__ 4324989b8cSPeter Brune #define __FUNCT__ "DMDATSGetContext" 44942e3340SBarry Smith static PetscErrorCode DMDATSGetContext(DM dm,DMTS sdm,DMTS_DA **dmdats) 4524989b8cSPeter Brune { 4624989b8cSPeter Brune PetscErrorCode ierr; 4724989b8cSPeter Brune 4824989b8cSPeter Brune PetscFunctionBegin; 4924989b8cSPeter Brune *dmdats = PETSC_NULL; 5024989b8cSPeter Brune if (!sdm->data) { 51942e3340SBarry Smith ierr = PetscNewLog(dm,DMTS_DA,&sdm->data);CHKERRQ(ierr); 52*d74926cbSBarry Smith sdm->ops->destroy = DMTSDestroy_DMDA; 53*d74926cbSBarry Smith sdm->ops->duplicate = DMTSDuplicate_DMDA; 5424989b8cSPeter Brune } 55942e3340SBarry Smith *dmdats = (DMTS_DA*)sdm->data; 5624989b8cSPeter Brune PetscFunctionReturn(0); 5724989b8cSPeter Brune } 5824989b8cSPeter Brune 5924989b8cSPeter Brune #undef __FUNCT__ 6024989b8cSPeter Brune #define __FUNCT__ "TSComputeIFunction_DMDA" 6124989b8cSPeter Brune /* 6224989b8cSPeter Brune This function should eventually replace: 6324989b8cSPeter Brune DMDAComputeFunction() and DMDAComputeFunction1() 6424989b8cSPeter Brune */ 6524989b8cSPeter Brune static PetscErrorCode TSComputeIFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,Vec F,void *ctx) 6624989b8cSPeter Brune { 6724989b8cSPeter Brune PetscErrorCode ierr; 6824989b8cSPeter Brune DM dm; 69942e3340SBarry Smith DMTS_DA *dmdats = (DMTS_DA*)ctx; 7024989b8cSPeter Brune DMDALocalInfo info; 7124989b8cSPeter Brune Vec Xloc; 7224989b8cSPeter Brune void *x,*f,*xdot; 7324989b8cSPeter Brune 7424989b8cSPeter Brune PetscFunctionBegin; 7524989b8cSPeter Brune PetscValidHeaderSpecific(ts,TS_CLASSID,1); 7624989b8cSPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 7724989b8cSPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 7824989b8cSPeter Brune if (!dmdats->ifunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context"); 7924989b8cSPeter Brune ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 8024989b8cSPeter Brune ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 8124989b8cSPeter Brune ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 8224989b8cSPeter Brune ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 8324989b8cSPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 8424989b8cSPeter Brune ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 8524989b8cSPeter Brune ierr = DMDAVecGetArray(dm,Xdot,&xdot);CHKERRQ(ierr); 8624989b8cSPeter Brune switch (dmdats->ifunctionlocalimode) { 8724989b8cSPeter Brune case INSERT_VALUES: { 8824989b8cSPeter Brune ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr); 8924989b8cSPeter Brune CHKMEMQ; 9024989b8cSPeter Brune ierr = (*dmdats->ifunctionlocal)(&info,ptime,x,xdot,f,dmdats->ifunctionlocalctx);CHKERRQ(ierr); 9124989b8cSPeter Brune CHKMEMQ; 9224989b8cSPeter Brune ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr); 9324989b8cSPeter Brune } break; 9424989b8cSPeter Brune case ADD_VALUES: { 9524989b8cSPeter Brune Vec Floc; 9624989b8cSPeter Brune ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 9724989b8cSPeter Brune ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 9824989b8cSPeter Brune ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr); 9924989b8cSPeter Brune CHKMEMQ; 10024989b8cSPeter Brune ierr = (*dmdats->ifunctionlocal)(&info,ptime,x,xdot,f,dmdats->ifunctionlocalctx);CHKERRQ(ierr); 10124989b8cSPeter Brune CHKMEMQ; 10224989b8cSPeter Brune ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr); 10324989b8cSPeter Brune ierr = VecZeroEntries(F);CHKERRQ(ierr); 10424989b8cSPeter Brune ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 10524989b8cSPeter Brune ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 10624989b8cSPeter Brune ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 10724989b8cSPeter Brune } break; 10824989b8cSPeter Brune default: SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->ifunctionlocalimode); 10924989b8cSPeter Brune } 11024989b8cSPeter Brune ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 11124989b8cSPeter Brune ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 11211c2e1cfSJed Brown ierr = DMDAVecRestoreArray(dm,Xdot,&xdot);CHKERRQ(ierr); 11324989b8cSPeter Brune PetscFunctionReturn(0); 11424989b8cSPeter Brune } 11524989b8cSPeter Brune 11624989b8cSPeter Brune #undef __FUNCT__ 11724989b8cSPeter Brune #define __FUNCT__ "TSComputeIJacobian_DMDA" 11824989b8cSPeter Brune static PetscErrorCode TSComputeIJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *mstr,void *ctx) 11924989b8cSPeter Brune { 12024989b8cSPeter Brune PetscErrorCode ierr; 12124989b8cSPeter Brune DM dm; 122942e3340SBarry Smith DMTS_DA *dmdats = (DMTS_DA*)ctx; 12324989b8cSPeter Brune DMDALocalInfo info; 12424989b8cSPeter Brune Vec Xloc; 12524989b8cSPeter Brune void *x,*xdot; 12624989b8cSPeter Brune 12724989b8cSPeter Brune PetscFunctionBegin; 12824989b8cSPeter Brune if (!dmdats->ifunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context"); 12924989b8cSPeter Brune ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 13024989b8cSPeter Brune 13124989b8cSPeter Brune if (dmdats->ijacobianlocal) { 13224989b8cSPeter Brune ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 13324989b8cSPeter Brune ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 13424989b8cSPeter Brune ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 13524989b8cSPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 13624989b8cSPeter Brune ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 13724989b8cSPeter Brune ierr = DMDAVecGetArray(dm,Xdot,&xdot);CHKERRQ(ierr); 13824989b8cSPeter Brune CHKMEMQ; 13924989b8cSPeter Brune ierr = (*dmdats->ijacobianlocal)(&info,ptime,x,xdot,shift,*A,*B,mstr,dmdats->ijacobianlocalctx);CHKERRQ(ierr); 14024989b8cSPeter Brune CHKMEMQ; 14124989b8cSPeter Brune ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 14224989b8cSPeter Brune ierr = DMDAVecRestoreArray(dm,Xdot,&xdot);CHKERRQ(ierr); 14324989b8cSPeter Brune ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 14424989b8cSPeter Brune } else { 14524989b8cSPeter Brune SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"TSComputeIJacobian_DMDA() called without calling DMDATSSetIJacobian()"); 14624989b8cSPeter Brune } 14724989b8cSPeter Brune /* This will be redundant if the user called both, but it's too common to forget. */ 14824989b8cSPeter Brune if (*A != *B) { 14924989b8cSPeter Brune ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15024989b8cSPeter Brune ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15124989b8cSPeter Brune } 15224989b8cSPeter Brune PetscFunctionReturn(0); 15324989b8cSPeter Brune } 15424989b8cSPeter Brune 15524989b8cSPeter Brune #undef __FUNCT__ 15624989b8cSPeter Brune #define __FUNCT__ "TSComputeRHSFunction_DMDA" 15724989b8cSPeter Brune /* 15824989b8cSPeter Brune This function should eventually replace: 15924989b8cSPeter Brune DMDAComputeFunction() and DMDAComputeFunction1() 16024989b8cSPeter Brune */ 16124989b8cSPeter Brune static PetscErrorCode TSComputeRHSFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec F,void *ctx) 16224989b8cSPeter Brune { 16324989b8cSPeter Brune PetscErrorCode ierr; 16424989b8cSPeter Brune DM dm; 165942e3340SBarry Smith DMTS_DA *dmdats = (DMTS_DA*)ctx; 16624989b8cSPeter Brune DMDALocalInfo info; 16724989b8cSPeter Brune Vec Xloc; 16824989b8cSPeter Brune void *x,*f; 16924989b8cSPeter Brune 17024989b8cSPeter Brune PetscFunctionBegin; 17124989b8cSPeter Brune PetscValidHeaderSpecific(ts,TS_CLASSID,1); 17224989b8cSPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 17324989b8cSPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 17424989b8cSPeter Brune if (!dmdats->rhsfunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context"); 17524989b8cSPeter Brune ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 17624989b8cSPeter Brune ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 17724989b8cSPeter Brune ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 17824989b8cSPeter Brune ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 17924989b8cSPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 18024989b8cSPeter Brune ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 18124989b8cSPeter Brune switch (dmdats->rhsfunctionlocalimode) { 18224989b8cSPeter Brune case INSERT_VALUES: { 18324989b8cSPeter Brune ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr); 18424989b8cSPeter Brune CHKMEMQ; 18524989b8cSPeter Brune ierr = (*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx);CHKERRQ(ierr); 18624989b8cSPeter Brune CHKMEMQ; 18724989b8cSPeter Brune ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr); 18824989b8cSPeter Brune } break; 18924989b8cSPeter Brune case ADD_VALUES: { 19024989b8cSPeter Brune Vec Floc; 19124989b8cSPeter Brune ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 19224989b8cSPeter Brune ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 19324989b8cSPeter Brune ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr); 19424989b8cSPeter Brune CHKMEMQ; 19524989b8cSPeter Brune ierr = (*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx);CHKERRQ(ierr); 19624989b8cSPeter Brune CHKMEMQ; 19724989b8cSPeter Brune ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr); 19824989b8cSPeter Brune ierr = VecZeroEntries(F);CHKERRQ(ierr); 19924989b8cSPeter Brune ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 20024989b8cSPeter Brune ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 20124989b8cSPeter Brune ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 20224989b8cSPeter Brune } break; 20324989b8cSPeter Brune default: SETERRQ1(((PetscObject)ts)->comm,PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->rhsfunctionlocalimode); 20424989b8cSPeter Brune } 20524989b8cSPeter Brune ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 20624989b8cSPeter Brune ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 20724989b8cSPeter Brune PetscFunctionReturn(0); 20824989b8cSPeter Brune } 20924989b8cSPeter Brune 21024989b8cSPeter Brune #undef __FUNCT__ 21124989b8cSPeter Brune #define __FUNCT__ "TSComputeRHSJacobian_DMDA" 21224989b8cSPeter Brune static PetscErrorCode TSComputeRHSJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx) 21324989b8cSPeter Brune { 21424989b8cSPeter Brune PetscErrorCode ierr; 21524989b8cSPeter Brune DM dm; 216942e3340SBarry Smith DMTS_DA *dmdats = (DMTS_DA*)ctx; 21724989b8cSPeter Brune DMDALocalInfo info; 21824989b8cSPeter Brune Vec Xloc; 21924989b8cSPeter Brune void *x; 22024989b8cSPeter Brune 22124989b8cSPeter Brune PetscFunctionBegin; 22224989b8cSPeter Brune if (!dmdats->rhsfunctionlocal) SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"Corrupt context"); 22324989b8cSPeter Brune ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 22424989b8cSPeter Brune 22524989b8cSPeter Brune if (dmdats->rhsjacobianlocal) { 22624989b8cSPeter Brune ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 22724989b8cSPeter Brune ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 22824989b8cSPeter Brune ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 22924989b8cSPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 23024989b8cSPeter Brune ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 23124989b8cSPeter Brune CHKMEMQ; 23224989b8cSPeter Brune ierr = (*dmdats->rhsjacobianlocal)(&info,ptime,x,*A,*B,mstr,dmdats->rhsjacobianlocalctx);CHKERRQ(ierr); 23324989b8cSPeter Brune CHKMEMQ; 23424989b8cSPeter Brune ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 23524989b8cSPeter Brune ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 23624989b8cSPeter Brune } else { 23724989b8cSPeter Brune SETERRQ(((PetscObject)ts)->comm,PETSC_ERR_PLIB,"TSComputeRHSJacobian_DMDA() called without calling DMDATSSetRHSJacobian()"); 23824989b8cSPeter Brune } 23924989b8cSPeter Brune /* This will be redundant if the user called both, but it's too common to forget. */ 24024989b8cSPeter Brune if (*A != *B) { 24124989b8cSPeter Brune ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 24224989b8cSPeter Brune ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 24324989b8cSPeter Brune } 24424989b8cSPeter Brune PetscFunctionReturn(0); 24524989b8cSPeter Brune } 24624989b8cSPeter Brune 24724989b8cSPeter Brune 24824989b8cSPeter Brune #undef __FUNCT__ 24924989b8cSPeter Brune #define __FUNCT__ "DMDATSSetRHSFunctionLocal" 25024989b8cSPeter Brune /*@C 25124989b8cSPeter Brune DMDATSSetRHSFunctionLocal - set a local residual evaluation function 25224989b8cSPeter Brune 25324989b8cSPeter Brune Logically Collective 25424989b8cSPeter Brune 25524989b8cSPeter Brune Input Arguments: 25624989b8cSPeter Brune + dm - DM to associate callback with 257b398ffe5SPeter Brune . imode - insert mode for the residual 25824989b8cSPeter Brune . func - local residual evaluation 25924989b8cSPeter Brune - ctx - optional context for local residual evaluation 26024989b8cSPeter Brune 26124989b8cSPeter Brune Calling sequence for func: 262b398ffe5SPeter Brune 263b398ffe5SPeter Brune $ func(DMDALocalInfo info,PetscReal t,void *x,void *f,void *ctx) 264b398ffe5SPeter Brune 26524989b8cSPeter Brune + info - DMDALocalInfo defining the subdomain to evaluate the residual on 266b398ffe5SPeter Brune . t - time at which to evaluate residual 267b398ffe5SPeter Brune . x - array of local state information 268b398ffe5SPeter Brune . f - output array of local residual information 269b398ffe5SPeter Brune - ctx - optional user context 27024989b8cSPeter Brune 27124989b8cSPeter Brune Level: beginner 27224989b8cSPeter Brune 2738e625dcfSPeter Brune .seealso: DMTSSetRHSFunction(), DMDATSSetRHSJacobianLocal(), DMDASNESSetFunctionLocal() 27424989b8cSPeter Brune @*/ 2756c6b9e74SPeter Brune PetscErrorCode DMDATSSetRHSFunctionLocal(DM dm,InsertMode imode,DMDATSRHSFunctionLocal func,void *ctx) 27624989b8cSPeter Brune { 27724989b8cSPeter Brune PetscErrorCode ierr; 278942e3340SBarry Smith DMTS sdm; 279942e3340SBarry Smith DMTS_DA *dmdats; 28024989b8cSPeter Brune 28124989b8cSPeter Brune PetscFunctionBegin; 28224989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 283942e3340SBarry Smith ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr); 28424989b8cSPeter Brune ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr); 28524989b8cSPeter Brune dmdats->rhsfunctionlocalimode = imode; 28624989b8cSPeter Brune dmdats->rhsfunctionlocal = func; 28724989b8cSPeter Brune dmdats->rhsfunctionlocalctx = ctx; 28824989b8cSPeter Brune ierr = DMTSSetRHSFunction(dm,TSComputeRHSFunction_DMDA,dmdats);CHKERRQ(ierr); 28924989b8cSPeter Brune PetscFunctionReturn(0); 29024989b8cSPeter Brune } 29124989b8cSPeter Brune 29224989b8cSPeter Brune #undef __FUNCT__ 29324989b8cSPeter Brune #define __FUNCT__ "DMDATSSetRHSJacobianLocal" 29424989b8cSPeter Brune /*@C 29524989b8cSPeter Brune DMDATSSetRHSJacobianLocal - set a local residual evaluation function 29624989b8cSPeter Brune 29724989b8cSPeter Brune Logically Collective 29824989b8cSPeter Brune 29924989b8cSPeter Brune Input Arguments: 30024989b8cSPeter Brune + dm - DM to associate callback with 301b398ffe5SPeter Brune . func - local RHS Jacobian evaluation routine 302b398ffe5SPeter Brune - ctx - optional context for local jacobian evaluation 30324989b8cSPeter Brune 30424989b8cSPeter Brune Calling sequence for func: 305b398ffe5SPeter Brune 306b398ffe5SPeter Brune $ func(DMDALocalInfo* info,PetscReal t,void* x,Mat J,Mat B,MatStructure *flg,void *ctx); 307b398ffe5SPeter Brune 30824989b8cSPeter Brune + info - DMDALocalInfo defining the subdomain to evaluate the residual on 309b398ffe5SPeter Brune . t - time at which to evaluate residual 310b398ffe5SPeter Brune . x - array of local state information 311b398ffe5SPeter Brune . J - Jacobian matrix 312b398ffe5SPeter Brune . B - preconditioner matrix; often same as J 313b398ffe5SPeter Brune . flg - flag indicating information about the preconditioner matrix structure (same as flag in KSPSetOperators()) 31424989b8cSPeter Brune - ctx - optional context passed above 31524989b8cSPeter Brune 31624989b8cSPeter Brune Level: beginner 31724989b8cSPeter Brune 3188e625dcfSPeter Brune .seealso: DMTSSetRHSJacobian(), DMDATSSetRHSFunctionLocal(), DMDASNESSetJacobianLocal() 31924989b8cSPeter Brune @*/ 3206c6b9e74SPeter Brune PetscErrorCode DMDATSSetRHSJacobianLocal(DM dm,DMDATSRHSJacobianLocal func,void *ctx) 32124989b8cSPeter Brune { 32224989b8cSPeter Brune PetscErrorCode ierr; 323942e3340SBarry Smith DMTS sdm; 324942e3340SBarry Smith DMTS_DA *dmdats; 32524989b8cSPeter Brune 32624989b8cSPeter Brune PetscFunctionBegin; 32724989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 328942e3340SBarry Smith ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr); 32924989b8cSPeter Brune ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr); 33024989b8cSPeter Brune dmdats->rhsjacobianlocal = func; 33124989b8cSPeter Brune dmdats->rhsjacobianlocalctx = ctx; 33224989b8cSPeter Brune ierr = DMTSSetRHSJacobian(dm,TSComputeRHSJacobian_DMDA,dmdats);CHKERRQ(ierr); 33324989b8cSPeter Brune PetscFunctionReturn(0); 33424989b8cSPeter Brune } 33524989b8cSPeter Brune 33624989b8cSPeter Brune 33724989b8cSPeter Brune #undef __FUNCT__ 33824989b8cSPeter Brune #define __FUNCT__ "DMDATSSetIFunctionLocal" 33924989b8cSPeter Brune /*@C 34024989b8cSPeter Brune DMDATSSetIFunctionLocal - set a local residual evaluation function 34124989b8cSPeter Brune 34224989b8cSPeter Brune Logically Collective 34324989b8cSPeter Brune 34424989b8cSPeter Brune Input Arguments: 34524989b8cSPeter Brune + dm - DM to associate callback with 34624989b8cSPeter Brune . func - local residual evaluation 34724989b8cSPeter Brune - ctx - optional context for local residual evaluation 34824989b8cSPeter Brune 34924989b8cSPeter Brune Calling sequence for func: 35024989b8cSPeter Brune + info - DMDALocalInfo defining the subdomain to evaluate the residual on 351b398ffe5SPeter Brune . t - time at which to evaluate residual 352b398ffe5SPeter Brune . x - array of local state information 353b398ffe5SPeter Brune . xdot - array of local time derivative information 354b398ffe5SPeter Brune . f - output array of local function evaluation information 35524989b8cSPeter Brune - ctx - optional context passed above 35624989b8cSPeter Brune 35724989b8cSPeter Brune Level: beginner 35824989b8cSPeter Brune 3598e625dcfSPeter Brune .seealso: DMTSSetIFunction(), DMDATSSetIJacobianLocal(), DMDASNESSetFunctionLocal() 36024989b8cSPeter Brune @*/ 3616c6b9e74SPeter Brune PetscErrorCode DMDATSSetIFunctionLocal(DM dm,InsertMode imode,DMDATSIFunctionLocal func,void *ctx) 36224989b8cSPeter Brune { 36324989b8cSPeter Brune PetscErrorCode ierr; 364942e3340SBarry Smith DMTS sdm; 365942e3340SBarry Smith DMTS_DA *dmdats; 36624989b8cSPeter Brune 36724989b8cSPeter Brune PetscFunctionBegin; 36824989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 369942e3340SBarry Smith ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr); 37024989b8cSPeter Brune ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr); 37124989b8cSPeter Brune dmdats->ifunctionlocalimode = imode; 37224989b8cSPeter Brune dmdats->ifunctionlocal = func; 37324989b8cSPeter Brune dmdats->ifunctionlocalctx = ctx; 37424989b8cSPeter Brune ierr = DMTSSetIFunction(dm,TSComputeIFunction_DMDA,dmdats);CHKERRQ(ierr); 37524989b8cSPeter Brune PetscFunctionReturn(0); 37624989b8cSPeter Brune } 37724989b8cSPeter Brune 37824989b8cSPeter Brune #undef __FUNCT__ 37924989b8cSPeter Brune #define __FUNCT__ "DMDATSSetIJacobianLocal" 38024989b8cSPeter Brune /*@C 38124989b8cSPeter Brune DMDATSSetIJacobianLocal - set a local residual evaluation function 38224989b8cSPeter Brune 38324989b8cSPeter Brune Logically Collective 38424989b8cSPeter Brune 38524989b8cSPeter Brune Input Arguments: 38624989b8cSPeter Brune + dm - DM to associate callback with 38724989b8cSPeter Brune . func - local residual evaluation 38824989b8cSPeter Brune - ctx - optional context for local residual evaluation 38924989b8cSPeter Brune 39024989b8cSPeter Brune Calling sequence for func: 391b398ffe5SPeter Brune 392b398ffe5SPeter Brune $ func(DMDALocalInfo* info,PetscReal t,void* x,void *xdot,Mat J,Mat B,MatStructure *flg,void *ctx); 393b398ffe5SPeter Brune 39424989b8cSPeter Brune + info - DMDALocalInfo defining the subdomain to evaluate the residual on 395b398ffe5SPeter Brune . t - time at which to evaluate the jacobian 396b398ffe5SPeter Brune . x - array of local state information 397b398ffe5SPeter Brune . xdot - time derivative at this state 398b398ffe5SPeter Brune . J - Jacobian matrix 399b398ffe5SPeter Brune . B - preconditioner matrix; often same as J 400b398ffe5SPeter Brune . flg - flag indicating information about the preconditioner matrix structure (same as flag in KSPSetOperators()) 40124989b8cSPeter Brune - ctx - optional context passed above 40224989b8cSPeter Brune 40324989b8cSPeter Brune Level: beginner 40424989b8cSPeter Brune 4058e625dcfSPeter Brune .seealso: DMTSSetJacobian(), DMDATSSetIFunctionLocal(), DMDASNESSetJacobianLocal() 40624989b8cSPeter Brune @*/ 4076c6b9e74SPeter Brune PetscErrorCode DMDATSSetIJacobianLocal(DM dm,DMDATSIJacobianLocal func,void *ctx) 40824989b8cSPeter Brune { 40924989b8cSPeter Brune PetscErrorCode ierr; 410942e3340SBarry Smith DMTS sdm; 411942e3340SBarry Smith DMTS_DA *dmdats; 41224989b8cSPeter Brune 41324989b8cSPeter Brune PetscFunctionBegin; 41424989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 415942e3340SBarry Smith ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr); 41624989b8cSPeter Brune ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr); 41724989b8cSPeter Brune dmdats->ijacobianlocal = func; 41824989b8cSPeter Brune dmdats->ijacobianlocalctx = ctx; 41924989b8cSPeter Brune ierr = DMTSSetIJacobian(dm,TSComputeIJacobian_DMDA,dmdats);CHKERRQ(ierr); 42024989b8cSPeter Brune PetscFunctionReturn(0); 42124989b8cSPeter Brune } 422