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 47d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeIFunction_DMLocal(TS ts, PetscReal time, Vec X, Vec X_t, Vec F, void *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)); 793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 80d433e6cbSMatthew G. Knepley } 81d433e6cbSMatthew G. Knepley 82d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeRHSFunction_DMLocal(TS ts, PetscReal time, Vec X, Vec F, void *ctx) 83d71ae5a4SJacob Faibussowitsch { 84d433e6cbSMatthew G. Knepley DM dm; 85b4937a87SMatthew G. Knepley Vec locX, locF; 86d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts = (DMTS_Local *)ctx; 87d433e6cbSMatthew G. Knepley 88d433e6cbSMatthew G. Knepley PetscFunctionBegin; 89d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 90d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(X, VEC_CLASSID, 3); 91064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(F, VEC_CLASSID, 4); 929566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 939566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &locX)); 949566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &locF)); 959566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(locX)); 969566063dSJacob Faibussowitsch if (dmlocalts->boundarylocal) PetscCall((*dmlocalts->boundarylocal)(dm, time, locX, NULL, dmlocalts->boundarylocalctx)); 979566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX)); 989566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX)); 999566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(locF)); 100b4937a87SMatthew G. Knepley CHKMEMQ; 1019566063dSJacob Faibussowitsch PetscCall((*dmlocalts->rhsfunctionlocal)(dm, time, locX, locF, dmlocalts->rhsfunctionlocalctx)); 102b4937a87SMatthew G. Knepley CHKMEMQ; 1039566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(F)); 1049566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalBegin(dm, locF, ADD_VALUES, F)); 1059566063dSJacob Faibussowitsch PetscCall(DMLocalToGlobalEnd(dm, locF, ADD_VALUES, F)); 106b4937a87SMatthew G. Knepley if (dmlocalts->lumpedmassinv) { 1079566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(F, dmlocalts->lumpedmassinv, F)); 108b4937a87SMatthew G. Knepley } else if (dmlocalts->kspmass) { 109b4937a87SMatthew G. Knepley Vec tmp; 110b4937a87SMatthew G. Knepley 1119566063dSJacob Faibussowitsch PetscCall(DMGetGlobalVector(dm, &tmp)); 1129566063dSJacob Faibussowitsch PetscCall(KSPSolve(dmlocalts->kspmass, F, tmp)); 1139566063dSJacob Faibussowitsch PetscCall(VecCopy(tmp, F)); 1149566063dSJacob Faibussowitsch PetscCall(DMRestoreGlobalVector(dm, &tmp)); 115b4937a87SMatthew G. Knepley } 1169566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &locX)); 1179566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &locF)); 1183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 119d433e6cbSMatthew G. Knepley } 120d433e6cbSMatthew G. Knepley 121d71ae5a4SJacob Faibussowitsch static PetscErrorCode TSComputeIJacobian_DMLocal(TS ts, PetscReal time, Vec X, Vec X_t, PetscReal a, Mat A, Mat B, void *ctx) 122d71ae5a4SJacob Faibussowitsch { 123d433e6cbSMatthew G. Knepley DM dm; 124d433e6cbSMatthew G. Knepley Vec locX, locX_t; 125d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts = (DMTS_Local *)ctx; 126d433e6cbSMatthew G. Knepley 127d433e6cbSMatthew G. Knepley PetscFunctionBegin; 1289566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 129d433e6cbSMatthew G. Knepley if (dmlocalts->ijacobianlocal) { 1309566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &locX)); 1319566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(dm, &locX_t)); 1329566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(locX)); 1339566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(locX_t)); 1349566063dSJacob Faibussowitsch if (dmlocalts->boundarylocal) PetscCall((*dmlocalts->boundarylocal)(dm, time, locX, locX_t, dmlocalts->boundarylocalctx)); 1359566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX)); 1369566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX)); 1379566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(dm, X_t, INSERT_VALUES, locX_t)); 1389566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(dm, X_t, INSERT_VALUES, locX_t)); 139d433e6cbSMatthew G. Knepley CHKMEMQ; 1409566063dSJacob Faibussowitsch PetscCall((*dmlocalts->ijacobianlocal)(dm, time, locX, locX_t, a, A, B, dmlocalts->ijacobianlocalctx)); 141d433e6cbSMatthew G. Knepley CHKMEMQ; 1429566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &locX)); 1439566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(dm, &locX_t)); 144d433e6cbSMatthew G. Knepley } else { 145d433e6cbSMatthew G. Knepley MatFDColoring fdcoloring; 1469566063dSJacob Faibussowitsch PetscCall(PetscObjectQuery((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject *)&fdcoloring)); 147d433e6cbSMatthew G. Knepley if (!fdcoloring) { 148d433e6cbSMatthew G. Knepley ISColoring coloring; 149d433e6cbSMatthew G. Knepley 1509566063dSJacob Faibussowitsch PetscCall(DMCreateColoring(dm, dm->coloringtype, &coloring)); 1519566063dSJacob Faibussowitsch PetscCall(MatFDColoringCreate(B, coloring, &fdcoloring)); 1529566063dSJacob Faibussowitsch PetscCall(ISColoringDestroy(&coloring)); 153d433e6cbSMatthew G. Knepley switch (dm->coloringtype) { 154d71ae5a4SJacob Faibussowitsch case IS_COLORING_GLOBAL: 155d71ae5a4SJacob Faibussowitsch PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))TSComputeIFunction_DMLocal, dmlocalts)); 156d71ae5a4SJacob Faibussowitsch break; 157d71ae5a4SJacob Faibussowitsch default: 158d71ae5a4SJacob Faibussowitsch SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "No support for coloring type '%s'", ISColoringTypes[dm->coloringtype]); 159d433e6cbSMatthew G. Knepley } 1609566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring, ((PetscObject)dm)->prefix)); 1619566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetFromOptions(fdcoloring)); 1629566063dSJacob Faibussowitsch PetscCall(MatFDColoringSetUp(B, coloring, fdcoloring)); 1639566063dSJacob Faibussowitsch PetscCall(PetscObjectCompose((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject)fdcoloring)); 1649566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)fdcoloring)); 165d433e6cbSMatthew G. Knepley 166d433e6cbSMatthew G. Knepley /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 167d433e6cbSMatthew G. Knepley * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 168d433e6cbSMatthew G. Knepley * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 169d433e6cbSMatthew G. Knepley * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 170d433e6cbSMatthew G. Knepley * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed. 171d433e6cbSMatthew G. Knepley */ 1729566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)dm)); 173d433e6cbSMatthew G. Knepley } 1749566063dSJacob Faibussowitsch PetscCall(MatFDColoringApply(B, fdcoloring, X, ts)); 175d433e6cbSMatthew G. Knepley } 176d433e6cbSMatthew G. Knepley /* This will be redundant if the user called both, but it's too common to forget. */ 177d433e6cbSMatthew G. Knepley if (A != B) { 1789566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1799566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 180d433e6cbSMatthew G. Knepley } 1813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 182d433e6cbSMatthew G. Knepley } 183d433e6cbSMatthew G. Knepley 1847eb36faeSMatthew G. Knepley /*@C 185f2b97b97SMatthew G. Knepley DMTSSetBoundaryLocal - set the function for essential boundary data for a local implicit function evaluation. 186a40652d4SToby Isaac 1877eb36faeSMatthew G. Knepley Logically Collective 1887eb36faeSMatthew G. Knepley 1894165533cSJose E. Roman Input Parameters: 190bcf0153eSBarry Smith + dm - `DM` to associate callback with 191a40652d4SToby Isaac . func - local function evaluation 192a40652d4SToby Isaac - ctx - context for function evaluation 1937eb36faeSMatthew G. Knepley 1947eb36faeSMatthew G. Knepley Level: intermediate 1957eb36faeSMatthew G. Knepley 19614d0ab18SJacob Faibussowitsch Notes: 19714d0ab18SJacob Faibussowitsch `func` should set the essential boundary data for the local portion of the solution, as 19814d0ab18SJacob Faibussowitsch well its time derivative (if it is not `NULL`). 19914d0ab18SJacob Faibussowitsch 20014d0ab18SJacob Faibussowitsch Vectors are initialized to zero before this function, so it is only needed for non 20114d0ab18SJacob Faibussowitsch homogeneous data. 20214d0ab18SJacob Faibussowitsch 20314d0ab18SJacob Faibussowitsch This function is somewhat optional: boundary data could potentially be inserted by a function 20414d0ab18SJacob Faibussowitsch passed to `DMTSSetIFunctionLocal()`. The use case for this function is for discretizations 20514d0ab18SJacob Faibussowitsch with constraints (see `DMGetDefaultConstraints()`): this function inserts boundary values 20614d0ab18SJacob Faibussowitsch before constraint interpolation. 2072fe279fdSBarry Smith 2081cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `TS`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 2097eb36faeSMatthew G. Knepley @*/ 210d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetBoundaryLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, void *), void *ctx) 211d71ae5a4SJacob Faibussowitsch { 2127eb36faeSMatthew G. Knepley DMTS tdm; 2137eb36faeSMatthew G. Knepley DMTS_Local *dmlocalts; 2147eb36faeSMatthew G. Knepley 2157eb36faeSMatthew G. Knepley PetscFunctionBegin; 2167eb36faeSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2179566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 2189566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 2197eb36faeSMatthew G. Knepley 2207eb36faeSMatthew G. Knepley dmlocalts->boundarylocal = func; 2217eb36faeSMatthew G. Knepley dmlocalts->boundarylocalctx = ctx; 2223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2237eb36faeSMatthew G. Knepley } 2247eb36faeSMatthew G. Knepley 225d433e6cbSMatthew G. Knepley /*@C 2260c9409e4SMatthew G. Knepley DMTSGetIFunctionLocal - get the local implicit function evaluation function. This function is called with local vector 2270c9409e4SMatthew G. Knepley containing the local vector information PLUS ghost point information. It should compute a result for all local 228bcf0153eSBarry Smith elements and `DM` will automatically accumulate the overlapping values. 2290c9409e4SMatthew G. Knepley 2300c9409e4SMatthew G. Knepley Logically Collective 2310c9409e4SMatthew G. Knepley 2320c9409e4SMatthew G. Knepley Input Parameter: 233bcf0153eSBarry Smith . dm - `DM` to associate callback with 2340c9409e4SMatthew G. Knepley 2350c9409e4SMatthew G. Knepley Output Parameters: 2360c9409e4SMatthew G. Knepley + func - local function evaluation 2370c9409e4SMatthew G. Knepley - ctx - context for function evaluation 2380c9409e4SMatthew G. Knepley 2390c9409e4SMatthew G. Knepley Level: beginner 2400c9409e4SMatthew G. Knepley 2411cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSSetIFunctionLocal()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 2420c9409e4SMatthew G. Knepley @*/ 243d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetIFunctionLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, Vec, void *), void **ctx) 244d71ae5a4SJacob Faibussowitsch { 2450c9409e4SMatthew G. Knepley DMTS tdm; 2460c9409e4SMatthew G. Knepley DMTS_Local *dmlocalts; 2470c9409e4SMatthew G. Knepley 2480c9409e4SMatthew G. Knepley PetscFunctionBegin; 2490c9409e4SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 250a86d55f5SMatthew G. Knepley PetscCall(DMGetDMTS(dm, &tdm)); 251a86d55f5SMatthew G. Knepley PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 2529371c9d4SSatish Balay if (func) { 2534f572ea9SToby Isaac PetscAssertPointer(func, 2); 2549371c9d4SSatish Balay *func = dmlocalts->ifunctionlocal; 2559371c9d4SSatish Balay } 2569371c9d4SSatish Balay if (ctx) { 2574f572ea9SToby Isaac PetscAssertPointer(ctx, 3); 2589371c9d4SSatish Balay *ctx = dmlocalts->ifunctionlocalctx; 2599371c9d4SSatish Balay } 2603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2610c9409e4SMatthew G. Knepley } 2620c9409e4SMatthew G. Knepley 2630c9409e4SMatthew G. Knepley /*@C 264d433e6cbSMatthew G. Knepley DMTSSetIFunctionLocal - set a local implicit function evaluation function. This function is called with local vector 265d433e6cbSMatthew G. Knepley containing the local vector information PLUS ghost point information. It should compute a result for all local 266bcf0153eSBarry Smith elements and `DM` will automatically accumulate the overlapping values. 267d433e6cbSMatthew G. Knepley 268d433e6cbSMatthew G. Knepley Logically Collective 269d433e6cbSMatthew G. Knepley 2704165533cSJose E. Roman Input Parameters: 271bcf0153eSBarry Smith + dm - `DM` to associate callback with 272d433e6cbSMatthew G. Knepley . func - local function evaluation 273d433e6cbSMatthew G. Knepley - ctx - context for function evaluation 274d433e6cbSMatthew G. Knepley 275d433e6cbSMatthew G. Knepley Level: beginner 276d433e6cbSMatthew G. Knepley 2771cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSGetIFunctionLocal()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 278d433e6cbSMatthew G. Knepley @*/ 279d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, Vec, void *), void *ctx) 280d71ae5a4SJacob Faibussowitsch { 281d433e6cbSMatthew G. Knepley DMTS tdm; 282d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts; 283d433e6cbSMatthew G. Knepley 284d433e6cbSMatthew G. Knepley PetscFunctionBegin; 285d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2869566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 2879566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 288d433e6cbSMatthew G. Knepley 289d433e6cbSMatthew G. Knepley dmlocalts->ifunctionlocal = func; 290d433e6cbSMatthew G. Knepley dmlocalts->ifunctionlocalctx = ctx; 291d433e6cbSMatthew G. Knepley 2929566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunction(dm, TSComputeIFunction_DMLocal, dmlocalts)); 293d433e6cbSMatthew G. Knepley if (!tdm->ops->ijacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 2949566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts)); 295d433e6cbSMatthew G. Knepley } 2963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 297d433e6cbSMatthew G. Knepley } 298d433e6cbSMatthew G. Knepley 299d433e6cbSMatthew G. Knepley /*@C 3000c9409e4SMatthew G. Knepley DMTSGetIJacobianLocal - get a local Jacobian evaluation function 3010c9409e4SMatthew G. Knepley 3020c9409e4SMatthew G. Knepley Logically Collective 3030c9409e4SMatthew G. Knepley 3040c9409e4SMatthew G. Knepley Input Parameter: 305bcf0153eSBarry Smith . dm - `DM` to associate callback with 3060c9409e4SMatthew G. Knepley 3070c9409e4SMatthew G. Knepley Output Parameters: 3080c9409e4SMatthew G. Knepley + func - local Jacobian evaluation 3090c9409e4SMatthew G. Knepley - ctx - optional context for local Jacobian evaluation 3100c9409e4SMatthew G. Knepley 3110c9409e4SMatthew G. Knepley Level: beginner 3120c9409e4SMatthew G. Knepley 3131cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSSetIJacobianLocal()`, `DMTSSetIFunctionLocal()`, `DMTSSetIJacobian()`, `DMTSSetIFunction()` 3140c9409e4SMatthew G. Knepley @*/ 315d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetIJacobianLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void **ctx) 316d71ae5a4SJacob Faibussowitsch { 3170c9409e4SMatthew G. Knepley DMTS tdm; 3180c9409e4SMatthew G. Knepley DMTS_Local *dmlocalts; 3190c9409e4SMatthew G. Knepley 3200c9409e4SMatthew G. Knepley PetscFunctionBegin; 3210c9409e4SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 322a86d55f5SMatthew G. Knepley PetscCall(DMGetDMTS(dm, &tdm)); 323a86d55f5SMatthew G. Knepley PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 3249371c9d4SSatish Balay if (func) { 3254f572ea9SToby Isaac PetscAssertPointer(func, 2); 3269371c9d4SSatish Balay *func = dmlocalts->ijacobianlocal; 3279371c9d4SSatish Balay } 3289371c9d4SSatish Balay if (ctx) { 3294f572ea9SToby Isaac PetscAssertPointer(ctx, 3); 3309371c9d4SSatish Balay *ctx = dmlocalts->ijacobianlocalctx; 3319371c9d4SSatish Balay } 3323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3330c9409e4SMatthew G. Knepley } 3340c9409e4SMatthew G. Knepley 3350c9409e4SMatthew G. Knepley /*@C 336d433e6cbSMatthew G. Knepley DMTSSetIJacobianLocal - set a local Jacobian evaluation function 337d433e6cbSMatthew G. Knepley 338d433e6cbSMatthew G. Knepley Logically Collective 339d433e6cbSMatthew G. Knepley 3404165533cSJose E. Roman Input Parameters: 341bcf0153eSBarry Smith + dm - `DM` to associate callback with 342d433e6cbSMatthew G. Knepley . func - local Jacobian evaluation 343d433e6cbSMatthew G. Knepley - ctx - optional context for local Jacobian evaluation 344d433e6cbSMatthew G. Knepley 345d433e6cbSMatthew G. Knepley Level: beginner 346d433e6cbSMatthew G. Knepley 3471cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSGetIJacobianLocal()`, `DMTSSetIFunctionLocal()`, `DMTSSetIJacobian()`, `DMTSSetIFunction()` 348d433e6cbSMatthew G. Knepley @*/ 349d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIJacobianLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *ctx) 350d71ae5a4SJacob Faibussowitsch { 351d433e6cbSMatthew G. Knepley DMTS tdm; 352d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts; 353d433e6cbSMatthew G. Knepley 354d433e6cbSMatthew G. Knepley PetscFunctionBegin; 355d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3569566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 3579566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 358d433e6cbSMatthew G. Knepley 359d433e6cbSMatthew G. Knepley dmlocalts->ijacobianlocal = func; 360d433e6cbSMatthew G. Knepley dmlocalts->ijacobianlocalctx = ctx; 361d433e6cbSMatthew G. Knepley 3629566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts)); 3633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 364d433e6cbSMatthew G. Knepley } 365d433e6cbSMatthew G. Knepley 366d433e6cbSMatthew G. Knepley /*@C 3670c9409e4SMatthew G. Knepley DMTSGetRHSFunctionLocal - get a local rhs function evaluation function. This function is called with local vector 3680c9409e4SMatthew G. Knepley containing the local vector information PLUS ghost point information. It should compute a result for all local 369bcf0153eSBarry Smith elements and `DM` will automatically accumulate the overlapping values. 3700c9409e4SMatthew G. Knepley 3710c9409e4SMatthew G. Knepley Logically Collective 3720c9409e4SMatthew G. Knepley 3730c9409e4SMatthew G. Knepley Input Parameter: 374bcf0153eSBarry Smith . dm - `DM` to associate callback with 3750c9409e4SMatthew G. Knepley 3760c9409e4SMatthew G. Knepley Output Parameters: 3770c9409e4SMatthew G. Knepley + func - local function evaluation 3780c9409e4SMatthew G. Knepley - ctx - context for function evaluation 3790c9409e4SMatthew G. Knepley 3800c9409e4SMatthew G. Knepley Level: beginner 3810c9409e4SMatthew G. Knepley 3821cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSSetRHSFunctionLocal()`, `DMTSSetRHSFunction()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 3830c9409e4SMatthew G. Knepley @*/ 384d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetRHSFunctionLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, void *), void **ctx) 385d71ae5a4SJacob Faibussowitsch { 3860c9409e4SMatthew G. Knepley DMTS tdm; 3870c9409e4SMatthew G. Knepley DMTS_Local *dmlocalts; 3880c9409e4SMatthew G. Knepley 3890c9409e4SMatthew G. Knepley PetscFunctionBegin; 3900c9409e4SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 391a86d55f5SMatthew G. Knepley PetscCall(DMGetDMTS(dm, &tdm)); 392a86d55f5SMatthew G. Knepley PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 3939371c9d4SSatish Balay if (func) { 3944f572ea9SToby Isaac PetscAssertPointer(func, 2); 3959371c9d4SSatish Balay *func = dmlocalts->rhsfunctionlocal; 3969371c9d4SSatish Balay } 3979371c9d4SSatish Balay if (ctx) { 3984f572ea9SToby Isaac PetscAssertPointer(ctx, 3); 3999371c9d4SSatish Balay *ctx = dmlocalts->rhsfunctionlocalctx; 4009371c9d4SSatish Balay } 4013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4020c9409e4SMatthew G. Knepley } 4030c9409e4SMatthew G. Knepley 4040c9409e4SMatthew G. Knepley /*@C 405d433e6cbSMatthew G. Knepley DMTSSetRHSFunctionLocal - set a local rhs function evaluation function. This function is called with local vector 406d433e6cbSMatthew G. Knepley containing the local vector information PLUS ghost point information. It should compute a result for all local 407bcf0153eSBarry Smith elements and `DM` will automatically accumulate the overlapping values. 408d433e6cbSMatthew G. Knepley 409d433e6cbSMatthew G. Knepley Logically Collective 410d433e6cbSMatthew G. Knepley 4114165533cSJose E. Roman Input Parameters: 412bcf0153eSBarry Smith + dm - `DM` to associate callback with 413d433e6cbSMatthew G. Knepley . func - local function evaluation 414d433e6cbSMatthew G. Knepley - ctx - context for function evaluation 415d433e6cbSMatthew G. Knepley 416d433e6cbSMatthew G. Knepley Level: beginner 417d433e6cbSMatthew G. Knepley 4181cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSGetRHSFunctionLocal()`, `DMTSSetRHSFunction()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 419d433e6cbSMatthew G. Knepley @*/ 420d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetRHSFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, void *), void *ctx) 421d71ae5a4SJacob Faibussowitsch { 422d433e6cbSMatthew G. Knepley DMTS tdm; 423d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts; 424d433e6cbSMatthew G. Knepley 425d433e6cbSMatthew G. Knepley PetscFunctionBegin; 426d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4279566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 4289566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 429d433e6cbSMatthew G. Knepley 430d433e6cbSMatthew G. Knepley dmlocalts->rhsfunctionlocal = func; 431d433e6cbSMatthew G. Knepley dmlocalts->rhsfunctionlocalctx = ctx; 432d433e6cbSMatthew G. Knepley 4339566063dSJacob Faibussowitsch PetscCall(DMTSSetRHSFunction(dm, TSComputeRHSFunction_DMLocal, dmlocalts)); 4343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 435d433e6cbSMatthew G. Knepley } 436b4937a87SMatthew G. Knepley 437cc4c1da9SBarry Smith /*@ 438bcf0153eSBarry 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. 439b4937a87SMatthew G. Knepley 440c3339decSBarry Smith Collective 441b4937a87SMatthew G. Knepley 4422fe279fdSBarry Smith Input Parameter: 443bcf0153eSBarry Smith . dm - `DM` providing the mass matrix 444b4937a87SMatthew G. Knepley 445b4937a87SMatthew G. Knepley Level: developer 446b4937a87SMatthew G. Knepley 447bcf0153eSBarry Smith Note: 448bcf0153eSBarry 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. 449bcf0153eSBarry Smith 4501cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSCreateRHSMassMatrixLumped()`, `DMTSDestroyRHSMassMatrix()`, `DMCreateMassMatrix()`, `DMTS` 451b4937a87SMatthew G. Knepley @*/ 452d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSCreateRHSMassMatrix(DM dm) 453d71ae5a4SJacob Faibussowitsch { 454b4937a87SMatthew G. Knepley DMTS tdm; 455b4937a87SMatthew G. Knepley DMTS_Local *dmlocalts; 456b4937a87SMatthew G. Knepley const char *prefix; 457b4937a87SMatthew G. Knepley 458b4937a87SMatthew G. Knepley PetscFunctionBegin; 459b4937a87SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4609566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 4619566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 4629566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(dm, dm, &dmlocalts->mass)); 4639566063dSJacob Faibussowitsch PetscCall(KSPCreate(PetscObjectComm((PetscObject)dm), &dmlocalts->kspmass)); 4649566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 4659566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(dmlocalts->kspmass, prefix)); 4669566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(dmlocalts->kspmass, "mass_")); 4679566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(dmlocalts->kspmass)); 4689566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(dmlocalts->kspmass, dmlocalts->mass, dmlocalts->mass)); 4693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 470b4937a87SMatthew G. Knepley } 471b4937a87SMatthew G. Knepley 472cc4c1da9SBarry Smith /*@ 473bcf0153eSBarry 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. 474b4937a87SMatthew G. Knepley 475c3339decSBarry Smith Collective 476b4937a87SMatthew G. Knepley 4772fe279fdSBarry Smith Input Parameter: 478bcf0153eSBarry Smith . dm - `DM` providing the mass matrix 479b4937a87SMatthew G. Knepley 480b4937a87SMatthew G. Knepley Level: developer 481b4937a87SMatthew G. Knepley 482bcf0153eSBarry Smith Note: 483bcf0153eSBarry 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. 484bcf0153eSBarry Smith Since the matrix is lumped, inversion is trivial. 485bcf0153eSBarry Smith 4861cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSCreateRHSMassMatrix()`, `DMTSDestroyRHSMassMatrix()`, `DMCreateMassMatrix()`, `DMTS` 487b4937a87SMatthew G. Knepley @*/ 488d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM dm) 489d71ae5a4SJacob Faibussowitsch { 490b4937a87SMatthew G. Knepley DMTS tdm; 491b4937a87SMatthew G. Knepley DMTS_Local *dmlocalts; 492b4937a87SMatthew G. Knepley 493b4937a87SMatthew G. Knepley PetscFunctionBegin; 494b4937a87SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4959566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 4969566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 497*8e9849d2SStefano Zampini PetscCall(DMCreateMassMatrixLumped(dm, NULL, &dmlocalts->lumpedmassinv)); 4989566063dSJacob Faibussowitsch PetscCall(VecReciprocal(dmlocalts->lumpedmassinv)); 4999566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(dmlocalts->lumpedmassinv, NULL, "-lumped_mass_inv_view")); 5003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 501b4937a87SMatthew G. Knepley } 502b4937a87SMatthew G. Knepley 503cc4c1da9SBarry Smith /*@ 504bcf0153eSBarry Smith DMTSDestroyRHSMassMatrix - Destroys the mass matrix and solver stored in the `DM` context, if they exist. 505b4937a87SMatthew G. Knepley 506b4937a87SMatthew G. Knepley Logically Collective 507b4937a87SMatthew G. Knepley 5082fe279fdSBarry Smith Input Parameter: 509bcf0153eSBarry Smith . dm - `DM` providing the mass matrix 510b4937a87SMatthew G. Knepley 511b4937a87SMatthew G. Knepley Level: developer 512b4937a87SMatthew G. Knepley 513b43aa488SJacob Faibussowitsch .seealso: [](ch_ts), `DM`, `DMTSCreateRHSMassMatrixLumped()`, `DMCreateMassMatrix()`, `DMTS` 514b4937a87SMatthew G. Knepley @*/ 515d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSDestroyRHSMassMatrix(DM dm) 516d71ae5a4SJacob Faibussowitsch { 517b4937a87SMatthew G. Knepley DMTS tdm; 518b4937a87SMatthew G. Knepley DMTS_Local *dmlocalts; 519b4937a87SMatthew G. Knepley 520b4937a87SMatthew G. Knepley PetscFunctionBegin; 521b4937a87SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5229566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 5239566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 5249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dmlocalts->lumpedmassinv)); 5259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&dmlocalts->mass)); 5269566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&dmlocalts->kspmass)); 5273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 528b4937a87SMatthew G. Knepley } 529