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