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