1af0996ceSBarry Smith #include <petsc/private/dmimpl.h> 2af0996ceSBarry Smith #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 3d433e6cbSMatthew G. Knepley 4d433e6cbSMatthew G. Knepley typedef struct { 5f2b97b97SMatthew G. Knepley PetscErrorCode (*boundarylocal)(DM, PetscReal, Vec, Vec, void *); 6d433e6cbSMatthew G. Knepley PetscErrorCode (*ifunctionlocal)(DM, PetscReal, Vec, Vec, Vec, void *); 7d433e6cbSMatthew G. Knepley PetscErrorCode (*ijacobianlocal)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 8d433e6cbSMatthew G. Knepley PetscErrorCode (*rhsfunctionlocal)(DM, PetscReal, Vec, Vec, void *); 9f2b97b97SMatthew G. Knepley void *boundarylocalctx; 10d433e6cbSMatthew G. Knepley void *ifunctionlocalctx; 11d433e6cbSMatthew G. Knepley void *ijacobianlocalctx; 12d433e6cbSMatthew G. Knepley void *rhsfunctionlocalctx; 13b4937a87SMatthew G. Knepley Vec lumpedmassinv; 14b4937a87SMatthew G. Knepley Mat mass; 15b4937a87SMatthew G. Knepley KSP kspmass; 16d433e6cbSMatthew G. Knepley } DMTS_Local; 17d433e6cbSMatthew G. Knepley 18d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMTSDestroy_DMLocal(DMTS tdm) 19d71ae5a4SJacob Faibussowitsch { 20d433e6cbSMatthew G. Knepley PetscFunctionBegin; 219566063dSJacob Faibussowitsch PetscCall(PetscFree(tdm->data)); 223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 23d433e6cbSMatthew G. Knepley } 24d433e6cbSMatthew G. Knepley 25d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMTSDuplicate_DMLocal(DMTS oldtdm, DMTS tdm) 26d71ae5a4SJacob Faibussowitsch { 27d433e6cbSMatthew G. Knepley PetscFunctionBegin; 284dfa11a4SJacob Faibussowitsch PetscCall(PetscNew((DMTS_Local **)&tdm->data)); 299566063dSJacob Faibussowitsch if (oldtdm->data) PetscCall(PetscMemcpy(tdm->data, oldtdm->data, sizeof(DMTS_Local))); 303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 31d433e6cbSMatthew G. Knepley } 32d433e6cbSMatthew G. Knepley 33d71ae5a4SJacob Faibussowitsch static PetscErrorCode DMLocalTSGetContext(DM dm, DMTS tdm, DMTS_Local **dmlocalts) 34d71ae5a4SJacob Faibussowitsch { 35d433e6cbSMatthew G. Knepley PetscFunctionBegin; 36d433e6cbSMatthew G. Knepley *dmlocalts = NULL; 37d433e6cbSMatthew G. Knepley if (!tdm->data) { 384dfa11a4SJacob Faibussowitsch PetscCall(PetscNew((DMTS_Local **)&tdm->data)); 39d433e6cbSMatthew G. Knepley 40d433e6cbSMatthew G. Knepley tdm->ops->destroy = DMTSDestroy_DMLocal; 41d433e6cbSMatthew G. Knepley tdm->ops->duplicate = DMTSDuplicate_DMLocal; 42d433e6cbSMatthew G. Knepley } 43d433e6cbSMatthew G. Knepley *dmlocalts = (DMTS_Local *)tdm->data; 443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 45d433e6cbSMatthew G. Knepley } 46d433e6cbSMatthew G. Knepley 47*2a8381b2SBarry Smith static PetscErrorCode TSComputeIFunction_DMLocal(TS ts, PetscReal time, Vec X, Vec X_t, Vec F, PetscCtx ctx) 48d71ae5a4SJacob Faibussowitsch { 49d433e6cbSMatthew G. Knepley DM dm; 50d433e6cbSMatthew G. Knepley Vec locX, locX_t, locF; 51d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts = (DMTS_Local *)ctx; 52d433e6cbSMatthew G. Knepley 53d433e6cbSMatthew G. Knepley PetscFunctionBegin; 54d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 55d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(X, VEC_CLASSID, 3); 56d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(X_t, VEC_CLASSID, 4); 57d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(F, VEC_CLASSID, 5); 589566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 599566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &locX)); 609566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &locX_t)); 619566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &locF)); 629566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(locX)); 639566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(locX_t)); 649566063dSJacob Faibussowitsch if (dmlocalts->boundarylocal) PetscCall((*dmlocalts->boundarylocal)(dm, time, locX, locX_t, dmlocalts->boundarylocalctx)); 659566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX)); 669566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX)); 679566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, X_t, INSERT_VALUES, locX_t)); 689566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, X_t, INSERT_VALUES, locX_t)); 699566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(locF)); 70d433e6cbSMatthew G. Knepley CHKMEMQ; 719566063dSJacob Faibussowitsch PetscCall((*dmlocalts->ifunctionlocal)(dm, time, locX, locX_t, locF, dmlocalts->ifunctionlocalctx)); 72d433e6cbSMatthew G. Knepley CHKMEMQ; 739566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 749566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm, locF, ADD_VALUES, F)); 759566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm, locF, ADD_VALUES, F)); 769566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &locX)); 779566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &locX_t)); 789566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &locF)); 798d78ee70SStefano Zampini 808d78ee70SStefano Zampini /* remove nullspace from residual */ 818d78ee70SStefano Zampini { 828d78ee70SStefano Zampini MatNullSpace nullsp; 838d78ee70SStefano Zampini PetscCall(PetscObjectQuery((PetscObject)dm, "__dmtsnullspace", (PetscObject *)&nullsp)); 848d78ee70SStefano Zampini if (nullsp) PetscCall(MatNullSpaceRemove(nullsp, F)); 858d78ee70SStefano Zampini } 863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 87d433e6cbSMatthew G. Knepley } 88d433e6cbSMatthew G. Knepley 89*2a8381b2SBarry Smith static PetscErrorCode TSComputeRHSFunction_DMLocal(TS ts, PetscReal time, Vec X, Vec F, PetscCtx ctx) 90d71ae5a4SJacob Faibussowitsch { 91d433e6cbSMatthew G. Knepley DM dm; 92b4937a87SMatthew G. Knepley Vec locX, locF; 93d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts = (DMTS_Local *)ctx; 94d433e6cbSMatthew G. Knepley 95d433e6cbSMatthew G. Knepley PetscFunctionBegin; 96d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 97d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(X, VEC_CLASSID, 3); 98064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(F, VEC_CLASSID, 4); 999566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 1009566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &locX)); 1019566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &locF)); 1029566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(locX)); 1039566063dSJacob Faibussowitsch if (dmlocalts->boundarylocal) PetscCall((*dmlocalts->boundarylocal)(dm, time, locX, NULL, dmlocalts->boundarylocalctx)); 1049566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX)); 1059566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX)); 1069566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(locF)); 107b4937a87SMatthew G. Knepley CHKMEMQ; 1089566063dSJacob Faibussowitsch PetscCall((*dmlocalts->rhsfunctionlocal)(dm, time, locX, locF, dmlocalts->rhsfunctionlocalctx)); 109b4937a87SMatthew G. Knepley CHKMEMQ; 1109566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 1119566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm, locF, ADD_VALUES, F)); 1129566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm, locF, ADD_VALUES, F)); 113b4937a87SMatthew G. Knepley if (dmlocalts->lumpedmassinv) { 1149566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(F, dmlocalts->lumpedmassinv, F)); 115b4937a87SMatthew G. Knepley } else if (dmlocalts->kspmass) { 116b4937a87SMatthew G. Knepley Vec tmp; 117b4937a87SMatthew G. Knepley 1189566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &tmp)); 1199566063dSJacob Faibussowitsch PetscCall(KSPSolve(dmlocalts->kspmass, F, tmp)); 1209566063dSJacob Faibussowitsch PetscCall(VecCopy(tmp, F)); 1219566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &tmp)); 122b4937a87SMatthew G. Knepley } 1239566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &locX)); 1249566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &locF)); 1253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 126d433e6cbSMatthew G. Knepley } 127d433e6cbSMatthew G. Knepley 128*2a8381b2SBarry Smith static PetscErrorCode TSComputeIJacobian_DMLocal(TS ts, PetscReal time, Vec X, Vec X_t, PetscReal a, Mat A, Mat B, PetscCtx ctx) 129d71ae5a4SJacob Faibussowitsch { 130d433e6cbSMatthew G. Knepley DM dm; 131d433e6cbSMatthew G. Knepley Vec locX, locX_t; 132d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts = (DMTS_Local *)ctx; 133d433e6cbSMatthew G. Knepley 134d433e6cbSMatthew G. Knepley PetscFunctionBegin; 1359566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 136d433e6cbSMatthew G. Knepley if (dmlocalts->ijacobianlocal) { 1379566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &locX)); 1389566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &locX_t)); 1399566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(locX)); 1409566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(locX_t)); 1419566063dSJacob Faibussowitsch if (dmlocalts->boundarylocal) PetscCall((*dmlocalts->boundarylocal)(dm, time, locX, locX_t, dmlocalts->boundarylocalctx)); 1429566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX)); 1439566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX)); 1449566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, X_t, INSERT_VALUES, locX_t)); 1459566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, X_t, INSERT_VALUES, locX_t)); 146d433e6cbSMatthew G. Knepley CHKMEMQ; 1479566063dSJacob Faibussowitsch PetscCall((*dmlocalts->ijacobianlocal)(dm, time, locX, locX_t, a, A, B, dmlocalts->ijacobianlocalctx)); 148d433e6cbSMatthew G. Knepley CHKMEMQ; 1499566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &locX)); 1509566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &locX_t)); 151d433e6cbSMatthew G. Knepley } else { 152d433e6cbSMatthew G. Knepley MatFDColoring fdcoloring; 1539566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject *)&fdcoloring)); 154d433e6cbSMatthew G. Knepley if (!fdcoloring) { 155d433e6cbSMatthew G. Knepley ISColoring coloring; 156d433e6cbSMatthew G. Knepley 1579566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(dm, dm->coloringtype, &coloring)); 1589566063dSJacob Faibussowitsch PetscCall(MatFDColoringCreate(B, coloring, &fdcoloring)); 1599566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&coloring)); 160d433e6cbSMatthew G. Knepley switch (dm->coloringtype) { 161d71ae5a4SJacob Faibussowitsch case IS_COLORING_GLOBAL: 162453a69bbSBarry Smith PetscCall(MatFDColoringSetFunction(fdcoloring, (MatFDColoringFn *)(PetscVoidFn *)TSComputeIFunction_DMLocal, dmlocalts)); 163d71ae5a4SJacob Faibussowitsch break; 164d71ae5a4SJacob Faibussowitsch default: 165d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "No support for coloring type '%s'", ISColoringTypes[dm->coloringtype]); 166d433e6cbSMatthew G. Knepley } 1679566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring, ((PetscObject)dm)->prefix)); 1689566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFromOptions(fdcoloring)); 1699566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetUp(B, coloring, fdcoloring)); 1709566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject)fdcoloring)); 1719566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)fdcoloring)); 172d433e6cbSMatthew G. Knepley 173d433e6cbSMatthew G. Knepley /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 174d433e6cbSMatthew G. Knepley * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 175d433e6cbSMatthew G. Knepley * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 176d433e6cbSMatthew G. Knepley * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 177d433e6cbSMatthew G. Knepley * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed. 178d433e6cbSMatthew G. Knepley */ 1799566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dm)); 180d433e6cbSMatthew G. Knepley } 1819566063dSJacob Faibussowitsch PetscCall(MatFDColoringApply(B, fdcoloring, X, ts)); 182d433e6cbSMatthew G. Knepley } 183d433e6cbSMatthew G. Knepley /* This will be redundant if the user called both, but it's too common to forget. */ 184d433e6cbSMatthew G. Knepley if (A != B) { 1859566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1869566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 187d433e6cbSMatthew G. Knepley } 1883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 189d433e6cbSMatthew G. Knepley } 190d433e6cbSMatthew G. Knepley 1917eb36faeSMatthew G. Knepley /*@C 192f2b97b97SMatthew G. Knepley DMTSSetBoundaryLocal - set the function for essential boundary data for a local implicit function evaluation. 193a40652d4SToby Isaac 1947eb36faeSMatthew G. Knepley Logically Collective 1957eb36faeSMatthew G. Knepley 1964165533cSJose E. Roman Input Parameters: 197bcf0153eSBarry Smith + dm - `DM` to associate callback with 198a40652d4SToby Isaac . func - local function evaluation 199a40652d4SToby Isaac - ctx - context for function evaluation 2007eb36faeSMatthew G. Knepley 2017eb36faeSMatthew G. Knepley Level: intermediate 2027eb36faeSMatthew G. Knepley 20314d0ab18SJacob Faibussowitsch Notes: 20414d0ab18SJacob Faibussowitsch `func` should set the essential boundary data for the local portion of the solution, as 20514d0ab18SJacob Faibussowitsch well its time derivative (if it is not `NULL`). 20614d0ab18SJacob Faibussowitsch 20714d0ab18SJacob Faibussowitsch Vectors are initialized to zero before this function, so it is only needed for non 20814d0ab18SJacob Faibussowitsch homogeneous data. 20914d0ab18SJacob Faibussowitsch 21014d0ab18SJacob Faibussowitsch This function is somewhat optional: boundary data could potentially be inserted by a function 21114d0ab18SJacob Faibussowitsch passed to `DMTSSetIFunctionLocal()`. The use case for this function is for discretizations 21214d0ab18SJacob Faibussowitsch with constraints (see `DMGetDefaultConstraints()`): this function inserts boundary values 21314d0ab18SJacob Faibussowitsch before constraint interpolation. 2142fe279fdSBarry Smith 2151cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 2167eb36faeSMatthew G. Knepley @*/ 217*2a8381b2SBarry Smith PetscErrorCode DMTSSetBoundaryLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, void *), PetscCtx ctx) 218d71ae5a4SJacob Faibussowitsch { 2197eb36faeSMatthew G. Knepley DMTS tdm; 2207eb36faeSMatthew G. Knepley DMTS_Local *dmlocalts; 2217eb36faeSMatthew G. Knepley 2227eb36faeSMatthew G. Knepley PetscFunctionBegin; 2237eb36faeSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2249566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 2259566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 2267eb36faeSMatthew G. Knepley 2277eb36faeSMatthew G. Knepley dmlocalts->boundarylocal = func; 2287eb36faeSMatthew G. Knepley dmlocalts->boundarylocalctx = ctx; 2293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2307eb36faeSMatthew G. Knepley } 2317eb36faeSMatthew G. Knepley 232d433e6cbSMatthew G. Knepley /*@C 2330c9409e4SMatthew G. Knepley DMTSGetIFunctionLocal - get the local implicit function evaluation function. This function is called with local vector 2340c9409e4SMatthew G. Knepley containing the local vector information PLUS ghost point information. It should compute a result for all local 235bcf0153eSBarry Smith elements and `DM` will automatically accumulate the overlapping values. 2360c9409e4SMatthew G. Knepley 2370c9409e4SMatthew G. Knepley Logically Collective 2380c9409e4SMatthew G. Knepley 2390c9409e4SMatthew G. Knepley Input Parameter: 240bcf0153eSBarry Smith . dm - `DM` to associate callback with 2410c9409e4SMatthew G. Knepley 2420c9409e4SMatthew G. Knepley Output Parameters: 2430c9409e4SMatthew G. Knepley + func - local function evaluation 2440c9409e4SMatthew G. Knepley - ctx - context for function evaluation 2450c9409e4SMatthew G. Knepley 2460c9409e4SMatthew G. Knepley Level: beginner 2470c9409e4SMatthew G. Knepley 2481cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSSetIFunctionLocal()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 2490c9409e4SMatthew G. Knepley @*/ 250*2a8381b2SBarry Smith PetscErrorCode DMTSGetIFunctionLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, Vec, void *), PetscCtxRt ctx) 251d71ae5a4SJacob Faibussowitsch { 2520c9409e4SMatthew G. Knepley DMTS tdm; 2530c9409e4SMatthew G. Knepley DMTS_Local *dmlocalts; 2540c9409e4SMatthew G. Knepley 2550c9409e4SMatthew G. Knepley PetscFunctionBegin; 2560c9409e4SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 257a86d55f5SMatthew G. Knepley PetscCall(DMGetDMTS(dm, &tdm)); 258a86d55f5SMatthew G. Knepley PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 2599371c9d4SSatish Balay if (func) { 2604f572ea9SToby Isaac PetscAssertPointer(func, 2); 2619371c9d4SSatish Balay *func = dmlocalts->ifunctionlocal; 2629371c9d4SSatish Balay } 2639371c9d4SSatish Balay if (ctx) { 2644f572ea9SToby Isaac PetscAssertPointer(ctx, 3); 265*2a8381b2SBarry Smith *(void **)ctx = dmlocalts->ifunctionlocalctx; 2669371c9d4SSatish Balay } 2673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2680c9409e4SMatthew G. Knepley } 2690c9409e4SMatthew G. Knepley 2700c9409e4SMatthew G. Knepley /*@C 271d433e6cbSMatthew G. Knepley DMTSSetIFunctionLocal - set a local implicit function evaluation function. This function is called with local vector 272d433e6cbSMatthew G. Knepley containing the local vector information PLUS ghost point information. It should compute a result for all local 273bcf0153eSBarry Smith elements and `DM` will automatically accumulate the overlapping values. 274d433e6cbSMatthew G. Knepley 275d433e6cbSMatthew G. Knepley Logically Collective 276d433e6cbSMatthew G. Knepley 2774165533cSJose E. Roman Input Parameters: 278bcf0153eSBarry Smith + dm - `DM` to associate callback with 279d433e6cbSMatthew G. Knepley . func - local function evaluation 280d433e6cbSMatthew G. Knepley - ctx - context for function evaluation 281d433e6cbSMatthew G. Knepley 282d433e6cbSMatthew G. Knepley Level: beginner 283d433e6cbSMatthew G. Knepley 2841cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSGetIFunctionLocal()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 285d433e6cbSMatthew G. Knepley @*/ 286*2a8381b2SBarry Smith PetscErrorCode DMTSSetIFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, Vec, void *), PetscCtx ctx) 287d71ae5a4SJacob Faibussowitsch { 288d433e6cbSMatthew G. Knepley DMTS tdm; 289d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts; 290d433e6cbSMatthew G. Knepley 291d433e6cbSMatthew G. Knepley PetscFunctionBegin; 292d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2939566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 2949566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 295d433e6cbSMatthew G. Knepley 296d433e6cbSMatthew G. Knepley dmlocalts->ifunctionlocal = func; 297d433e6cbSMatthew G. Knepley dmlocalts->ifunctionlocalctx = ctx; 298d433e6cbSMatthew G. Knepley 2999566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunction(dm, TSComputeIFunction_DMLocal, dmlocalts)); 300d433e6cbSMatthew G. Knepley if (!tdm->ops->ijacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 3019566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts)); 302d433e6cbSMatthew G. Knepley } 3033ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 304d433e6cbSMatthew G. Knepley } 305d433e6cbSMatthew G. Knepley 306d433e6cbSMatthew G. Knepley /*@C 3070c9409e4SMatthew G. Knepley DMTSGetIJacobianLocal - get a local Jacobian evaluation function 3080c9409e4SMatthew G. Knepley 3090c9409e4SMatthew G. Knepley Logically Collective 3100c9409e4SMatthew G. Knepley 3110c9409e4SMatthew G. Knepley Input Parameter: 312bcf0153eSBarry Smith . dm - `DM` to associate callback with 3130c9409e4SMatthew G. Knepley 3140c9409e4SMatthew G. Knepley Output Parameters: 3150c9409e4SMatthew G. Knepley + func - local Jacobian evaluation 3160c9409e4SMatthew G. Knepley - ctx - optional context for local Jacobian evaluation 3170c9409e4SMatthew G. Knepley 3180c9409e4SMatthew G. Knepley Level: beginner 3190c9409e4SMatthew G. Knepley 3201cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSSetIJacobianLocal()`, `DMTSSetIFunctionLocal()`, `DMTSSetIJacobian()`, `DMTSSetIFunction()` 3210c9409e4SMatthew G. Knepley @*/ 322*2a8381b2SBarry Smith PetscErrorCode DMTSGetIJacobianLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), PetscCtxRt ctx) 323d71ae5a4SJacob Faibussowitsch { 3240c9409e4SMatthew G. Knepley DMTS tdm; 3250c9409e4SMatthew G. Knepley DMTS_Local *dmlocalts; 3260c9409e4SMatthew G. Knepley 3270c9409e4SMatthew G. Knepley PetscFunctionBegin; 3280c9409e4SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 329a86d55f5SMatthew G. Knepley PetscCall(DMGetDMTS(dm, &tdm)); 330a86d55f5SMatthew G. Knepley PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 3319371c9d4SSatish Balay if (func) { 3324f572ea9SToby Isaac PetscAssertPointer(func, 2); 3339371c9d4SSatish Balay *func = dmlocalts->ijacobianlocal; 3349371c9d4SSatish Balay } 3359371c9d4SSatish Balay if (ctx) { 3364f572ea9SToby Isaac PetscAssertPointer(ctx, 3); 337*2a8381b2SBarry Smith *(void **)ctx = dmlocalts->ijacobianlocalctx; 3389371c9d4SSatish Balay } 3393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3400c9409e4SMatthew G. Knepley } 3410c9409e4SMatthew G. Knepley 3420c9409e4SMatthew G. Knepley /*@C 343d433e6cbSMatthew G. Knepley DMTSSetIJacobianLocal - set a local Jacobian evaluation function 344d433e6cbSMatthew G. Knepley 345d433e6cbSMatthew G. Knepley Logically Collective 346d433e6cbSMatthew G. Knepley 3474165533cSJose E. Roman Input Parameters: 348bcf0153eSBarry Smith + dm - `DM` to associate callback with 349d433e6cbSMatthew G. Knepley . func - local Jacobian evaluation 350d433e6cbSMatthew G. Knepley - ctx - optional context for local Jacobian evaluation 351d433e6cbSMatthew G. Knepley 352d433e6cbSMatthew G. Knepley Level: beginner 353d433e6cbSMatthew G. Knepley 3541cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSGetIJacobianLocal()`, `DMTSSetIFunctionLocal()`, `DMTSSetIJacobian()`, `DMTSSetIFunction()` 355d433e6cbSMatthew G. Knepley @*/ 356*2a8381b2SBarry Smith PetscErrorCode DMTSSetIJacobianLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), PetscCtx ctx) 357d71ae5a4SJacob Faibussowitsch { 358d433e6cbSMatthew G. Knepley DMTS tdm; 359d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts; 360d433e6cbSMatthew G. Knepley 361d433e6cbSMatthew G. Knepley PetscFunctionBegin; 362d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3639566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 3649566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 365d433e6cbSMatthew G. Knepley 366d433e6cbSMatthew G. Knepley dmlocalts->ijacobianlocal = func; 367d433e6cbSMatthew G. Knepley dmlocalts->ijacobianlocalctx = ctx; 368d433e6cbSMatthew G. Knepley 3699566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts)); 3703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 371d433e6cbSMatthew G. Knepley } 372d433e6cbSMatthew G. Knepley 373d433e6cbSMatthew G. Knepley /*@C 3740c9409e4SMatthew G. Knepley DMTSGetRHSFunctionLocal - get a local rhs function evaluation function. This function is called with local vector 3750c9409e4SMatthew G. Knepley containing the local vector information PLUS ghost point information. It should compute a result for all local 376bcf0153eSBarry Smith elements and `DM` will automatically accumulate the overlapping values. 3770c9409e4SMatthew G. Knepley 3780c9409e4SMatthew G. Knepley Logically Collective 3790c9409e4SMatthew G. Knepley 3800c9409e4SMatthew G. Knepley Input Parameter: 381bcf0153eSBarry Smith . dm - `DM` to associate callback with 3820c9409e4SMatthew G. Knepley 3830c9409e4SMatthew G. Knepley Output Parameters: 3840c9409e4SMatthew G. Knepley + func - local function evaluation 3850c9409e4SMatthew G. Knepley - ctx - context for function evaluation 3860c9409e4SMatthew G. Knepley 3870c9409e4SMatthew G. Knepley Level: beginner 3880c9409e4SMatthew G. Knepley 3891cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSSetRHSFunctionLocal()`, `DMTSSetRHSFunction()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 3900c9409e4SMatthew G. Knepley @*/ 391*2a8381b2SBarry Smith PetscErrorCode DMTSGetRHSFunctionLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, void *), PetscCtxRt ctx) 392d71ae5a4SJacob Faibussowitsch { 3930c9409e4SMatthew G. Knepley DMTS tdm; 3940c9409e4SMatthew G. Knepley DMTS_Local *dmlocalts; 3950c9409e4SMatthew G. Knepley 3960c9409e4SMatthew G. Knepley PetscFunctionBegin; 3970c9409e4SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 398a86d55f5SMatthew G. Knepley PetscCall(DMGetDMTS(dm, &tdm)); 399a86d55f5SMatthew G. Knepley PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 4009371c9d4SSatish Balay if (func) { 4014f572ea9SToby Isaac PetscAssertPointer(func, 2); 4029371c9d4SSatish Balay *func = dmlocalts->rhsfunctionlocal; 4039371c9d4SSatish Balay } 4049371c9d4SSatish Balay if (ctx) { 4054f572ea9SToby Isaac PetscAssertPointer(ctx, 3); 406*2a8381b2SBarry Smith *(void **)ctx = dmlocalts->rhsfunctionlocalctx; 4079371c9d4SSatish Balay } 4083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4090c9409e4SMatthew G. Knepley } 4100c9409e4SMatthew G. Knepley 4110c9409e4SMatthew G. Knepley /*@C 412d433e6cbSMatthew G. Knepley DMTSSetRHSFunctionLocal - set a local rhs function evaluation function. This function is called with local vector 413d433e6cbSMatthew G. Knepley containing the local vector information PLUS ghost point information. It should compute a result for all local 414bcf0153eSBarry Smith elements and `DM` will automatically accumulate the overlapping values. 415d433e6cbSMatthew G. Knepley 416d433e6cbSMatthew G. Knepley Logically Collective 417d433e6cbSMatthew G. Knepley 4184165533cSJose E. Roman Input Parameters: 419bcf0153eSBarry Smith + dm - `DM` to associate callback with 420d433e6cbSMatthew G. Knepley . func - local function evaluation 421d433e6cbSMatthew G. Knepley - ctx - context for function evaluation 422d433e6cbSMatthew G. Knepley 423d433e6cbSMatthew G. Knepley Level: beginner 424d433e6cbSMatthew G. Knepley 4251cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSGetRHSFunctionLocal()`, `DMTSSetRHSFunction()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 426d433e6cbSMatthew G. Knepley @*/ 427*2a8381b2SBarry Smith PetscErrorCode DMTSSetRHSFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, void *), PetscCtx ctx) 428d71ae5a4SJacob Faibussowitsch { 429d433e6cbSMatthew G. Knepley DMTS tdm; 430d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts; 431d433e6cbSMatthew G. Knepley 432d433e6cbSMatthew G. Knepley PetscFunctionBegin; 433d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4349566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 4359566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 436d433e6cbSMatthew G. Knepley 437d433e6cbSMatthew G. Knepley dmlocalts->rhsfunctionlocal = func; 438d433e6cbSMatthew G. Knepley dmlocalts->rhsfunctionlocalctx = ctx; 439d433e6cbSMatthew G. Knepley 4409566063dSJacob Faibussowitsch PetscCall(DMTSSetRHSFunction(dm, TSComputeRHSFunction_DMLocal, dmlocalts)); 4413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 442d433e6cbSMatthew G. Knepley } 443b4937a87SMatthew G. Knepley 444cc4c1da9SBarry Smith /*@ 445bcf0153eSBarry Smith DMTSCreateRHSMassMatrix - This creates the mass matrix associated with the given `DM`, and a solver to invert it, and stores them in the `DM` context. 446b4937a87SMatthew G. Knepley 447c3339decSBarry Smith Collective 448b4937a87SMatthew G. Knepley 4492fe279fdSBarry Smith Input Parameter: 450bcf0153eSBarry Smith . dm - `DM` providing the mass matrix 451b4937a87SMatthew G. Knepley 452b4937a87SMatthew G. Knepley Level: developer 453b4937a87SMatthew G. Knepley 454bcf0153eSBarry Smith Note: 455bcf0153eSBarry Smith The idea here is that an explicit system can be given a mass matrix, based on the `DM`, which is inverted on the RHS at each step. 456bcf0153eSBarry Smith 4571cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSCreateRHSMassMatrixLumped()`, `DMTSDestroyRHSMassMatrix()`, `DMCreateMassMatrix()`, `DMTS` 458b4937a87SMatthew G. Knepley @*/ 459d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSCreateRHSMassMatrix(DM dm) 460d71ae5a4SJacob Faibussowitsch { 461b4937a87SMatthew G. Knepley DMTS tdm; 462b4937a87SMatthew G. Knepley DMTS_Local *dmlocalts; 463b4937a87SMatthew G. Knepley const char *prefix; 464b4937a87SMatthew G. Knepley 465b4937a87SMatthew G. Knepley PetscFunctionBegin; 466b4937a87SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4679566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 4689566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 4699566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(dm, dm, &dmlocalts->mass)); 4709566063dSJacob Faibussowitsch PetscCall(KSPCreate(PetscObjectComm((PetscObject)dm), &dmlocalts->kspmass)); 4719566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 4729566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(dmlocalts->kspmass, prefix)); 4739566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(dmlocalts->kspmass, "mass_")); 4749566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(dmlocalts->kspmass)); 4759566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(dmlocalts->kspmass, dmlocalts->mass, dmlocalts->mass)); 4763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 477b4937a87SMatthew G. Knepley } 478b4937a87SMatthew G. Knepley 479cc4c1da9SBarry Smith /*@ 480bcf0153eSBarry Smith DMTSCreateRHSMassMatrixLumped - This creates the lumped mass matrix associated with the given `DM`, and a solver to invert it, and stores them in the `DM` context. 481b4937a87SMatthew G. Knepley 482c3339decSBarry Smith Collective 483b4937a87SMatthew G. Knepley 4842fe279fdSBarry Smith Input Parameter: 485bcf0153eSBarry Smith . dm - `DM` providing the mass matrix 486b4937a87SMatthew G. Knepley 487b4937a87SMatthew G. Knepley Level: developer 488b4937a87SMatthew G. Knepley 489bcf0153eSBarry Smith Note: 490bcf0153eSBarry Smith The idea here is that an explicit system can be given a mass matrix, based on the `DM`, which is inverted on the RHS at each step. 491bcf0153eSBarry Smith Since the matrix is lumped, inversion is trivial. 492bcf0153eSBarry Smith 4931cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSCreateRHSMassMatrix()`, `DMTSDestroyRHSMassMatrix()`, `DMCreateMassMatrix()`, `DMTS` 494b4937a87SMatthew G. Knepley @*/ 495d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM dm) 496d71ae5a4SJacob Faibussowitsch { 497b4937a87SMatthew G. Knepley DMTS tdm; 498b4937a87SMatthew G. Knepley DMTS_Local *dmlocalts; 499b4937a87SMatthew G. Knepley 500b4937a87SMatthew G. Knepley PetscFunctionBegin; 501b4937a87SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5029566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 5039566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 5048e9849d2SStefano Zampini PetscCall(DMCreateMassMatrixLumped(dm, NULL, &dmlocalts->lumpedmassinv)); 5059566063dSJacob Faibussowitsch PetscCall(VecReciprocal(dmlocalts->lumpedmassinv)); 5069566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(dmlocalts->lumpedmassinv, NULL, "-lumped_mass_inv_view")); 5073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 508b4937a87SMatthew G. Knepley } 509b4937a87SMatthew G. Knepley 510cc4c1da9SBarry Smith /*@ 511bcf0153eSBarry Smith DMTSDestroyRHSMassMatrix - Destroys the mass matrix and solver stored in the `DM` context, if they exist. 512b4937a87SMatthew G. Knepley 513b4937a87SMatthew G. Knepley Logically Collective 514b4937a87SMatthew G. Knepley 5152fe279fdSBarry Smith Input Parameter: 516bcf0153eSBarry Smith . dm - `DM` providing the mass matrix 517b4937a87SMatthew G. Knepley 518b4937a87SMatthew G. Knepley Level: developer 519b4937a87SMatthew G. Knepley 520b43aa488SJacob Faibussowitsch .seealso: [](ch_ts), `DM`, `DMTSCreateRHSMassMatrixLumped()`, `DMCreateMassMatrix()`, `DMTS` 521b4937a87SMatthew G. Knepley @*/ 522d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSDestroyRHSMassMatrix(DM dm) 523d71ae5a4SJacob Faibussowitsch { 524b4937a87SMatthew G. Knepley DMTS tdm; 525b4937a87SMatthew G. Knepley DMTS_Local *dmlocalts; 526b4937a87SMatthew G. Knepley 527b4937a87SMatthew G. Knepley PetscFunctionBegin; 528b4937a87SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5299566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 5309566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 5319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dmlocalts->lumpedmassinv)); 5329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&dmlocalts->mass)); 5339566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&dmlocalts->kspmass)); 5343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 535b4937a87SMatthew G. Knepley } 536