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*/ 451b4a12fSMatthew G. Knepley #include <petscdraw.h> 524989b8cSPeter Brune 624989b8cSPeter Brune /* This structure holds the user-provided DMDA callbacks */ 724989b8cSPeter Brune typedef struct { 824989b8cSPeter Brune PetscErrorCode (*ifunctionlocal)(DMDALocalInfo*,PetscReal,void*,void*,void*,void*); 924989b8cSPeter Brune PetscErrorCode (*rhsfunctionlocal)(DMDALocalInfo*,PetscReal,void*,void*,void*); 1024989b8cSPeter Brune PetscErrorCode (*ijacobianlocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,MatStructure*,void*); 1124989b8cSPeter Brune PetscErrorCode (*rhsjacobianlocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,MatStructure*,void*); 1224989b8cSPeter Brune void *ifunctionlocalctx; 1324989b8cSPeter Brune void *ijacobianlocalctx; 1424989b8cSPeter Brune void *rhsfunctionlocalctx; 1524989b8cSPeter Brune void *rhsjacobianlocalctx; 1624989b8cSPeter Brune InsertMode ifunctionlocalimode; 1724989b8cSPeter Brune InsertMode rhsfunctionlocalimode; 18942e3340SBarry Smith } DMTS_DA; 1924989b8cSPeter Brune 2024989b8cSPeter Brune #undef __FUNCT__ 21942e3340SBarry Smith #define __FUNCT__ "DMTSDestroy_DMDA" 22942e3340SBarry Smith static PetscErrorCode DMTSDestroy_DMDA(DMTS sdm) 2324989b8cSPeter Brune { 2424989b8cSPeter Brune PetscErrorCode ierr; 2524989b8cSPeter Brune 2624989b8cSPeter Brune PetscFunctionBegin; 2724989b8cSPeter Brune ierr = PetscFree(sdm->data);CHKERRQ(ierr); 2824989b8cSPeter Brune PetscFunctionReturn(0); 2924989b8cSPeter Brune } 3024989b8cSPeter Brune 3124989b8cSPeter Brune #undef __FUNCT__ 32942e3340SBarry Smith #define __FUNCT__ "DMTSDuplicate_DMDA" 33d74926cbSBarry Smith static PetscErrorCode DMTSDuplicate_DMDA(DMTS oldsdm,DMTS sdm) 340313d0e0SJed Brown { 350313d0e0SJed Brown PetscErrorCode ierr; 360313d0e0SJed Brown 370313d0e0SJed Brown PetscFunctionBegin; 38*b00a9115SJed Brown ierr = PetscNewLog(sdm,(DMTS_DA**)&sdm->data);CHKERRQ(ierr); 39942e3340SBarry Smith if (oldsdm->data) {ierr = PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMTS_DA));CHKERRQ(ierr);} 400313d0e0SJed Brown PetscFunctionReturn(0); 410313d0e0SJed Brown } 420313d0e0SJed Brown 430313d0e0SJed Brown #undef __FUNCT__ 4424989b8cSPeter Brune #define __FUNCT__ "DMDATSGetContext" 45942e3340SBarry Smith static PetscErrorCode DMDATSGetContext(DM dm,DMTS sdm,DMTS_DA **dmdats) 4624989b8cSPeter Brune { 4724989b8cSPeter Brune PetscErrorCode ierr; 4824989b8cSPeter Brune 4924989b8cSPeter Brune PetscFunctionBegin; 500298fd71SBarry Smith *dmdats = NULL; 5124989b8cSPeter Brune if (!sdm->data) { 52*b00a9115SJed Brown ierr = PetscNewLog(dm,(DMTS_DA**)&sdm->data);CHKERRQ(ierr); 53d74926cbSBarry Smith sdm->ops->destroy = DMTSDestroy_DMDA; 54d74926cbSBarry Smith sdm->ops->duplicate = DMTSDuplicate_DMDA; 5524989b8cSPeter Brune } 56942e3340SBarry Smith *dmdats = (DMTS_DA*)sdm->data; 5724989b8cSPeter Brune PetscFunctionReturn(0); 5824989b8cSPeter Brune } 5924989b8cSPeter Brune 6024989b8cSPeter Brune #undef __FUNCT__ 6124989b8cSPeter Brune #define __FUNCT__ "TSComputeIFunction_DMDA" 6224989b8cSPeter Brune /* 6324989b8cSPeter Brune This function should eventually replace: 6424989b8cSPeter Brune DMDAComputeFunction() and DMDAComputeFunction1() 6524989b8cSPeter Brune */ 6624989b8cSPeter Brune static PetscErrorCode TSComputeIFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,Vec F,void *ctx) 6724989b8cSPeter Brune { 6824989b8cSPeter Brune PetscErrorCode ierr; 6924989b8cSPeter Brune DM dm; 70942e3340SBarry Smith DMTS_DA *dmdats = (DMTS_DA*)ctx; 7124989b8cSPeter Brune DMDALocalInfo info; 7224989b8cSPeter Brune Vec Xloc; 7324989b8cSPeter Brune void *x,*f,*xdot; 7424989b8cSPeter Brune 7524989b8cSPeter Brune PetscFunctionBegin; 7624989b8cSPeter Brune PetscValidHeaderSpecific(ts,TS_CLASSID,1); 7724989b8cSPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 7824989b8cSPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 79ce94432eSBarry Smith if (!dmdats->ifunctionlocal) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context"); 8024989b8cSPeter Brune ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 8124989b8cSPeter Brune ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 8224989b8cSPeter Brune ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 8324989b8cSPeter Brune ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 8424989b8cSPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 8524989b8cSPeter Brune ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 8624989b8cSPeter Brune ierr = DMDAVecGetArray(dm,Xdot,&xdot);CHKERRQ(ierr); 8724989b8cSPeter Brune switch (dmdats->ifunctionlocalimode) { 8824989b8cSPeter Brune case INSERT_VALUES: { 8924989b8cSPeter Brune ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr); 9024989b8cSPeter Brune CHKMEMQ; 9124989b8cSPeter Brune ierr = (*dmdats->ifunctionlocal)(&info,ptime,x,xdot,f,dmdats->ifunctionlocalctx);CHKERRQ(ierr); 9224989b8cSPeter Brune CHKMEMQ; 9324989b8cSPeter Brune ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr); 9424989b8cSPeter Brune } break; 9524989b8cSPeter Brune case ADD_VALUES: { 9624989b8cSPeter Brune Vec Floc; 9724989b8cSPeter Brune ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 9824989b8cSPeter Brune ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 9924989b8cSPeter Brune ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr); 10024989b8cSPeter Brune CHKMEMQ; 10124989b8cSPeter Brune ierr = (*dmdats->ifunctionlocal)(&info,ptime,x,xdot,f,dmdats->ifunctionlocalctx);CHKERRQ(ierr); 10224989b8cSPeter Brune CHKMEMQ; 10324989b8cSPeter Brune ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr); 10424989b8cSPeter Brune ierr = VecZeroEntries(F);CHKERRQ(ierr); 10524989b8cSPeter Brune ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 10624989b8cSPeter Brune ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 10724989b8cSPeter Brune ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 10824989b8cSPeter Brune } break; 109ce94432eSBarry Smith default: SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->ifunctionlocalimode); 11024989b8cSPeter Brune } 11124989b8cSPeter Brune ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 11224989b8cSPeter Brune ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 11311c2e1cfSJed Brown ierr = DMDAVecRestoreArray(dm,Xdot,&xdot);CHKERRQ(ierr); 11424989b8cSPeter Brune PetscFunctionReturn(0); 11524989b8cSPeter Brune } 11624989b8cSPeter Brune 11724989b8cSPeter Brune #undef __FUNCT__ 11824989b8cSPeter Brune #define __FUNCT__ "TSComputeIJacobian_DMDA" 11924989b8cSPeter Brune static PetscErrorCode TSComputeIJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,PetscReal shift,Mat *A,Mat *B,MatStructure *mstr,void *ctx) 12024989b8cSPeter Brune { 12124989b8cSPeter Brune PetscErrorCode ierr; 12224989b8cSPeter Brune DM dm; 123942e3340SBarry Smith DMTS_DA *dmdats = (DMTS_DA*)ctx; 12424989b8cSPeter Brune DMDALocalInfo info; 12524989b8cSPeter Brune Vec Xloc; 12624989b8cSPeter Brune void *x,*xdot; 12724989b8cSPeter Brune 12824989b8cSPeter Brune PetscFunctionBegin; 129ce94432eSBarry Smith if (!dmdats->ifunctionlocal) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context"); 13024989b8cSPeter Brune ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 13124989b8cSPeter Brune 13224989b8cSPeter Brune if (dmdats->ijacobianlocal) { 13324989b8cSPeter Brune ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 13424989b8cSPeter Brune ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 13524989b8cSPeter Brune ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 13624989b8cSPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 13724989b8cSPeter Brune ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 13824989b8cSPeter Brune ierr = DMDAVecGetArray(dm,Xdot,&xdot);CHKERRQ(ierr); 13924989b8cSPeter Brune CHKMEMQ; 14024989b8cSPeter Brune ierr = (*dmdats->ijacobianlocal)(&info,ptime,x,xdot,shift,*A,*B,mstr,dmdats->ijacobianlocalctx);CHKERRQ(ierr); 14124989b8cSPeter Brune CHKMEMQ; 14224989b8cSPeter Brune ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 14324989b8cSPeter Brune ierr = DMDAVecRestoreArray(dm,Xdot,&xdot);CHKERRQ(ierr); 14424989b8cSPeter Brune ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 145ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"TSComputeIJacobian_DMDA() called without calling DMDATSSetIJacobian()"); 14624989b8cSPeter Brune /* This will be redundant if the user called both, but it's too common to forget. */ 14724989b8cSPeter Brune if (*A != *B) { 14824989b8cSPeter Brune ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14924989b8cSPeter Brune ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 15024989b8cSPeter Brune } 15124989b8cSPeter Brune PetscFunctionReturn(0); 15224989b8cSPeter Brune } 15324989b8cSPeter Brune 15424989b8cSPeter Brune #undef __FUNCT__ 15524989b8cSPeter Brune #define __FUNCT__ "TSComputeRHSFunction_DMDA" 15624989b8cSPeter Brune /* 15724989b8cSPeter Brune This function should eventually replace: 15824989b8cSPeter Brune DMDAComputeFunction() and DMDAComputeFunction1() 15924989b8cSPeter Brune */ 16024989b8cSPeter Brune static PetscErrorCode TSComputeRHSFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec F,void *ctx) 16124989b8cSPeter Brune { 16224989b8cSPeter Brune PetscErrorCode ierr; 16324989b8cSPeter Brune DM dm; 164942e3340SBarry Smith DMTS_DA *dmdats = (DMTS_DA*)ctx; 16524989b8cSPeter Brune DMDALocalInfo info; 16624989b8cSPeter Brune Vec Xloc; 16724989b8cSPeter Brune void *x,*f; 16824989b8cSPeter Brune 16924989b8cSPeter Brune PetscFunctionBegin; 17024989b8cSPeter Brune PetscValidHeaderSpecific(ts,TS_CLASSID,1); 17124989b8cSPeter Brune PetscValidHeaderSpecific(X,VEC_CLASSID,2); 17224989b8cSPeter Brune PetscValidHeaderSpecific(F,VEC_CLASSID,3); 173ce94432eSBarry Smith if (!dmdats->rhsfunctionlocal) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context"); 17424989b8cSPeter Brune ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 17524989b8cSPeter Brune ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 17624989b8cSPeter Brune ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 17724989b8cSPeter Brune ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 17824989b8cSPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 17924989b8cSPeter Brune ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 18024989b8cSPeter Brune switch (dmdats->rhsfunctionlocalimode) { 18124989b8cSPeter Brune case INSERT_VALUES: { 18224989b8cSPeter Brune ierr = DMDAVecGetArray(dm,F,&f);CHKERRQ(ierr); 18324989b8cSPeter Brune CHKMEMQ; 18424989b8cSPeter Brune ierr = (*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx);CHKERRQ(ierr); 18524989b8cSPeter Brune CHKMEMQ; 18624989b8cSPeter Brune ierr = DMDAVecRestoreArray(dm,F,&f);CHKERRQ(ierr); 18724989b8cSPeter Brune } break; 18824989b8cSPeter Brune case ADD_VALUES: { 18924989b8cSPeter Brune Vec Floc; 19024989b8cSPeter Brune ierr = DMGetLocalVector(dm,&Floc);CHKERRQ(ierr); 19124989b8cSPeter Brune ierr = VecZeroEntries(Floc);CHKERRQ(ierr); 19224989b8cSPeter Brune ierr = DMDAVecGetArray(dm,Floc,&f);CHKERRQ(ierr); 19324989b8cSPeter Brune CHKMEMQ; 19424989b8cSPeter Brune ierr = (*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx);CHKERRQ(ierr); 19524989b8cSPeter Brune CHKMEMQ; 19624989b8cSPeter Brune ierr = DMDAVecRestoreArray(dm,Floc,&f);CHKERRQ(ierr); 19724989b8cSPeter Brune ierr = VecZeroEntries(F);CHKERRQ(ierr); 19824989b8cSPeter Brune ierr = DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 19924989b8cSPeter Brune ierr = DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);CHKERRQ(ierr); 20024989b8cSPeter Brune ierr = DMRestoreLocalVector(dm,&Floc);CHKERRQ(ierr); 20124989b8cSPeter Brune } break; 202ce94432eSBarry Smith default: SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->rhsfunctionlocalimode); 20324989b8cSPeter Brune } 20424989b8cSPeter Brune ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 20524989b8cSPeter Brune ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 20624989b8cSPeter Brune PetscFunctionReturn(0); 20724989b8cSPeter Brune } 20824989b8cSPeter Brune 20924989b8cSPeter Brune #undef __FUNCT__ 21024989b8cSPeter Brune #define __FUNCT__ "TSComputeRHSJacobian_DMDA" 21124989b8cSPeter Brune static PetscErrorCode TSComputeRHSJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Mat *A,Mat *B,MatStructure *mstr,void *ctx) 21224989b8cSPeter Brune { 21324989b8cSPeter Brune PetscErrorCode ierr; 21424989b8cSPeter Brune DM dm; 215942e3340SBarry Smith DMTS_DA *dmdats = (DMTS_DA*)ctx; 21624989b8cSPeter Brune DMDALocalInfo info; 21724989b8cSPeter Brune Vec Xloc; 21824989b8cSPeter Brune void *x; 21924989b8cSPeter Brune 22024989b8cSPeter Brune PetscFunctionBegin; 221ce94432eSBarry Smith if (!dmdats->rhsfunctionlocal) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context"); 22224989b8cSPeter Brune ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 22324989b8cSPeter Brune 22424989b8cSPeter Brune if (dmdats->rhsjacobianlocal) { 22524989b8cSPeter Brune ierr = DMGetLocalVector(dm,&Xloc);CHKERRQ(ierr); 22624989b8cSPeter Brune ierr = DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 22724989b8cSPeter Brune ierr = DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);CHKERRQ(ierr); 22824989b8cSPeter Brune ierr = DMDAGetLocalInfo(dm,&info);CHKERRQ(ierr); 22924989b8cSPeter Brune ierr = DMDAVecGetArray(dm,Xloc,&x);CHKERRQ(ierr); 23024989b8cSPeter Brune CHKMEMQ; 23124989b8cSPeter Brune ierr = (*dmdats->rhsjacobianlocal)(&info,ptime,x,*A,*B,mstr,dmdats->rhsjacobianlocalctx);CHKERRQ(ierr); 23224989b8cSPeter Brune CHKMEMQ; 23324989b8cSPeter Brune ierr = DMDAVecRestoreArray(dm,Xloc,&x);CHKERRQ(ierr); 23424989b8cSPeter Brune ierr = DMRestoreLocalVector(dm,&Xloc);CHKERRQ(ierr); 235ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"TSComputeRHSJacobian_DMDA() called without calling DMDATSSetRHSJacobian()"); 23624989b8cSPeter Brune /* This will be redundant if the user called both, but it's too common to forget. */ 23724989b8cSPeter Brune if (*A != *B) { 23824989b8cSPeter Brune ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 23924989b8cSPeter Brune ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 24024989b8cSPeter Brune } 24124989b8cSPeter Brune PetscFunctionReturn(0); 24224989b8cSPeter Brune } 24324989b8cSPeter Brune 24424989b8cSPeter Brune 24524989b8cSPeter Brune #undef __FUNCT__ 24624989b8cSPeter Brune #define __FUNCT__ "DMDATSSetRHSFunctionLocal" 24724989b8cSPeter Brune /*@C 24824989b8cSPeter Brune DMDATSSetRHSFunctionLocal - set a local residual evaluation function 24924989b8cSPeter Brune 25024989b8cSPeter Brune Logically Collective 25124989b8cSPeter Brune 25224989b8cSPeter Brune Input Arguments: 25324989b8cSPeter Brune + dm - DM to associate callback with 254b398ffe5SPeter Brune . imode - insert mode for the residual 25524989b8cSPeter Brune . func - local residual evaluation 25624989b8cSPeter Brune - ctx - optional context for local residual evaluation 25724989b8cSPeter Brune 25824989b8cSPeter Brune Calling sequence for func: 259b398ffe5SPeter Brune 260b398ffe5SPeter Brune $ func(DMDALocalInfo info,PetscReal t,void *x,void *f,void *ctx) 261b398ffe5SPeter Brune 26224989b8cSPeter Brune + info - DMDALocalInfo defining the subdomain to evaluate the residual on 263b398ffe5SPeter Brune . t - time at which to evaluate residual 264b398ffe5SPeter Brune . x - array of local state information 265b398ffe5SPeter Brune . f - output array of local residual information 266b398ffe5SPeter Brune - ctx - optional user context 26724989b8cSPeter Brune 26824989b8cSPeter Brune Level: beginner 26924989b8cSPeter Brune 2708e625dcfSPeter Brune .seealso: DMTSSetRHSFunction(), DMDATSSetRHSJacobianLocal(), DMDASNESSetFunctionLocal() 27124989b8cSPeter Brune @*/ 2726c6b9e74SPeter Brune PetscErrorCode DMDATSSetRHSFunctionLocal(DM dm,InsertMode imode,DMDATSRHSFunctionLocal func,void *ctx) 27324989b8cSPeter Brune { 27424989b8cSPeter Brune PetscErrorCode ierr; 275942e3340SBarry Smith DMTS sdm; 276942e3340SBarry Smith DMTS_DA *dmdats; 27724989b8cSPeter Brune 27824989b8cSPeter Brune PetscFunctionBegin; 27924989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 280942e3340SBarry Smith ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr); 28124989b8cSPeter Brune ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr); 28224989b8cSPeter Brune dmdats->rhsfunctionlocalimode = imode; 28324989b8cSPeter Brune dmdats->rhsfunctionlocal = func; 28424989b8cSPeter Brune dmdats->rhsfunctionlocalctx = ctx; 28524989b8cSPeter Brune ierr = DMTSSetRHSFunction(dm,TSComputeRHSFunction_DMDA,dmdats);CHKERRQ(ierr); 28624989b8cSPeter Brune PetscFunctionReturn(0); 28724989b8cSPeter Brune } 28824989b8cSPeter Brune 28924989b8cSPeter Brune #undef __FUNCT__ 29024989b8cSPeter Brune #define __FUNCT__ "DMDATSSetRHSJacobianLocal" 29124989b8cSPeter Brune /*@C 29224989b8cSPeter Brune DMDATSSetRHSJacobianLocal - set a local residual evaluation function 29324989b8cSPeter Brune 29424989b8cSPeter Brune Logically Collective 29524989b8cSPeter Brune 29624989b8cSPeter Brune Input Arguments: 29724989b8cSPeter Brune + dm - DM to associate callback with 298b398ffe5SPeter Brune . func - local RHS Jacobian evaluation routine 299b398ffe5SPeter Brune - ctx - optional context for local jacobian evaluation 30024989b8cSPeter Brune 30124989b8cSPeter Brune Calling sequence for func: 302b398ffe5SPeter Brune 303b398ffe5SPeter Brune $ func(DMDALocalInfo* info,PetscReal t,void* x,Mat J,Mat B,MatStructure *flg,void *ctx); 304b398ffe5SPeter Brune 30524989b8cSPeter Brune + info - DMDALocalInfo defining the subdomain to evaluate the residual on 306b398ffe5SPeter Brune . t - time at which to evaluate residual 307b398ffe5SPeter Brune . x - array of local state information 308b398ffe5SPeter Brune . J - Jacobian matrix 309b398ffe5SPeter Brune . B - preconditioner matrix; often same as J 310b398ffe5SPeter Brune . flg - flag indicating information about the preconditioner matrix structure (same as flag in KSPSetOperators()) 31124989b8cSPeter Brune - ctx - optional context passed above 31224989b8cSPeter Brune 31324989b8cSPeter Brune Level: beginner 31424989b8cSPeter Brune 3158e625dcfSPeter Brune .seealso: DMTSSetRHSJacobian(), DMDATSSetRHSFunctionLocal(), DMDASNESSetJacobianLocal() 31624989b8cSPeter Brune @*/ 3176c6b9e74SPeter Brune PetscErrorCode DMDATSSetRHSJacobianLocal(DM dm,DMDATSRHSJacobianLocal func,void *ctx) 31824989b8cSPeter Brune { 31924989b8cSPeter Brune PetscErrorCode ierr; 320942e3340SBarry Smith DMTS sdm; 321942e3340SBarry Smith DMTS_DA *dmdats; 32224989b8cSPeter Brune 32324989b8cSPeter Brune PetscFunctionBegin; 32424989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 325942e3340SBarry Smith ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr); 32624989b8cSPeter Brune ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr); 32724989b8cSPeter Brune dmdats->rhsjacobianlocal = func; 32824989b8cSPeter Brune dmdats->rhsjacobianlocalctx = ctx; 32924989b8cSPeter Brune ierr = DMTSSetRHSJacobian(dm,TSComputeRHSJacobian_DMDA,dmdats);CHKERRQ(ierr); 33024989b8cSPeter Brune PetscFunctionReturn(0); 33124989b8cSPeter Brune } 33224989b8cSPeter Brune 33324989b8cSPeter Brune 33424989b8cSPeter Brune #undef __FUNCT__ 33524989b8cSPeter Brune #define __FUNCT__ "DMDATSSetIFunctionLocal" 33624989b8cSPeter Brune /*@C 33724989b8cSPeter Brune DMDATSSetIFunctionLocal - set a local residual evaluation function 33824989b8cSPeter Brune 33924989b8cSPeter Brune Logically Collective 34024989b8cSPeter Brune 34124989b8cSPeter Brune Input Arguments: 34224989b8cSPeter Brune + dm - DM to associate callback with 34324989b8cSPeter Brune . func - local residual evaluation 34424989b8cSPeter Brune - ctx - optional context for local residual evaluation 34524989b8cSPeter Brune 34624989b8cSPeter Brune Calling sequence for func: 34724989b8cSPeter Brune + info - DMDALocalInfo defining the subdomain to evaluate the residual on 348b398ffe5SPeter Brune . t - time at which to evaluate residual 349b398ffe5SPeter Brune . x - array of local state information 350b398ffe5SPeter Brune . xdot - array of local time derivative information 351b398ffe5SPeter Brune . f - output array of local function evaluation information 35224989b8cSPeter Brune - ctx - optional context passed above 35324989b8cSPeter Brune 35424989b8cSPeter Brune Level: beginner 35524989b8cSPeter Brune 3568e625dcfSPeter Brune .seealso: DMTSSetIFunction(), DMDATSSetIJacobianLocal(), DMDASNESSetFunctionLocal() 35724989b8cSPeter Brune @*/ 3586c6b9e74SPeter Brune PetscErrorCode DMDATSSetIFunctionLocal(DM dm,InsertMode imode,DMDATSIFunctionLocal func,void *ctx) 35924989b8cSPeter Brune { 36024989b8cSPeter Brune PetscErrorCode ierr; 361942e3340SBarry Smith DMTS sdm; 362942e3340SBarry Smith DMTS_DA *dmdats; 36324989b8cSPeter Brune 36424989b8cSPeter Brune PetscFunctionBegin; 36524989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 366942e3340SBarry Smith ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr); 36724989b8cSPeter Brune ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr); 36824989b8cSPeter Brune dmdats->ifunctionlocalimode = imode; 36924989b8cSPeter Brune dmdats->ifunctionlocal = func; 37024989b8cSPeter Brune dmdats->ifunctionlocalctx = ctx; 37124989b8cSPeter Brune ierr = DMTSSetIFunction(dm,TSComputeIFunction_DMDA,dmdats);CHKERRQ(ierr); 37224989b8cSPeter Brune PetscFunctionReturn(0); 37324989b8cSPeter Brune } 37424989b8cSPeter Brune 37524989b8cSPeter Brune #undef __FUNCT__ 37624989b8cSPeter Brune #define __FUNCT__ "DMDATSSetIJacobianLocal" 37724989b8cSPeter Brune /*@C 37824989b8cSPeter Brune DMDATSSetIJacobianLocal - set a local residual evaluation function 37924989b8cSPeter Brune 38024989b8cSPeter Brune Logically Collective 38124989b8cSPeter Brune 38224989b8cSPeter Brune Input Arguments: 38324989b8cSPeter Brune + dm - DM to associate callback with 38424989b8cSPeter Brune . func - local residual evaluation 38524989b8cSPeter Brune - ctx - optional context for local residual evaluation 38624989b8cSPeter Brune 38724989b8cSPeter Brune Calling sequence for func: 388b398ffe5SPeter Brune 389b398ffe5SPeter Brune $ func(DMDALocalInfo* info,PetscReal t,void* x,void *xdot,Mat J,Mat B,MatStructure *flg,void *ctx); 390b398ffe5SPeter Brune 39124989b8cSPeter Brune + info - DMDALocalInfo defining the subdomain to evaluate the residual on 392b398ffe5SPeter Brune . t - time at which to evaluate the jacobian 393b398ffe5SPeter Brune . x - array of local state information 394b398ffe5SPeter Brune . xdot - time derivative at this state 395b398ffe5SPeter Brune . J - Jacobian matrix 396b398ffe5SPeter Brune . B - preconditioner matrix; often same as J 397b398ffe5SPeter Brune . flg - flag indicating information about the preconditioner matrix structure (same as flag in KSPSetOperators()) 39824989b8cSPeter Brune - ctx - optional context passed above 39924989b8cSPeter Brune 40024989b8cSPeter Brune Level: beginner 40124989b8cSPeter Brune 4028e625dcfSPeter Brune .seealso: DMTSSetJacobian(), DMDATSSetIFunctionLocal(), DMDASNESSetJacobianLocal() 40324989b8cSPeter Brune @*/ 4046c6b9e74SPeter Brune PetscErrorCode DMDATSSetIJacobianLocal(DM dm,DMDATSIJacobianLocal func,void *ctx) 40524989b8cSPeter Brune { 40624989b8cSPeter Brune PetscErrorCode ierr; 407942e3340SBarry Smith DMTS sdm; 408942e3340SBarry Smith DMTS_DA *dmdats; 40924989b8cSPeter Brune 41024989b8cSPeter Brune PetscFunctionBegin; 41124989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 412942e3340SBarry Smith ierr = DMGetDMTSWrite(dm,&sdm);CHKERRQ(ierr); 41324989b8cSPeter Brune ierr = DMDATSGetContext(dm,sdm,&dmdats);CHKERRQ(ierr); 41424989b8cSPeter Brune dmdats->ijacobianlocal = func; 41524989b8cSPeter Brune dmdats->ijacobianlocalctx = ctx; 41624989b8cSPeter Brune ierr = DMTSSetIJacobian(dm,TSComputeIJacobian_DMDA,dmdats);CHKERRQ(ierr); 41724989b8cSPeter Brune PetscFunctionReturn(0); 41824989b8cSPeter Brune } 419d1212d36SBarry Smith 420d1212d36SBarry Smith #undef __FUNCT__ 421d1212d36SBarry Smith #define __FUNCT__ "TSMonitorDMDARayDestroy" 422d1212d36SBarry Smith PetscErrorCode TSMonitorDMDARayDestroy(void **mctx) 423d1212d36SBarry Smith { 424d1212d36SBarry Smith TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *) *mctx; 425d1212d36SBarry Smith PetscErrorCode ierr; 426d1212d36SBarry Smith 427d1212d36SBarry Smith PetscFunctionBegin; 42851b4a12fSMatthew G. Knepley if (rayctx->lgctx) {ierr = TSMonitorLGCtxDestroy(&rayctx->lgctx);CHKERRQ(ierr);} 429d1212d36SBarry Smith ierr = VecDestroy(&rayctx->ray);CHKERRQ(ierr); 430d1212d36SBarry Smith ierr = VecScatterDestroy(&rayctx->scatter);CHKERRQ(ierr); 431d1212d36SBarry Smith ierr = PetscViewerDestroy(&rayctx->viewer);CHKERRQ(ierr); 432d1212d36SBarry Smith ierr = PetscFree(rayctx);CHKERRQ(ierr); 433d1212d36SBarry Smith PetscFunctionReturn(0); 434d1212d36SBarry Smith } 435d1212d36SBarry Smith 436d1212d36SBarry Smith #undef __FUNCT__ 437d1212d36SBarry Smith #define __FUNCT__ "TSMonitorDMDARay" 438d1212d36SBarry Smith PetscErrorCode TSMonitorDMDARay(TS ts,PetscInt steps,PetscReal time,Vec u,void *mctx) 439d1212d36SBarry Smith { 440d1212d36SBarry Smith TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx*)mctx; 441d1212d36SBarry Smith Vec solution; 442d1212d36SBarry Smith PetscErrorCode ierr; 443d1212d36SBarry Smith 444d1212d36SBarry Smith PetscFunctionBegin; 445d1212d36SBarry Smith ierr = TSGetSolution(ts,&solution);CHKERRQ(ierr); 446d1212d36SBarry Smith ierr = VecScatterBegin(rayctx->scatter,solution,rayctx->ray,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 447d1212d36SBarry Smith ierr = VecScatterEnd(rayctx->scatter,solution,rayctx->ray,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 448d1212d36SBarry Smith if (rayctx->viewer) { 449d1212d36SBarry Smith ierr = VecView(rayctx->ray,rayctx->viewer);CHKERRQ(ierr); 450d1212d36SBarry Smith } 451d1212d36SBarry Smith PetscFunctionReturn(0); 452d1212d36SBarry Smith } 45351b4a12fSMatthew G. Knepley 45451b4a12fSMatthew G. Knepley #undef __FUNCT__ 45551b4a12fSMatthew G. Knepley #define __FUNCT__ "TSMonitorLGDMDARay" 45651b4a12fSMatthew G. Knepley PetscErrorCode TSMonitorLGDMDARay(TS ts, PetscInt step, PetscReal ptime, Vec u, void *ctx) 45751b4a12fSMatthew G. Knepley { 45851b4a12fSMatthew G. Knepley TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *) ctx; 45951b4a12fSMatthew G. Knepley TSMonitorLGCtx lgctx = (TSMonitorLGCtx) rayctx->lgctx; 46051b4a12fSMatthew G. Knepley Vec v = rayctx->ray; 46151b4a12fSMatthew G. Knepley const PetscScalar *a; 46251b4a12fSMatthew G. Knepley PetscInt dim; 46351b4a12fSMatthew G. Knepley PetscErrorCode ierr; 46451b4a12fSMatthew G. Knepley 46551b4a12fSMatthew G. Knepley PetscFunctionBegin; 46651b4a12fSMatthew G. Knepley ierr = VecScatterBegin(rayctx->scatter, u, v, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 46751b4a12fSMatthew G. Knepley ierr = VecScatterEnd(rayctx->scatter, u, v, INSERT_VALUES, SCATTER_FORWARD);CHKERRQ(ierr); 46851b4a12fSMatthew G. Knepley if (!step) { 46951b4a12fSMatthew G. Knepley PetscDrawAxis axis; 47051b4a12fSMatthew G. Knepley 47151b4a12fSMatthew G. Knepley ierr = PetscDrawLGGetAxis(lgctx->lg, &axis);CHKERRQ(ierr); 47251b4a12fSMatthew G. Knepley ierr = PetscDrawAxisSetLabels(axis, "Solution Ray as function of time", "Time", "Solution");CHKERRQ(ierr); 47351b4a12fSMatthew G. Knepley ierr = VecGetLocalSize(rayctx->ray, &dim);CHKERRQ(ierr); 47451b4a12fSMatthew G. Knepley ierr = PetscDrawLGSetDimension(lgctx->lg, dim);CHKERRQ(ierr); 47551b4a12fSMatthew G. Knepley ierr = PetscDrawLGReset(lgctx->lg);CHKERRQ(ierr); 47651b4a12fSMatthew G. Knepley } 47751b4a12fSMatthew G. Knepley ierr = VecGetArrayRead(v, &a);CHKERRQ(ierr); 47851b4a12fSMatthew G. Knepley #if defined(PETSC_USE_COMPLEX) 47951b4a12fSMatthew G. Knepley { 48051b4a12fSMatthew G. Knepley PetscReal *areal; 48151b4a12fSMatthew G. Knepley PetscInt i,n; 48251b4a12fSMatthew G. Knepley ierr = VecGetLocalSize(v, &n);CHKERRQ(ierr); 483785e854fSJed Brown ierr = PetscMalloc1(n, &areal);CHKERRQ(ierr); 48451b4a12fSMatthew G. Knepley for (i = 0; i < n; ++i) areal[i] = PetscRealPart(a[i]); 48551b4a12fSMatthew G. Knepley ierr = PetscDrawLGAddCommonPoint(lgctx->lg, ptime, areal);CHKERRQ(ierr); 48651b4a12fSMatthew G. Knepley ierr = PetscFree(areal);CHKERRQ(ierr); 48751b4a12fSMatthew G. Knepley } 48851b4a12fSMatthew G. Knepley #else 48951b4a12fSMatthew G. Knepley ierr = PetscDrawLGAddCommonPoint(lgctx->lg, ptime, a);CHKERRQ(ierr); 49051b4a12fSMatthew G. Knepley #endif 49151b4a12fSMatthew G. Knepley ierr = VecRestoreArrayRead(v, &a);CHKERRQ(ierr); 49251b4a12fSMatthew G. Knepley if (((lgctx->howoften > 0) && (!(step % lgctx->howoften))) || ((lgctx->howoften == -1) && ts->reason)) { 49351b4a12fSMatthew G. Knepley ierr = PetscDrawLGDraw(lgctx->lg);CHKERRQ(ierr); 49451b4a12fSMatthew G. Knepley } 49551b4a12fSMatthew G. Knepley PetscFunctionReturn(0); 49651b4a12fSMatthew G. Knepley } 497