124989b8cSPeter Brune #include <petscdmda.h> /*I "petscdmda.h" I*/ 2af0996ceSBarry Smith #include <petsc/private/dmimpl.h> 3af0996ceSBarry Smith #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*); 10d1e9a80fSBarry Smith PetscErrorCode (*ijacobianlocal)(DMDALocalInfo*,PetscReal,void*,void*,PetscReal,Mat,Mat,void*); 11d1e9a80fSBarry Smith PetscErrorCode (*rhsjacobianlocal)(DMDALocalInfo*,PetscReal,void*,Mat,Mat,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 20942e3340SBarry Smith static PetscErrorCode DMTSDestroy_DMDA(DMTS sdm) 2124989b8cSPeter Brune { 2224989b8cSPeter Brune PetscFunctionBegin; 239566063dSJacob Faibussowitsch PetscCall(PetscFree(sdm->data)); 2424989b8cSPeter Brune PetscFunctionReturn(0); 2524989b8cSPeter Brune } 2624989b8cSPeter Brune 27d74926cbSBarry Smith static PetscErrorCode DMTSDuplicate_DMDA(DMTS oldsdm,DMTS sdm) 280313d0e0SJed Brown { 290313d0e0SJed Brown PetscFunctionBegin; 309566063dSJacob Faibussowitsch PetscCall(PetscNewLog(sdm,(DMTS_DA**)&sdm->data)); 319566063dSJacob Faibussowitsch if (oldsdm->data) PetscCall(PetscMemcpy(sdm->data,oldsdm->data,sizeof(DMTS_DA))); 320313d0e0SJed Brown PetscFunctionReturn(0); 330313d0e0SJed Brown } 340313d0e0SJed Brown 35942e3340SBarry Smith static PetscErrorCode DMDATSGetContext(DM dm,DMTS sdm,DMTS_DA **dmdats) 3624989b8cSPeter Brune { 3724989b8cSPeter Brune PetscFunctionBegin; 380298fd71SBarry Smith *dmdats = NULL; 3924989b8cSPeter Brune if (!sdm->data) { 409566063dSJacob Faibussowitsch PetscCall(PetscNewLog(dm,(DMTS_DA**)&sdm->data)); 41d74926cbSBarry Smith sdm->ops->destroy = DMTSDestroy_DMDA; 42d74926cbSBarry Smith sdm->ops->duplicate = DMTSDuplicate_DMDA; 4324989b8cSPeter Brune } 44942e3340SBarry Smith *dmdats = (DMTS_DA*)sdm->data; 4524989b8cSPeter Brune PetscFunctionReturn(0); 4624989b8cSPeter Brune } 4724989b8cSPeter Brune 4824989b8cSPeter Brune static PetscErrorCode TSComputeIFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,Vec F,void *ctx) 4924989b8cSPeter Brune { 5024989b8cSPeter Brune DM dm; 51942e3340SBarry Smith DMTS_DA *dmdats = (DMTS_DA*)ctx; 5224989b8cSPeter Brune DMDALocalInfo info; 53d2fc88d6Ssepideh Vec Xloc,Xdotloc; 5424989b8cSPeter Brune void *x,*f,*xdot; 5524989b8cSPeter Brune 5624989b8cSPeter Brune PetscFunctionBegin; 5724989b8cSPeter Brune PetscValidHeaderSpecific(ts,TS_CLASSID,1); 58064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(X,VEC_CLASSID,3); 59064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(F,VEC_CLASSID,5); 603c633725SBarry Smith PetscCheck(dmdats->ifunctionlocal,PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context"); 619566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 629566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Xdotloc)); 639566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm,Xdot,INSERT_VALUES,Xdotloc)); 649566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm,Xdot,INSERT_VALUES,Xdotloc)); 659566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Xloc)); 669566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc)); 679566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc)); 689566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm,&info)); 699566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm,Xloc,&x)); 709566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm,Xdotloc,&xdot)); 7124989b8cSPeter Brune switch (dmdats->ifunctionlocalimode) { 7224989b8cSPeter Brune case INSERT_VALUES: { 739566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm,F,&f)); 7424989b8cSPeter Brune CHKMEMQ; 759566063dSJacob Faibussowitsch PetscCall((*dmdats->ifunctionlocal)(&info,ptime,x,xdot,f,dmdats->ifunctionlocalctx)); 7624989b8cSPeter Brune CHKMEMQ; 779566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm,F,&f)); 7824989b8cSPeter Brune } break; 7924989b8cSPeter Brune case ADD_VALUES: { 8024989b8cSPeter Brune Vec Floc; 819566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Floc)); 829566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Floc)); 839566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm,Floc,&f)); 8424989b8cSPeter Brune CHKMEMQ; 859566063dSJacob Faibussowitsch PetscCall((*dmdats->ifunctionlocal)(&info,ptime,x,xdot,f,dmdats->ifunctionlocalctx)); 8624989b8cSPeter Brune CHKMEMQ; 879566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm,Floc,&f)); 889566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 899566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F)); 909566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F)); 919566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&Floc)); 9224989b8cSPeter Brune } break; 9398921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->ifunctionlocalimode); 9424989b8cSPeter Brune } 959566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm,Xloc,&x)); 969566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&Xloc)); 979566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm,Xdotloc,&xdot)); 989566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&Xdotloc)); 9924989b8cSPeter Brune PetscFunctionReturn(0); 10024989b8cSPeter Brune } 10124989b8cSPeter Brune 102d1e9a80fSBarry Smith static PetscErrorCode TSComputeIJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Vec Xdot,PetscReal shift,Mat A,Mat B,void *ctx) 10324989b8cSPeter Brune { 10424989b8cSPeter Brune DM dm; 105942e3340SBarry Smith DMTS_DA *dmdats = (DMTS_DA*)ctx; 10624989b8cSPeter Brune DMDALocalInfo info; 10724989b8cSPeter Brune Vec Xloc; 10824989b8cSPeter Brune void *x,*xdot; 10924989b8cSPeter Brune 11024989b8cSPeter Brune PetscFunctionBegin; 1113c633725SBarry Smith PetscCheck(dmdats->ifunctionlocal,PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context"); 1129566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 11324989b8cSPeter Brune 11424989b8cSPeter Brune if (dmdats->ijacobianlocal) { 1159566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Xloc)); 1169566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc)); 1179566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc)); 1189566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm,&info)); 1199566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm,Xloc,&x)); 1209566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm,Xdot,&xdot)); 12124989b8cSPeter Brune CHKMEMQ; 1229566063dSJacob Faibussowitsch PetscCall((*dmdats->ijacobianlocal)(&info,ptime,x,xdot,shift,A,B,dmdats->ijacobianlocalctx)); 12324989b8cSPeter Brune CHKMEMQ; 1249566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm,Xloc,&x)); 1259566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm,Xdot,&xdot)); 1269566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&Xloc)); 127ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"TSComputeIJacobian_DMDA() called without calling DMDATSSetIJacobian()"); 12824989b8cSPeter Brune /* This will be redundant if the user called both, but it's too common to forget. */ 12994ab13aaSBarry Smith if (A != B) { 1309566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1319566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 13224989b8cSPeter Brune } 13324989b8cSPeter Brune PetscFunctionReturn(0); 13424989b8cSPeter Brune } 13524989b8cSPeter Brune 13624989b8cSPeter Brune static PetscErrorCode TSComputeRHSFunction_DMDA(TS ts,PetscReal ptime,Vec X,Vec F,void *ctx) 13724989b8cSPeter Brune { 13824989b8cSPeter Brune DM dm; 139942e3340SBarry Smith DMTS_DA *dmdats = (DMTS_DA*)ctx; 14024989b8cSPeter Brune DMDALocalInfo info; 14124989b8cSPeter Brune Vec Xloc; 14224989b8cSPeter Brune void *x,*f; 14324989b8cSPeter Brune 14424989b8cSPeter Brune PetscFunctionBegin; 14524989b8cSPeter Brune PetscValidHeaderSpecific(ts,TS_CLASSID,1); 146064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(X,VEC_CLASSID,3); 147064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(F,VEC_CLASSID,4); 1483c633725SBarry Smith PetscCheck(dmdats->rhsfunctionlocal,PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context"); 1499566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 1509566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Xloc)); 1519566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc)); 1529566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc)); 1539566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm,&info)); 1549566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm,Xloc,&x)); 15524989b8cSPeter Brune switch (dmdats->rhsfunctionlocalimode) { 15624989b8cSPeter Brune case INSERT_VALUES: { 1579566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm,F,&f)); 15824989b8cSPeter Brune CHKMEMQ; 1599566063dSJacob Faibussowitsch PetscCall((*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx)); 16024989b8cSPeter Brune CHKMEMQ; 1619566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm,F,&f)); 16224989b8cSPeter Brune } break; 16324989b8cSPeter Brune case ADD_VALUES: { 16424989b8cSPeter Brune Vec Floc; 1659566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Floc)); 1669566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Floc)); 1679566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm,Floc,&f)); 16824989b8cSPeter Brune CHKMEMQ; 1699566063dSJacob Faibussowitsch PetscCall((*dmdats->rhsfunctionlocal)(&info,ptime,x,f,dmdats->rhsfunctionlocalctx)); 17024989b8cSPeter Brune CHKMEMQ; 1719566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm,Floc,&f)); 1729566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 1739566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F)); 1749566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F)); 1759566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&Floc)); 17624989b8cSPeter Brune } break; 17798921bdaSJacob Faibussowitsch default: SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_INCOMP,"Cannot use imode=%d",(int)dmdats->rhsfunctionlocalimode); 17824989b8cSPeter Brune } 1799566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm,Xloc,&x)); 1809566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&Xloc)); 18124989b8cSPeter Brune PetscFunctionReturn(0); 18224989b8cSPeter Brune } 18324989b8cSPeter Brune 184d1e9a80fSBarry Smith static PetscErrorCode TSComputeRHSJacobian_DMDA(TS ts,PetscReal ptime,Vec X,Mat A,Mat B,void *ctx) 18524989b8cSPeter Brune { 18624989b8cSPeter Brune DM dm; 187942e3340SBarry Smith DMTS_DA *dmdats = (DMTS_DA*)ctx; 18824989b8cSPeter Brune DMDALocalInfo info; 18924989b8cSPeter Brune Vec Xloc; 19024989b8cSPeter Brune void *x; 19124989b8cSPeter Brune 19224989b8cSPeter Brune PetscFunctionBegin; 1933c633725SBarry Smith PetscCheck(dmdats->rhsfunctionlocal,PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"Corrupt context"); 1949566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 19524989b8cSPeter Brune 19624989b8cSPeter Brune if (dmdats->rhsjacobianlocal) { 1979566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm,&Xloc)); 1989566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc)); 1999566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc)); 2009566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm,&info)); 2019566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm,Xloc,&x)); 20224989b8cSPeter Brune CHKMEMQ; 2039566063dSJacob Faibussowitsch PetscCall((*dmdats->rhsjacobianlocal)(&info,ptime,x,A,B,dmdats->rhsjacobianlocalctx)); 20424989b8cSPeter Brune CHKMEMQ; 2059566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm,Xloc,&x)); 2069566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm,&Xloc)); 207ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_PLIB,"TSComputeRHSJacobian_DMDA() called without calling DMDATSSetRHSJacobian()"); 20824989b8cSPeter Brune /* This will be redundant if the user called both, but it's too common to forget. */ 20994ab13aaSBarry Smith if (A != B) { 2109566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 2119566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 21224989b8cSPeter Brune } 21324989b8cSPeter Brune PetscFunctionReturn(0); 21424989b8cSPeter Brune } 21524989b8cSPeter Brune 21624989b8cSPeter Brune /*@C 21724989b8cSPeter Brune DMDATSSetRHSFunctionLocal - set a local residual evaluation function 21824989b8cSPeter Brune 21924989b8cSPeter Brune Logically Collective 22024989b8cSPeter Brune 2214165533cSJose E. Roman Input Parameters: 22224989b8cSPeter Brune + dm - DM to associate callback with 223b398ffe5SPeter Brune . imode - insert mode for the residual 22424989b8cSPeter Brune . func - local residual evaluation 22524989b8cSPeter Brune - ctx - optional context for local residual evaluation 22624989b8cSPeter Brune 22724989b8cSPeter Brune Calling sequence for func: 228b398ffe5SPeter Brune 229b398ffe5SPeter Brune $ func(DMDALocalInfo info,PetscReal t,void *x,void *f,void *ctx) 230b398ffe5SPeter Brune 23124989b8cSPeter Brune + info - DMDALocalInfo defining the subdomain to evaluate the residual on 232b398ffe5SPeter Brune . t - time at which to evaluate residual 233b398ffe5SPeter Brune . x - array of local state information 234b398ffe5SPeter Brune . f - output array of local residual information 235b398ffe5SPeter Brune - ctx - optional user context 23624989b8cSPeter Brune 23724989b8cSPeter Brune Level: beginner 23824989b8cSPeter Brune 239db781477SPatrick Sanan .seealso: `DMTSSetRHSFunction()`, `DMDATSSetRHSJacobianLocal()`, `DMDASNESSetFunctionLocal()` 24024989b8cSPeter Brune @*/ 2416c6b9e74SPeter Brune PetscErrorCode DMDATSSetRHSFunctionLocal(DM dm,InsertMode imode,DMDATSRHSFunctionLocal func,void *ctx) 24224989b8cSPeter Brune { 243942e3340SBarry Smith DMTS sdm; 244942e3340SBarry Smith DMTS_DA *dmdats; 24524989b8cSPeter Brune 24624989b8cSPeter Brune PetscFunctionBegin; 24724989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2489566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm,&sdm)); 2499566063dSJacob Faibussowitsch PetscCall(DMDATSGetContext(dm,sdm,&dmdats)); 25024989b8cSPeter Brune dmdats->rhsfunctionlocalimode = imode; 25124989b8cSPeter Brune dmdats->rhsfunctionlocal = func; 25224989b8cSPeter Brune dmdats->rhsfunctionlocalctx = ctx; 2539566063dSJacob Faibussowitsch PetscCall(DMTSSetRHSFunction(dm,TSComputeRHSFunction_DMDA,dmdats)); 25424989b8cSPeter Brune PetscFunctionReturn(0); 25524989b8cSPeter Brune } 25624989b8cSPeter Brune 25724989b8cSPeter Brune /*@C 25824989b8cSPeter Brune DMDATSSetRHSJacobianLocal - set a local residual evaluation function 25924989b8cSPeter Brune 26024989b8cSPeter Brune Logically Collective 26124989b8cSPeter Brune 2624165533cSJose E. Roman Input Parameters: 26324989b8cSPeter Brune + dm - DM to associate callback with 264b398ffe5SPeter Brune . func - local RHS Jacobian evaluation routine 265b398ffe5SPeter Brune - ctx - optional context for local jacobian evaluation 26624989b8cSPeter Brune 26724989b8cSPeter Brune Calling sequence for func: 268b398ffe5SPeter Brune 2697f6be190SBarry Smith $ func(DMDALocalInfo* info,PetscReal t,void* x,Mat J,Mat B,void *ctx); 270b398ffe5SPeter Brune 27124989b8cSPeter Brune + info - DMDALocalInfo defining the subdomain to evaluate the residual on 272b398ffe5SPeter Brune . t - time at which to evaluate residual 273b398ffe5SPeter Brune . x - array of local state information 274b398ffe5SPeter Brune . J - Jacobian matrix 275b398ffe5SPeter Brune . B - preconditioner matrix; often same as J 27624989b8cSPeter Brune - ctx - optional context passed above 27724989b8cSPeter Brune 27824989b8cSPeter Brune Level: beginner 27924989b8cSPeter Brune 280db781477SPatrick Sanan .seealso: `DMTSSetRHSJacobian()`, `DMDATSSetRHSFunctionLocal()`, `DMDASNESSetJacobianLocal()` 28124989b8cSPeter Brune @*/ 2826c6b9e74SPeter Brune PetscErrorCode DMDATSSetRHSJacobianLocal(DM dm,DMDATSRHSJacobianLocal func,void *ctx) 28324989b8cSPeter Brune { 284942e3340SBarry Smith DMTS sdm; 285942e3340SBarry Smith DMTS_DA *dmdats; 28624989b8cSPeter Brune 28724989b8cSPeter Brune PetscFunctionBegin; 28824989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 2899566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm,&sdm)); 2909566063dSJacob Faibussowitsch PetscCall(DMDATSGetContext(dm,sdm,&dmdats)); 29124989b8cSPeter Brune dmdats->rhsjacobianlocal = func; 29224989b8cSPeter Brune dmdats->rhsjacobianlocalctx = ctx; 2939566063dSJacob Faibussowitsch PetscCall(DMTSSetRHSJacobian(dm,TSComputeRHSJacobian_DMDA,dmdats)); 29424989b8cSPeter Brune PetscFunctionReturn(0); 29524989b8cSPeter Brune } 29624989b8cSPeter Brune 29724989b8cSPeter Brune /*@C 29824989b8cSPeter Brune DMDATSSetIFunctionLocal - set a local residual evaluation function 29924989b8cSPeter Brune 30024989b8cSPeter Brune Logically Collective 30124989b8cSPeter Brune 3024165533cSJose E. Roman Input Parameters: 30324989b8cSPeter Brune + dm - DM to associate callback with 30424989b8cSPeter Brune . func - local residual evaluation 30524989b8cSPeter Brune - ctx - optional context for local residual evaluation 30624989b8cSPeter Brune 30724989b8cSPeter Brune Calling sequence for func: 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 . xdot - array of local time derivative information 312b398ffe5SPeter Brune . f - output array of local function evaluation information 31324989b8cSPeter Brune - ctx - optional context passed above 31424989b8cSPeter Brune 31524989b8cSPeter Brune Level: beginner 31624989b8cSPeter Brune 317db781477SPatrick Sanan .seealso: `DMTSSetIFunction()`, `DMDATSSetIJacobianLocal()`, `DMDASNESSetFunctionLocal()` 31824989b8cSPeter Brune @*/ 3196c6b9e74SPeter Brune PetscErrorCode DMDATSSetIFunctionLocal(DM dm,InsertMode imode,DMDATSIFunctionLocal func,void *ctx) 32024989b8cSPeter Brune { 321942e3340SBarry Smith DMTS sdm; 322942e3340SBarry Smith DMTS_DA *dmdats; 32324989b8cSPeter Brune 32424989b8cSPeter Brune PetscFunctionBegin; 32524989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3269566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm,&sdm)); 3279566063dSJacob Faibussowitsch PetscCall(DMDATSGetContext(dm,sdm,&dmdats)); 32824989b8cSPeter Brune dmdats->ifunctionlocalimode = imode; 32924989b8cSPeter Brune dmdats->ifunctionlocal = func; 33024989b8cSPeter Brune dmdats->ifunctionlocalctx = ctx; 3319566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunction(dm,TSComputeIFunction_DMDA,dmdats)); 33224989b8cSPeter Brune PetscFunctionReturn(0); 33324989b8cSPeter Brune } 33424989b8cSPeter Brune 33524989b8cSPeter Brune /*@C 33624989b8cSPeter Brune DMDATSSetIJacobianLocal - set a local residual evaluation function 33724989b8cSPeter Brune 33824989b8cSPeter Brune Logically Collective 33924989b8cSPeter Brune 3404165533cSJose E. Roman Input Parameters: 34124989b8cSPeter Brune + dm - DM to associate callback with 34224989b8cSPeter Brune . func - local residual evaluation 34324989b8cSPeter Brune - ctx - optional context for local residual evaluation 34424989b8cSPeter Brune 34524989b8cSPeter Brune Calling sequence for func: 346b398ffe5SPeter Brune 347dab4423fSBarry Smith $ func(DMDALocalInfo* info,PetscReal t,void* x,void *xdot,PetscScalar shift,Mat J,Mat B,void *ctx); 348b398ffe5SPeter Brune 34924989b8cSPeter Brune + info - DMDALocalInfo defining the subdomain to evaluate the residual on 350b398ffe5SPeter Brune . t - time at which to evaluate the jacobian 351b398ffe5SPeter Brune . x - array of local state information 352b398ffe5SPeter Brune . xdot - time derivative at this state 353dab4423fSBarry Smith . shift - see TSSetIJacobian() for the meaning of this parameter 354b398ffe5SPeter Brune . J - Jacobian matrix 355b398ffe5SPeter Brune . B - preconditioner matrix; often same as J 35624989b8cSPeter Brune - ctx - optional context passed above 35724989b8cSPeter Brune 35824989b8cSPeter Brune Level: beginner 35924989b8cSPeter Brune 360db781477SPatrick Sanan .seealso: `DMTSSetJacobian()`, `DMDATSSetIFunctionLocal()`, `DMDASNESSetJacobianLocal()` 36124989b8cSPeter Brune @*/ 3626c6b9e74SPeter Brune PetscErrorCode DMDATSSetIJacobianLocal(DM dm,DMDATSIJacobianLocal func,void *ctx) 36324989b8cSPeter Brune { 364942e3340SBarry Smith DMTS sdm; 365942e3340SBarry Smith DMTS_DA *dmdats; 36624989b8cSPeter Brune 36724989b8cSPeter Brune PetscFunctionBegin; 36824989b8cSPeter Brune PetscValidHeaderSpecific(dm,DM_CLASSID,1); 3699566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm,&sdm)); 3709566063dSJacob Faibussowitsch PetscCall(DMDATSGetContext(dm,sdm,&dmdats)); 37124989b8cSPeter Brune dmdats->ijacobianlocal = func; 37224989b8cSPeter Brune dmdats->ijacobianlocalctx = ctx; 3739566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobian(dm,TSComputeIJacobian_DMDA,dmdats)); 37424989b8cSPeter Brune PetscFunctionReturn(0); 37524989b8cSPeter Brune } 376d1212d36SBarry Smith 377d1212d36SBarry Smith PetscErrorCode TSMonitorDMDARayDestroy(void **mctx) 378d1212d36SBarry Smith { 379d1212d36SBarry Smith TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *) *mctx; 380d1212d36SBarry Smith 381d1212d36SBarry Smith PetscFunctionBegin; 3829566063dSJacob Faibussowitsch if (rayctx->lgctx) PetscCall(TSMonitorLGCtxDestroy(&rayctx->lgctx)); 3839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rayctx->ray)); 3849566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&rayctx->scatter)); 3859566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&rayctx->viewer)); 3869566063dSJacob Faibussowitsch PetscCall(PetscFree(rayctx)); 387d1212d36SBarry Smith PetscFunctionReturn(0); 388d1212d36SBarry Smith } 389d1212d36SBarry Smith 390d1212d36SBarry Smith PetscErrorCode TSMonitorDMDARay(TS ts,PetscInt steps,PetscReal time,Vec u,void *mctx) 391d1212d36SBarry Smith { 392d1212d36SBarry Smith TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx*)mctx; 393d1212d36SBarry Smith Vec solution; 394d1212d36SBarry Smith 395d1212d36SBarry Smith PetscFunctionBegin; 3969566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts,&solution)); 3979566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(rayctx->scatter,solution,rayctx->ray,INSERT_VALUES,SCATTER_FORWARD)); 3989566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(rayctx->scatter,solution,rayctx->ray,INSERT_VALUES,SCATTER_FORWARD)); 399*1baa6e33SBarry Smith if (rayctx->viewer) PetscCall(VecView(rayctx->ray,rayctx->viewer)); 400d1212d36SBarry Smith PetscFunctionReturn(0); 401d1212d36SBarry Smith } 40251b4a12fSMatthew G. Knepley 40351b4a12fSMatthew G. Knepley PetscErrorCode TSMonitorLGDMDARay(TS ts, PetscInt step, PetscReal ptime, Vec u, void *ctx) 40451b4a12fSMatthew G. Knepley { 40551b4a12fSMatthew G. Knepley TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *) ctx; 40651b4a12fSMatthew G. Knepley TSMonitorLGCtx lgctx = (TSMonitorLGCtx) rayctx->lgctx; 40751b4a12fSMatthew G. Knepley Vec v = rayctx->ray; 40851b4a12fSMatthew G. Knepley const PetscScalar *a; 40951b4a12fSMatthew G. Knepley PetscInt dim; 41051b4a12fSMatthew G. Knepley 41151b4a12fSMatthew G. Knepley PetscFunctionBegin; 4129566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(rayctx->scatter, u, v, INSERT_VALUES, SCATTER_FORWARD)); 4139566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(rayctx->scatter, u, v, INSERT_VALUES, SCATTER_FORWARD)); 41451b4a12fSMatthew G. Knepley if (!step) { 41551b4a12fSMatthew G. Knepley PetscDrawAxis axis; 41651b4a12fSMatthew G. Knepley 4179566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(lgctx->lg, &axis)); 4189566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Solution Ray as function of time", "Time", "Solution")); 4199566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(rayctx->ray, &dim)); 4209566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(lgctx->lg, dim)); 4219566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(lgctx->lg)); 42251b4a12fSMatthew G. Knepley } 4239566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(v, &a)); 42451b4a12fSMatthew G. Knepley #if defined(PETSC_USE_COMPLEX) 42551b4a12fSMatthew G. Knepley { 42651b4a12fSMatthew G. Knepley PetscReal *areal; 42751b4a12fSMatthew G. Knepley PetscInt i,n; 4289566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(v, &n)); 4299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &areal)); 43051b4a12fSMatthew G. Knepley for (i = 0; i < n; ++i) areal[i] = PetscRealPart(a[i]); 4319566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(lgctx->lg, ptime, areal)); 4329566063dSJacob Faibussowitsch PetscCall(PetscFree(areal)); 43351b4a12fSMatthew G. Knepley } 43451b4a12fSMatthew G. Knepley #else 4359566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(lgctx->lg, ptime, a)); 43651b4a12fSMatthew G. Knepley #endif 4379566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(v, &a)); 43851b4a12fSMatthew G. Knepley if (((lgctx->howoften > 0) && (!(step % lgctx->howoften))) || ((lgctx->howoften == -1) && ts->reason)) { 4399566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(lgctx->lg)); 4409566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(lgctx->lg)); 44151b4a12fSMatthew G. Knepley } 44251b4a12fSMatthew G. Knepley PetscFunctionReturn(0); 44351b4a12fSMatthew G. Knepley } 444