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