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 20d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMTSDestroy_DMDA(DMTS sdm) 21d71ae5a4SJacob Faibussowitsch { 2224989b8cSPeter Brune PetscFunctionBegin; 239566063dSJacob Faibussowitsch PetscCall(PetscFree(sdm->data)); 243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2524989b8cSPeter Brune } 2624989b8cSPeter Brune 27d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMTSDuplicate_DMDA(DMTS oldsdm, DMTS sdm) 28d71ae5a4SJacob Faibussowitsch { 290313d0e0SJed Brown PetscFunctionBegin; 304dfa11a4SJacob Faibussowitsch PetscCall(PetscNew((DMTS_DA **)&sdm->data)); 319566063dSJacob Faibussowitsch if (oldsdm->data) PetscCall(PetscMemcpy(sdm->data, oldsdm->data, sizeof(DMTS_DA))); 323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 330313d0e0SJed Brown } 340313d0e0SJed Brown 35d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMDATSGetContext(DM dm, DMTS sdm, DMTS_DA **dmdats) 36d71ae5a4SJacob Faibussowitsch { 3724989b8cSPeter Brune PetscFunctionBegin; 380298fd71SBarry Smith *dmdats = NULL; 3924989b8cSPeter Brune if (!sdm->data) { 404dfa11a4SJacob Faibussowitsch PetscCall(PetscNew((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; 453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4624989b8cSPeter Brune } 4724989b8cSPeter Brune 48*2a8381b2SBarry Smith static PetscErrorCode TSComputeIFunction_DMDA(TS ts, PetscReal ptime, Vec X, Vec Xdot, Vec F, PetscCtx ctx) 49d71ae5a4SJacob Faibussowitsch { 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; 93d71ae5a4SJacob Faibussowitsch default: 94d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_INCOMP, "Cannot use imode=%d", (int)dmdats->ifunctionlocalimode); 9524989b8cSPeter Brune } 969566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, Xloc, &x)); 979566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &Xloc)); 989566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, Xdotloc, &xdot)); 999566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &Xdotloc)); 1003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10124989b8cSPeter Brune } 10224989b8cSPeter Brune 103*2a8381b2SBarry Smith static PetscErrorCode TSComputeIJacobian_DMDA(TS ts, PetscReal ptime, Vec X, Vec Xdot, PetscReal shift, Mat A, Mat B, PetscCtx ctx) 104d71ae5a4SJacob Faibussowitsch { 10524989b8cSPeter Brune DM dm; 106942e3340SBarry Smith DMTS_DA *dmdats = (DMTS_DA *)ctx; 10724989b8cSPeter Brune DMDALocalInfo info; 108e0ec5be5SStefano Zampini Vec Xloc, Xdotloc; 10924989b8cSPeter Brune void *x, *xdot; 11024989b8cSPeter Brune 11124989b8cSPeter Brune PetscFunctionBegin; 1123c633725SBarry Smith PetscCheck(dmdats->ifunctionlocal, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Corrupt context"); 1139566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 11424989b8cSPeter Brune 115966bd95aSPierre Jolivet PetscCheck(dmdats->ijacobianlocal, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "TSComputeIJacobian_DMDA() called without calling DMDATSSetIJacobian()"); 1169566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &Xloc)); 1179566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc)); 1189566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc)); 119e0ec5be5SStefano Zampini PetscCall(DMGetLocalVector(dm, &Xdotloc)); 120e0ec5be5SStefano Zampini PetscCall(DMGlobalToLocalBegin(dm, Xdot, INSERT_VALUES, Xdotloc)); 121e0ec5be5SStefano Zampini PetscCall(DMGlobalToLocalEnd(dm, Xdot, INSERT_VALUES, Xdotloc)); 1229566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 1239566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm, Xloc, &x)); 124e0ec5be5SStefano Zampini PetscCall(DMDAVecGetArray(dm, Xdotloc, &xdot)); 12524989b8cSPeter Brune CHKMEMQ; 1269566063dSJacob Faibussowitsch PetscCall((*dmdats->ijacobianlocal)(&info, ptime, x, xdot, shift, A, B, dmdats->ijacobianlocalctx)); 12724989b8cSPeter Brune CHKMEMQ; 1289566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, Xloc, &x)); 129e0ec5be5SStefano Zampini PetscCall(DMDAVecRestoreArray(dm, Xdotloc, &xdot)); 1309566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &Xloc)); 131e0ec5be5SStefano Zampini PetscCall(DMRestoreLocalVector(dm, &Xdotloc)); 13224989b8cSPeter Brune /* This will be redundant if the user called both, but it's too common to forget. */ 13394ab13aaSBarry Smith if (A != B) { 1349566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1359566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 13624989b8cSPeter Brune } 1373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13824989b8cSPeter Brune } 13924989b8cSPeter Brune 140*2a8381b2SBarry Smith static PetscErrorCode TSComputeRHSFunction_DMDA(TS ts, PetscReal ptime, Vec X, Vec F, PetscCtx ctx) 141d71ae5a4SJacob Faibussowitsch { 14224989b8cSPeter Brune DM dm; 143942e3340SBarry Smith DMTS_DA *dmdats = (DMTS_DA *)ctx; 14424989b8cSPeter Brune DMDALocalInfo info; 14524989b8cSPeter Brune Vec Xloc; 14624989b8cSPeter Brune void *x, *f; 14724989b8cSPeter Brune 14824989b8cSPeter Brune PetscFunctionBegin; 14924989b8cSPeter Brune PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 150064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(X, VEC_CLASSID, 3); 151064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(F, VEC_CLASSID, 4); 1523c633725SBarry Smith PetscCheck(dmdats->rhsfunctionlocal, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Corrupt context"); 1539566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 1549566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &Xloc)); 1559566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc)); 1569566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc)); 1579566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 1589566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm, Xloc, &x)); 15924989b8cSPeter Brune switch (dmdats->rhsfunctionlocalimode) { 16024989b8cSPeter Brune case INSERT_VALUES: { 1619566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm, F, &f)); 16224989b8cSPeter Brune CHKMEMQ; 1639566063dSJacob Faibussowitsch PetscCall((*dmdats->rhsfunctionlocal)(&info, ptime, x, f, dmdats->rhsfunctionlocalctx)); 16424989b8cSPeter Brune CHKMEMQ; 1659566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, F, &f)); 16624989b8cSPeter Brune } break; 16724989b8cSPeter Brune case ADD_VALUES: { 16824989b8cSPeter Brune Vec Floc; 1699566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &Floc)); 1709566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Floc)); 1719566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm, Floc, &f)); 17224989b8cSPeter Brune CHKMEMQ; 1739566063dSJacob Faibussowitsch PetscCall((*dmdats->rhsfunctionlocal)(&info, ptime, x, f, dmdats->rhsfunctionlocalctx)); 17424989b8cSPeter Brune CHKMEMQ; 1759566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, Floc, &f)); 1769566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 1779566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm, Floc, ADD_VALUES, F)); 1789566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm, Floc, ADD_VALUES, F)); 1799566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &Floc)); 18024989b8cSPeter Brune } break; 181d71ae5a4SJacob Faibussowitsch default: 182d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_INCOMP, "Cannot use imode=%d", (int)dmdats->rhsfunctionlocalimode); 18324989b8cSPeter Brune } 1849566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, Xloc, &x)); 1859566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &Xloc)); 1863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18724989b8cSPeter Brune } 18824989b8cSPeter Brune 189*2a8381b2SBarry Smith static PetscErrorCode TSComputeRHSJacobian_DMDA(TS ts, PetscReal ptime, Vec X, Mat A, Mat B, PetscCtx ctx) 190d71ae5a4SJacob Faibussowitsch { 19124989b8cSPeter Brune DM dm; 192942e3340SBarry Smith DMTS_DA *dmdats = (DMTS_DA *)ctx; 19324989b8cSPeter Brune DMDALocalInfo info; 19424989b8cSPeter Brune Vec Xloc; 19524989b8cSPeter Brune void *x; 19624989b8cSPeter Brune 19724989b8cSPeter Brune PetscFunctionBegin; 1983c633725SBarry Smith PetscCheck(dmdats->rhsfunctionlocal, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Corrupt context"); 1999566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 20024989b8cSPeter Brune 201966bd95aSPierre Jolivet PetscCheck(dmdats->rhsjacobianlocal, PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "TSComputeRHSJacobian_DMDA() called without calling DMDATSSetRHSJacobian()"); 2029566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &Xloc)); 2039566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, Xloc)); 2049566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, Xloc)); 2059566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(dm, &info)); 2069566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(dm, Xloc, &x)); 20724989b8cSPeter Brune CHKMEMQ; 2089566063dSJacob Faibussowitsch PetscCall((*dmdats->rhsjacobianlocal)(&info, ptime, x, A, B, dmdats->rhsjacobianlocalctx)); 20924989b8cSPeter Brune CHKMEMQ; 2109566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(dm, Xloc, &x)); 2119566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &Xloc)); 21224989b8cSPeter Brune /* This will be redundant if the user called both, but it's too common to forget. */ 21394ab13aaSBarry Smith if (A != B) { 2149566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 2159566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 21624989b8cSPeter Brune } 2173ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 21824989b8cSPeter Brune } 21924989b8cSPeter Brune 22024989b8cSPeter Brune /*@C 221bcf0153eSBarry Smith DMDATSSetRHSFunctionLocal - set a local residual evaluation function for use with `DMDA` 22224989b8cSPeter Brune 22324989b8cSPeter Brune Logically Collective 22424989b8cSPeter Brune 2254165533cSJose E. Roman Input Parameters: 226bcf0153eSBarry Smith + dm - `DM` to associate callback with 227b398ffe5SPeter Brune . imode - insert mode for the residual 2288434afd1SBarry Smith . func - local residual evaluation, see `DMDATSRHSFunctionLocalFn` for the calling sequence 22924989b8cSPeter Brune - ctx - optional context for local residual evaluation 23024989b8cSPeter Brune 23124989b8cSPeter Brune Level: beginner 23224989b8cSPeter Brune 2338434afd1SBarry Smith .seealso: [](ch_ts), `DMDA`, `DMDATSRHSFunctionLocalFn`, `TS`, `TSSetRHSFunction()`, `DMTSSetRHSFunction()`, `DMDATSSetRHSJacobianLocal()`, `DMDASNESSetFunctionLocal()` 23424989b8cSPeter Brune @*/ 235*2a8381b2SBarry Smith PetscErrorCode DMDATSSetRHSFunctionLocal(DM dm, InsertMode imode, DMDATSRHSFunctionLocalFn *func, PetscCtx ctx) 236d71ae5a4SJacob Faibussowitsch { 237942e3340SBarry Smith DMTS sdm; 238942e3340SBarry Smith DMTS_DA *dmdats; 23924989b8cSPeter Brune 24024989b8cSPeter Brune PetscFunctionBegin; 24124989b8cSPeter Brune PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2429566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &sdm)); 2439566063dSJacob Faibussowitsch PetscCall(DMDATSGetContext(dm, sdm, &dmdats)); 24424989b8cSPeter Brune dmdats->rhsfunctionlocalimode = imode; 24524989b8cSPeter Brune dmdats->rhsfunctionlocal = func; 24624989b8cSPeter Brune dmdats->rhsfunctionlocalctx = ctx; 2479566063dSJacob Faibussowitsch PetscCall(DMTSSetRHSFunction(dm, TSComputeRHSFunction_DMDA, dmdats)); 2483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 24924989b8cSPeter Brune } 25024989b8cSPeter Brune 25124989b8cSPeter Brune /*@C 252bcf0153eSBarry Smith DMDATSSetRHSJacobianLocal - set a local residual evaluation function for use with `DMDA` 25324989b8cSPeter Brune 25424989b8cSPeter Brune Logically Collective 25524989b8cSPeter Brune 2564165533cSJose E. Roman Input Parameters: 257bcf0153eSBarry Smith + dm - `DM` to associate callback with 2588434afd1SBarry Smith . func - local RHS Jacobian evaluation routine, see `DMDATSRHSJacobianLocalFn` for the calling sequence 259b398ffe5SPeter Brune - ctx - optional context for local jacobian evaluation 26024989b8cSPeter Brune 26124989b8cSPeter Brune Level: beginner 26224989b8cSPeter Brune 2638434afd1SBarry Smith .seealso: [](ch_ts), `DMDA`, `DMDATSRHSJacobianLocalFn`, `DMTSSetRHSJacobian()`, 26414d0ab18SJacob Faibussowitsch `DMDATSSetRHSFunctionLocal()`, `DMDASNESSetJacobianLocal()` 26524989b8cSPeter Brune @*/ 266*2a8381b2SBarry Smith PetscErrorCode DMDATSSetRHSJacobianLocal(DM dm, DMDATSRHSJacobianLocalFn *func, PetscCtx ctx) 267d71ae5a4SJacob Faibussowitsch { 268942e3340SBarry Smith DMTS sdm; 269942e3340SBarry Smith DMTS_DA *dmdats; 27024989b8cSPeter Brune 27124989b8cSPeter Brune PetscFunctionBegin; 27224989b8cSPeter Brune PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2739566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &sdm)); 2749566063dSJacob Faibussowitsch PetscCall(DMDATSGetContext(dm, sdm, &dmdats)); 27524989b8cSPeter Brune dmdats->rhsjacobianlocal = func; 27624989b8cSPeter Brune dmdats->rhsjacobianlocalctx = ctx; 2779566063dSJacob Faibussowitsch PetscCall(DMTSSetRHSJacobian(dm, TSComputeRHSJacobian_DMDA, dmdats)); 2783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 27924989b8cSPeter Brune } 28024989b8cSPeter Brune 28124989b8cSPeter Brune /*@C 282bcf0153eSBarry Smith DMDATSSetIFunctionLocal - set a local residual evaluation function for use with `DMDA` 28324989b8cSPeter Brune 28424989b8cSPeter Brune Logically Collective 28524989b8cSPeter Brune 2864165533cSJose E. Roman Input Parameters: 287bcf0153eSBarry Smith + dm - `DM` to associate callback with 28814d0ab18SJacob Faibussowitsch . imode - the insert mode of the function 2898434afd1SBarry Smith . func - local residual evaluation, see `DMDATSIFunctionLocalFn` for the calling sequence 29024989b8cSPeter Brune - ctx - optional context for local residual evaluation 29124989b8cSPeter Brune 29224989b8cSPeter Brune Level: beginner 29324989b8cSPeter Brune 2948434afd1SBarry Smith .seealso: [](ch_ts), `DMDA`, `DMDATSIFunctionLocalFn`, `DMTSSetIFunction()`, 29514d0ab18SJacob Faibussowitsch `DMDATSSetIJacobianLocal()`, `DMDASNESSetFunctionLocal()` 29624989b8cSPeter Brune @*/ 297*2a8381b2SBarry Smith PetscErrorCode DMDATSSetIFunctionLocal(DM dm, InsertMode imode, DMDATSIFunctionLocalFn *func, PetscCtx ctx) 298d71ae5a4SJacob Faibussowitsch { 299942e3340SBarry Smith DMTS sdm; 300942e3340SBarry Smith DMTS_DA *dmdats; 30124989b8cSPeter Brune 30224989b8cSPeter Brune PetscFunctionBegin; 30324989b8cSPeter Brune PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3049566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &sdm)); 3059566063dSJacob Faibussowitsch PetscCall(DMDATSGetContext(dm, sdm, &dmdats)); 30624989b8cSPeter Brune dmdats->ifunctionlocalimode = imode; 30724989b8cSPeter Brune dmdats->ifunctionlocal = func; 30824989b8cSPeter Brune dmdats->ifunctionlocalctx = ctx; 3099566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunction(dm, TSComputeIFunction_DMDA, dmdats)); 3103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 31124989b8cSPeter Brune } 31224989b8cSPeter Brune 31324989b8cSPeter Brune /*@C 314bcf0153eSBarry Smith DMDATSSetIJacobianLocal - set a local residual evaluation function for use with `DMDA` 31524989b8cSPeter Brune 31624989b8cSPeter Brune Logically Collective 31724989b8cSPeter Brune 3184165533cSJose E. Roman Input Parameters: 319bcf0153eSBarry Smith + dm - `DM` to associate callback with 3208434afd1SBarry Smith . func - local residual evaluation, see `DMDATSIJacobianLocalFn` for the calling sequence 32124989b8cSPeter Brune - ctx - optional context for local residual evaluation 32224989b8cSPeter Brune 32324989b8cSPeter Brune Level: beginner 32424989b8cSPeter Brune 32590907b5bSBarry Smith .seealso: [](ch_ts), `DMDA`, `DMDATSIJacobianLocalFn`, `DMTSSetIJacobian()`, 32614d0ab18SJacob Faibussowitsch `DMDATSSetIFunctionLocal()`, `DMDASNESSetJacobianLocal()` 32724989b8cSPeter Brune @*/ 328*2a8381b2SBarry Smith PetscErrorCode DMDATSSetIJacobianLocal(DM dm, DMDATSIJacobianLocalFn *func, PetscCtx ctx) 329d71ae5a4SJacob Faibussowitsch { 330942e3340SBarry Smith DMTS sdm; 331942e3340SBarry Smith DMTS_DA *dmdats; 33224989b8cSPeter Brune 33324989b8cSPeter Brune PetscFunctionBegin; 33424989b8cSPeter Brune PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3359566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &sdm)); 3369566063dSJacob Faibussowitsch PetscCall(DMDATSGetContext(dm, sdm, &dmdats)); 33724989b8cSPeter Brune dmdats->ijacobianlocal = func; 33824989b8cSPeter Brune dmdats->ijacobianlocalctx = ctx; 3399566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMDA, dmdats)); 3403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 34124989b8cSPeter Brune } 342d1212d36SBarry Smith 343*2a8381b2SBarry Smith PetscErrorCode TSMonitorDMDARayDestroy(PetscCtxRt mctx) 344d71ae5a4SJacob Faibussowitsch { 345*2a8381b2SBarry Smith TSMonitorDMDARayCtx *rayctx = *(TSMonitorDMDARayCtx **)mctx; 346d1212d36SBarry Smith 347d1212d36SBarry Smith PetscFunctionBegin; 3489566063dSJacob Faibussowitsch if (rayctx->lgctx) PetscCall(TSMonitorLGCtxDestroy(&rayctx->lgctx)); 3499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&rayctx->ray)); 3509566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&rayctx->scatter)); 3519566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&rayctx->viewer)); 3529566063dSJacob Faibussowitsch PetscCall(PetscFree(rayctx)); 3533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 354d1212d36SBarry Smith } 355d1212d36SBarry Smith 356d71ae5a4SJacob Faibussowitsch PetscErrorCode TSMonitorDMDARay(TS ts, PetscInt steps, PetscReal time, Vec u, void *mctx) 357d71ae5a4SJacob Faibussowitsch { 358d1212d36SBarry Smith TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *)mctx; 359d1212d36SBarry Smith Vec solution; 360d1212d36SBarry Smith 361d1212d36SBarry Smith PetscFunctionBegin; 3629566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &solution)); 3639566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(rayctx->scatter, solution, rayctx->ray, INSERT_VALUES, SCATTER_FORWARD)); 3649566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(rayctx->scatter, solution, rayctx->ray, INSERT_VALUES, SCATTER_FORWARD)); 3651baa6e33SBarry Smith if (rayctx->viewer) PetscCall(VecView(rayctx->ray, rayctx->viewer)); 3663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 367d1212d36SBarry Smith } 36851b4a12fSMatthew G. Knepley 369*2a8381b2SBarry Smith PetscErrorCode TSMonitorLGDMDARay(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscCtx ctx) 370d71ae5a4SJacob Faibussowitsch { 37151b4a12fSMatthew G. Knepley TSMonitorDMDARayCtx *rayctx = (TSMonitorDMDARayCtx *)ctx; 372835f2295SStefano Zampini TSMonitorLGCtx lgctx = rayctx->lgctx; 37351b4a12fSMatthew G. Knepley Vec v = rayctx->ray; 37451b4a12fSMatthew G. Knepley const PetscScalar *a; 37551b4a12fSMatthew G. Knepley PetscInt dim; 37651b4a12fSMatthew G. Knepley 37751b4a12fSMatthew G. Knepley PetscFunctionBegin; 3789566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(rayctx->scatter, u, v, INSERT_VALUES, SCATTER_FORWARD)); 3799566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(rayctx->scatter, u, v, INSERT_VALUES, SCATTER_FORWARD)); 38051b4a12fSMatthew G. Knepley if (!step) { 38151b4a12fSMatthew G. Knepley PetscDrawAxis axis; 38251b4a12fSMatthew G. Knepley 3839566063dSJacob Faibussowitsch PetscCall(PetscDrawLGGetAxis(lgctx->lg, &axis)); 3849566063dSJacob Faibussowitsch PetscCall(PetscDrawAxisSetLabels(axis, "Solution Ray as function of time", "Time", "Solution")); 3859566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(rayctx->ray, &dim)); 3869566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSetDimension(lgctx->lg, dim)); 3879566063dSJacob Faibussowitsch PetscCall(PetscDrawLGReset(lgctx->lg)); 38851b4a12fSMatthew G. Knepley } 3899566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(v, &a)); 39051b4a12fSMatthew G. Knepley #if defined(PETSC_USE_COMPLEX) 39151b4a12fSMatthew G. Knepley { 39251b4a12fSMatthew G. Knepley PetscReal *areal; 39351b4a12fSMatthew G. Knepley PetscInt i, n; 3949566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(v, &n)); 3959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &areal)); 39651b4a12fSMatthew G. Knepley for (i = 0; i < n; ++i) areal[i] = PetscRealPart(a[i]); 3979566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(lgctx->lg, ptime, areal)); 3989566063dSJacob Faibussowitsch PetscCall(PetscFree(areal)); 39951b4a12fSMatthew G. Knepley } 40051b4a12fSMatthew G. Knepley #else 4019566063dSJacob Faibussowitsch PetscCall(PetscDrawLGAddCommonPoint(lgctx->lg, ptime, a)); 40251b4a12fSMatthew G. Knepley #endif 4039566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(v, &a)); 40451b4a12fSMatthew G. Knepley if (((lgctx->howoften > 0) && (!(step % lgctx->howoften))) || ((lgctx->howoften == -1) && ts->reason)) { 4059566063dSJacob Faibussowitsch PetscCall(PetscDrawLGDraw(lgctx->lg)); 4069566063dSJacob Faibussowitsch PetscCall(PetscDrawLGSave(lgctx->lg)); 40751b4a12fSMatthew G. Knepley } 4083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 40951b4a12fSMatthew G. Knepley } 410