1 #include <petsc/private/dmimpl.h> 2 #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 3 4 typedef struct { 5 PetscErrorCode (*boundarylocal)(DM, PetscReal, Vec, Vec, void *); 6 PetscErrorCode (*ifunctionlocal)(DM, PetscReal, Vec, Vec, Vec, void *); 7 PetscErrorCode (*ijacobianlocal)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *); 8 PetscErrorCode (*rhsfunctionlocal)(DM, PetscReal, Vec, Vec, void *); 9 void *boundarylocalctx; 10 void *ifunctionlocalctx; 11 void *ijacobianlocalctx; 12 void *rhsfunctionlocalctx; 13 Vec lumpedmassinv; 14 Mat mass; 15 KSP kspmass; 16 } DMTS_Local; 17 18 static PetscErrorCode DMTSDestroy_DMLocal(DMTS tdm) 19 { 20 PetscFunctionBegin; 21 PetscCall(PetscFree(tdm->data)); 22 PetscFunctionReturn(PETSC_SUCCESS); 23 } 24 25 static PetscErrorCode DMTSDuplicate_DMLocal(DMTS oldtdm, DMTS tdm) 26 { 27 PetscFunctionBegin; 28 PetscCall(PetscNew((DMTS_Local **)&tdm->data)); 29 if (oldtdm->data) PetscCall(PetscMemcpy(tdm->data, oldtdm->data, sizeof(DMTS_Local))); 30 PetscFunctionReturn(PETSC_SUCCESS); 31 } 32 33 static PetscErrorCode DMLocalTSGetContext(DM dm, DMTS tdm, DMTS_Local **dmlocalts) 34 { 35 PetscFunctionBegin; 36 *dmlocalts = NULL; 37 if (!tdm->data) { 38 PetscCall(PetscNew((DMTS_Local **)&tdm->data)); 39 40 tdm->ops->destroy = DMTSDestroy_DMLocal; 41 tdm->ops->duplicate = DMTSDuplicate_DMLocal; 42 } 43 *dmlocalts = (DMTS_Local *)tdm->data; 44 PetscFunctionReturn(PETSC_SUCCESS); 45 } 46 47 static PetscErrorCode TSComputeIFunction_DMLocal(TS ts, PetscReal time, Vec X, Vec X_t, Vec F, void *ctx) 48 { 49 DM dm; 50 Vec locX, locX_t, locF; 51 DMTS_Local *dmlocalts = (DMTS_Local *)ctx; 52 53 PetscFunctionBegin; 54 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 55 PetscValidHeaderSpecific(X, VEC_CLASSID, 3); 56 PetscValidHeaderSpecific(X_t, VEC_CLASSID, 4); 57 PetscValidHeaderSpecific(F, VEC_CLASSID, 5); 58 PetscCall(TSGetDM(ts, &dm)); 59 PetscCall(DMGetLocalVector(dm, &locX)); 60 PetscCall(DMGetLocalVector(dm, &locX_t)); 61 PetscCall(DMGetLocalVector(dm, &locF)); 62 PetscCall(VecZeroEntries(locX)); 63 PetscCall(VecZeroEntries(locX_t)); 64 if (dmlocalts->boundarylocal) PetscCall((*dmlocalts->boundarylocal)(dm, time, locX, locX_t, dmlocalts->boundarylocalctx)); 65 PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX)); 66 PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX)); 67 PetscCall(DMGlobalToLocalBegin(dm, X_t, INSERT_VALUES, locX_t)); 68 PetscCall(DMGlobalToLocalEnd(dm, X_t, INSERT_VALUES, locX_t)); 69 PetscCall(VecZeroEntries(locF)); 70 CHKMEMQ; 71 PetscCall((*dmlocalts->ifunctionlocal)(dm, time, locX, locX_t, locF, dmlocalts->ifunctionlocalctx)); 72 CHKMEMQ; 73 PetscCall(VecZeroEntries(F)); 74 PetscCall(DMLocalToGlobalBegin(dm, locF, ADD_VALUES, F)); 75 PetscCall(DMLocalToGlobalEnd(dm, locF, ADD_VALUES, F)); 76 PetscCall(DMRestoreLocalVector(dm, &locX)); 77 PetscCall(DMRestoreLocalVector(dm, &locX_t)); 78 PetscCall(DMRestoreLocalVector(dm, &locF)); 79 PetscFunctionReturn(PETSC_SUCCESS); 80 } 81 82 static PetscErrorCode TSComputeRHSFunction_DMLocal(TS ts, PetscReal time, Vec X, Vec F, void *ctx) 83 { 84 DM dm; 85 Vec locX, locF; 86 DMTS_Local *dmlocalts = (DMTS_Local *)ctx; 87 88 PetscFunctionBegin; 89 PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 90 PetscValidHeaderSpecific(X, VEC_CLASSID, 3); 91 PetscValidHeaderSpecific(F, VEC_CLASSID, 4); 92 PetscCall(TSGetDM(ts, &dm)); 93 PetscCall(DMGetLocalVector(dm, &locX)); 94 PetscCall(DMGetLocalVector(dm, &locF)); 95 PetscCall(VecZeroEntries(locX)); 96 if (dmlocalts->boundarylocal) PetscCall((*dmlocalts->boundarylocal)(dm, time, locX, NULL, dmlocalts->boundarylocalctx)); 97 PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX)); 98 PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX)); 99 PetscCall(VecZeroEntries(locF)); 100 CHKMEMQ; 101 PetscCall((*dmlocalts->rhsfunctionlocal)(dm, time, locX, locF, dmlocalts->rhsfunctionlocalctx)); 102 CHKMEMQ; 103 PetscCall(VecZeroEntries(F)); 104 PetscCall(DMLocalToGlobalBegin(dm, locF, ADD_VALUES, F)); 105 PetscCall(DMLocalToGlobalEnd(dm, locF, ADD_VALUES, F)); 106 if (dmlocalts->lumpedmassinv) { 107 PetscCall(VecPointwiseMult(F, dmlocalts->lumpedmassinv, F)); 108 } else if (dmlocalts->kspmass) { 109 Vec tmp; 110 111 PetscCall(DMGetGlobalVector(dm, &tmp)); 112 PetscCall(KSPSolve(dmlocalts->kspmass, F, tmp)); 113 PetscCall(VecCopy(tmp, F)); 114 PetscCall(DMRestoreGlobalVector(dm, &tmp)); 115 } 116 PetscCall(DMRestoreLocalVector(dm, &locX)); 117 PetscCall(DMRestoreLocalVector(dm, &locF)); 118 PetscFunctionReturn(PETSC_SUCCESS); 119 } 120 121 static PetscErrorCode TSComputeIJacobian_DMLocal(TS ts, PetscReal time, Vec X, Vec X_t, PetscReal a, Mat A, Mat B, void *ctx) 122 { 123 DM dm; 124 Vec locX, locX_t; 125 DMTS_Local *dmlocalts = (DMTS_Local *)ctx; 126 127 PetscFunctionBegin; 128 PetscCall(TSGetDM(ts, &dm)); 129 if (dmlocalts->ijacobianlocal) { 130 PetscCall(DMGetLocalVector(dm, &locX)); 131 PetscCall(DMGetLocalVector(dm, &locX_t)); 132 PetscCall(VecZeroEntries(locX)); 133 PetscCall(VecZeroEntries(locX_t)); 134 if (dmlocalts->boundarylocal) PetscCall((*dmlocalts->boundarylocal)(dm, time, locX, locX_t, dmlocalts->boundarylocalctx)); 135 PetscCall(DMGlobalToLocalBegin(dm, X, INSERT_VALUES, locX)); 136 PetscCall(DMGlobalToLocalEnd(dm, X, INSERT_VALUES, locX)); 137 PetscCall(DMGlobalToLocalBegin(dm, X_t, INSERT_VALUES, locX_t)); 138 PetscCall(DMGlobalToLocalEnd(dm, X_t, INSERT_VALUES, locX_t)); 139 CHKMEMQ; 140 PetscCall((*dmlocalts->ijacobianlocal)(dm, time, locX, locX_t, a, A, B, dmlocalts->ijacobianlocalctx)); 141 CHKMEMQ; 142 PetscCall(DMRestoreLocalVector(dm, &locX)); 143 PetscCall(DMRestoreLocalVector(dm, &locX_t)); 144 } else { 145 MatFDColoring fdcoloring; 146 PetscCall(PetscObjectQuery((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject *)&fdcoloring)); 147 if (!fdcoloring) { 148 ISColoring coloring; 149 150 PetscCall(DMCreateColoring(dm, dm->coloringtype, &coloring)); 151 PetscCall(MatFDColoringCreate(B, coloring, &fdcoloring)); 152 PetscCall(ISColoringDestroy(&coloring)); 153 switch (dm->coloringtype) { 154 case IS_COLORING_GLOBAL: 155 PetscCall(MatFDColoringSetFunction(fdcoloring, (PetscErrorCode(*)(void))TSComputeIFunction_DMLocal, dmlocalts)); 156 break; 157 default: 158 SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "No support for coloring type '%s'", ISColoringTypes[dm->coloringtype]); 159 } 160 PetscCall(PetscObjectSetOptionsPrefix((PetscObject)fdcoloring, ((PetscObject)dm)->prefix)); 161 PetscCall(MatFDColoringSetFromOptions(fdcoloring)); 162 PetscCall(MatFDColoringSetUp(B, coloring, fdcoloring)); 163 PetscCall(PetscObjectCompose((PetscObject)dm, "DMDASNES_FDCOLORING", (PetscObject)fdcoloring)); 164 PetscCall(PetscObjectDereference((PetscObject)fdcoloring)); 165 166 /* The following breaks an ugly reference counting loop that deserves a paragraph. MatFDColoringApply() will call 167 * VecDuplicate() with the state Vec and store inside the MatFDColoring. This Vec will duplicate the Vec, but the 168 * MatFDColoring is composed with the DM. We dereference the DM here so that the reference count will eventually 169 * drop to 0. Note the code in DMDestroy() that exits early for a negative reference count. That code path will be 170 * taken when the PetscObjectList for the Vec inside MatFDColoring is destroyed. 171 */ 172 PetscCall(PetscObjectDereference((PetscObject)dm)); 173 } 174 PetscCall(MatFDColoringApply(B, fdcoloring, X, ts)); 175 } 176 /* This will be redundant if the user called both, but it's too common to forget. */ 177 if (A != B) { 178 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 179 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 180 } 181 PetscFunctionReturn(PETSC_SUCCESS); 182 } 183 184 /*@C 185 DMTSSetBoundaryLocal - set the function for essential boundary data for a local implicit function evaluation. 186 It should set the essential boundary data for the local portion of the solution X, as well its time derivative X_t (if it is not NULL). 187 Vectors are initialized to zero before this function, so it is only needed for non homogeneous data. 188 189 Note that this function is somewhat optional: boundary data could potentially be inserted by a function passed to 190 `DMTSSetIFunctionLocal()`. The use case for this function is for discretizations with constraints (see 191 `DMGetDefaultConstraints()`): this function inserts boundary values before constraint interpolation. 192 193 Logically Collective 194 195 Input Parameters: 196 + dm - `DM` to associate callback with 197 . func - local function evaluation 198 - ctx - context for function evaluation 199 200 Level: intermediate 201 202 .seealso: [](chapter_ts), `DM`, `TS`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 203 @*/ 204 PetscErrorCode DMTSSetBoundaryLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, void *), void *ctx) 205 { 206 DMTS tdm; 207 DMTS_Local *dmlocalts; 208 209 PetscFunctionBegin; 210 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 211 PetscCall(DMGetDMTSWrite(dm, &tdm)); 212 PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 213 214 dmlocalts->boundarylocal = func; 215 dmlocalts->boundarylocalctx = ctx; 216 217 PetscFunctionReturn(PETSC_SUCCESS); 218 } 219 220 /*@C 221 DMTSGetIFunctionLocal - get the local implicit function evaluation function. This function is called with local vector 222 containing the local vector information PLUS ghost point information. It should compute a result for all local 223 elements and `DM` will automatically accumulate the overlapping values. 224 225 Logically Collective 226 227 Input Parameter: 228 . dm - `DM` to associate callback with 229 230 Output Parameters: 231 + func - local function evaluation 232 - ctx - context for function evaluation 233 234 Level: beginner 235 236 .seealso: [](chapter_ts), `DM`, `DMTSSetIFunctionLocal()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 237 @*/ 238 PetscErrorCode DMTSGetIFunctionLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, Vec, void *), void **ctx) 239 { 240 DMTS tdm; 241 DMTS_Local *dmlocalts; 242 PetscErrorCode ierr; 243 244 PetscFunctionBegin; 245 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 246 ierr = DMGetDMTS(dm, &tdm); 247 CHKERRQ(ierr); 248 ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts); 249 CHKERRQ(ierr); 250 if (func) { 251 PetscValidPointer(func, 2); 252 *func = dmlocalts->ifunctionlocal; 253 } 254 if (ctx) { 255 PetscValidPointer(ctx, 3); 256 *ctx = dmlocalts->ifunctionlocalctx; 257 } 258 PetscFunctionReturn(PETSC_SUCCESS); 259 } 260 261 /*@C 262 DMTSSetIFunctionLocal - set a local implicit function evaluation function. This function is called with local vector 263 containing the local vector information PLUS ghost point information. It should compute a result for all local 264 elements and `DM` will automatically accumulate the overlapping values. 265 266 Logically Collective 267 268 Input Parameters: 269 + dm - `DM` to associate callback with 270 . func - local function evaluation 271 - ctx - context for function evaluation 272 273 Level: beginner 274 275 .seealso: [](chapter_ts), `DM`, `DMTSGetIFunctionLocal()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 276 @*/ 277 PetscErrorCode DMTSSetIFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, Vec, void *), void *ctx) 278 { 279 DMTS tdm; 280 DMTS_Local *dmlocalts; 281 282 PetscFunctionBegin; 283 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 284 PetscCall(DMGetDMTSWrite(dm, &tdm)); 285 PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 286 287 dmlocalts->ifunctionlocal = func; 288 dmlocalts->ifunctionlocalctx = ctx; 289 290 PetscCall(DMTSSetIFunction(dm, TSComputeIFunction_DMLocal, dmlocalts)); 291 if (!tdm->ops->ijacobian) { /* Call us for the Jacobian too, can be overridden by the user. */ 292 PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts)); 293 } 294 PetscFunctionReturn(PETSC_SUCCESS); 295 } 296 297 /*@C 298 DMTSGetIJacobianLocal - get a local Jacobian evaluation function 299 300 Logically Collective 301 302 Input Parameter: 303 . dm - `DM` to associate callback with 304 305 Output Parameters: 306 + func - local Jacobian evaluation 307 - ctx - optional context for local Jacobian evaluation 308 309 Level: beginner 310 311 .seealso: [](chapter_ts), `DM`, `DMTSSetIJacobianLocal()`, `DMTSSetIFunctionLocal()`, `DMTSSetIJacobian()`, `DMTSSetIFunction()` 312 @*/ 313 PetscErrorCode DMTSGetIJacobianLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void **ctx) 314 { 315 DMTS tdm; 316 DMTS_Local *dmlocalts; 317 PetscErrorCode ierr; 318 319 PetscFunctionBegin; 320 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 321 ierr = DMGetDMTS(dm, &tdm); 322 CHKERRQ(ierr); 323 ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts); 324 CHKERRQ(ierr); 325 if (func) { 326 PetscValidPointer(func, 2); 327 *func = dmlocalts->ijacobianlocal; 328 } 329 if (ctx) { 330 PetscValidPointer(ctx, 3); 331 *ctx = dmlocalts->ijacobianlocalctx; 332 } 333 PetscFunctionReturn(PETSC_SUCCESS); 334 } 335 336 /*@C 337 DMTSSetIJacobianLocal - set a local Jacobian evaluation function 338 339 Logically Collective 340 341 Input Parameters: 342 + dm - `DM` to associate callback with 343 . func - local Jacobian evaluation 344 - ctx - optional context for local Jacobian evaluation 345 346 Level: beginner 347 348 .seealso: [](chapter_ts), `DM`, `DMTSGetIJacobianLocal()`, `DMTSSetIFunctionLocal()`, `DMTSSetIJacobian()`, `DMTSSetIFunction()` 349 @*/ 350 PetscErrorCode DMTSSetIJacobianLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *), void *ctx) 351 { 352 DMTS tdm; 353 DMTS_Local *dmlocalts; 354 355 PetscFunctionBegin; 356 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 357 PetscCall(DMGetDMTSWrite(dm, &tdm)); 358 PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 359 360 dmlocalts->ijacobianlocal = func; 361 dmlocalts->ijacobianlocalctx = ctx; 362 363 PetscCall(DMTSSetIJacobian(dm, TSComputeIJacobian_DMLocal, dmlocalts)); 364 PetscFunctionReturn(PETSC_SUCCESS); 365 } 366 367 /*@C 368 DMTSGetRHSFunctionLocal - get a local rhs function evaluation function. This function is called with local vector 369 containing the local vector information PLUS ghost point information. It should compute a result for all local 370 elements and `DM` will automatically accumulate the overlapping values. 371 372 Logically Collective 373 374 Input Parameter: 375 . dm - `DM` to associate callback with 376 377 Output Parameters: 378 + func - local function evaluation 379 - ctx - context for function evaluation 380 381 Level: beginner 382 383 .seealso: [](chapter_ts), `DM`, `DMTSSetRHSFunctionLocal()`, `DMTSSetRHSFunction()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 384 @*/ 385 PetscErrorCode DMTSGetRHSFunctionLocal(DM dm, PetscErrorCode (**func)(DM, PetscReal, Vec, Vec, void *), void **ctx) 386 { 387 DMTS tdm; 388 DMTS_Local *dmlocalts; 389 PetscErrorCode ierr; 390 391 PetscFunctionBegin; 392 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 393 ierr = DMGetDMTS(dm, &tdm); 394 CHKERRQ(ierr); 395 ierr = DMLocalTSGetContext(dm, tdm, &dmlocalts); 396 CHKERRQ(ierr); 397 if (func) { 398 PetscValidPointer(func, 2); 399 *func = dmlocalts->rhsfunctionlocal; 400 } 401 if (ctx) { 402 PetscValidPointer(ctx, 3); 403 *ctx = dmlocalts->rhsfunctionlocalctx; 404 } 405 PetscFunctionReturn(PETSC_SUCCESS); 406 } 407 408 /*@C 409 DMTSSetRHSFunctionLocal - set a local rhs function evaluation function. This function is called with local vector 410 containing the local vector information PLUS ghost point information. It should compute a result for all local 411 elements and `DM` will automatically accumulate the overlapping values. 412 413 Logically Collective 414 415 Input Parameters: 416 + dm - `DM` to associate callback with 417 . func - local function evaluation 418 - ctx - context for function evaluation 419 420 Level: beginner 421 422 .seealso: [](chapter_ts), `DM`, `DMTSGetRHSFunctionLocal()`, `DMTSSetRHSFunction()`, `DMTSSetIFunction()`, `DMTSSetIJacobianLocal()` 423 @*/ 424 PetscErrorCode DMTSSetRHSFunctionLocal(DM dm, PetscErrorCode (*func)(DM, PetscReal, Vec, Vec, void *), void *ctx) 425 { 426 DMTS tdm; 427 DMTS_Local *dmlocalts; 428 429 PetscFunctionBegin; 430 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 431 PetscCall(DMGetDMTSWrite(dm, &tdm)); 432 PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 433 434 dmlocalts->rhsfunctionlocal = func; 435 dmlocalts->rhsfunctionlocalctx = ctx; 436 437 PetscCall(DMTSSetRHSFunction(dm, TSComputeRHSFunction_DMLocal, dmlocalts)); 438 PetscFunctionReturn(PETSC_SUCCESS); 439 } 440 441 /*@C 442 DMTSCreateRHSMassMatrix - This creates the mass matrix associated with the given `DM`, and a solver to invert it, and stores them in the `DM` context. 443 444 Collective 445 446 Input Parameters: 447 . dm - `DM` providing the mass matrix 448 449 Level: developer 450 451 Note: 452 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. 453 454 .seealso: [](chapter_ts), `DM`, `DMTSCreateRHSMassMatrixLumped()`, `DMTSDestroyRHSMassMatrix()`, `DMCreateMassMatrix()`, `DMTS` 455 @*/ 456 PetscErrorCode DMTSCreateRHSMassMatrix(DM dm) 457 { 458 DMTS tdm; 459 DMTS_Local *dmlocalts; 460 const char *prefix; 461 462 PetscFunctionBegin; 463 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 464 PetscCall(DMGetDMTSWrite(dm, &tdm)); 465 PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 466 PetscCall(DMCreateMassMatrix(dm, dm, &dmlocalts->mass)); 467 PetscCall(KSPCreate(PetscObjectComm((PetscObject)dm), &dmlocalts->kspmass)); 468 PetscCall(PetscObjectGetOptionsPrefix((PetscObject)dm, &prefix)); 469 PetscCall(KSPSetOptionsPrefix(dmlocalts->kspmass, prefix)); 470 PetscCall(KSPAppendOptionsPrefix(dmlocalts->kspmass, "mass_")); 471 PetscCall(KSPSetFromOptions(dmlocalts->kspmass)); 472 PetscCall(KSPSetOperators(dmlocalts->kspmass, dmlocalts->mass, dmlocalts->mass)); 473 PetscFunctionReturn(PETSC_SUCCESS); 474 } 475 476 /*@C 477 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. 478 479 Collective 480 481 Input Parameters: 482 . dm - `DM` providing the mass matrix 483 484 Level: developer 485 486 Note: 487 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. 488 Since the matrix is lumped, inversion is trivial. 489 490 .seealso: [](chapter_ts), `DM`, `DMTSCreateRHSMassMatrix()`, `DMTSDestroyRHSMassMatrix()`, `DMCreateMassMatrix()`, `DMTS` 491 @*/ 492 PetscErrorCode DMTSCreateRHSMassMatrixLumped(DM dm) 493 { 494 DMTS tdm; 495 DMTS_Local *dmlocalts; 496 497 PetscFunctionBegin; 498 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 499 PetscCall(DMGetDMTSWrite(dm, &tdm)); 500 PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 501 PetscCall(DMCreateMassMatrixLumped(dm, &dmlocalts->lumpedmassinv)); 502 PetscCall(VecReciprocal(dmlocalts->lumpedmassinv)); 503 PetscCall(VecViewFromOptions(dmlocalts->lumpedmassinv, NULL, "-lumped_mass_inv_view")); 504 PetscFunctionReturn(PETSC_SUCCESS); 505 } 506 507 /*@C 508 DMTSDestroyRHSMassMatrix - Destroys the mass matrix and solver stored in the `DM` context, if they exist. 509 510 Logically Collective 511 512 Input Parameters: 513 . dm - `DM` providing the mass matrix 514 515 Level: developer 516 517 .seealso: [](chapter_ts), `DM`, `DMTSCreateRHSMassMatrixLumped()`, `DMCreateMassMatrix()`, `DMCreateMassMatrix()`, `DMTS` 518 @*/ 519 PetscErrorCode DMTSDestroyRHSMassMatrix(DM dm) 520 { 521 DMTS tdm; 522 DMTS_Local *dmlocalts; 523 524 PetscFunctionBegin; 525 PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 526 PetscCall(DMGetDMTSWrite(dm, &tdm)); 527 PetscCall(DMLocalTSGetContext(dm, tdm, &dmlocalts)); 528 PetscCall(VecDestroy(&dmlocalts->lumpedmassinv)); 529 PetscCall(MatDestroy(&dmlocalts->mass)); 530 PetscCall(KSPDestroy(&dmlocalts->kspmass)); 531 PetscFunctionReturn(PETSC_SUCCESS); 532 } 533