1 #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 2 #include <petscdmplex.h> 3 4 PetscClassId PETSCDUALSPACE_CLASSID = 0; 5 6 PetscFunctionList PetscDualSpaceList = NULL; 7 PetscBool PetscDualSpaceRegisterAllCalled = PETSC_FALSE; 8 9 const char *const PetscDualSpaceReferenceCells[] = {"SIMPLEX", "TENSOR", "PetscDualSpaceReferenceCell", "PETSCDUALSPACE_REFCELL_",0}; 10 11 /*@C 12 PetscDualSpaceRegister - Adds a new PetscDualSpace implementation 13 14 Not Collective 15 16 Input Parameters: 17 + name - The name of a new user-defined creation routine 18 - create_func - The creation routine itself 19 20 Notes: 21 PetscDualSpaceRegister() may be called multiple times to add several user-defined PetscDualSpaces 22 23 Sample usage: 24 .vb 25 PetscDualSpaceRegister("my_space", MyPetscDualSpaceCreate); 26 .ve 27 28 Then, your PetscDualSpace type can be chosen with the procedural interface via 29 .vb 30 PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *); 31 PetscDualSpaceSetType(PetscDualSpace, "my_dual_space"); 32 .ve 33 or at runtime via the option 34 .vb 35 -petscdualspace_type my_dual_space 36 .ve 37 38 Level: advanced 39 40 .keywords: PetscDualSpace, register 41 .seealso: PetscDualSpaceRegisterAll(), PetscDualSpaceRegisterDestroy() 42 43 @*/ 44 PetscErrorCode PetscDualSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscDualSpace)) 45 { 46 PetscErrorCode ierr; 47 48 PetscFunctionBegin; 49 ierr = PetscFunctionListAdd(&PetscDualSpaceList, sname, function);CHKERRQ(ierr); 50 PetscFunctionReturn(0); 51 } 52 53 /*@C 54 PetscDualSpaceSetType - Builds a particular PetscDualSpace 55 56 Collective on PetscDualSpace 57 58 Input Parameters: 59 + sp - The PetscDualSpace object 60 - name - The kind of space 61 62 Options Database Key: 63 . -petscdualspace_type <type> - Sets the PetscDualSpace type; use -help for a list of available types 64 65 Level: intermediate 66 67 .keywords: PetscDualSpace, set, type 68 .seealso: PetscDualSpaceGetType(), PetscDualSpaceCreate() 69 @*/ 70 PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name) 71 { 72 PetscErrorCode (*r)(PetscDualSpace); 73 PetscBool match; 74 PetscErrorCode ierr; 75 76 PetscFunctionBegin; 77 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 78 ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr); 79 if (match) PetscFunctionReturn(0); 80 81 if (!PetscDualSpaceRegisterAllCalled) {ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);} 82 ierr = PetscFunctionListFind(PetscDualSpaceList, name, &r);CHKERRQ(ierr); 83 if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name); 84 85 if (sp->ops->destroy) { 86 ierr = (*sp->ops->destroy)(sp);CHKERRQ(ierr); 87 sp->ops->destroy = NULL; 88 } 89 ierr = (*r)(sp);CHKERRQ(ierr); 90 ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr); 91 PetscFunctionReturn(0); 92 } 93 94 /*@C 95 PetscDualSpaceGetType - Gets the PetscDualSpace type name (as a string) from the object. 96 97 Not Collective 98 99 Input Parameter: 100 . sp - The PetscDualSpace 101 102 Output Parameter: 103 . name - The PetscDualSpace type name 104 105 Level: intermediate 106 107 .keywords: PetscDualSpace, get, type, name 108 .seealso: PetscDualSpaceSetType(), PetscDualSpaceCreate() 109 @*/ 110 PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name) 111 { 112 PetscErrorCode ierr; 113 114 PetscFunctionBegin; 115 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 116 PetscValidPointer(name, 2); 117 if (!PetscDualSpaceRegisterAllCalled) { 118 ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr); 119 } 120 *name = ((PetscObject) sp)->type_name; 121 PetscFunctionReturn(0); 122 } 123 124 static PetscErrorCode PetscDualSpaceView_ASCII(PetscDualSpace sp, PetscViewer v) 125 { 126 PetscViewerFormat format; 127 PetscInt pdim, f; 128 PetscErrorCode ierr; 129 130 PetscFunctionBegin; 131 ierr = PetscDualSpaceGetDimension(sp, &pdim);CHKERRQ(ierr); 132 ierr = PetscObjectPrintClassNamePrefixType((PetscObject) sp, v);CHKERRQ(ierr); 133 ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 134 ierr = PetscViewerASCIIPrintf(v, "Dual space with %D components, size %D\n", sp->Nc, pdim);CHKERRQ(ierr); 135 if (sp->ops->view) {ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr);} 136 ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr); 137 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 138 ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 139 for (f = 0; f < pdim; ++f) { 140 ierr = PetscViewerASCIIPrintf(v, "Dual basis vector %D\n", f);CHKERRQ(ierr); 141 ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 142 ierr = PetscQuadratureView(sp->functional[f], v);CHKERRQ(ierr); 143 ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 144 } 145 ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 146 } 147 ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 148 PetscFunctionReturn(0); 149 } 150 151 /*@ 152 PetscDualSpaceView - Views a PetscDualSpace 153 154 Collective on PetscDualSpace 155 156 Input Parameter: 157 + sp - the PetscDualSpace object to view 158 - v - the viewer 159 160 Level: developer 161 162 .seealso PetscDualSpaceDestroy() 163 @*/ 164 PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v) 165 { 166 PetscBool iascii; 167 PetscErrorCode ierr; 168 169 PetscFunctionBegin; 170 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 171 if (v) PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2); 172 if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr);} 173 ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 174 if (iascii) {ierr = PetscDualSpaceView_ASCII(sp, v);CHKERRQ(ierr);} 175 PetscFunctionReturn(0); 176 } 177 178 /*@ 179 PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database 180 181 Collective on PetscDualSpace 182 183 Input Parameter: 184 . sp - the PetscDualSpace object to set options for 185 186 Options Database: 187 . -petscspace_degree the approximation order of the space 188 189 Level: developer 190 191 .seealso PetscDualSpaceView() 192 @*/ 193 PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp) 194 { 195 PetscDualSpaceReferenceCell refCell = PETSCDUALSPACE_REFCELL_SIMPLEX; 196 PetscInt refDim = 0; 197 PetscBool flg; 198 const char *defaultType; 199 char name[256]; 200 PetscErrorCode ierr; 201 202 PetscFunctionBegin; 203 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 204 if (!((PetscObject) sp)->type_name) { 205 defaultType = PETSCDUALSPACELAGRANGE; 206 } else { 207 defaultType = ((PetscObject) sp)->type_name; 208 } 209 if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);} 210 211 ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr); 212 ierr = PetscOptionsFList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, name, 256, &flg);CHKERRQ(ierr); 213 if (flg) { 214 ierr = PetscDualSpaceSetType(sp, name);CHKERRQ(ierr); 215 } else if (!((PetscObject) sp)->type_name) { 216 ierr = PetscDualSpaceSetType(sp, defaultType);CHKERRQ(ierr); 217 } 218 ierr = PetscOptionsInt("-petscdualspace_degree", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL);CHKERRQ(ierr); 219 ierr = PetscOptionsInt("-petscdualspace_components", "The number of components", "PetscDualSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL);CHKERRQ(ierr); 220 if (sp->ops->setfromoptions) { 221 ierr = (*sp->ops->setfromoptions)(PetscOptionsObject,sp);CHKERRQ(ierr); 222 } 223 ierr = PetscOptionsInt("-petscdualspace_refdim", "The spatial dimension of the reference cell", "PetscDualSpaceSetReferenceCell", refDim, &refDim, NULL);CHKERRQ(ierr); 224 ierr = PetscOptionsEnum("-petscdualspace_refcell", "Reference cell", "PetscDualSpaceSetReferenceCell", PetscDualSpaceReferenceCells, (PetscEnum) refCell, (PetscEnum *) &refCell, &flg);CHKERRQ(ierr); 225 if (flg) { 226 DM K; 227 228 if (!refDim) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_INCOMP, "Reference cell specified without a dimension. Use -petscdualspace_refdim."); 229 ierr = PetscDualSpaceCreateReferenceCell(sp, refDim, refCell == PETSCDUALSPACE_REFCELL_SIMPLEX ? PETSC_TRUE : PETSC_FALSE, &K);CHKERRQ(ierr); 230 ierr = PetscDualSpaceSetDM(sp, K);CHKERRQ(ierr); 231 ierr = DMDestroy(&K);CHKERRQ(ierr); 232 } 233 234 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 235 ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) sp);CHKERRQ(ierr); 236 ierr = PetscOptionsEnd();CHKERRQ(ierr); 237 sp->setfromoptionscalled = PETSC_TRUE; 238 PetscFunctionReturn(0); 239 } 240 241 /*@ 242 PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace 243 244 Collective on PetscDualSpace 245 246 Input Parameter: 247 . sp - the PetscDualSpace object to setup 248 249 Level: developer 250 251 .seealso PetscDualSpaceView(), PetscDualSpaceDestroy() 252 @*/ 253 PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp) 254 { 255 PetscErrorCode ierr; 256 257 PetscFunctionBegin; 258 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 259 if (sp->setupcalled) PetscFunctionReturn(0); 260 sp->setupcalled = PETSC_TRUE; 261 if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);} 262 if (sp->setfromoptionscalled) {ierr = PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view");CHKERRQ(ierr);} 263 PetscFunctionReturn(0); 264 } 265 266 /*@ 267 PetscDualSpaceDestroy - Destroys a PetscDualSpace object 268 269 Collective on PetscDualSpace 270 271 Input Parameter: 272 . sp - the PetscDualSpace object to destroy 273 274 Level: developer 275 276 .seealso PetscDualSpaceView() 277 @*/ 278 PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp) 279 { 280 PetscInt dim, f; 281 PetscErrorCode ierr; 282 283 PetscFunctionBegin; 284 if (!*sp) PetscFunctionReturn(0); 285 PetscValidHeaderSpecific((*sp), PETSCDUALSPACE_CLASSID, 1); 286 287 if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);} 288 ((PetscObject) (*sp))->refct = 0; 289 290 ierr = PetscDualSpaceGetDimension(*sp, &dim);CHKERRQ(ierr); 291 for (f = 0; f < dim; ++f) { 292 ierr = PetscQuadratureDestroy(&(*sp)->functional[f]);CHKERRQ(ierr); 293 } 294 ierr = PetscFree((*sp)->functional);CHKERRQ(ierr); 295 ierr = PetscQuadratureDestroy(&(*sp)->allPoints);CHKERRQ(ierr); 296 ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr); 297 298 if ((*sp)->ops->destroy) {ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr);} 299 ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr); 300 PetscFunctionReturn(0); 301 } 302 303 /*@ 304 PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType(). 305 306 Collective on MPI_Comm 307 308 Input Parameter: 309 . comm - The communicator for the PetscDualSpace object 310 311 Output Parameter: 312 . sp - The PetscDualSpace object 313 314 Level: beginner 315 316 .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE 317 @*/ 318 PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp) 319 { 320 PetscDualSpace s; 321 PetscErrorCode ierr; 322 323 PetscFunctionBegin; 324 PetscValidPointer(sp, 2); 325 ierr = PetscCitationsRegister(FECitation,&FEcite);CHKERRQ(ierr); 326 *sp = NULL; 327 ierr = PetscFEInitializePackage();CHKERRQ(ierr); 328 329 ierr = PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView);CHKERRQ(ierr); 330 331 s->order = 0; 332 s->Nc = 1; 333 s->setupcalled = PETSC_FALSE; 334 335 *sp = s; 336 PetscFunctionReturn(0); 337 } 338 339 /*@ 340 PetscDualSpaceDuplicate - Creates a duplicate PetscDualSpace object, however it is not setup. 341 342 Collective on PetscDualSpace 343 344 Input Parameter: 345 . sp - The original PetscDualSpace 346 347 Output Parameter: 348 . spNew - The duplicate PetscDualSpace 349 350 Level: beginner 351 352 .seealso: PetscDualSpaceCreate(), PetscDualSpaceSetType() 353 @*/ 354 PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew) 355 { 356 PetscErrorCode ierr; 357 358 PetscFunctionBegin; 359 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 360 PetscValidPointer(spNew, 2); 361 ierr = (*sp->ops->duplicate)(sp, spNew);CHKERRQ(ierr); 362 PetscFunctionReturn(0); 363 } 364 365 /*@ 366 PetscDualSpaceGetDM - Get the DM representing the reference cell 367 368 Not collective 369 370 Input Parameter: 371 . sp - The PetscDualSpace 372 373 Output Parameter: 374 . dm - The reference cell 375 376 Level: intermediate 377 378 .seealso: PetscDualSpaceSetDM(), PetscDualSpaceCreate() 379 @*/ 380 PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm) 381 { 382 PetscFunctionBegin; 383 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 384 PetscValidPointer(dm, 2); 385 *dm = sp->dm; 386 PetscFunctionReturn(0); 387 } 388 389 /*@ 390 PetscDualSpaceSetDM - Get the DM representing the reference cell 391 392 Not collective 393 394 Input Parameters: 395 + sp - The PetscDualSpace 396 - dm - The reference cell 397 398 Level: intermediate 399 400 .seealso: PetscDualSpaceGetDM(), PetscDualSpaceCreate() 401 @*/ 402 PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm) 403 { 404 PetscErrorCode ierr; 405 406 PetscFunctionBegin; 407 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 408 PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 409 ierr = DMDestroy(&sp->dm);CHKERRQ(ierr); 410 ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 411 sp->dm = dm; 412 PetscFunctionReturn(0); 413 } 414 415 /*@ 416 PetscDualSpaceGetOrder - Get the order of the dual space 417 418 Not collective 419 420 Input Parameter: 421 . sp - The PetscDualSpace 422 423 Output Parameter: 424 . order - The order 425 426 Level: intermediate 427 428 .seealso: PetscDualSpaceSetOrder(), PetscDualSpaceCreate() 429 @*/ 430 PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order) 431 { 432 PetscFunctionBegin; 433 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 434 PetscValidPointer(order, 2); 435 *order = sp->order; 436 PetscFunctionReturn(0); 437 } 438 439 /*@ 440 PetscDualSpaceSetOrder - Set the order of the dual space 441 442 Not collective 443 444 Input Parameters: 445 + sp - The PetscDualSpace 446 - order - The order 447 448 Level: intermediate 449 450 .seealso: PetscDualSpaceGetOrder(), PetscDualSpaceCreate() 451 @*/ 452 PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order) 453 { 454 PetscFunctionBegin; 455 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 456 sp->order = order; 457 PetscFunctionReturn(0); 458 } 459 460 /*@ 461 PetscDualSpaceGetNumComponents - Return the number of components for this space 462 463 Input Parameter: 464 . sp - The PetscDualSpace 465 466 Output Parameter: 467 . Nc - The number of components 468 469 Note: A vector space, for example, will have d components, where d is the spatial dimension 470 471 Level: intermediate 472 473 .seealso: PetscDualSpaceSetNumComponents(), PetscDualSpaceGetDimension(), PetscDualSpaceCreate(), PetscDualSpace 474 @*/ 475 PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc) 476 { 477 PetscFunctionBegin; 478 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 479 PetscValidPointer(Nc, 2); 480 *Nc = sp->Nc; 481 PetscFunctionReturn(0); 482 } 483 484 /*@ 485 PetscDualSpaceSetNumComponents - Set the number of components for this space 486 487 Input Parameters: 488 + sp - The PetscDualSpace 489 - order - The number of components 490 491 Level: intermediate 492 493 .seealso: PetscDualSpaceGetNumComponents(), PetscDualSpaceCreate(), PetscDualSpace 494 @*/ 495 PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc) 496 { 497 PetscFunctionBegin; 498 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 499 sp->Nc = Nc; 500 PetscFunctionReturn(0); 501 } 502 503 /*@ 504 PetscDualSpaceLagrangeGetTensor - Get the tensor nature of the dual space 505 506 Not collective 507 508 Input Parameter: 509 . sp - The PetscDualSpace 510 511 Output Parameter: 512 . tensor - Whether the dual space has tensor layout (vs. simplicial) 513 514 Level: intermediate 515 516 .seealso: PetscDualSpaceLagrangeSetTensor(), PetscDualSpaceCreate() 517 @*/ 518 PetscErrorCode PetscDualSpaceLagrangeGetTensor(PetscDualSpace sp, PetscBool *tensor) 519 { 520 PetscErrorCode ierr; 521 522 PetscFunctionBegin; 523 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 524 PetscValidPointer(tensor, 2); 525 ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeGetTensor_C",(PetscDualSpace,PetscBool *),(sp,tensor));CHKERRQ(ierr); 526 PetscFunctionReturn(0); 527 } 528 529 /*@ 530 PetscDualSpaceLagrangeSetTensor - Set the tensor nature of the dual space 531 532 Not collective 533 534 Input Parameters: 535 + sp - The PetscDualSpace 536 - tensor - Whether the dual space has tensor layout (vs. simplicial) 537 538 Level: intermediate 539 540 .seealso: PetscDualSpaceLagrangeGetTensor(), PetscDualSpaceCreate() 541 @*/ 542 PetscErrorCode PetscDualSpaceLagrangeSetTensor(PetscDualSpace sp, PetscBool tensor) 543 { 544 PetscErrorCode ierr; 545 546 PetscFunctionBegin; 547 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 548 ierr = PetscTryMethod(sp,"PetscDualSpaceLagrangeSetTensor_C",(PetscDualSpace,PetscBool),(sp,tensor));CHKERRQ(ierr); 549 PetscFunctionReturn(0); 550 } 551 552 /*@ 553 PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space 554 555 Not collective 556 557 Input Parameters: 558 + sp - The PetscDualSpace 559 - i - The basis number 560 561 Output Parameter: 562 . functional - The basis functional 563 564 Level: intermediate 565 566 .seealso: PetscDualSpaceGetDimension(), PetscDualSpaceCreate() 567 @*/ 568 PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional) 569 { 570 PetscInt dim; 571 PetscErrorCode ierr; 572 573 PetscFunctionBegin; 574 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 575 PetscValidPointer(functional, 3); 576 ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr); 577 if ((i < 0) || (i >= dim)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %d must be in [0, %d)", i, dim); 578 *functional = sp->functional[i]; 579 PetscFunctionReturn(0); 580 } 581 582 /*@ 583 PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals 584 585 Not collective 586 587 Input Parameter: 588 . sp - The PetscDualSpace 589 590 Output Parameter: 591 . dim - The dimension 592 593 Level: intermediate 594 595 .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate() 596 @*/ 597 PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim) 598 { 599 PetscErrorCode ierr; 600 601 PetscFunctionBegin; 602 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 603 PetscValidPointer(dim, 2); 604 *dim = 0; 605 if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);} 606 PetscFunctionReturn(0); 607 } 608 609 /*@C 610 PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension 611 612 Not collective 613 614 Input Parameter: 615 . sp - The PetscDualSpace 616 617 Output Parameter: 618 . numDof - An array of length dim+1 which holds the number of dofs for each dimension 619 620 Level: intermediate 621 622 .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate() 623 @*/ 624 PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt **numDof) 625 { 626 PetscErrorCode ierr; 627 628 PetscFunctionBegin; 629 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 630 PetscValidPointer(numDof, 2); 631 ierr = (*sp->ops->getnumdof)(sp, numDof);CHKERRQ(ierr); 632 if (!*numDof) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation"); 633 PetscFunctionReturn(0); 634 } 635 636 PetscErrorCode PetscDualSpaceCreateSection(PetscDualSpace sp, PetscSection *section) 637 { 638 DM dm; 639 PetscInt pStart, pEnd, depth, h, offset; 640 const PetscInt *numDof; 641 PetscErrorCode ierr; 642 643 PetscFunctionBegin; 644 ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr); 645 ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 646 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)sp),section);CHKERRQ(ierr); 647 ierr = PetscSectionSetChart(*section,pStart,pEnd);CHKERRQ(ierr); 648 ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 649 ierr = PetscDualSpaceGetNumDof(sp,&numDof);CHKERRQ(ierr); 650 for (h = 0; h <= depth; h++) { 651 PetscInt hStart, hEnd, p, dof; 652 653 ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr); 654 dof = numDof[depth - h]; 655 for (p = hStart; p < hEnd; p++) { 656 ierr = PetscSectionSetDof(*section,p,dof);CHKERRQ(ierr); 657 } 658 } 659 ierr = PetscSectionSetUp(*section);CHKERRQ(ierr); 660 for (h = 0, offset = 0; h <= depth; h++) { 661 PetscInt hStart, hEnd, p, dof; 662 663 ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr); 664 dof = numDof[depth - h]; 665 for (p = hStart; p < hEnd; p++) { 666 ierr = PetscSectionGetDof(*section,p,&dof);CHKERRQ(ierr); 667 ierr = PetscSectionSetOffset(*section,p,offset);CHKERRQ(ierr); 668 offset += dof; 669 } 670 } 671 PetscFunctionReturn(0); 672 } 673 674 /*@ 675 PetscDualSpaceCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell 676 677 Collective on PetscDualSpace 678 679 Input Parameters: 680 + sp - The PetscDualSpace 681 . dim - The spatial dimension 682 - simplex - Flag for simplex, otherwise use a tensor-product cell 683 684 Output Parameter: 685 . refdm - The reference cell 686 687 Level: advanced 688 689 .keywords: PetscDualSpace, reference cell 690 .seealso: PetscDualSpaceCreate(), DMPLEX 691 @*/ 692 PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace sp, PetscInt dim, PetscBool simplex, DM *refdm) 693 { 694 PetscErrorCode ierr; 695 696 PetscFunctionBeginUser; 697 ierr = DMPlexCreateReferenceCell(PetscObjectComm((PetscObject) sp), dim, simplex, refdm);CHKERRQ(ierr); 698 PetscFunctionReturn(0); 699 } 700 701 /*@C 702 PetscDualSpaceApply - Apply a functional from the dual space basis to an input function 703 704 Input Parameters: 705 + sp - The PetscDualSpace object 706 . f - The basis functional index 707 . time - The time 708 . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian) (or evaluated at the coordinates of the functional) 709 . numComp - The number of components for the function 710 . func - The input function 711 - ctx - A context for the function 712 713 Output Parameter: 714 . value - numComp output values 715 716 Note: The calling sequence for the callback func is given by: 717 718 $ func(PetscInt dim, PetscReal time, const PetscReal x[], 719 $ PetscInt numComponents, PetscScalar values[], void *ctx) 720 721 Level: developer 722 723 .seealso: PetscDualSpaceCreate() 724 @*/ 725 PetscErrorCode PetscDualSpaceApply(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFEGeom *cgeom, PetscInt numComp, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx, PetscScalar *value) 726 { 727 PetscErrorCode ierr; 728 729 PetscFunctionBegin; 730 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 731 PetscValidPointer(cgeom, 4); 732 PetscValidPointer(value, 8); 733 ierr = (*sp->ops->apply)(sp, f, time, cgeom, numComp, func, ctx, value);CHKERRQ(ierr); 734 PetscFunctionReturn(0); 735 } 736 737 /*@C 738 PetscDualSpaceApplyAll - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints() 739 740 Input Parameters: 741 + sp - The PetscDualSpace object 742 - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints() 743 744 Output Parameter: 745 . spValue - The values of all dual space functionals 746 747 Level: developer 748 749 .seealso: PetscDualSpaceCreate() 750 @*/ 751 PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue) 752 { 753 PetscErrorCode ierr; 754 755 PetscFunctionBegin; 756 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 757 ierr = (*sp->ops->applyall)(sp, pointEval, spValue);CHKERRQ(ierr); 758 PetscFunctionReturn(0); 759 } 760 761 /*@C 762 PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional. 763 764 Input Parameters: 765 + sp - The PetscDualSpace object 766 . f - The basis functional index 767 . time - The time 768 . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian) 769 . Nc - The number of components for the function 770 . func - The input function 771 - ctx - A context for the function 772 773 Output Parameter: 774 . value - The output value 775 776 Note: The calling sequence for the callback func is given by: 777 778 $ func(PetscInt dim, PetscReal time, const PetscReal x[], 779 $ PetscInt numComponents, PetscScalar values[], void *ctx) 780 781 and the idea is to evaluate the functional as an integral 782 783 $ n(f) = int dx n(x) . f(x) 784 785 where both n and f have Nc components. 786 787 Level: developer 788 789 .seealso: PetscDualSpaceCreate() 790 @*/ 791 PetscErrorCode PetscDualSpaceApplyDefault(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFEGeom *cgeom, PetscInt Nc, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx, PetscScalar *value) 792 { 793 DM dm; 794 PetscQuadrature n; 795 const PetscReal *points, *weights; 796 PetscReal x[3]; 797 PetscScalar *val; 798 PetscInt dim, dE, qNc, c, Nq, q; 799 PetscBool isAffine; 800 PetscErrorCode ierr; 801 802 PetscFunctionBegin; 803 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 804 PetscValidPointer(value, 5); 805 ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 806 ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr); 807 ierr = PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights);CHKERRQ(ierr); 808 if (dim != cgeom->dim) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature spatial dimension %D != cell geometry dimension %D", dim, cgeom->dim); 809 if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc); 810 ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr); 811 *value = 0.0; 812 isAffine = cgeom->isAffine; 813 dE = cgeom->dimEmbed; 814 for (q = 0; q < Nq; ++q) { 815 if (isAffine) { 816 CoordinatesRefToReal(dE, cgeom->dim, cgeom->xi, cgeom->v, cgeom->J, &points[q*dim], x); 817 ierr = (*func)(dE, time, x, Nc, val, ctx);CHKERRQ(ierr); 818 } else { 819 ierr = (*func)(dE, time, &cgeom->v[dE*q], Nc, val, ctx);CHKERRQ(ierr); 820 } 821 for (c = 0; c < Nc; ++c) { 822 *value += val[c]*weights[q*Nc+c]; 823 } 824 } 825 ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr); 826 PetscFunctionReturn(0); 827 } 828 829 /*@C 830 PetscDualSpaceApplyAllDefault - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints() 831 832 Input Parameters: 833 + sp - The PetscDualSpace object 834 - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints() 835 836 Output Parameter: 837 . spValue - The values of all dual space functionals 838 839 Level: developer 840 841 .seealso: PetscDualSpaceCreate() 842 @*/ 843 PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue) 844 { 845 PetscQuadrature n; 846 const PetscReal *points, *weights; 847 PetscInt qNc, c, Nq, q, f, spdim, Nc; 848 PetscInt offset; 849 PetscErrorCode ierr; 850 851 PetscFunctionBegin; 852 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 853 PetscValidScalarPointer(pointEval, 2); 854 PetscValidScalarPointer(spValue, 5); 855 ierr = PetscDualSpaceGetDimension(sp, &spdim);CHKERRQ(ierr); 856 ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr); 857 for (f = 0, offset = 0; f < spdim; f++) { 858 ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr); 859 ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr); 860 if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc); 861 spValue[f] = 0.0; 862 for (q = 0; q < Nq; ++q) { 863 for (c = 0; c < Nc; ++c) { 864 spValue[f] += pointEval[offset++]*weights[q*Nc+c]; 865 } 866 } 867 } 868 PetscFunctionReturn(0); 869 } 870 871 PetscErrorCode PetscDualSpaceGetAllPoints(PetscDualSpace sp, PetscQuadrature *allPoints) 872 { 873 PetscErrorCode ierr; 874 875 PetscFunctionBegin; 876 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 877 PetscValidPointer(allPoints,2); 878 if (!sp->allPoints && sp->ops->createallpoints) { 879 ierr = (*sp->ops->createallpoints)(sp,&sp->allPoints);CHKERRQ(ierr); 880 } 881 *allPoints = sp->allPoints; 882 PetscFunctionReturn(0); 883 } 884 885 PetscErrorCode PetscDualSpaceCreateAllPointsDefault(PetscDualSpace sp, PetscQuadrature *allPoints) 886 { 887 PetscInt spdim; 888 PetscInt numPoints, offset; 889 PetscReal *points; 890 PetscInt f, dim; 891 PetscQuadrature q; 892 PetscErrorCode ierr; 893 894 PetscFunctionBegin; 895 ierr = PetscDualSpaceGetDimension(sp,&spdim);CHKERRQ(ierr); 896 if (!spdim) { 897 ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr); 898 ierr = PetscQuadratureSetData(*allPoints,0,0,0,NULL,NULL);CHKERRQ(ierr); 899 } 900 ierr = PetscDualSpaceGetFunctional(sp,0,&q);CHKERRQ(ierr); 901 ierr = PetscQuadratureGetData(q,&dim,NULL,&numPoints,NULL,NULL);CHKERRQ(ierr); 902 for (f = 1; f < spdim; f++) { 903 PetscInt Np; 904 905 ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr); 906 ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL);CHKERRQ(ierr); 907 numPoints += Np; 908 } 909 ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr); 910 for (f = 0, offset = 0; f < spdim; f++) { 911 const PetscReal *p; 912 PetscInt Np, i; 913 914 ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr); 915 ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,&p,NULL);CHKERRQ(ierr); 916 for (i = 0; i < Np * dim; i++) { 917 points[offset + i] = p[i]; 918 } 919 offset += Np * dim; 920 } 921 ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr); 922 ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr); 923 PetscFunctionReturn(0); 924 } 925 926 /*@C 927 PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid. 928 929 Input Parameters: 930 + sp - The PetscDualSpace object 931 . f - The basis functional index 932 . time - The time 933 . cgeom - A context with geometric information for this cell, we currently just use the centroid 934 . Nc - The number of components for the function 935 . func - The input function 936 - ctx - A context for the function 937 938 Output Parameter: 939 . value - The output value (scalar) 940 941 Note: The calling sequence for the callback func is given by: 942 943 $ func(PetscInt dim, PetscReal time, const PetscReal x[], 944 $ PetscInt numComponents, PetscScalar values[], void *ctx) 945 946 and the idea is to evaluate the functional as an integral 947 948 $ n(f) = int dx n(x) . f(x) 949 950 where both n and f have Nc components. 951 952 Level: developer 953 954 .seealso: PetscDualSpaceCreate() 955 @*/ 956 PetscErrorCode PetscDualSpaceApplyFVM(PetscDualSpace sp, PetscInt f, PetscReal time, PetscFVCellGeom *cgeom, PetscInt Nc, PetscErrorCode (*func)(PetscInt, PetscReal, const PetscReal [], PetscInt, PetscScalar *, void *), void *ctx, PetscScalar *value) 957 { 958 DM dm; 959 PetscQuadrature n; 960 const PetscReal *points, *weights; 961 PetscScalar *val; 962 PetscInt dimEmbed, qNc, c, Nq, q; 963 PetscErrorCode ierr; 964 965 PetscFunctionBegin; 966 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 967 PetscValidPointer(value, 5); 968 ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 969 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 970 ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr); 971 ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr); 972 if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc); 973 ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr); 974 *value = 0.; 975 for (q = 0; q < Nq; ++q) { 976 ierr = (*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx);CHKERRQ(ierr); 977 for (c = 0; c < Nc; ++c) { 978 *value += val[c]*weights[q*Nc+c]; 979 } 980 } 981 ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr); 982 PetscFunctionReturn(0); 983 } 984 985 /*@ 986 PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a 987 given height. This assumes that the reference cell is symmetric over points of this height. 988 989 If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and 990 pointwise values are not defined on the element boundaries), or if the implementation of PetscDualSpace does not 991 support extracting subspaces, then NULL is returned. 992 993 This does not increment the reference count on the returned dual space, and the user should not destroy it. 994 995 Not collective 996 997 Input Parameters: 998 + sp - the PetscDualSpace object 999 - height - the height of the mesh point for which the subspace is desired 1000 1001 Output Parameter: 1002 . subsp - the subspace. Note that the functionals in the subspace are with respect to the intrinsic geometry of the 1003 point, which will be of lesser dimension if height > 0. 1004 1005 Level: advanced 1006 1007 .seealso: PetscSpaceGetHeightSubspace(), PetscDualSpace 1008 @*/ 1009 PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp) 1010 { 1011 PetscErrorCode ierr; 1012 1013 PetscFunctionBegin; 1014 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1015 PetscValidPointer(subsp, 3); 1016 *subsp = NULL; 1017 if (sp->ops->getheightsubspace) { 1018 ierr = (*sp->ops->getheightsubspace)(sp, height, subsp);CHKERRQ(ierr); 1019 } 1020 PetscFunctionReturn(0); 1021 } 1022 1023 /*@ 1024 PetscDualSpaceGetPointSubspace - Get the subset of the dual space basis that is supported on a particular mesh point. 1025 1026 If the dual space is not defined on the mesh point (e.g. if the space is discontinuous and pointwise values are not 1027 defined on the element boundaries), or if the implementation of PetscDualSpace does not support extracting 1028 subspaces, then NULL is returned. 1029 1030 This does not increment the reference count on the returned dual space, and the user should not destroy it. 1031 1032 Not collective 1033 1034 Input Parameters: 1035 + sp - the PetscDualSpace object 1036 - point - the point (in the dual space's DM) for which the subspace is desired 1037 1038 Output Parameters: 1039 bdsp - the subspace. Note that the functionals in the subspace are with respect to the intrinsic geometry of the 1040 point, which will be of lesser dimension if height > 0. 1041 1042 Level: advanced 1043 1044 .seealso: PetscDualSpace 1045 @*/ 1046 PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace sp, PetscInt point, PetscDualSpace *bdsp) 1047 { 1048 PetscErrorCode ierr; 1049 1050 PetscFunctionBegin; 1051 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1052 PetscValidPointer(bdsp,2); 1053 *bdsp = NULL; 1054 if (sp->ops->getpointsubspace) { 1055 ierr = (*sp->ops->getpointsubspace)(sp,point,bdsp);CHKERRQ(ierr); 1056 } else if (sp->ops->getheightsubspace) { 1057 DM dm; 1058 DMLabel label; 1059 PetscInt dim, depth, height; 1060 1061 ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr); 1062 ierr = DMPlexGetDepth(dm,&dim);CHKERRQ(ierr); 1063 ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr); 1064 ierr = DMLabelGetValue(label,point,&depth);CHKERRQ(ierr); 1065 height = dim - depth; 1066 ierr = (*sp->ops->getheightsubspace)(sp,height,bdsp);CHKERRQ(ierr); 1067 } 1068 PetscFunctionReturn(0); 1069 } 1070 1071