1 #include <petsc/private/petscfvimpl.h> /*I "petscfv.h" I*/ 2 #include <petscdmplex.h> 3 #include <petscdmplextransform.h> 4 #include <petscds.h> 5 6 PetscClassId PETSCLIMITER_CLASSID = 0; 7 8 PetscFunctionList PetscLimiterList = NULL; 9 PetscBool PetscLimiterRegisterAllCalled = PETSC_FALSE; 10 11 PetscBool Limitercite = PETSC_FALSE; 12 const char LimiterCitation[] = "@article{BergerAftosmisMurman2005,\n" 13 " title = {Analysis of slope limiters on irregular grids},\n" 14 " journal = {AIAA paper},\n" 15 " author = {Marsha Berger and Michael J. Aftosmis and Scott M. Murman},\n" 16 " volume = {490},\n" 17 " year = {2005}\n}\n"; 18 19 /*@C 20 PetscLimiterRegister - Adds a new PetscLimiter implementation 21 22 Not Collective 23 24 Input Parameters: 25 + name - The name of a new user-defined creation routine 26 - create_func - The creation routine itself 27 28 Notes: 29 PetscLimiterRegister() may be called multiple times to add several user-defined PetscLimiters 30 31 Sample usage: 32 .vb 33 PetscLimiterRegister("my_lim", MyPetscLimiterCreate); 34 .ve 35 36 Then, your PetscLimiter type can be chosen with the procedural interface via 37 .vb 38 PetscLimiterCreate(MPI_Comm, PetscLimiter *); 39 PetscLimiterSetType(PetscLimiter, "my_lim"); 40 .ve 41 or at runtime via the option 42 .vb 43 -petsclimiter_type my_lim 44 .ve 45 46 Level: advanced 47 48 .seealso: PetscLimiterRegisterAll(), PetscLimiterRegisterDestroy() 49 50 @*/ 51 PetscErrorCode PetscLimiterRegister(const char sname[], PetscErrorCode (*function)(PetscLimiter)) 52 { 53 PetscFunctionBegin; 54 PetscCall(PetscFunctionListAdd(&PetscLimiterList, sname, function)); 55 PetscFunctionReturn(0); 56 } 57 58 /*@C 59 PetscLimiterSetType - Builds a particular PetscLimiter 60 61 Collective on lim 62 63 Input Parameters: 64 + lim - The PetscLimiter object 65 - name - The kind of limiter 66 67 Options Database Key: 68 . -petsclimiter_type <type> - Sets the PetscLimiter type; use -help for a list of available types 69 70 Level: intermediate 71 72 .seealso: PetscLimiterGetType(), PetscLimiterCreate() 73 @*/ 74 PetscErrorCode PetscLimiterSetType(PetscLimiter lim, PetscLimiterType name) 75 { 76 PetscErrorCode (*r)(PetscLimiter); 77 PetscBool match; 78 79 PetscFunctionBegin; 80 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 81 PetscCall(PetscObjectTypeCompare((PetscObject) lim, name, &match)); 82 if (match) PetscFunctionReturn(0); 83 84 PetscCall(PetscLimiterRegisterAll()); 85 PetscCall(PetscFunctionListFind(PetscLimiterList, name, &r)); 86 PetscCheck(r,PetscObjectComm((PetscObject) lim), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscLimiter type: %s", name); 87 88 if (lim->ops->destroy) { 89 PetscCall((*lim->ops->destroy)(lim)); 90 lim->ops->destroy = NULL; 91 } 92 PetscCall((*r)(lim)); 93 PetscCall(PetscObjectChangeTypeName((PetscObject) lim, name)); 94 PetscFunctionReturn(0); 95 } 96 97 /*@C 98 PetscLimiterGetType - Gets the PetscLimiter type name (as a string) from the object. 99 100 Not Collective 101 102 Input Parameter: 103 . lim - The PetscLimiter 104 105 Output Parameter: 106 . name - The PetscLimiter type name 107 108 Level: intermediate 109 110 .seealso: PetscLimiterSetType(), PetscLimiterCreate() 111 @*/ 112 PetscErrorCode PetscLimiterGetType(PetscLimiter lim, PetscLimiterType *name) 113 { 114 PetscFunctionBegin; 115 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 116 PetscValidPointer(name, 2); 117 PetscCall(PetscLimiterRegisterAll()); 118 *name = ((PetscObject) lim)->type_name; 119 PetscFunctionReturn(0); 120 } 121 122 /*@C 123 PetscLimiterViewFromOptions - View from Options 124 125 Collective on PetscLimiter 126 127 Input Parameters: 128 + A - the PetscLimiter object to view 129 . obj - Optional object 130 - name - command line option 131 132 Level: intermediate 133 .seealso: PetscLimiter, PetscLimiterView, PetscObjectViewFromOptions(), PetscLimiterCreate() 134 @*/ 135 PetscErrorCode PetscLimiterViewFromOptions(PetscLimiter A,PetscObject obj,const char name[]) 136 { 137 PetscFunctionBegin; 138 PetscValidHeaderSpecific(A,PETSCLIMITER_CLASSID,1); 139 PetscCall(PetscObjectViewFromOptions((PetscObject)A,obj,name)); 140 PetscFunctionReturn(0); 141 } 142 143 /*@C 144 PetscLimiterView - Views a PetscLimiter 145 146 Collective on lim 147 148 Input Parameters: 149 + lim - the PetscLimiter object to view 150 - v - the viewer 151 152 Level: beginner 153 154 .seealso: PetscLimiterDestroy() 155 @*/ 156 PetscErrorCode PetscLimiterView(PetscLimiter lim, PetscViewer v) 157 { 158 PetscFunctionBegin; 159 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 160 if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) lim), &v)); 161 if (lim->ops->view) PetscCall((*lim->ops->view)(lim, v)); 162 PetscFunctionReturn(0); 163 } 164 165 /*@ 166 PetscLimiterSetFromOptions - sets parameters in a PetscLimiter from the options database 167 168 Collective on lim 169 170 Input Parameter: 171 . lim - the PetscLimiter object to set options for 172 173 Level: intermediate 174 175 .seealso: PetscLimiterView() 176 @*/ 177 PetscErrorCode PetscLimiterSetFromOptions(PetscLimiter lim) 178 { 179 const char *defaultType; 180 char name[256]; 181 PetscBool flg; 182 183 PetscFunctionBegin; 184 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 185 if (!((PetscObject) lim)->type_name) defaultType = PETSCLIMITERSIN; 186 else defaultType = ((PetscObject) lim)->type_name; 187 PetscCall(PetscLimiterRegisterAll()); 188 189 PetscObjectOptionsBegin((PetscObject) lim); 190 PetscCall(PetscOptionsFList("-petsclimiter_type", "Finite volume slope limiter", "PetscLimiterSetType", PetscLimiterList, defaultType, name, 256, &flg)); 191 if (flg) { 192 PetscCall(PetscLimiterSetType(lim, name)); 193 } else if (!((PetscObject) lim)->type_name) { 194 PetscCall(PetscLimiterSetType(lim, defaultType)); 195 } 196 if (lim->ops->setfromoptions) PetscCall((*lim->ops->setfromoptions)(lim)); 197 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 198 PetscCall(PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) lim)); 199 PetscOptionsEnd(); 200 PetscCall(PetscLimiterViewFromOptions(lim, NULL, "-petsclimiter_view")); 201 PetscFunctionReturn(0); 202 } 203 204 /*@C 205 PetscLimiterSetUp - Construct data structures for the PetscLimiter 206 207 Collective on lim 208 209 Input Parameter: 210 . lim - the PetscLimiter object to setup 211 212 Level: intermediate 213 214 .seealso: PetscLimiterView(), PetscLimiterDestroy() 215 @*/ 216 PetscErrorCode PetscLimiterSetUp(PetscLimiter lim) 217 { 218 PetscFunctionBegin; 219 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 220 if (lim->ops->setup) PetscCall((*lim->ops->setup)(lim)); 221 PetscFunctionReturn(0); 222 } 223 224 /*@ 225 PetscLimiterDestroy - Destroys a PetscLimiter object 226 227 Collective on lim 228 229 Input Parameter: 230 . lim - the PetscLimiter object to destroy 231 232 Level: beginner 233 234 .seealso: PetscLimiterView() 235 @*/ 236 PetscErrorCode PetscLimiterDestroy(PetscLimiter *lim) 237 { 238 PetscFunctionBegin; 239 if (!*lim) PetscFunctionReturn(0); 240 PetscValidHeaderSpecific((*lim), PETSCLIMITER_CLASSID, 1); 241 242 if (--((PetscObject)(*lim))->refct > 0) {*lim = NULL; PetscFunctionReturn(0);} 243 ((PetscObject) (*lim))->refct = 0; 244 245 if ((*lim)->ops->destroy) PetscCall((*(*lim)->ops->destroy)(*lim)); 246 PetscCall(PetscHeaderDestroy(lim)); 247 PetscFunctionReturn(0); 248 } 249 250 /*@ 251 PetscLimiterCreate - Creates an empty PetscLimiter object. The type can then be set with PetscLimiterSetType(). 252 253 Collective 254 255 Input Parameter: 256 . comm - The communicator for the PetscLimiter object 257 258 Output Parameter: 259 . lim - The PetscLimiter object 260 261 Level: beginner 262 263 .seealso: PetscLimiterSetType(), PETSCLIMITERSIN 264 @*/ 265 PetscErrorCode PetscLimiterCreate(MPI_Comm comm, PetscLimiter *lim) 266 { 267 PetscLimiter l; 268 269 PetscFunctionBegin; 270 PetscValidPointer(lim, 2); 271 PetscCall(PetscCitationsRegister(LimiterCitation,&Limitercite)); 272 *lim = NULL; 273 PetscCall(PetscFVInitializePackage()); 274 275 PetscCall(PetscHeaderCreate(l, PETSCLIMITER_CLASSID, "PetscLimiter", "Finite Volume Slope Limiter", "PetscLimiter", comm, PetscLimiterDestroy, PetscLimiterView)); 276 277 *lim = l; 278 PetscFunctionReturn(0); 279 } 280 281 /*@ 282 PetscLimiterLimit - Limit the flux 283 284 Input Parameters: 285 + lim - The PetscLimiter 286 - flim - The input field 287 288 Output Parameter: 289 . phi - The limited field 290 291 Note: Limiters given in symmetric form following Berger, Aftosmis, and Murman 2005 292 $ The classical flux-limited formulation is psi(r) where 293 $ 294 $ r = (u[0] - u[-1]) / (u[1] - u[0]) 295 $ 296 $ The second order TVD region is bounded by 297 $ 298 $ psi_minmod(r) = min(r,1) and psi_superbee(r) = min(2, 2r, max(1,r)) 299 $ 300 $ where all limiters are implicitly clipped to be non-negative. A more convenient slope-limited form is psi(r) = 301 $ phi(r)(r+1)/2 in which the reconstructed interface values are 302 $ 303 $ u(v) = u[0] + phi(r) (grad u)[0] v 304 $ 305 $ where v is the vector from centroid to quadrature point. In these variables, the usual limiters become 306 $ 307 $ phi_minmod(r) = 2 min(1/(1+r),r/(1+r)) phi_superbee(r) = 2 min(2/(1+r), 2r/(1+r), max(1,r)/(1+r)) 308 $ 309 $ For a nicer symmetric formulation, rewrite in terms of 310 $ 311 $ f = (u[0] - u[-1]) / (u[1] - u[-1]) 312 $ 313 $ where r(f) = f/(1-f). Not that r(1-f) = (1-f)/f = 1/r(f) so the symmetry condition 314 $ 315 $ phi(r) = phi(1/r) 316 $ 317 $ becomes 318 $ 319 $ w(f) = w(1-f). 320 $ 321 $ The limiters below implement this final form w(f). The reference methods are 322 $ 323 $ w_minmod(f) = 2 min(f,(1-f)) w_superbee(r) = 4 min((1-f), f) 324 325 Level: beginner 326 327 .seealso: PetscLimiterSetType(), PetscLimiterCreate() 328 @*/ 329 PetscErrorCode PetscLimiterLimit(PetscLimiter lim, PetscReal flim, PetscReal *phi) 330 { 331 PetscFunctionBegin; 332 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 333 PetscValidRealPointer(phi, 3); 334 PetscCall((*lim->ops->limit)(lim, flim, phi)); 335 PetscFunctionReturn(0); 336 } 337 338 static PetscErrorCode PetscLimiterDestroy_Sin(PetscLimiter lim) 339 { 340 PetscLimiter_Sin *l = (PetscLimiter_Sin *) lim->data; 341 342 PetscFunctionBegin; 343 PetscCall(PetscFree(l)); 344 PetscFunctionReturn(0); 345 } 346 347 static PetscErrorCode PetscLimiterView_Sin_Ascii(PetscLimiter lim, PetscViewer viewer) 348 { 349 PetscViewerFormat format; 350 351 PetscFunctionBegin; 352 PetscCall(PetscViewerGetFormat(viewer, &format)); 353 PetscCall(PetscViewerASCIIPrintf(viewer, "Sin Slope Limiter:\n")); 354 PetscFunctionReturn(0); 355 } 356 357 static PetscErrorCode PetscLimiterView_Sin(PetscLimiter lim, PetscViewer viewer) 358 { 359 PetscBool iascii; 360 361 PetscFunctionBegin; 362 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 363 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 364 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 365 if (iascii) PetscCall(PetscLimiterView_Sin_Ascii(lim, viewer)); 366 PetscFunctionReturn(0); 367 } 368 369 static PetscErrorCode PetscLimiterLimit_Sin(PetscLimiter lim, PetscReal f, PetscReal *phi) 370 { 371 PetscFunctionBegin; 372 *phi = PetscSinReal(PETSC_PI*PetscMax(0, PetscMin(f, 1))); 373 PetscFunctionReturn(0); 374 } 375 376 static PetscErrorCode PetscLimiterInitialize_Sin(PetscLimiter lim) 377 { 378 PetscFunctionBegin; 379 lim->ops->view = PetscLimiterView_Sin; 380 lim->ops->destroy = PetscLimiterDestroy_Sin; 381 lim->ops->limit = PetscLimiterLimit_Sin; 382 PetscFunctionReturn(0); 383 } 384 385 /*MC 386 PETSCLIMITERSIN = "sin" - A PetscLimiter object 387 388 Level: intermediate 389 390 .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType() 391 M*/ 392 393 PETSC_EXTERN PetscErrorCode PetscLimiterCreate_Sin(PetscLimiter lim) 394 { 395 PetscLimiter_Sin *l; 396 397 PetscFunctionBegin; 398 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 399 PetscCall(PetscNewLog(lim, &l)); 400 lim->data = l; 401 402 PetscCall(PetscLimiterInitialize_Sin(lim)); 403 PetscFunctionReturn(0); 404 } 405 406 static PetscErrorCode PetscLimiterDestroy_Zero(PetscLimiter lim) 407 { 408 PetscLimiter_Zero *l = (PetscLimiter_Zero *) lim->data; 409 410 PetscFunctionBegin; 411 PetscCall(PetscFree(l)); 412 PetscFunctionReturn(0); 413 } 414 415 static PetscErrorCode PetscLimiterView_Zero_Ascii(PetscLimiter lim, PetscViewer viewer) 416 { 417 PetscViewerFormat format; 418 419 PetscFunctionBegin; 420 PetscCall(PetscViewerGetFormat(viewer, &format)); 421 PetscCall(PetscViewerASCIIPrintf(viewer, "Zero Slope Limiter:\n")); 422 PetscFunctionReturn(0); 423 } 424 425 static PetscErrorCode PetscLimiterView_Zero(PetscLimiter lim, PetscViewer viewer) 426 { 427 PetscBool iascii; 428 429 PetscFunctionBegin; 430 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 431 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 432 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 433 if (iascii) PetscCall(PetscLimiterView_Zero_Ascii(lim, viewer)); 434 PetscFunctionReturn(0); 435 } 436 437 static PetscErrorCode PetscLimiterLimit_Zero(PetscLimiter lim, PetscReal f, PetscReal *phi) 438 { 439 PetscFunctionBegin; 440 *phi = 0.0; 441 PetscFunctionReturn(0); 442 } 443 444 static PetscErrorCode PetscLimiterInitialize_Zero(PetscLimiter lim) 445 { 446 PetscFunctionBegin; 447 lim->ops->view = PetscLimiterView_Zero; 448 lim->ops->destroy = PetscLimiterDestroy_Zero; 449 lim->ops->limit = PetscLimiterLimit_Zero; 450 PetscFunctionReturn(0); 451 } 452 453 /*MC 454 PETSCLIMITERZERO = "zero" - A PetscLimiter object 455 456 Level: intermediate 457 458 .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType() 459 M*/ 460 461 PETSC_EXTERN PetscErrorCode PetscLimiterCreate_Zero(PetscLimiter lim) 462 { 463 PetscLimiter_Zero *l; 464 465 PetscFunctionBegin; 466 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 467 PetscCall(PetscNewLog(lim, &l)); 468 lim->data = l; 469 470 PetscCall(PetscLimiterInitialize_Zero(lim)); 471 PetscFunctionReturn(0); 472 } 473 474 static PetscErrorCode PetscLimiterDestroy_None(PetscLimiter lim) 475 { 476 PetscLimiter_None *l = (PetscLimiter_None *) lim->data; 477 478 PetscFunctionBegin; 479 PetscCall(PetscFree(l)); 480 PetscFunctionReturn(0); 481 } 482 483 static PetscErrorCode PetscLimiterView_None_Ascii(PetscLimiter lim, PetscViewer viewer) 484 { 485 PetscViewerFormat format; 486 487 PetscFunctionBegin; 488 PetscCall(PetscViewerGetFormat(viewer, &format)); 489 PetscCall(PetscViewerASCIIPrintf(viewer, "None Slope Limiter:\n")); 490 PetscFunctionReturn(0); 491 } 492 493 static PetscErrorCode PetscLimiterView_None(PetscLimiter lim, PetscViewer viewer) 494 { 495 PetscBool iascii; 496 497 PetscFunctionBegin; 498 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 499 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 500 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 501 if (iascii) PetscCall(PetscLimiterView_None_Ascii(lim, viewer)); 502 PetscFunctionReturn(0); 503 } 504 505 static PetscErrorCode PetscLimiterLimit_None(PetscLimiter lim, PetscReal f, PetscReal *phi) 506 { 507 PetscFunctionBegin; 508 *phi = 1.0; 509 PetscFunctionReturn(0); 510 } 511 512 static PetscErrorCode PetscLimiterInitialize_None(PetscLimiter lim) 513 { 514 PetscFunctionBegin; 515 lim->ops->view = PetscLimiterView_None; 516 lim->ops->destroy = PetscLimiterDestroy_None; 517 lim->ops->limit = PetscLimiterLimit_None; 518 PetscFunctionReturn(0); 519 } 520 521 /*MC 522 PETSCLIMITERNONE = "none" - A PetscLimiter object 523 524 Level: intermediate 525 526 .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType() 527 M*/ 528 529 PETSC_EXTERN PetscErrorCode PetscLimiterCreate_None(PetscLimiter lim) 530 { 531 PetscLimiter_None *l; 532 533 PetscFunctionBegin; 534 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 535 PetscCall(PetscNewLog(lim, &l)); 536 lim->data = l; 537 538 PetscCall(PetscLimiterInitialize_None(lim)); 539 PetscFunctionReturn(0); 540 } 541 542 static PetscErrorCode PetscLimiterDestroy_Minmod(PetscLimiter lim) 543 { 544 PetscLimiter_Minmod *l = (PetscLimiter_Minmod *) lim->data; 545 546 PetscFunctionBegin; 547 PetscCall(PetscFree(l)); 548 PetscFunctionReturn(0); 549 } 550 551 static PetscErrorCode PetscLimiterView_Minmod_Ascii(PetscLimiter lim, PetscViewer viewer) 552 { 553 PetscViewerFormat format; 554 555 PetscFunctionBegin; 556 PetscCall(PetscViewerGetFormat(viewer, &format)); 557 PetscCall(PetscViewerASCIIPrintf(viewer, "Minmod Slope Limiter:\n")); 558 PetscFunctionReturn(0); 559 } 560 561 static PetscErrorCode PetscLimiterView_Minmod(PetscLimiter lim, PetscViewer viewer) 562 { 563 PetscBool iascii; 564 565 PetscFunctionBegin; 566 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 567 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 568 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 569 if (iascii) PetscCall(PetscLimiterView_Minmod_Ascii(lim, viewer)); 570 PetscFunctionReturn(0); 571 } 572 573 static PetscErrorCode PetscLimiterLimit_Minmod(PetscLimiter lim, PetscReal f, PetscReal *phi) 574 { 575 PetscFunctionBegin; 576 *phi = 2*PetscMax(0, PetscMin(f, 1-f)); 577 PetscFunctionReturn(0); 578 } 579 580 static PetscErrorCode PetscLimiterInitialize_Minmod(PetscLimiter lim) 581 { 582 PetscFunctionBegin; 583 lim->ops->view = PetscLimiterView_Minmod; 584 lim->ops->destroy = PetscLimiterDestroy_Minmod; 585 lim->ops->limit = PetscLimiterLimit_Minmod; 586 PetscFunctionReturn(0); 587 } 588 589 /*MC 590 PETSCLIMITERMINMOD = "minmod" - A PetscLimiter object 591 592 Level: intermediate 593 594 .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType() 595 M*/ 596 597 PETSC_EXTERN PetscErrorCode PetscLimiterCreate_Minmod(PetscLimiter lim) 598 { 599 PetscLimiter_Minmod *l; 600 601 PetscFunctionBegin; 602 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 603 PetscCall(PetscNewLog(lim, &l)); 604 lim->data = l; 605 606 PetscCall(PetscLimiterInitialize_Minmod(lim)); 607 PetscFunctionReturn(0); 608 } 609 610 static PetscErrorCode PetscLimiterDestroy_VanLeer(PetscLimiter lim) 611 { 612 PetscLimiter_VanLeer *l = (PetscLimiter_VanLeer *) lim->data; 613 614 PetscFunctionBegin; 615 PetscCall(PetscFree(l)); 616 PetscFunctionReturn(0); 617 } 618 619 static PetscErrorCode PetscLimiterView_VanLeer_Ascii(PetscLimiter lim, PetscViewer viewer) 620 { 621 PetscViewerFormat format; 622 623 PetscFunctionBegin; 624 PetscCall(PetscViewerGetFormat(viewer, &format)); 625 PetscCall(PetscViewerASCIIPrintf(viewer, "Van Leer Slope Limiter:\n")); 626 PetscFunctionReturn(0); 627 } 628 629 static PetscErrorCode PetscLimiterView_VanLeer(PetscLimiter lim, PetscViewer viewer) 630 { 631 PetscBool iascii; 632 633 PetscFunctionBegin; 634 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 635 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 636 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 637 if (iascii) PetscCall(PetscLimiterView_VanLeer_Ascii(lim, viewer)); 638 PetscFunctionReturn(0); 639 } 640 641 static PetscErrorCode PetscLimiterLimit_VanLeer(PetscLimiter lim, PetscReal f, PetscReal *phi) 642 { 643 PetscFunctionBegin; 644 *phi = PetscMax(0, 4*f*(1-f)); 645 PetscFunctionReturn(0); 646 } 647 648 static PetscErrorCode PetscLimiterInitialize_VanLeer(PetscLimiter lim) 649 { 650 PetscFunctionBegin; 651 lim->ops->view = PetscLimiterView_VanLeer; 652 lim->ops->destroy = PetscLimiterDestroy_VanLeer; 653 lim->ops->limit = PetscLimiterLimit_VanLeer; 654 PetscFunctionReturn(0); 655 } 656 657 /*MC 658 PETSCLIMITERVANLEER = "vanleer" - A PetscLimiter object 659 660 Level: intermediate 661 662 .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType() 663 M*/ 664 665 PETSC_EXTERN PetscErrorCode PetscLimiterCreate_VanLeer(PetscLimiter lim) 666 { 667 PetscLimiter_VanLeer *l; 668 669 PetscFunctionBegin; 670 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 671 PetscCall(PetscNewLog(lim, &l)); 672 lim->data = l; 673 674 PetscCall(PetscLimiterInitialize_VanLeer(lim)); 675 PetscFunctionReturn(0); 676 } 677 678 static PetscErrorCode PetscLimiterDestroy_VanAlbada(PetscLimiter lim) 679 { 680 PetscLimiter_VanAlbada *l = (PetscLimiter_VanAlbada *) lim->data; 681 682 PetscFunctionBegin; 683 PetscCall(PetscFree(l)); 684 PetscFunctionReturn(0); 685 } 686 687 static PetscErrorCode PetscLimiterView_VanAlbada_Ascii(PetscLimiter lim, PetscViewer viewer) 688 { 689 PetscViewerFormat format; 690 691 PetscFunctionBegin; 692 PetscCall(PetscViewerGetFormat(viewer, &format)); 693 PetscCall(PetscViewerASCIIPrintf(viewer, "Van Albada Slope Limiter:\n")); 694 PetscFunctionReturn(0); 695 } 696 697 static PetscErrorCode PetscLimiterView_VanAlbada(PetscLimiter lim, PetscViewer viewer) 698 { 699 PetscBool iascii; 700 701 PetscFunctionBegin; 702 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 703 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 704 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 705 if (iascii) PetscCall(PetscLimiterView_VanAlbada_Ascii(lim, viewer)); 706 PetscFunctionReturn(0); 707 } 708 709 static PetscErrorCode PetscLimiterLimit_VanAlbada(PetscLimiter lim, PetscReal f, PetscReal *phi) 710 { 711 PetscFunctionBegin; 712 *phi = PetscMax(0, 2*f*(1-f) / (PetscSqr(f) + PetscSqr(1-f))); 713 PetscFunctionReturn(0); 714 } 715 716 static PetscErrorCode PetscLimiterInitialize_VanAlbada(PetscLimiter lim) 717 { 718 PetscFunctionBegin; 719 lim->ops->view = PetscLimiterView_VanAlbada; 720 lim->ops->destroy = PetscLimiterDestroy_VanAlbada; 721 lim->ops->limit = PetscLimiterLimit_VanAlbada; 722 PetscFunctionReturn(0); 723 } 724 725 /*MC 726 PETSCLIMITERVANALBADA = "vanalbada" - A PetscLimiter object 727 728 Level: intermediate 729 730 .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType() 731 M*/ 732 733 PETSC_EXTERN PetscErrorCode PetscLimiterCreate_VanAlbada(PetscLimiter lim) 734 { 735 PetscLimiter_VanAlbada *l; 736 737 PetscFunctionBegin; 738 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 739 PetscCall(PetscNewLog(lim, &l)); 740 lim->data = l; 741 742 PetscCall(PetscLimiterInitialize_VanAlbada(lim)); 743 PetscFunctionReturn(0); 744 } 745 746 static PetscErrorCode PetscLimiterDestroy_Superbee(PetscLimiter lim) 747 { 748 PetscLimiter_Superbee *l = (PetscLimiter_Superbee *) lim->data; 749 750 PetscFunctionBegin; 751 PetscCall(PetscFree(l)); 752 PetscFunctionReturn(0); 753 } 754 755 static PetscErrorCode PetscLimiterView_Superbee_Ascii(PetscLimiter lim, PetscViewer viewer) 756 { 757 PetscViewerFormat format; 758 759 PetscFunctionBegin; 760 PetscCall(PetscViewerGetFormat(viewer, &format)); 761 PetscCall(PetscViewerASCIIPrintf(viewer, "Superbee Slope Limiter:\n")); 762 PetscFunctionReturn(0); 763 } 764 765 static PetscErrorCode PetscLimiterView_Superbee(PetscLimiter lim, PetscViewer viewer) 766 { 767 PetscBool iascii; 768 769 PetscFunctionBegin; 770 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 771 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 772 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 773 if (iascii) PetscCall(PetscLimiterView_Superbee_Ascii(lim, viewer)); 774 PetscFunctionReturn(0); 775 } 776 777 static PetscErrorCode PetscLimiterLimit_Superbee(PetscLimiter lim, PetscReal f, PetscReal *phi) 778 { 779 PetscFunctionBegin; 780 *phi = 4*PetscMax(0, PetscMin(f, 1-f)); 781 PetscFunctionReturn(0); 782 } 783 784 static PetscErrorCode PetscLimiterInitialize_Superbee(PetscLimiter lim) 785 { 786 PetscFunctionBegin; 787 lim->ops->view = PetscLimiterView_Superbee; 788 lim->ops->destroy = PetscLimiterDestroy_Superbee; 789 lim->ops->limit = PetscLimiterLimit_Superbee; 790 PetscFunctionReturn(0); 791 } 792 793 /*MC 794 PETSCLIMITERSUPERBEE = "superbee" - A PetscLimiter object 795 796 Level: intermediate 797 798 .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType() 799 M*/ 800 801 PETSC_EXTERN PetscErrorCode PetscLimiterCreate_Superbee(PetscLimiter lim) 802 { 803 PetscLimiter_Superbee *l; 804 805 PetscFunctionBegin; 806 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 807 PetscCall(PetscNewLog(lim, &l)); 808 lim->data = l; 809 810 PetscCall(PetscLimiterInitialize_Superbee(lim)); 811 PetscFunctionReturn(0); 812 } 813 814 static PetscErrorCode PetscLimiterDestroy_MC(PetscLimiter lim) 815 { 816 PetscLimiter_MC *l = (PetscLimiter_MC *) lim->data; 817 818 PetscFunctionBegin; 819 PetscCall(PetscFree(l)); 820 PetscFunctionReturn(0); 821 } 822 823 static PetscErrorCode PetscLimiterView_MC_Ascii(PetscLimiter lim, PetscViewer viewer) 824 { 825 PetscViewerFormat format; 826 827 PetscFunctionBegin; 828 PetscCall(PetscViewerGetFormat(viewer, &format)); 829 PetscCall(PetscViewerASCIIPrintf(viewer, "MC Slope Limiter:\n")); 830 PetscFunctionReturn(0); 831 } 832 833 static PetscErrorCode PetscLimiterView_MC(PetscLimiter lim, PetscViewer viewer) 834 { 835 PetscBool iascii; 836 837 PetscFunctionBegin; 838 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 839 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 840 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 841 if (iascii) PetscCall(PetscLimiterView_MC_Ascii(lim, viewer)); 842 PetscFunctionReturn(0); 843 } 844 845 /* aka Barth-Jespersen */ 846 static PetscErrorCode PetscLimiterLimit_MC(PetscLimiter lim, PetscReal f, PetscReal *phi) 847 { 848 PetscFunctionBegin; 849 *phi = PetscMin(1, 4*PetscMax(0, PetscMin(f, 1-f))); 850 PetscFunctionReturn(0); 851 } 852 853 static PetscErrorCode PetscLimiterInitialize_MC(PetscLimiter lim) 854 { 855 PetscFunctionBegin; 856 lim->ops->view = PetscLimiterView_MC; 857 lim->ops->destroy = PetscLimiterDestroy_MC; 858 lim->ops->limit = PetscLimiterLimit_MC; 859 PetscFunctionReturn(0); 860 } 861 862 /*MC 863 PETSCLIMITERMC = "mc" - A PetscLimiter object 864 865 Level: intermediate 866 867 .seealso: PetscLimiterType, PetscLimiterCreate(), PetscLimiterSetType() 868 M*/ 869 870 PETSC_EXTERN PetscErrorCode PetscLimiterCreate_MC(PetscLimiter lim) 871 { 872 PetscLimiter_MC *l; 873 874 PetscFunctionBegin; 875 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 1); 876 PetscCall(PetscNewLog(lim, &l)); 877 lim->data = l; 878 879 PetscCall(PetscLimiterInitialize_MC(lim)); 880 PetscFunctionReturn(0); 881 } 882 883 PetscClassId PETSCFV_CLASSID = 0; 884 885 PetscFunctionList PetscFVList = NULL; 886 PetscBool PetscFVRegisterAllCalled = PETSC_FALSE; 887 888 /*@C 889 PetscFVRegister - Adds a new PetscFV implementation 890 891 Not Collective 892 893 Input Parameters: 894 + name - The name of a new user-defined creation routine 895 - create_func - The creation routine itself 896 897 Notes: 898 PetscFVRegister() may be called multiple times to add several user-defined PetscFVs 899 900 Sample usage: 901 .vb 902 PetscFVRegister("my_fv", MyPetscFVCreate); 903 .ve 904 905 Then, your PetscFV type can be chosen with the procedural interface via 906 .vb 907 PetscFVCreate(MPI_Comm, PetscFV *); 908 PetscFVSetType(PetscFV, "my_fv"); 909 .ve 910 or at runtime via the option 911 .vb 912 -petscfv_type my_fv 913 .ve 914 915 Level: advanced 916 917 .seealso: PetscFVRegisterAll(), PetscFVRegisterDestroy() 918 919 @*/ 920 PetscErrorCode PetscFVRegister(const char sname[], PetscErrorCode (*function)(PetscFV)) 921 { 922 PetscFunctionBegin; 923 PetscCall(PetscFunctionListAdd(&PetscFVList, sname, function)); 924 PetscFunctionReturn(0); 925 } 926 927 /*@C 928 PetscFVSetType - Builds a particular PetscFV 929 930 Collective on fvm 931 932 Input Parameters: 933 + fvm - The PetscFV object 934 - name - The kind of FVM space 935 936 Options Database Key: 937 . -petscfv_type <type> - Sets the PetscFV type; use -help for a list of available types 938 939 Level: intermediate 940 941 .seealso: PetscFVGetType(), PetscFVCreate() 942 @*/ 943 PetscErrorCode PetscFVSetType(PetscFV fvm, PetscFVType name) 944 { 945 PetscErrorCode (*r)(PetscFV); 946 PetscBool match; 947 948 PetscFunctionBegin; 949 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 950 PetscCall(PetscObjectTypeCompare((PetscObject) fvm, name, &match)); 951 if (match) PetscFunctionReturn(0); 952 953 PetscCall(PetscFVRegisterAll()); 954 PetscCall(PetscFunctionListFind(PetscFVList, name, &r)); 955 PetscCheck(r,PetscObjectComm((PetscObject) fvm), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscFV type: %s", name); 956 957 if (fvm->ops->destroy) { 958 PetscCall((*fvm->ops->destroy)(fvm)); 959 fvm->ops->destroy = NULL; 960 } 961 PetscCall((*r)(fvm)); 962 PetscCall(PetscObjectChangeTypeName((PetscObject) fvm, name)); 963 PetscFunctionReturn(0); 964 } 965 966 /*@C 967 PetscFVGetType - Gets the PetscFV type name (as a string) from the object. 968 969 Not Collective 970 971 Input Parameter: 972 . fvm - The PetscFV 973 974 Output Parameter: 975 . name - The PetscFV type name 976 977 Level: intermediate 978 979 .seealso: PetscFVSetType(), PetscFVCreate() 980 @*/ 981 PetscErrorCode PetscFVGetType(PetscFV fvm, PetscFVType *name) 982 { 983 PetscFunctionBegin; 984 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 985 PetscValidPointer(name, 2); 986 PetscCall(PetscFVRegisterAll()); 987 *name = ((PetscObject) fvm)->type_name; 988 PetscFunctionReturn(0); 989 } 990 991 /*@C 992 PetscFVViewFromOptions - View from Options 993 994 Collective on PetscFV 995 996 Input Parameters: 997 + A - the PetscFV object 998 . obj - Optional object 999 - name - command line option 1000 1001 Level: intermediate 1002 .seealso: PetscFV, PetscFVView, PetscObjectViewFromOptions(), PetscFVCreate() 1003 @*/ 1004 PetscErrorCode PetscFVViewFromOptions(PetscFV A,PetscObject obj,const char name[]) 1005 { 1006 PetscFunctionBegin; 1007 PetscValidHeaderSpecific(A,PETSCFV_CLASSID,1); 1008 PetscCall(PetscObjectViewFromOptions((PetscObject)A,obj,name)); 1009 PetscFunctionReturn(0); 1010 } 1011 1012 /*@C 1013 PetscFVView - Views a PetscFV 1014 1015 Collective on fvm 1016 1017 Input Parameters: 1018 + fvm - the PetscFV object to view 1019 - v - the viewer 1020 1021 Level: beginner 1022 1023 .seealso: PetscFVDestroy() 1024 @*/ 1025 PetscErrorCode PetscFVView(PetscFV fvm, PetscViewer v) 1026 { 1027 PetscFunctionBegin; 1028 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1029 if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) fvm), &v)); 1030 if (fvm->ops->view) PetscCall((*fvm->ops->view)(fvm, v)); 1031 PetscFunctionReturn(0); 1032 } 1033 1034 /*@ 1035 PetscFVSetFromOptions - sets parameters in a PetscFV from the options database 1036 1037 Collective on fvm 1038 1039 Input Parameter: 1040 . fvm - the PetscFV object to set options for 1041 1042 Options Database Key: 1043 . -petscfv_compute_gradients <bool> - Determines whether cell gradients are calculated 1044 1045 Level: intermediate 1046 1047 .seealso: PetscFVView() 1048 @*/ 1049 PetscErrorCode PetscFVSetFromOptions(PetscFV fvm) 1050 { 1051 const char *defaultType; 1052 char name[256]; 1053 PetscBool flg; 1054 1055 PetscFunctionBegin; 1056 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1057 if (!((PetscObject) fvm)->type_name) defaultType = PETSCFVUPWIND; 1058 else defaultType = ((PetscObject) fvm)->type_name; 1059 PetscCall(PetscFVRegisterAll()); 1060 1061 PetscObjectOptionsBegin((PetscObject) fvm); 1062 PetscCall(PetscOptionsFList("-petscfv_type", "Finite volume discretization", "PetscFVSetType", PetscFVList, defaultType, name, 256, &flg)); 1063 if (flg) { 1064 PetscCall(PetscFVSetType(fvm, name)); 1065 } else if (!((PetscObject) fvm)->type_name) { 1066 PetscCall(PetscFVSetType(fvm, defaultType)); 1067 1068 } 1069 PetscCall(PetscOptionsBool("-petscfv_compute_gradients", "Compute cell gradients", "PetscFVSetComputeGradients", fvm->computeGradients, &fvm->computeGradients, NULL)); 1070 if (fvm->ops->setfromoptions) PetscCall((*fvm->ops->setfromoptions)(fvm)); 1071 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 1072 PetscCall(PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) fvm)); 1073 PetscCall(PetscLimiterSetFromOptions(fvm->limiter)); 1074 PetscOptionsEnd(); 1075 PetscCall(PetscFVViewFromOptions(fvm, NULL, "-petscfv_view")); 1076 PetscFunctionReturn(0); 1077 } 1078 1079 /*@ 1080 PetscFVSetUp - Construct data structures for the PetscFV 1081 1082 Collective on fvm 1083 1084 Input Parameter: 1085 . fvm - the PetscFV object to setup 1086 1087 Level: intermediate 1088 1089 .seealso: PetscFVView(), PetscFVDestroy() 1090 @*/ 1091 PetscErrorCode PetscFVSetUp(PetscFV fvm) 1092 { 1093 PetscFunctionBegin; 1094 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1095 PetscCall(PetscLimiterSetUp(fvm->limiter)); 1096 if (fvm->ops->setup) PetscCall((*fvm->ops->setup)(fvm)); 1097 PetscFunctionReturn(0); 1098 } 1099 1100 /*@ 1101 PetscFVDestroy - Destroys a PetscFV object 1102 1103 Collective on fvm 1104 1105 Input Parameter: 1106 . fvm - the PetscFV object to destroy 1107 1108 Level: beginner 1109 1110 .seealso: PetscFVView() 1111 @*/ 1112 PetscErrorCode PetscFVDestroy(PetscFV *fvm) 1113 { 1114 PetscInt i; 1115 1116 PetscFunctionBegin; 1117 if (!*fvm) PetscFunctionReturn(0); 1118 PetscValidHeaderSpecific((*fvm), PETSCFV_CLASSID, 1); 1119 1120 if (--((PetscObject)(*fvm))->refct > 0) {*fvm = NULL; PetscFunctionReturn(0);} 1121 ((PetscObject) (*fvm))->refct = 0; 1122 1123 for (i = 0; i < (*fvm)->numComponents; i++) { 1124 PetscCall(PetscFree((*fvm)->componentNames[i])); 1125 } 1126 PetscCall(PetscFree((*fvm)->componentNames)); 1127 PetscCall(PetscLimiterDestroy(&(*fvm)->limiter)); 1128 PetscCall(PetscDualSpaceDestroy(&(*fvm)->dualSpace)); 1129 PetscCall(PetscFree((*fvm)->fluxWork)); 1130 PetscCall(PetscQuadratureDestroy(&(*fvm)->quadrature)); 1131 PetscCall(PetscTabulationDestroy(&(*fvm)->T)); 1132 1133 if ((*fvm)->ops->destroy) PetscCall((*(*fvm)->ops->destroy)(*fvm)); 1134 PetscCall(PetscHeaderDestroy(fvm)); 1135 PetscFunctionReturn(0); 1136 } 1137 1138 /*@ 1139 PetscFVCreate - Creates an empty PetscFV object. The type can then be set with PetscFVSetType(). 1140 1141 Collective 1142 1143 Input Parameter: 1144 . comm - The communicator for the PetscFV object 1145 1146 Output Parameter: 1147 . fvm - The PetscFV object 1148 1149 Level: beginner 1150 1151 .seealso: PetscFVSetType(), PETSCFVUPWIND 1152 @*/ 1153 PetscErrorCode PetscFVCreate(MPI_Comm comm, PetscFV *fvm) 1154 { 1155 PetscFV f; 1156 1157 PetscFunctionBegin; 1158 PetscValidPointer(fvm, 2); 1159 *fvm = NULL; 1160 PetscCall(PetscFVInitializePackage()); 1161 1162 PetscCall(PetscHeaderCreate(f, PETSCFV_CLASSID, "PetscFV", "Finite Volume", "PetscFV", comm, PetscFVDestroy, PetscFVView)); 1163 PetscCall(PetscMemzero(f->ops, sizeof(struct _PetscFVOps))); 1164 1165 PetscCall(PetscLimiterCreate(comm, &f->limiter)); 1166 f->numComponents = 1; 1167 f->dim = 0; 1168 f->computeGradients = PETSC_FALSE; 1169 f->fluxWork = NULL; 1170 PetscCall(PetscCalloc1(f->numComponents,&f->componentNames)); 1171 1172 *fvm = f; 1173 PetscFunctionReturn(0); 1174 } 1175 1176 /*@ 1177 PetscFVSetLimiter - Set the limiter object 1178 1179 Logically collective on fvm 1180 1181 Input Parameters: 1182 + fvm - the PetscFV object 1183 - lim - The PetscLimiter 1184 1185 Level: intermediate 1186 1187 .seealso: PetscFVGetLimiter() 1188 @*/ 1189 PetscErrorCode PetscFVSetLimiter(PetscFV fvm, PetscLimiter lim) 1190 { 1191 PetscFunctionBegin; 1192 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1193 PetscValidHeaderSpecific(lim, PETSCLIMITER_CLASSID, 2); 1194 PetscCall(PetscLimiterDestroy(&fvm->limiter)); 1195 PetscCall(PetscObjectReference((PetscObject) lim)); 1196 fvm->limiter = lim; 1197 PetscFunctionReturn(0); 1198 } 1199 1200 /*@ 1201 PetscFVGetLimiter - Get the limiter object 1202 1203 Not collective 1204 1205 Input Parameter: 1206 . fvm - the PetscFV object 1207 1208 Output Parameter: 1209 . lim - The PetscLimiter 1210 1211 Level: intermediate 1212 1213 .seealso: PetscFVSetLimiter() 1214 @*/ 1215 PetscErrorCode PetscFVGetLimiter(PetscFV fvm, PetscLimiter *lim) 1216 { 1217 PetscFunctionBegin; 1218 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1219 PetscValidPointer(lim, 2); 1220 *lim = fvm->limiter; 1221 PetscFunctionReturn(0); 1222 } 1223 1224 /*@ 1225 PetscFVSetNumComponents - Set the number of field components 1226 1227 Logically collective on fvm 1228 1229 Input Parameters: 1230 + fvm - the PetscFV object 1231 - comp - The number of components 1232 1233 Level: intermediate 1234 1235 .seealso: PetscFVGetNumComponents() 1236 @*/ 1237 PetscErrorCode PetscFVSetNumComponents(PetscFV fvm, PetscInt comp) 1238 { 1239 PetscFunctionBegin; 1240 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1241 if (fvm->numComponents != comp) { 1242 PetscInt i; 1243 1244 for (i = 0; i < fvm->numComponents; i++) { 1245 PetscCall(PetscFree(fvm->componentNames[i])); 1246 } 1247 PetscCall(PetscFree(fvm->componentNames)); 1248 PetscCall(PetscCalloc1(comp,&fvm->componentNames)); 1249 } 1250 fvm->numComponents = comp; 1251 PetscCall(PetscFree(fvm->fluxWork)); 1252 PetscCall(PetscMalloc1(comp, &fvm->fluxWork)); 1253 PetscFunctionReturn(0); 1254 } 1255 1256 /*@ 1257 PetscFVGetNumComponents - Get the number of field components 1258 1259 Not collective 1260 1261 Input Parameter: 1262 . fvm - the PetscFV object 1263 1264 Output Parameter: 1265 , comp - The number of components 1266 1267 Level: intermediate 1268 1269 .seealso: PetscFVSetNumComponents() 1270 @*/ 1271 PetscErrorCode PetscFVGetNumComponents(PetscFV fvm, PetscInt *comp) 1272 { 1273 PetscFunctionBegin; 1274 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1275 PetscValidIntPointer(comp, 2); 1276 *comp = fvm->numComponents; 1277 PetscFunctionReturn(0); 1278 } 1279 1280 /*@C 1281 PetscFVSetComponentName - Set the name of a component (used in output and viewing) 1282 1283 Logically collective on fvm 1284 Input Parameters: 1285 + fvm - the PetscFV object 1286 . comp - the component number 1287 - name - the component name 1288 1289 Level: intermediate 1290 1291 .seealso: PetscFVGetComponentName() 1292 @*/ 1293 PetscErrorCode PetscFVSetComponentName(PetscFV fvm, PetscInt comp, const char *name) 1294 { 1295 PetscFunctionBegin; 1296 PetscCall(PetscFree(fvm->componentNames[comp])); 1297 PetscCall(PetscStrallocpy(name,&fvm->componentNames[comp])); 1298 PetscFunctionReturn(0); 1299 } 1300 1301 /*@C 1302 PetscFVGetComponentName - Get the name of a component (used in output and viewing) 1303 1304 Logically collective on fvm 1305 Input Parameters: 1306 + fvm - the PetscFV object 1307 - comp - the component number 1308 1309 Output Parameter: 1310 . name - the component name 1311 1312 Level: intermediate 1313 1314 .seealso: PetscFVSetComponentName() 1315 @*/ 1316 PetscErrorCode PetscFVGetComponentName(PetscFV fvm, PetscInt comp, const char **name) 1317 { 1318 PetscFunctionBegin; 1319 *name = fvm->componentNames[comp]; 1320 PetscFunctionReturn(0); 1321 } 1322 1323 /*@ 1324 PetscFVSetSpatialDimension - Set the spatial dimension 1325 1326 Logically collective on fvm 1327 1328 Input Parameters: 1329 + fvm - the PetscFV object 1330 - dim - The spatial dimension 1331 1332 Level: intermediate 1333 1334 .seealso: PetscFVGetSpatialDimension() 1335 @*/ 1336 PetscErrorCode PetscFVSetSpatialDimension(PetscFV fvm, PetscInt dim) 1337 { 1338 PetscFunctionBegin; 1339 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1340 fvm->dim = dim; 1341 PetscFunctionReturn(0); 1342 } 1343 1344 /*@ 1345 PetscFVGetSpatialDimension - Get the spatial dimension 1346 1347 Logically collective on fvm 1348 1349 Input Parameter: 1350 . fvm - the PetscFV object 1351 1352 Output Parameter: 1353 . dim - The spatial dimension 1354 1355 Level: intermediate 1356 1357 .seealso: PetscFVSetSpatialDimension() 1358 @*/ 1359 PetscErrorCode PetscFVGetSpatialDimension(PetscFV fvm, PetscInt *dim) 1360 { 1361 PetscFunctionBegin; 1362 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1363 PetscValidIntPointer(dim, 2); 1364 *dim = fvm->dim; 1365 PetscFunctionReturn(0); 1366 } 1367 1368 /*@ 1369 PetscFVSetComputeGradients - Toggle computation of cell gradients 1370 1371 Logically collective on fvm 1372 1373 Input Parameters: 1374 + fvm - the PetscFV object 1375 - computeGradients - Flag to compute cell gradients 1376 1377 Level: intermediate 1378 1379 .seealso: PetscFVGetComputeGradients() 1380 @*/ 1381 PetscErrorCode PetscFVSetComputeGradients(PetscFV fvm, PetscBool computeGradients) 1382 { 1383 PetscFunctionBegin; 1384 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1385 fvm->computeGradients = computeGradients; 1386 PetscFunctionReturn(0); 1387 } 1388 1389 /*@ 1390 PetscFVGetComputeGradients - Return flag for computation of cell gradients 1391 1392 Not collective 1393 1394 Input Parameter: 1395 . fvm - the PetscFV object 1396 1397 Output Parameter: 1398 . computeGradients - Flag to compute cell gradients 1399 1400 Level: intermediate 1401 1402 .seealso: PetscFVSetComputeGradients() 1403 @*/ 1404 PetscErrorCode PetscFVGetComputeGradients(PetscFV fvm, PetscBool *computeGradients) 1405 { 1406 PetscFunctionBegin; 1407 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1408 PetscValidBoolPointer(computeGradients, 2); 1409 *computeGradients = fvm->computeGradients; 1410 PetscFunctionReturn(0); 1411 } 1412 1413 /*@ 1414 PetscFVSetQuadrature - Set the quadrature object 1415 1416 Logically collective on fvm 1417 1418 Input Parameters: 1419 + fvm - the PetscFV object 1420 - q - The PetscQuadrature 1421 1422 Level: intermediate 1423 1424 .seealso: PetscFVGetQuadrature() 1425 @*/ 1426 PetscErrorCode PetscFVSetQuadrature(PetscFV fvm, PetscQuadrature q) 1427 { 1428 PetscFunctionBegin; 1429 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1430 PetscCall(PetscQuadratureDestroy(&fvm->quadrature)); 1431 PetscCall(PetscObjectReference((PetscObject) q)); 1432 fvm->quadrature = q; 1433 PetscFunctionReturn(0); 1434 } 1435 1436 /*@ 1437 PetscFVGetQuadrature - Get the quadrature object 1438 1439 Not collective 1440 1441 Input Parameter: 1442 . fvm - the PetscFV object 1443 1444 Output Parameter: 1445 . lim - The PetscQuadrature 1446 1447 Level: intermediate 1448 1449 .seealso: PetscFVSetQuadrature() 1450 @*/ 1451 PetscErrorCode PetscFVGetQuadrature(PetscFV fvm, PetscQuadrature *q) 1452 { 1453 PetscFunctionBegin; 1454 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1455 PetscValidPointer(q, 2); 1456 if (!fvm->quadrature) { 1457 /* Create default 1-point quadrature */ 1458 PetscReal *points, *weights; 1459 1460 PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &fvm->quadrature)); 1461 PetscCall(PetscCalloc1(fvm->dim, &points)); 1462 PetscCall(PetscMalloc1(1, &weights)); 1463 weights[0] = 1.0; 1464 PetscCall(PetscQuadratureSetData(fvm->quadrature, fvm->dim, 1, 1, points, weights)); 1465 } 1466 *q = fvm->quadrature; 1467 PetscFunctionReturn(0); 1468 } 1469 1470 /*@ 1471 PetscFVGetDualSpace - Returns the PetscDualSpace used to define the inner product 1472 1473 Not collective 1474 1475 Input Parameter: 1476 . fvm - The PetscFV object 1477 1478 Output Parameter: 1479 . sp - The PetscDualSpace object 1480 1481 Note: A simple dual space is provided automatically, and the user typically will not need to override it. 1482 1483 Level: intermediate 1484 1485 .seealso: PetscFVCreate() 1486 @*/ 1487 PetscErrorCode PetscFVGetDualSpace(PetscFV fvm, PetscDualSpace *sp) 1488 { 1489 PetscFunctionBegin; 1490 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1491 PetscValidPointer(sp, 2); 1492 if (!fvm->dualSpace) { 1493 DM K; 1494 PetscInt dim, Nc, c; 1495 1496 PetscCall(PetscFVGetSpatialDimension(fvm, &dim)); 1497 PetscCall(PetscFVGetNumComponents(fvm, &Nc)); 1498 PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject) fvm), &fvm->dualSpace)); 1499 PetscCall(PetscDualSpaceSetType(fvm->dualSpace, PETSCDUALSPACESIMPLE)); 1500 PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, DMPolytopeTypeSimpleShape(dim, PETSC_FALSE), &K)); 1501 PetscCall(PetscDualSpaceSetNumComponents(fvm->dualSpace, Nc)); 1502 PetscCall(PetscDualSpaceSetDM(fvm->dualSpace, K)); 1503 PetscCall(DMDestroy(&K)); 1504 PetscCall(PetscDualSpaceSimpleSetDimension(fvm->dualSpace, Nc)); 1505 /* Should we be using PetscFVGetQuadrature() here? */ 1506 for (c = 0; c < Nc; ++c) { 1507 PetscQuadrature qc; 1508 PetscReal *points, *weights; 1509 1510 PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &qc)); 1511 PetscCall(PetscCalloc1(dim, &points)); 1512 PetscCall(PetscCalloc1(Nc, &weights)); 1513 weights[c] = 1.0; 1514 PetscCall(PetscQuadratureSetData(qc, dim, Nc, 1, points, weights)); 1515 PetscCall(PetscDualSpaceSimpleSetFunctional(fvm->dualSpace, c, qc)); 1516 PetscCall(PetscQuadratureDestroy(&qc)); 1517 } 1518 PetscCall(PetscDualSpaceSetUp(fvm->dualSpace)); 1519 } 1520 *sp = fvm->dualSpace; 1521 PetscFunctionReturn(0); 1522 } 1523 1524 /*@ 1525 PetscFVSetDualSpace - Sets the PetscDualSpace used to define the inner product 1526 1527 Not collective 1528 1529 Input Parameters: 1530 + fvm - The PetscFV object 1531 - sp - The PetscDualSpace object 1532 1533 Level: intermediate 1534 1535 Note: A simple dual space is provided automatically, and the user typically will not need to override it. 1536 1537 .seealso: PetscFVCreate() 1538 @*/ 1539 PetscErrorCode PetscFVSetDualSpace(PetscFV fvm, PetscDualSpace sp) 1540 { 1541 PetscFunctionBegin; 1542 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1543 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 2); 1544 PetscCall(PetscDualSpaceDestroy(&fvm->dualSpace)); 1545 fvm->dualSpace = sp; 1546 PetscCall(PetscObjectReference((PetscObject) fvm->dualSpace)); 1547 PetscFunctionReturn(0); 1548 } 1549 1550 /*@C 1551 PetscFVGetCellTabulation - Returns the tabulation of the basis functions at the quadrature points 1552 1553 Not collective 1554 1555 Input Parameter: 1556 . fvm - The PetscFV object 1557 1558 Output Parameter: 1559 . T - The basis function values and derivatives at quadrature points 1560 1561 Note: 1562 $ T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 1563 $ T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d 1564 $ T->T[2] = H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis function i, component c, in directions d and e 1565 1566 Level: intermediate 1567 1568 .seealso: PetscFEGetCellTabulation(), PetscFVCreateTabulation(), PetscFVGetQuadrature(), PetscQuadratureGetData() 1569 @*/ 1570 PetscErrorCode PetscFVGetCellTabulation(PetscFV fvm, PetscTabulation *T) 1571 { 1572 PetscInt npoints; 1573 const PetscReal *points; 1574 1575 PetscFunctionBegin; 1576 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1577 PetscValidPointer(T, 2); 1578 PetscCall(PetscQuadratureGetData(fvm->quadrature, NULL, NULL, &npoints, &points, NULL)); 1579 if (!fvm->T) PetscCall(PetscFVCreateTabulation(fvm, 1, npoints, points, 1, &fvm->T)); 1580 *T = fvm->T; 1581 PetscFunctionReturn(0); 1582 } 1583 1584 /*@C 1585 PetscFVCreateTabulation - Tabulates the basis functions, and perhaps derivatives, at the points provided. 1586 1587 Not collective 1588 1589 Input Parameters: 1590 + fvm - The PetscFV object 1591 . nrepl - The number of replicas 1592 . npoints - The number of tabulation points in a replica 1593 . points - The tabulation point coordinates 1594 - K - The order of derivative to tabulate 1595 1596 Output Parameter: 1597 . T - The basis function values and derivative at tabulation points 1598 1599 Note: 1600 $ T->T[0] = B[(p*pdim + i)*Nc + c] is the value at point p for basis function i and component c 1601 $ T->T[1] = D[((p*pdim + i)*Nc + c)*dim + d] is the derivative value at point p for basis function i, component c, in direction d 1602 $ T->T[2] = H[(((p*pdim + i)*Nc + c)*dim + d)*dim + e] is the value at point p for basis function i, component c, in directions d and e 1603 1604 Level: intermediate 1605 1606 .seealso: PetscFECreateTabulation(), PetscTabulationDestroy(), PetscFEGetCellTabulation() 1607 @*/ 1608 PetscErrorCode PetscFVCreateTabulation(PetscFV fvm, PetscInt nrepl, PetscInt npoints, const PetscReal points[], PetscInt K, PetscTabulation *T) 1609 { 1610 PetscInt pdim = 1; /* Dimension of approximation space P */ 1611 PetscInt cdim; /* Spatial dimension */ 1612 PetscInt Nc; /* Field components */ 1613 PetscInt k, p, d, c, e; 1614 1615 PetscFunctionBegin; 1616 if (!npoints || K < 0) { 1617 *T = NULL; 1618 PetscFunctionReturn(0); 1619 } 1620 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1621 PetscValidRealPointer(points, 4); 1622 PetscValidPointer(T, 6); 1623 PetscCall(PetscFVGetSpatialDimension(fvm, &cdim)); 1624 PetscCall(PetscFVGetNumComponents(fvm, &Nc)); 1625 PetscCall(PetscMalloc1(1, T)); 1626 (*T)->K = !cdim ? 0 : K; 1627 (*T)->Nr = nrepl; 1628 (*T)->Np = npoints; 1629 (*T)->Nb = pdim; 1630 (*T)->Nc = Nc; 1631 (*T)->cdim = cdim; 1632 PetscCall(PetscMalloc1((*T)->K+1, &(*T)->T)); 1633 for (k = 0; k <= (*T)->K; ++k) { 1634 PetscCall(PetscMalloc1(nrepl*npoints*pdim*Nc*PetscPowInt(cdim, k), &(*T)->T[k])); 1635 } 1636 if (K >= 0) {for (p = 0; p < nrepl*npoints; ++p) for (d = 0; d < pdim; ++d) for (c = 0; c < Nc; ++c) (*T)->T[0][(p*pdim + d)*Nc + c] = 1.0;} 1637 if (K >= 1) {for (p = 0; p < nrepl*npoints; ++p) for (d = 0; d < pdim; ++d) for (c = 0; c < Nc; ++c) for (e = 0; e < cdim; ++e) (*T)->T[1][((p*pdim + d)*Nc + c)*cdim + e] = 0.0;} 1638 if (K >= 2) {for (p = 0; p < nrepl*npoints; ++p) for (d = 0; d < pdim; ++d) for (c = 0; c < Nc; ++c) for (e = 0; e < cdim*cdim; ++e) (*T)->T[2][((p*pdim + d)*Nc + c)*cdim*cdim + e] = 0.0;} 1639 PetscFunctionReturn(0); 1640 } 1641 1642 /*@C 1643 PetscFVComputeGradient - Compute the gradient reconstruction matrix for a given cell 1644 1645 Input Parameters: 1646 + fvm - The PetscFV object 1647 . numFaces - The number of cell faces which are not constrained 1648 - dx - The vector from the cell centroid to the neighboring cell centroid for each face 1649 1650 Level: advanced 1651 1652 .seealso: PetscFVCreate() 1653 @*/ 1654 PetscErrorCode PetscFVComputeGradient(PetscFV fvm, PetscInt numFaces, PetscScalar dx[], PetscScalar grad[]) 1655 { 1656 PetscFunctionBegin; 1657 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1658 if (fvm->ops->computegradient) PetscCall((*fvm->ops->computegradient)(fvm, numFaces, dx, grad)); 1659 PetscFunctionReturn(0); 1660 } 1661 1662 /*@C 1663 PetscFVIntegrateRHSFunction - Produce the cell residual vector for a chunk of elements by quadrature integration 1664 1665 Not collective 1666 1667 Input Parameters: 1668 + fvm - The PetscFV object for the field being integrated 1669 . prob - The PetscDS specifing the discretizations and continuum functions 1670 . field - The field being integrated 1671 . Nf - The number of faces in the chunk 1672 . fgeom - The face geometry for each face in the chunk 1673 . neighborVol - The volume for each pair of cells in the chunk 1674 . uL - The state from the cell on the left 1675 - uR - The state from the cell on the right 1676 1677 Output Parameters: 1678 + fluxL - the left fluxes for each face 1679 - fluxR - the right fluxes for each face 1680 1681 Level: developer 1682 1683 .seealso: PetscFVCreate() 1684 @*/ 1685 PetscErrorCode PetscFVIntegrateRHSFunction(PetscFV fvm, PetscDS prob, PetscInt field, PetscInt Nf, PetscFVFaceGeom *fgeom, PetscReal *neighborVol, 1686 PetscScalar uL[], PetscScalar uR[], PetscScalar fluxL[], PetscScalar fluxR[]) 1687 { 1688 PetscFunctionBegin; 1689 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1690 if (fvm->ops->integraterhsfunction) PetscCall((*fvm->ops->integraterhsfunction)(fvm, prob, field, Nf, fgeom, neighborVol, uL, uR, fluxL, fluxR)); 1691 PetscFunctionReturn(0); 1692 } 1693 1694 /*@ 1695 PetscFVRefine - Create a "refined" PetscFV object that refines the reference cell into smaller copies. This is typically used 1696 to precondition a higher order method with a lower order method on a refined mesh having the same number of dofs (but more 1697 sparsity). It is also used to create an interpolation between regularly refined meshes. 1698 1699 Input Parameter: 1700 . fv - The initial PetscFV 1701 1702 Output Parameter: 1703 . fvRef - The refined PetscFV 1704 1705 Level: advanced 1706 1707 .seealso: PetscFVType, PetscFVCreate(), PetscFVSetType() 1708 @*/ 1709 PetscErrorCode PetscFVRefine(PetscFV fv, PetscFV *fvRef) 1710 { 1711 PetscDualSpace Q, Qref; 1712 DM K, Kref; 1713 PetscQuadrature q, qref; 1714 DMPolytopeType ct; 1715 DMPlexTransform tr; 1716 PetscReal *v0; 1717 PetscReal *jac, *invjac; 1718 PetscInt numComp, numSubelements, s; 1719 1720 PetscFunctionBegin; 1721 PetscCall(PetscFVGetDualSpace(fv, &Q)); 1722 PetscCall(PetscFVGetQuadrature(fv, &q)); 1723 PetscCall(PetscDualSpaceGetDM(Q, &K)); 1724 /* Create dual space */ 1725 PetscCall(PetscDualSpaceDuplicate(Q, &Qref)); 1726 PetscCall(DMRefine(K, PetscObjectComm((PetscObject) fv), &Kref)); 1727 PetscCall(PetscDualSpaceSetDM(Qref, Kref)); 1728 PetscCall(DMDestroy(&Kref)); 1729 PetscCall(PetscDualSpaceSetUp(Qref)); 1730 /* Create volume */ 1731 PetscCall(PetscFVCreate(PetscObjectComm((PetscObject) fv), fvRef)); 1732 PetscCall(PetscFVSetDualSpace(*fvRef, Qref)); 1733 PetscCall(PetscFVGetNumComponents(fv, &numComp)); 1734 PetscCall(PetscFVSetNumComponents(*fvRef, numComp)); 1735 PetscCall(PetscFVSetUp(*fvRef)); 1736 /* Create quadrature */ 1737 PetscCall(DMPlexGetCellType(K, 0, &ct)); 1738 PetscCall(DMPlexTransformCreate(PETSC_COMM_SELF, &tr)); 1739 PetscCall(DMPlexTransformSetType(tr, DMPLEXREFINEREGULAR)); 1740 PetscCall(DMPlexRefineRegularGetAffineTransforms(tr, ct, &numSubelements, &v0, &jac, &invjac)); 1741 PetscCall(PetscQuadratureExpandComposite(q, numSubelements, v0, jac, &qref)); 1742 PetscCall(PetscDualSpaceSimpleSetDimension(Qref, numSubelements)); 1743 for (s = 0; s < numSubelements; ++s) { 1744 PetscQuadrature qs; 1745 const PetscReal *points, *weights; 1746 PetscReal *p, *w; 1747 PetscInt dim, Nc, npoints, np; 1748 1749 PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF, &qs)); 1750 PetscCall(PetscQuadratureGetData(q, &dim, &Nc, &npoints, &points, &weights)); 1751 np = npoints/numSubelements; 1752 PetscCall(PetscMalloc1(np*dim,&p)); 1753 PetscCall(PetscMalloc1(np*Nc,&w)); 1754 PetscCall(PetscArraycpy(p, &points[s*np*dim], np*dim)); 1755 PetscCall(PetscArraycpy(w, &weights[s*np*Nc], np*Nc)); 1756 PetscCall(PetscQuadratureSetData(qs, dim, Nc, np, p, w)); 1757 PetscCall(PetscDualSpaceSimpleSetFunctional(Qref, s, qs)); 1758 PetscCall(PetscQuadratureDestroy(&qs)); 1759 } 1760 PetscCall(PetscFVSetQuadrature(*fvRef, qref)); 1761 PetscCall(DMPlexTransformDestroy(&tr)); 1762 PetscCall(PetscQuadratureDestroy(&qref)); 1763 PetscCall(PetscDualSpaceDestroy(&Qref)); 1764 PetscFunctionReturn(0); 1765 } 1766 1767 static PetscErrorCode PetscFVDestroy_Upwind(PetscFV fvm) 1768 { 1769 PetscFV_Upwind *b = (PetscFV_Upwind *) fvm->data; 1770 1771 PetscFunctionBegin; 1772 PetscCall(PetscFree(b)); 1773 PetscFunctionReturn(0); 1774 } 1775 1776 static PetscErrorCode PetscFVView_Upwind_Ascii(PetscFV fv, PetscViewer viewer) 1777 { 1778 PetscInt Nc, c; 1779 PetscViewerFormat format; 1780 1781 PetscFunctionBegin; 1782 PetscCall(PetscFVGetNumComponents(fv, &Nc)); 1783 PetscCall(PetscViewerGetFormat(viewer, &format)); 1784 PetscCall(PetscViewerASCIIPrintf(viewer, "Upwind Finite Volume:\n")); 1785 PetscCall(PetscViewerASCIIPrintf(viewer, " num components: %" PetscInt_FMT "\n", Nc)); 1786 for (c = 0; c < Nc; c++) { 1787 if (fv->componentNames[c]) { 1788 PetscCall(PetscViewerASCIIPrintf(viewer, " component %" PetscInt_FMT ": %s\n", c, fv->componentNames[c])); 1789 } 1790 } 1791 PetscFunctionReturn(0); 1792 } 1793 1794 static PetscErrorCode PetscFVView_Upwind(PetscFV fv, PetscViewer viewer) 1795 { 1796 PetscBool iascii; 1797 1798 PetscFunctionBegin; 1799 PetscValidHeaderSpecific(fv, PETSCFV_CLASSID, 1); 1800 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1801 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 1802 if (iascii) PetscCall(PetscFVView_Upwind_Ascii(fv, viewer)); 1803 PetscFunctionReturn(0); 1804 } 1805 1806 static PetscErrorCode PetscFVSetUp_Upwind(PetscFV fvm) 1807 { 1808 PetscFunctionBegin; 1809 PetscFunctionReturn(0); 1810 } 1811 1812 /* 1813 neighborVol[f*2+0] contains the left geom 1814 neighborVol[f*2+1] contains the right geom 1815 */ 1816 static PetscErrorCode PetscFVIntegrateRHSFunction_Upwind(PetscFV fvm, PetscDS prob, PetscInt field, PetscInt Nf, PetscFVFaceGeom *fgeom, PetscReal *neighborVol, 1817 PetscScalar uL[], PetscScalar uR[], PetscScalar fluxL[], PetscScalar fluxR[]) 1818 { 1819 void (*riemann)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *); 1820 void *rctx; 1821 PetscScalar *flux = fvm->fluxWork; 1822 const PetscScalar *constants; 1823 PetscInt dim, numConstants, pdim, totDim, Nc, off, f, d; 1824 1825 PetscFunctionBegin; 1826 PetscCall(PetscDSGetTotalComponents(prob, &Nc)); 1827 PetscCall(PetscDSGetTotalDimension(prob, &totDim)); 1828 PetscCall(PetscDSGetFieldOffset(prob, field, &off)); 1829 PetscCall(PetscDSGetRiemannSolver(prob, field, &riemann)); 1830 PetscCall(PetscDSGetContext(prob, field, &rctx)); 1831 PetscCall(PetscDSGetConstants(prob, &numConstants, &constants)); 1832 PetscCall(PetscFVGetSpatialDimension(fvm, &dim)); 1833 PetscCall(PetscFVGetNumComponents(fvm, &pdim)); 1834 for (f = 0; f < Nf; ++f) { 1835 (*riemann)(dim, pdim, fgeom[f].centroid, fgeom[f].normal, &uL[f*Nc], &uR[f*Nc], numConstants, constants, flux, rctx); 1836 for (d = 0; d < pdim; ++d) { 1837 fluxL[f*totDim+off+d] = flux[d] / neighborVol[f*2+0]; 1838 fluxR[f*totDim+off+d] = flux[d] / neighborVol[f*2+1]; 1839 } 1840 } 1841 PetscFunctionReturn(0); 1842 } 1843 1844 static PetscErrorCode PetscFVInitialize_Upwind(PetscFV fvm) 1845 { 1846 PetscFunctionBegin; 1847 fvm->ops->setfromoptions = NULL; 1848 fvm->ops->setup = PetscFVSetUp_Upwind; 1849 fvm->ops->view = PetscFVView_Upwind; 1850 fvm->ops->destroy = PetscFVDestroy_Upwind; 1851 fvm->ops->integraterhsfunction = PetscFVIntegrateRHSFunction_Upwind; 1852 PetscFunctionReturn(0); 1853 } 1854 1855 /*MC 1856 PETSCFVUPWIND = "upwind" - A PetscFV object 1857 1858 Level: intermediate 1859 1860 .seealso: PetscFVType, PetscFVCreate(), PetscFVSetType() 1861 M*/ 1862 1863 PETSC_EXTERN PetscErrorCode PetscFVCreate_Upwind(PetscFV fvm) 1864 { 1865 PetscFV_Upwind *b; 1866 1867 PetscFunctionBegin; 1868 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 1869 PetscCall(PetscNewLog(fvm,&b)); 1870 fvm->data = b; 1871 1872 PetscCall(PetscFVInitialize_Upwind(fvm)); 1873 PetscFunctionReturn(0); 1874 } 1875 1876 #include <petscblaslapack.h> 1877 1878 static PetscErrorCode PetscFVDestroy_LeastSquares(PetscFV fvm) 1879 { 1880 PetscFV_LeastSquares *ls = (PetscFV_LeastSquares *) fvm->data; 1881 1882 PetscFunctionBegin; 1883 PetscCall(PetscObjectComposeFunction((PetscObject) fvm, "PetscFVLeastSquaresSetMaxFaces_C", NULL)); 1884 PetscCall(PetscFree4(ls->B, ls->Binv, ls->tau, ls->work)); 1885 PetscCall(PetscFree(ls)); 1886 PetscFunctionReturn(0); 1887 } 1888 1889 static PetscErrorCode PetscFVView_LeastSquares_Ascii(PetscFV fv, PetscViewer viewer) 1890 { 1891 PetscInt Nc, c; 1892 PetscViewerFormat format; 1893 1894 PetscFunctionBegin; 1895 PetscCall(PetscFVGetNumComponents(fv, &Nc)); 1896 PetscCall(PetscViewerGetFormat(viewer, &format)); 1897 PetscCall(PetscViewerASCIIPrintf(viewer, "Finite Volume with Least Squares Reconstruction:\n")); 1898 PetscCall(PetscViewerASCIIPrintf(viewer, " num components: %" PetscInt_FMT "\n", Nc)); 1899 for (c = 0; c < Nc; c++) { 1900 if (fv->componentNames[c]) { 1901 PetscCall(PetscViewerASCIIPrintf(viewer, " component %" PetscInt_FMT ": %s\n", c, fv->componentNames[c])); 1902 } 1903 } 1904 PetscFunctionReturn(0); 1905 } 1906 1907 static PetscErrorCode PetscFVView_LeastSquares(PetscFV fv, PetscViewer viewer) 1908 { 1909 PetscBool iascii; 1910 1911 PetscFunctionBegin; 1912 PetscValidHeaderSpecific(fv, PETSCFV_CLASSID, 1); 1913 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 1914 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii)); 1915 if (iascii) PetscCall(PetscFVView_LeastSquares_Ascii(fv, viewer)); 1916 PetscFunctionReturn(0); 1917 } 1918 1919 static PetscErrorCode PetscFVSetUp_LeastSquares(PetscFV fvm) 1920 { 1921 PetscFunctionBegin; 1922 PetscFunctionReturn(0); 1923 } 1924 1925 /* Overwrites A. Can only handle full-rank problems with m>=n */ 1926 static PetscErrorCode PetscFVLeastSquaresPseudoInverse_Static(PetscInt m,PetscInt mstride,PetscInt n,PetscScalar *A,PetscScalar *Ainv,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 1927 { 1928 PetscBool debug = PETSC_FALSE; 1929 PetscBLASInt M,N,K,lda,ldb,ldwork,info; 1930 PetscScalar *R,*Q,*Aback,Alpha; 1931 1932 PetscFunctionBegin; 1933 if (debug) { 1934 PetscCall(PetscMalloc1(m*n,&Aback)); 1935 PetscCall(PetscArraycpy(Aback,A,m*n)); 1936 } 1937 1938 PetscCall(PetscBLASIntCast(m,&M)); 1939 PetscCall(PetscBLASIntCast(n,&N)); 1940 PetscCall(PetscBLASIntCast(mstride,&lda)); 1941 PetscCall(PetscBLASIntCast(worksize,&ldwork)); 1942 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 1943 PetscStackCallBLAS("LAPACKgeqrf",LAPACKgeqrf_(&M,&N,A,&lda,tau,work,&ldwork,&info)); 1944 PetscCall(PetscFPTrapPop()); 1945 PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_LIB,"xGEQRF error"); 1946 R = A; /* Upper triangular part of A now contains R, the rest contains the elementary reflectors */ 1947 1948 /* Extract an explicit representation of Q */ 1949 Q = Ainv; 1950 PetscCall(PetscArraycpy(Q,A,mstride*n)); 1951 K = N; /* full rank */ 1952 PetscStackCallBLAS("LAPACKorgqr",LAPACKorgqr_(&M,&N,&K,Q,&lda,tau,work,&ldwork,&info)); 1953 PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_LIB,"xORGQR/xUNGQR error"); 1954 1955 /* Compute A^{-T} = (R^{-1} Q^T)^T = Q R^{-T} */ 1956 Alpha = 1.0; 1957 ldb = lda; 1958 BLAStrsm_("Right","Upper","ConjugateTranspose","NotUnitTriangular",&M,&N,&Alpha,R,&lda,Q,&ldb); 1959 /* Ainv is Q, overwritten with inverse */ 1960 1961 if (debug) { /* Check that pseudo-inverse worked */ 1962 PetscScalar Beta = 0.0; 1963 PetscBLASInt ldc; 1964 K = N; 1965 ldc = N; 1966 BLASgemm_("ConjugateTranspose","Normal",&N,&K,&M,&Alpha,Ainv,&lda,Aback,&ldb,&Beta,work,&ldc); 1967 PetscCall(PetscScalarView(n*n,work,PETSC_VIEWER_STDOUT_SELF)); 1968 PetscCall(PetscFree(Aback)); 1969 } 1970 PetscFunctionReturn(0); 1971 } 1972 1973 /* Overwrites A. Can handle degenerate problems and m<n. */ 1974 static PetscErrorCode PetscFVLeastSquaresPseudoInverseSVD_Static(PetscInt m,PetscInt mstride,PetscInt n,PetscScalar *A,PetscScalar *Ainv,PetscScalar *tau,PetscInt worksize,PetscScalar *work) 1975 { 1976 PetscBool debug = PETSC_FALSE; 1977 PetscScalar *Brhs, *Aback; 1978 PetscScalar *tmpwork; 1979 PetscReal rcond; 1980 #if defined (PETSC_USE_COMPLEX) 1981 PetscInt rworkSize; 1982 PetscReal *rwork; 1983 #endif 1984 PetscInt i, j, maxmn; 1985 PetscBLASInt M, N, lda, ldb, ldwork; 1986 PetscBLASInt nrhs, irank, info; 1987 1988 PetscFunctionBegin; 1989 if (debug) { 1990 PetscCall(PetscMalloc1(m*n,&Aback)); 1991 PetscCall(PetscArraycpy(Aback,A,m*n)); 1992 } 1993 1994 /* initialize to identity */ 1995 tmpwork = work; 1996 Brhs = Ainv; 1997 maxmn = PetscMax(m,n); 1998 for (j=0; j<maxmn; j++) { 1999 for (i=0; i<maxmn; i++) Brhs[i + j*maxmn] = 1.0*(i == j); 2000 } 2001 2002 PetscCall(PetscBLASIntCast(m,&M)); 2003 PetscCall(PetscBLASIntCast(n,&N)); 2004 PetscCall(PetscBLASIntCast(mstride,&lda)); 2005 PetscCall(PetscBLASIntCast(maxmn,&ldb)); 2006 PetscCall(PetscBLASIntCast(worksize,&ldwork)); 2007 rcond = -1; 2008 PetscCall(PetscFPTrapPush(PETSC_FP_TRAP_OFF)); 2009 nrhs = M; 2010 #if defined(PETSC_USE_COMPLEX) 2011 rworkSize = 5 * PetscMin(M,N); 2012 PetscCall(PetscMalloc1(rworkSize,&rwork)); 2013 PetscStackCallBLAS("LAPACKgelss",LAPACKgelss_(&M,&N,&nrhs,A,&lda,Brhs,&ldb, (PetscReal *) tau,&rcond,&irank,tmpwork,&ldwork,rwork,&info)); 2014 PetscCall(PetscFPTrapPop()); 2015 PetscCall(PetscFree(rwork)); 2016 #else 2017 nrhs = M; 2018 PetscStackCallBLAS("LAPACKgelss",LAPACKgelss_(&M,&N,&nrhs,A,&lda,Brhs,&ldb, (PetscReal *) tau,&rcond,&irank,tmpwork,&ldwork,&info)); 2019 PetscCall(PetscFPTrapPop()); 2020 #endif 2021 PetscCheck(!info,PETSC_COMM_SELF,PETSC_ERR_LIB,"xGELSS error"); 2022 /* The following check should be turned into a diagnostic as soon as someone wants to do this intentionally */ 2023 PetscCheck(irank >= PetscMin(M,N),PETSC_COMM_SELF,PETSC_ERR_USER,"Rank deficient least squares fit, indicates an isolated cell with two colinear points"); 2024 PetscFunctionReturn(0); 2025 } 2026 2027 #if 0 2028 static PetscErrorCode PetscFVLeastSquaresDebugCell_Static(PetscFV fvm, PetscInt cell, DM dm, DM dmFace, PetscScalar *fgeom, DM dmCell, PetscScalar *cgeom) 2029 { 2030 PetscReal grad[2] = {0, 0}; 2031 const PetscInt *faces; 2032 PetscInt numFaces, f; 2033 2034 PetscFunctionBegin; 2035 PetscCall(DMPlexGetConeSize(dm, cell, &numFaces)); 2036 PetscCall(DMPlexGetCone(dm, cell, &faces)); 2037 for (f = 0; f < numFaces; ++f) { 2038 const PetscInt *fcells; 2039 const CellGeom *cg1; 2040 const FaceGeom *fg; 2041 2042 PetscCall(DMPlexGetSupport(dm, faces[f], &fcells)); 2043 PetscCall(DMPlexPointLocalRead(dmFace, faces[f], fgeom, &fg)); 2044 for (i = 0; i < 2; ++i) { 2045 PetscScalar du; 2046 2047 if (fcells[i] == c) continue; 2048 PetscCall(DMPlexPointLocalRead(dmCell, fcells[i], cgeom, &cg1)); 2049 du = cg1->centroid[0] + 3*cg1->centroid[1] - (cg->centroid[0] + 3*cg->centroid[1]); 2050 grad[0] += fg->grad[!i][0] * du; 2051 grad[1] += fg->grad[!i][1] * du; 2052 } 2053 } 2054 PetscCall(PetscPrintf(PETSC_COMM_SELF, "cell[%d] grad (%g, %g)\n", cell, grad[0], grad[1])); 2055 PetscFunctionReturn(0); 2056 } 2057 #endif 2058 2059 /* 2060 PetscFVComputeGradient - Compute the gradient reconstruction matrix for a given cell 2061 2062 Input Parameters: 2063 + fvm - The PetscFV object 2064 . numFaces - The number of cell faces which are not constrained 2065 . dx - The vector from the cell centroid to the neighboring cell centroid for each face 2066 2067 Level: developer 2068 2069 .seealso: PetscFVCreate() 2070 */ 2071 static PetscErrorCode PetscFVComputeGradient_LeastSquares(PetscFV fvm, PetscInt numFaces, const PetscScalar dx[], PetscScalar grad[]) 2072 { 2073 PetscFV_LeastSquares *ls = (PetscFV_LeastSquares *) fvm->data; 2074 const PetscBool useSVD = PETSC_TRUE; 2075 const PetscInt maxFaces = ls->maxFaces; 2076 PetscInt dim, f, d; 2077 2078 PetscFunctionBegin; 2079 if (numFaces > maxFaces) { 2080 PetscCheck(maxFaces >= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "Reconstruction has not been initialized, call PetscFVLeastSquaresSetMaxFaces()"); 2081 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Number of input faces %" PetscInt_FMT " > %" PetscInt_FMT " maxfaces", numFaces, maxFaces); 2082 } 2083 PetscCall(PetscFVGetSpatialDimension(fvm, &dim)); 2084 for (f = 0; f < numFaces; ++f) { 2085 for (d = 0; d < dim; ++d) ls->B[d*maxFaces+f] = dx[f*dim+d]; 2086 } 2087 /* Overwrites B with garbage, returns Binv in row-major format */ 2088 if (useSVD) { 2089 PetscInt maxmn = PetscMax(numFaces, dim); 2090 PetscCall(PetscFVLeastSquaresPseudoInverseSVD_Static(numFaces, maxFaces, dim, ls->B, ls->Binv, ls->tau, ls->workSize, ls->work)); 2091 /* Binv shaped in column-major, coldim=maxmn.*/ 2092 for (f = 0; f < numFaces; ++f) { 2093 for (d = 0; d < dim; ++d) grad[f*dim+d] = ls->Binv[d + maxmn*f]; 2094 } 2095 } else { 2096 PetscCall(PetscFVLeastSquaresPseudoInverse_Static(numFaces, maxFaces, dim, ls->B, ls->Binv, ls->tau, ls->workSize, ls->work)); 2097 /* Binv shaped in row-major, rowdim=maxFaces.*/ 2098 for (f = 0; f < numFaces; ++f) { 2099 for (d = 0; d < dim; ++d) grad[f*dim+d] = ls->Binv[d*maxFaces + f]; 2100 } 2101 } 2102 PetscFunctionReturn(0); 2103 } 2104 2105 /* 2106 neighborVol[f*2+0] contains the left geom 2107 neighborVol[f*2+1] contains the right geom 2108 */ 2109 static PetscErrorCode PetscFVIntegrateRHSFunction_LeastSquares(PetscFV fvm, PetscDS prob, PetscInt field, PetscInt Nf, PetscFVFaceGeom *fgeom, PetscReal *neighborVol, 2110 PetscScalar uL[], PetscScalar uR[], PetscScalar fluxL[], PetscScalar fluxR[]) 2111 { 2112 void (*riemann)(PetscInt, PetscInt, const PetscReal[], const PetscReal[], const PetscScalar[], const PetscScalar[], PetscInt, const PetscScalar[], PetscScalar[], void *); 2113 void *rctx; 2114 PetscScalar *flux = fvm->fluxWork; 2115 const PetscScalar *constants; 2116 PetscInt dim, numConstants, pdim, Nc, totDim, off, f, d; 2117 2118 PetscFunctionBegin; 2119 PetscCall(PetscDSGetTotalComponents(prob, &Nc)); 2120 PetscCall(PetscDSGetTotalDimension(prob, &totDim)); 2121 PetscCall(PetscDSGetFieldOffset(prob, field, &off)); 2122 PetscCall(PetscDSGetRiemannSolver(prob, field, &riemann)); 2123 PetscCall(PetscDSGetContext(prob, field, &rctx)); 2124 PetscCall(PetscDSGetConstants(prob, &numConstants, &constants)); 2125 PetscCall(PetscFVGetSpatialDimension(fvm, &dim)); 2126 PetscCall(PetscFVGetNumComponents(fvm, &pdim)); 2127 for (f = 0; f < Nf; ++f) { 2128 (*riemann)(dim, pdim, fgeom[f].centroid, fgeom[f].normal, &uL[f*Nc], &uR[f*Nc], numConstants, constants, flux, rctx); 2129 for (d = 0; d < pdim; ++d) { 2130 fluxL[f*totDim+off+d] = flux[d] / neighborVol[f*2+0]; 2131 fluxR[f*totDim+off+d] = flux[d] / neighborVol[f*2+1]; 2132 } 2133 } 2134 PetscFunctionReturn(0); 2135 } 2136 2137 static PetscErrorCode PetscFVLeastSquaresSetMaxFaces_LS(PetscFV fvm, PetscInt maxFaces) 2138 { 2139 PetscFV_LeastSquares *ls = (PetscFV_LeastSquares *) fvm->data; 2140 PetscInt dim,m,n,nrhs,minmn,maxmn; 2141 2142 PetscFunctionBegin; 2143 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 2144 PetscCall(PetscFVGetSpatialDimension(fvm, &dim)); 2145 PetscCall(PetscFree4(ls->B, ls->Binv, ls->tau, ls->work)); 2146 ls->maxFaces = maxFaces; 2147 m = ls->maxFaces; 2148 n = dim; 2149 nrhs = ls->maxFaces; 2150 minmn = PetscMin(m,n); 2151 maxmn = PetscMax(m,n); 2152 ls->workSize = 3*minmn + PetscMax(2*minmn, PetscMax(maxmn, nrhs)); /* required by LAPACK */ 2153 PetscCall(PetscMalloc4(m*n,&ls->B,maxmn*maxmn,&ls->Binv,minmn,&ls->tau,ls->workSize,&ls->work)); 2154 PetscFunctionReturn(0); 2155 } 2156 2157 PetscErrorCode PetscFVInitialize_LeastSquares(PetscFV fvm) 2158 { 2159 PetscFunctionBegin; 2160 fvm->ops->setfromoptions = NULL; 2161 fvm->ops->setup = PetscFVSetUp_LeastSquares; 2162 fvm->ops->view = PetscFVView_LeastSquares; 2163 fvm->ops->destroy = PetscFVDestroy_LeastSquares; 2164 fvm->ops->computegradient = PetscFVComputeGradient_LeastSquares; 2165 fvm->ops->integraterhsfunction = PetscFVIntegrateRHSFunction_LeastSquares; 2166 PetscFunctionReturn(0); 2167 } 2168 2169 /*MC 2170 PETSCFVLEASTSQUARES = "leastsquares" - A PetscFV object 2171 2172 Level: intermediate 2173 2174 .seealso: PetscFVType, PetscFVCreate(), PetscFVSetType() 2175 M*/ 2176 2177 PETSC_EXTERN PetscErrorCode PetscFVCreate_LeastSquares(PetscFV fvm) 2178 { 2179 PetscFV_LeastSquares *ls; 2180 2181 PetscFunctionBegin; 2182 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 2183 PetscCall(PetscNewLog(fvm, &ls)); 2184 fvm->data = ls; 2185 2186 ls->maxFaces = -1; 2187 ls->workSize = -1; 2188 ls->B = NULL; 2189 ls->Binv = NULL; 2190 ls->tau = NULL; 2191 ls->work = NULL; 2192 2193 PetscCall(PetscFVSetComputeGradients(fvm, PETSC_TRUE)); 2194 PetscCall(PetscFVInitialize_LeastSquares(fvm)); 2195 PetscCall(PetscObjectComposeFunction((PetscObject) fvm, "PetscFVLeastSquaresSetMaxFaces_C", PetscFVLeastSquaresSetMaxFaces_LS)); 2196 PetscFunctionReturn(0); 2197 } 2198 2199 /*@ 2200 PetscFVLeastSquaresSetMaxFaces - Set the maximum number of cell faces for gradient reconstruction 2201 2202 Not collective 2203 2204 Input parameters: 2205 + fvm - The PetscFV object 2206 - maxFaces - The maximum number of cell faces 2207 2208 Level: intermediate 2209 2210 .seealso: PetscFVCreate(), PETSCFVLEASTSQUARES 2211 @*/ 2212 PetscErrorCode PetscFVLeastSquaresSetMaxFaces(PetscFV fvm, PetscInt maxFaces) 2213 { 2214 PetscFunctionBegin; 2215 PetscValidHeaderSpecific(fvm, PETSCFV_CLASSID, 1); 2216 PetscTryMethod(fvm, "PetscFVLeastSquaresSetMaxFaces_C", (PetscFV,PetscInt), (fvm,maxFaces)); 2217 PetscFunctionReturn(0); 2218 } 2219