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; 222a40652d4SToby Isaac 2233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2247eb36faeSMatthew G. Knepley } 2257eb36faeSMatthew G. Knepley 226d433e6cbSMatthew G. Knepley /*@C 2270c9409e4SMatthew G. Knepley DMTSGetIFunctionLocal - get the local implicit function evaluation function. This function is called with local vector 2280c9409e4SMatthew G. Knepley containing the local vector information PLUS ghost point information. It should compute a result for all local 229bcf0153eSBarry Smith elements and `DM` will automatically accumulate the overlapping values. 2300c9409e4SMatthew G. Knepley 2310c9409e4SMatthew G. Knepley Logically Collective 2320c9409e4SMatthew G. Knepley 2330c9409e4SMatthew G. Knepley Input Parameter: 234bcf0153eSBarry Smith . dm - `DM` to associate callback with 2350c9409e4SMatthew G. Knepley 2360c9409e4SMatthew G. Knepley Output Parameters: 2370c9409e4SMatthew G. Knepley + func - local function evaluation 2380c9409e4SMatthew G. Knepley - ctx - context for function evaluation 2390c9409e4SMatthew G. Knepley 2400c9409e4SMatthew G. Knepley Level: beginner 2410c9409e4SMatthew G. Knepley 2421cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSSetIFunctionLocal()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 2430c9409e4SMatthew G. Knepley @*/ 244d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetIFunctionLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, Vec, void *), void **ctx) 245d71ae5a4SJacob Faibussowitsch { 2460c9409e4SMatthew G. Knepley DMTS tdm; 2470c9409e4SMatthew G. Knepley DMTS_Local *dmlocalts; 2480c9409e4SMatthew G. Knepley 2490c9409e4SMatthew G. Knepley PetscFunctionBegin; 2500c9409e4SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 251*a86d55f5SMatthew G. Knepley PetscCall(DMGetDMTS(dm, &tdm)); 252*a86d55f5SMatthew G. Knepley PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 2539371c9d4SSatish Balay if (func) { 2544f572ea9SToby Isaac PetscAssertPointer(func, 2); 2559371c9d4SSatish Balay *func = dmlocalts->ifunctionlocal; 2569371c9d4SSatish Balay } 2579371c9d4SSatish Balay if (ctx) { 2584f572ea9SToby Isaac PetscAssertPointer(ctx, 3); 2599371c9d4SSatish Balay *ctx = dmlocalts->ifunctionlocalctx; 2609371c9d4SSatish Balay } 2613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 2620c9409e4SMatthew G. Knepley } 2630c9409e4SMatthew G. Knepley 2640c9409e4SMatthew G. Knepley /*@C 265d433e6cbSMatthew G. Knepley DMTSSetIFunctionLocal - set a local implicit function evaluation function. This function is called with local vector 266d433e6cbSMatthew G. Knepley containing the local vector information PLUS ghost point information. It should compute a result for all local 267bcf0153eSBarry Smith elements and `DM` will automatically accumulate the overlapping values. 268d433e6cbSMatthew G. Knepley 269d433e6cbSMatthew G. Knepley Logically Collective 270d433e6cbSMatthew G. Knepley 2714165533cSJose E. Roman Input Parameters: 272bcf0153eSBarry Smith + dm - `DM` to associate callback with 273d433e6cbSMatthew G. Knepley . func - local function evaluation 274d433e6cbSMatthew G. Knepley - ctx - context for function evaluation 275d433e6cbSMatthew G. Knepley 276d433e6cbSMatthew G. Knepley Level: beginner 277d433e6cbSMatthew G. Knepley 2781cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSGetIFunctionLocal()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 279d433e6cbSMatthew G. Knepley @*/ 280d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, Vec, void *), void *ctx) 281d71ae5a4SJacob Faibussowitsch { 282d433e6cbSMatthew G. Knepley DMTS tdm; 283d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts; 284d433e6cbSMatthew G. Knepley 285d433e6cbSMatthew G. Knepley PetscFunctionBegin; 286d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 2879566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 2889566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 289d433e6cbSMatthew G. Knepley 290d433e6cbSMatthew G. Knepley dmlocalts->ifunctionlocal = func; 291d433e6cbSMatthew G. Knepley dmlocalts->ifunctionlocalctx = ctx; 292d433e6cbSMatthew G. Knepley 2939566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunction(dm, TSComputeIFunction_DMLocal, dmlocalts)); 294d433e6cbSMatthew G. Knepley if (!tdm->ops->ijacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 2959566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts)); 296d433e6cbSMatthew G. Knepley } 2973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 298d433e6cbSMatthew G. Knepley } 299d433e6cbSMatthew G. Knepley 300d433e6cbSMatthew G. Knepley /*@C 3010c9409e4SMatthew G. Knepley DMTSGetIJacobianLocal - get a local Jacobian evaluation function 3020c9409e4SMatthew G. Knepley 3030c9409e4SMatthew G. Knepley Logically Collective 3040c9409e4SMatthew G. Knepley 3050c9409e4SMatthew G. Knepley Input Parameter: 306bcf0153eSBarry Smith . dm - `DM` to associate callback with 3070c9409e4SMatthew G. Knepley 3080c9409e4SMatthew G. Knepley Output Parameters: 3090c9409e4SMatthew G. Knepley + func - local Jacobian evaluation 3100c9409e4SMatthew G. Knepley - ctx - optional context for local Jacobian evaluation 3110c9409e4SMatthew G. Knepley 3120c9409e4SMatthew G. Knepley Level: beginner 3130c9409e4SMatthew G. Knepley 3141cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSSetIJacobianLocal()`, `DMTSSetIFunctionLocal()`, `DMTSSetIJacobian()`, `DMTSSetIFunction()` 3150c9409e4SMatthew G. Knepley @*/ 316d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetIJacobianLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void **ctx) 317d71ae5a4SJacob Faibussowitsch { 3180c9409e4SMatthew G. Knepley DMTS tdm; 3190c9409e4SMatthew G. Knepley DMTS_Local *dmlocalts; 3200c9409e4SMatthew G. Knepley 3210c9409e4SMatthew G. Knepley PetscFunctionBegin; 3220c9409e4SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 323*a86d55f5SMatthew G. Knepley PetscCall(DMGetDMTS(dm, &tdm)); 324*a86d55f5SMatthew G. Knepley PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 3259371c9d4SSatish Balay if (func) { 3264f572ea9SToby Isaac PetscAssertPointer(func, 2); 3279371c9d4SSatish Balay *func = dmlocalts->ijacobianlocal; 3289371c9d4SSatish Balay } 3299371c9d4SSatish Balay if (ctx) { 3304f572ea9SToby Isaac PetscAssertPointer(ctx, 3); 3319371c9d4SSatish Balay *ctx = dmlocalts->ijacobianlocalctx; 3329371c9d4SSatish Balay } 3333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 3340c9409e4SMatthew G. Knepley } 3350c9409e4SMatthew G. Knepley 3360c9409e4SMatthew G. Knepley /*@C 337d433e6cbSMatthew G. Knepley DMTSSetIJacobianLocal - set a local Jacobian evaluation function 338d433e6cbSMatthew G. Knepley 339d433e6cbSMatthew G. Knepley Logically Collective 340d433e6cbSMatthew G. Knepley 3414165533cSJose E. Roman Input Parameters: 342bcf0153eSBarry Smith + dm - `DM` to associate callback with 343d433e6cbSMatthew G. Knepley . func - local Jacobian evaluation 344d433e6cbSMatthew G. Knepley - ctx - optional context for local Jacobian evaluation 345d433e6cbSMatthew G. Knepley 346d433e6cbSMatthew G. Knepley Level: beginner 347d433e6cbSMatthew G. Knepley 3481cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSGetIJacobianLocal()`, `DMTSSetIFunctionLocal()`, `DMTSSetIJacobian()`, `DMTSSetIFunction()` 349d433e6cbSMatthew G. Knepley @*/ 350d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetIJacobianLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *ctx) 351d71ae5a4SJacob Faibussowitsch { 352d433e6cbSMatthew G. Knepley DMTS tdm; 353d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts; 354d433e6cbSMatthew G. Knepley 355d433e6cbSMatthew G. Knepley PetscFunctionBegin; 356d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 3579566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 3589566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 359d433e6cbSMatthew G. Knepley 360d433e6cbSMatthew G. Knepley dmlocalts->ijacobianlocal = func; 361d433e6cbSMatthew G. Knepley dmlocalts->ijacobianlocalctx = ctx; 362d433e6cbSMatthew G. Knepley 3639566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts)); 3643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 365d433e6cbSMatthew G. Knepley } 366d433e6cbSMatthew G. Knepley 367d433e6cbSMatthew G. Knepley /*@C 3680c9409e4SMatthew G. Knepley DMTSGetRHSFunctionLocal - get a local rhs function evaluation function. This function is called with local vector 3690c9409e4SMatthew G. Knepley containing the local vector information PLUS ghost point information. It should compute a result for all local 370bcf0153eSBarry Smith elements and `DM` will automatically accumulate the overlapping values. 3710c9409e4SMatthew G. Knepley 3720c9409e4SMatthew G. Knepley Logically Collective 3730c9409e4SMatthew G. Knepley 3740c9409e4SMatthew G. Knepley Input Parameter: 375bcf0153eSBarry Smith . dm - `DM` to associate callback with 3760c9409e4SMatthew G. Knepley 3770c9409e4SMatthew G. Knepley Output Parameters: 3780c9409e4SMatthew G. Knepley + func - local function evaluation 3790c9409e4SMatthew G. Knepley - ctx - context for function evaluation 3800c9409e4SMatthew G. Knepley 3810c9409e4SMatthew G. Knepley Level: beginner 3820c9409e4SMatthew G. Knepley 3831cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSSetRHSFunctionLocal()`, `DMTSSetRHSFunction()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 3840c9409e4SMatthew G. Knepley @*/ 385d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSGetRHSFunctionLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, void *), void **ctx) 386d71ae5a4SJacob Faibussowitsch { 3870c9409e4SMatthew G. Knepley DMTS tdm; 3880c9409e4SMatthew G. Knepley DMTS_Local *dmlocalts; 3890c9409e4SMatthew G. Knepley 3900c9409e4SMatthew G. Knepley PetscFunctionBegin; 3910c9409e4SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 392*a86d55f5SMatthew G. Knepley PetscCall(DMGetDMTS(dm, &tdm)); 393*a86d55f5SMatthew G. Knepley PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 3949371c9d4SSatish Balay if (func) { 3954f572ea9SToby Isaac PetscAssertPointer(func, 2); 3969371c9d4SSatish Balay *func = dmlocalts->rhsfunctionlocal; 3979371c9d4SSatish Balay } 3989371c9d4SSatish Balay if (ctx) { 3994f572ea9SToby Isaac PetscAssertPointer(ctx, 3); 4009371c9d4SSatish Balay *ctx = dmlocalts->rhsfunctionlocalctx; 4019371c9d4SSatish Balay } 4023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4030c9409e4SMatthew G. Knepley } 4040c9409e4SMatthew G. Knepley 4050c9409e4SMatthew G. Knepley /*@C 406d433e6cbSMatthew G. Knepley DMTSSetRHSFunctionLocal - set a local rhs function evaluation function. This function is called with local vector 407d433e6cbSMatthew G. Knepley containing the local vector information PLUS ghost point information. It should compute a result for all local 408bcf0153eSBarry Smith elements and `DM` will automatically accumulate the overlapping values. 409d433e6cbSMatthew G. Knepley 410d433e6cbSMatthew G. Knepley Logically Collective 411d433e6cbSMatthew G. Knepley 4124165533cSJose E. Roman Input Parameters: 413bcf0153eSBarry Smith + dm - `DM` to associate callback with 414d433e6cbSMatthew G. Knepley . func - local function evaluation 415d433e6cbSMatthew G. Knepley - ctx - context for function evaluation 416d433e6cbSMatthew G. Knepley 417d433e6cbSMatthew G. Knepley Level: beginner 418d433e6cbSMatthew G. Knepley 4191cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSGetRHSFunctionLocal()`, `DMTSSetRHSFunction()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 420d433e6cbSMatthew G. Knepley @*/ 421d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSSetRHSFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, void *), void *ctx) 422d71ae5a4SJacob Faibussowitsch { 423d433e6cbSMatthew G. Knepley DMTS tdm; 424d433e6cbSMatthew G. Knepley DMTS_Local *dmlocalts; 425d433e6cbSMatthew G. Knepley 426d433e6cbSMatthew G. Knepley PetscFunctionBegin; 427d433e6cbSMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4289566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 4299566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 430d433e6cbSMatthew G. Knepley 431d433e6cbSMatthew G. Knepley dmlocalts->rhsfunctionlocal = func; 432d433e6cbSMatthew G. Knepley dmlocalts->rhsfunctionlocalctx = ctx; 433d433e6cbSMatthew G. Knepley 4349566063dSJacob Faibussowitsch PetscCall(DMTSSetRHSFunction(dm, TSComputeRHSFunction_DMLocal, dmlocalts)); 4353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 436d433e6cbSMatthew G. Knepley } 437b4937a87SMatthew G. Knepley 438b4937a87SMatthew G. Knepley /*@C 439bcf0153eSBarry 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. 440b4937a87SMatthew G. Knepley 441c3339decSBarry Smith Collective 442b4937a87SMatthew G. Knepley 4432fe279fdSBarry Smith Input Parameter: 444bcf0153eSBarry Smith . dm - `DM` providing the mass matrix 445b4937a87SMatthew G. Knepley 446b4937a87SMatthew G. Knepley Level: developer 447b4937a87SMatthew G. Knepley 448bcf0153eSBarry Smith Note: 449bcf0153eSBarry 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. 450bcf0153eSBarry Smith 4511cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSCreateRHSMassMatrixLumped()`, `DMTSDestroyRHSMassMatrix()`, `DMCreateMassMatrix()`, `DMTS` 452b4937a87SMatthew G. Knepley @*/ 453d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSCreateRHSMassMatrix(DM dm) 454d71ae5a4SJacob Faibussowitsch { 455b4937a87SMatthew G. Knepley DMTS tdm; 456b4937a87SMatthew G. Knepley DMTS_Local *dmlocalts; 457b4937a87SMatthew G. Knepley const char *prefix; 458b4937a87SMatthew G. Knepley 459b4937a87SMatthew G. Knepley PetscFunctionBegin; 460b4937a87SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4619566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 4629566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 4639566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrix(dm, dm, &dmlocalts->mass)); 4649566063dSJacob Faibussowitsch PetscCall(KSPCreate(PetscObjectComm((PetscObject)dm), &dmlocalts->kspmass)); 4659566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 4669566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(dmlocalts->kspmass, prefix)); 4679566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(dmlocalts->kspmass, "mass_")); 4689566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(dmlocalts->kspmass)); 4699566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(dmlocalts->kspmass, dmlocalts->mass, dmlocalts->mass)); 4703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 471b4937a87SMatthew G. Knepley } 472b4937a87SMatthew G. Knepley 473b4937a87SMatthew G. Knepley /*@C 474bcf0153eSBarry 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. 475b4937a87SMatthew G. Knepley 476c3339decSBarry Smith Collective 477b4937a87SMatthew G. Knepley 4782fe279fdSBarry Smith Input Parameter: 479bcf0153eSBarry Smith . dm - `DM` providing the mass matrix 480b4937a87SMatthew G. Knepley 481b4937a87SMatthew G. Knepley Level: developer 482b4937a87SMatthew G. Knepley 483bcf0153eSBarry Smith Note: 484bcf0153eSBarry 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. 485bcf0153eSBarry Smith Since the matrix is lumped, inversion is trivial. 486bcf0153eSBarry Smith 4871cc06b55SBarry Smith .seealso: [](ch_ts), `DM`, `DMTSCreateRHSMassMatrix()`, `DMTSDestroyRHSMassMatrix()`, `DMCreateMassMatrix()`, `DMTS` 488b4937a87SMatthew G. Knepley @*/ 489d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM dm) 490d71ae5a4SJacob Faibussowitsch { 491b4937a87SMatthew G. Knepley DMTS tdm; 492b4937a87SMatthew G. Knepley DMTS_Local *dmlocalts; 493b4937a87SMatthew G. Knepley 494b4937a87SMatthew G. Knepley PetscFunctionBegin; 495b4937a87SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 4969566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 4979566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 4989566063dSJacob Faibussowitsch PetscCall(DMCreateMassMatrixLumped(dm, &dmlocalts->lumpedmassinv)); 4999566063dSJacob Faibussowitsch PetscCall(VecReciprocal(dmlocalts->lumpedmassinv)); 5009566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(dmlocalts->lumpedmassinv, NULL, "-lumped_mass_inv_view")); 5013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 502b4937a87SMatthew G. Knepley } 503b4937a87SMatthew G. Knepley 504b4937a87SMatthew G. Knepley /*@C 505bcf0153eSBarry Smith DMTSDestroyRHSMassMatrix - Destroys the mass matrix and solver stored in the `DM` context, if they exist. 506b4937a87SMatthew G. Knepley 507b4937a87SMatthew G. Knepley Logically Collective 508b4937a87SMatthew G. Knepley 5092fe279fdSBarry Smith Input Parameter: 510bcf0153eSBarry Smith . dm - `DM` providing the mass matrix 511b4937a87SMatthew G. Knepley 512b4937a87SMatthew G. Knepley Level: developer 513b4937a87SMatthew G. Knepley 514b43aa488SJacob Faibussowitsch .seealso: [](ch_ts), `DM`, `DMTSCreateRHSMassMatrixLumped()`, `DMCreateMassMatrix()`, `DMTS` 515b4937a87SMatthew G. Knepley @*/ 516d71ae5a4SJacob Faibussowitsch PetscErrorCode DMTSDestroyRHSMassMatrix(DM dm) 517d71ae5a4SJacob Faibussowitsch { 518b4937a87SMatthew G. Knepley DMTS tdm; 519b4937a87SMatthew G. Knepley DMTS_Local *dmlocalts; 520b4937a87SMatthew G. Knepley 521b4937a87SMatthew G. Knepley PetscFunctionBegin; 522b4937a87SMatthew G. Knepley PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 5239566063dSJacob Faibussowitsch PetscCall(DMGetDMTSWrite(dm, &tdm)); 5249566063dSJacob Faibussowitsch PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 5259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&dmlocalts->lumpedmassinv)); 5269566063dSJacob Faibussowitsch PetscCall(MatDestroy(&dmlocalts->mass)); 5279566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&dmlocalts->kspmass)); 5283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 529b4937a87SMatthew G. Knepley } 530