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 /* 12 PetscDualSpaceLatticePointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that sum up to at most 'max'. 13 Ordering is lexicographic with lowest index as least significant in ordering. 14 e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,0}. 15 16 Input Parameters: 17 + len - The length of the tuple 18 . max - The maximum sum 19 - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition 20 21 Output Parameter: 22 . tup - A tuple of len integers whos sum is at most 'max' 23 24 Level: developer 25 26 .seealso: PetscDualSpaceTensorPointLexicographic_Internal() 27 */ 28 PetscErrorCode PetscDualSpaceLatticePointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[]) 29 { 30 PetscFunctionBegin; 31 while (len--) { 32 max -= tup[len]; 33 if (!max) { 34 tup[len] = 0; 35 break; 36 } 37 } 38 tup[++len]++; 39 PetscFunctionReturn(0); 40 } 41 42 /* 43 PetscDualSpaceTensorPointLexicographic_Internal - Returns all tuples of size 'len' with nonnegative integers that are all less than or equal to 'max'. 44 Ordering is lexicographic with lowest index as least significant in ordering. 45 e.g. for len == 2 and max == 2, this will return, in order, {0,0}, {1,0}, {2,0}, {0,1}, {1,1}, {2,1}, {0,2}, {1,2}, {2,2}. 46 47 Input Parameters: 48 + len - The length of the tuple 49 . max - The maximum value 50 - tup - A tuple of length len+1: tup[len] > 0 indicates a stopping condition 51 52 Output Parameter: 53 . tup - A tuple of len integers whos sum is at most 'max' 54 55 Level: developer 56 57 .seealso: PetscDualSpaceLatticePointLexicographic_Internal() 58 */ 59 PetscErrorCode PetscDualSpaceTensorPointLexicographic_Internal(PetscInt len, PetscInt max, PetscInt tup[]) 60 { 61 PetscInt i; 62 63 PetscFunctionBegin; 64 for (i = 0; i < len; i++) { 65 if (tup[i] < max) { 66 break; 67 } else { 68 tup[i] = 0; 69 } 70 } 71 tup[i]++; 72 PetscFunctionReturn(0); 73 } 74 75 /*@C 76 PetscDualSpaceRegister - Adds a new PetscDualSpace implementation 77 78 Not Collective 79 80 Input Parameters: 81 + name - The name of a new user-defined creation routine 82 - create_func - The creation routine itself 83 84 Notes: 85 PetscDualSpaceRegister() may be called multiple times to add several user-defined PetscDualSpaces 86 87 Sample usage: 88 .vb 89 PetscDualSpaceRegister("my_space", MyPetscDualSpaceCreate); 90 .ve 91 92 Then, your PetscDualSpace type can be chosen with the procedural interface via 93 .vb 94 PetscDualSpaceCreate(MPI_Comm, PetscDualSpace *); 95 PetscDualSpaceSetType(PetscDualSpace, "my_dual_space"); 96 .ve 97 or at runtime via the option 98 .vb 99 -petscdualspace_type my_dual_space 100 .ve 101 102 Level: advanced 103 104 .seealso: PetscDualSpaceRegisterAll(), PetscDualSpaceRegisterDestroy() 105 106 @*/ 107 PetscErrorCode PetscDualSpaceRegister(const char sname[], PetscErrorCode (*function)(PetscDualSpace)) 108 { 109 PetscErrorCode ierr; 110 111 PetscFunctionBegin; 112 ierr = PetscFunctionListAdd(&PetscDualSpaceList, sname, function);CHKERRQ(ierr); 113 PetscFunctionReturn(0); 114 } 115 116 /*@C 117 PetscDualSpaceSetType - Builds a particular PetscDualSpace 118 119 Collective on sp 120 121 Input Parameters: 122 + sp - The PetscDualSpace object 123 - name - The kind of space 124 125 Options Database Key: 126 . -petscdualspace_type <type> - Sets the PetscDualSpace type; use -help for a list of available types 127 128 Level: intermediate 129 130 .seealso: PetscDualSpaceGetType(), PetscDualSpaceCreate() 131 @*/ 132 PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name) 133 { 134 PetscErrorCode (*r)(PetscDualSpace); 135 PetscBool match; 136 PetscErrorCode ierr; 137 138 PetscFunctionBegin; 139 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 140 ierr = PetscObjectTypeCompare((PetscObject) sp, name, &match);CHKERRQ(ierr); 141 if (match) PetscFunctionReturn(0); 142 143 if (!PetscDualSpaceRegisterAllCalled) {ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr);} 144 ierr = PetscFunctionListFind(PetscDualSpaceList, name, &r);CHKERRQ(ierr); 145 if (!r) SETERRQ1(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name); 146 147 if (sp->ops->destroy) { 148 ierr = (*sp->ops->destroy)(sp);CHKERRQ(ierr); 149 sp->ops->destroy = NULL; 150 } 151 ierr = (*r)(sp);CHKERRQ(ierr); 152 ierr = PetscObjectChangeTypeName((PetscObject) sp, name);CHKERRQ(ierr); 153 PetscFunctionReturn(0); 154 } 155 156 /*@C 157 PetscDualSpaceGetType - Gets the PetscDualSpace type name (as a string) from the object. 158 159 Not Collective 160 161 Input Parameter: 162 . sp - The PetscDualSpace 163 164 Output Parameter: 165 . name - The PetscDualSpace type name 166 167 Level: intermediate 168 169 .seealso: PetscDualSpaceSetType(), PetscDualSpaceCreate() 170 @*/ 171 PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name) 172 { 173 PetscErrorCode ierr; 174 175 PetscFunctionBegin; 176 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 177 PetscValidPointer(name, 2); 178 if (!PetscDualSpaceRegisterAllCalled) { 179 ierr = PetscDualSpaceRegisterAll();CHKERRQ(ierr); 180 } 181 *name = ((PetscObject) sp)->type_name; 182 PetscFunctionReturn(0); 183 } 184 185 static PetscErrorCode PetscDualSpaceView_ASCII(PetscDualSpace sp, PetscViewer v) 186 { 187 PetscViewerFormat format; 188 PetscInt pdim, f; 189 PetscErrorCode ierr; 190 191 PetscFunctionBegin; 192 ierr = PetscDualSpaceGetDimension(sp, &pdim);CHKERRQ(ierr); 193 ierr = PetscObjectPrintClassNamePrefixType((PetscObject) sp, v);CHKERRQ(ierr); 194 ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 195 ierr = PetscViewerASCIIPrintf(v, "Dual space with %D components, size %D\n", sp->Nc, pdim);CHKERRQ(ierr); 196 if (sp->ops->view) {ierr = (*sp->ops->view)(sp, v);CHKERRQ(ierr);} 197 ierr = PetscViewerGetFormat(v, &format);CHKERRQ(ierr); 198 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 199 ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 200 for (f = 0; f < pdim; ++f) { 201 ierr = PetscViewerASCIIPrintf(v, "Dual basis vector %D\n", f);CHKERRQ(ierr); 202 ierr = PetscViewerASCIIPushTab(v);CHKERRQ(ierr); 203 ierr = PetscQuadratureView(sp->functional[f], v);CHKERRQ(ierr); 204 ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 205 } 206 ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 207 } 208 ierr = PetscViewerASCIIPopTab(v);CHKERRQ(ierr); 209 PetscFunctionReturn(0); 210 } 211 212 /*@ 213 PetscDualSpaceView - Views a PetscDualSpace 214 215 Collective on sp 216 217 Input Parameter: 218 + sp - the PetscDualSpace object to view 219 - v - the viewer 220 221 Level: developer 222 223 .seealso PetscDualSpaceDestroy() 224 @*/ 225 PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v) 226 { 227 PetscBool iascii; 228 PetscErrorCode ierr; 229 230 PetscFunctionBegin; 231 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 232 if (v) PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2); 233 if (!v) {ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v);CHKERRQ(ierr);} 234 ierr = PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 235 if (iascii) {ierr = PetscDualSpaceView_ASCII(sp, v);CHKERRQ(ierr);} 236 PetscFunctionReturn(0); 237 } 238 239 /*@ 240 PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database 241 242 Collective on sp 243 244 Input Parameter: 245 . sp - the PetscDualSpace object to set options for 246 247 Options Database: 248 . -petscspace_degree the approximation order of the space 249 250 Level: developer 251 252 .seealso PetscDualSpaceView() 253 @*/ 254 PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp) 255 { 256 PetscDualSpaceReferenceCell refCell = PETSCDUALSPACE_REFCELL_SIMPLEX; 257 PetscInt refDim = 0; 258 PetscBool flg; 259 const char *defaultType; 260 char name[256]; 261 PetscErrorCode ierr; 262 263 PetscFunctionBegin; 264 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 265 if (!((PetscObject) sp)->type_name) { 266 defaultType = PETSCDUALSPACELAGRANGE; 267 } else { 268 defaultType = ((PetscObject) sp)->type_name; 269 } 270 if (!PetscSpaceRegisterAllCalled) {ierr = PetscSpaceRegisterAll();CHKERRQ(ierr);} 271 272 ierr = PetscObjectOptionsBegin((PetscObject) sp);CHKERRQ(ierr); 273 ierr = PetscOptionsFList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, name, 256, &flg);CHKERRQ(ierr); 274 if (flg) { 275 ierr = PetscDualSpaceSetType(sp, name);CHKERRQ(ierr); 276 } else if (!((PetscObject) sp)->type_name) { 277 ierr = PetscDualSpaceSetType(sp, defaultType);CHKERRQ(ierr); 278 } 279 ierr = PetscOptionsBoundedInt("-petscdualspace_degree", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL,0);CHKERRQ(ierr); 280 ierr = PetscOptionsBoundedInt("-petscdualspace_components", "The number of components", "PetscDualSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL,1);CHKERRQ(ierr); 281 if (sp->ops->setfromoptions) { 282 ierr = (*sp->ops->setfromoptions)(PetscOptionsObject,sp);CHKERRQ(ierr); 283 } 284 ierr = PetscOptionsBoundedInt("-petscdualspace_refdim", "The spatial dimension of the reference cell", "PetscDualSpaceSetReferenceCell", refDim, &refDim, NULL,0);CHKERRQ(ierr); 285 ierr = PetscOptionsEnum("-petscdualspace_refcell", "Reference cell", "PetscDualSpaceSetReferenceCell", PetscDualSpaceReferenceCells, (PetscEnum) refCell, (PetscEnum *) &refCell, &flg);CHKERRQ(ierr); 286 if (flg) { 287 DM K; 288 289 if (!refDim) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_INCOMP, "Reference cell specified without a dimension. Use -petscdualspace_refdim."); 290 ierr = PetscDualSpaceCreateReferenceCell(sp, refDim, refCell == PETSCDUALSPACE_REFCELL_SIMPLEX ? PETSC_TRUE : PETSC_FALSE, &K);CHKERRQ(ierr); 291 ierr = PetscDualSpaceSetDM(sp, K);CHKERRQ(ierr); 292 ierr = DMDestroy(&K);CHKERRQ(ierr); 293 } 294 295 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 296 ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) sp);CHKERRQ(ierr); 297 ierr = PetscOptionsEnd();CHKERRQ(ierr); 298 sp->setfromoptionscalled = PETSC_TRUE; 299 PetscFunctionReturn(0); 300 } 301 302 /*@ 303 PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace 304 305 Collective on sp 306 307 Input Parameter: 308 . sp - the PetscDualSpace object to setup 309 310 Level: developer 311 312 .seealso PetscDualSpaceView(), PetscDualSpaceDestroy() 313 @*/ 314 PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp) 315 { 316 PetscErrorCode ierr; 317 318 PetscFunctionBegin; 319 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 320 if (sp->setupcalled) PetscFunctionReturn(0); 321 sp->setupcalled = PETSC_TRUE; 322 if (sp->ops->setup) {ierr = (*sp->ops->setup)(sp);CHKERRQ(ierr);} 323 if (sp->setfromoptionscalled) {ierr = PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view");CHKERRQ(ierr);} 324 PetscFunctionReturn(0); 325 } 326 327 /*@ 328 PetscDualSpaceDestroy - Destroys a PetscDualSpace object 329 330 Collective on sp 331 332 Input Parameter: 333 . sp - the PetscDualSpace object to destroy 334 335 Level: developer 336 337 .seealso PetscDualSpaceView() 338 @*/ 339 PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp) 340 { 341 PetscInt dim, f; 342 PetscErrorCode ierr; 343 344 PetscFunctionBegin; 345 if (!*sp) PetscFunctionReturn(0); 346 PetscValidHeaderSpecific((*sp), PETSCDUALSPACE_CLASSID, 1); 347 348 if (--((PetscObject)(*sp))->refct > 0) {*sp = 0; PetscFunctionReturn(0);} 349 ((PetscObject) (*sp))->refct = 0; 350 351 ierr = PetscDualSpaceGetDimension(*sp, &dim);CHKERRQ(ierr); 352 for (f = 0; f < dim; ++f) { 353 ierr = PetscQuadratureDestroy(&(*sp)->functional[f]);CHKERRQ(ierr); 354 } 355 ierr = PetscFree((*sp)->functional);CHKERRQ(ierr); 356 ierr = PetscQuadratureDestroy(&(*sp)->allPoints);CHKERRQ(ierr); 357 ierr = DMDestroy(&(*sp)->dm);CHKERRQ(ierr); 358 359 if ((*sp)->ops->destroy) {ierr = (*(*sp)->ops->destroy)(*sp);CHKERRQ(ierr);} 360 ierr = PetscHeaderDestroy(sp);CHKERRQ(ierr); 361 PetscFunctionReturn(0); 362 } 363 364 /*@ 365 PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType(). 366 367 Collective 368 369 Input Parameter: 370 . comm - The communicator for the PetscDualSpace object 371 372 Output Parameter: 373 . sp - The PetscDualSpace object 374 375 Level: beginner 376 377 .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE 378 @*/ 379 PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp) 380 { 381 PetscDualSpace s; 382 PetscErrorCode ierr; 383 384 PetscFunctionBegin; 385 PetscValidPointer(sp, 2); 386 ierr = PetscCitationsRegister(FECitation,&FEcite);CHKERRQ(ierr); 387 *sp = NULL; 388 ierr = PetscFEInitializePackage();CHKERRQ(ierr); 389 390 ierr = PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView);CHKERRQ(ierr); 391 392 s->order = 0; 393 s->Nc = 1; 394 s->k = 0; 395 s->setupcalled = PETSC_FALSE; 396 397 *sp = s; 398 PetscFunctionReturn(0); 399 } 400 401 /*@ 402 PetscDualSpaceDuplicate - Creates a duplicate PetscDualSpace object, however it is not setup. 403 404 Collective on sp 405 406 Input Parameter: 407 . sp - The original PetscDualSpace 408 409 Output Parameter: 410 . spNew - The duplicate PetscDualSpace 411 412 Level: beginner 413 414 .seealso: PetscDualSpaceCreate(), PetscDualSpaceSetType() 415 @*/ 416 PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew) 417 { 418 PetscErrorCode ierr; 419 420 PetscFunctionBegin; 421 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 422 PetscValidPointer(spNew, 2); 423 ierr = (*sp->ops->duplicate)(sp, spNew);CHKERRQ(ierr); 424 PetscFunctionReturn(0); 425 } 426 427 /*@ 428 PetscDualSpaceGetDM - Get the DM representing the reference cell 429 430 Not collective 431 432 Input Parameter: 433 . sp - The PetscDualSpace 434 435 Output Parameter: 436 . dm - The reference cell 437 438 Level: intermediate 439 440 .seealso: PetscDualSpaceSetDM(), PetscDualSpaceCreate() 441 @*/ 442 PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm) 443 { 444 PetscFunctionBegin; 445 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 446 PetscValidPointer(dm, 2); 447 *dm = sp->dm; 448 PetscFunctionReturn(0); 449 } 450 451 /*@ 452 PetscDualSpaceSetDM - Get the DM representing the reference cell 453 454 Not collective 455 456 Input Parameters: 457 + sp - The PetscDualSpace 458 - dm - The reference cell 459 460 Level: intermediate 461 462 .seealso: PetscDualSpaceGetDM(), PetscDualSpaceCreate() 463 @*/ 464 PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm) 465 { 466 PetscErrorCode ierr; 467 468 PetscFunctionBegin; 469 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 470 PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 471 ierr = DMDestroy(&sp->dm);CHKERRQ(ierr); 472 ierr = PetscObjectReference((PetscObject) dm);CHKERRQ(ierr); 473 sp->dm = dm; 474 PetscFunctionReturn(0); 475 } 476 477 /*@ 478 PetscDualSpaceGetOrder - Get the order of the dual space 479 480 Not collective 481 482 Input Parameter: 483 . sp - The PetscDualSpace 484 485 Output Parameter: 486 . order - The order 487 488 Level: intermediate 489 490 .seealso: PetscDualSpaceSetOrder(), PetscDualSpaceCreate() 491 @*/ 492 PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order) 493 { 494 PetscFunctionBegin; 495 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 496 PetscValidPointer(order, 2); 497 *order = sp->order; 498 PetscFunctionReturn(0); 499 } 500 501 /*@ 502 PetscDualSpaceSetOrder - Set the order of the dual space 503 504 Not collective 505 506 Input Parameters: 507 + sp - The PetscDualSpace 508 - order - The order 509 510 Level: intermediate 511 512 .seealso: PetscDualSpaceGetOrder(), PetscDualSpaceCreate() 513 @*/ 514 PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order) 515 { 516 PetscFunctionBegin; 517 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 518 sp->order = order; 519 PetscFunctionReturn(0); 520 } 521 522 /*@ 523 PetscDualSpaceGetNumComponents - Return the number of components for this space 524 525 Input Parameter: 526 . sp - The PetscDualSpace 527 528 Output Parameter: 529 . Nc - The number of components 530 531 Note: A vector space, for example, will have d components, where d is the spatial dimension 532 533 Level: intermediate 534 535 .seealso: PetscDualSpaceSetNumComponents(), PetscDualSpaceGetDimension(), PetscDualSpaceCreate(), PetscDualSpace 536 @*/ 537 PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc) 538 { 539 PetscFunctionBegin; 540 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 541 PetscValidPointer(Nc, 2); 542 *Nc = sp->Nc; 543 PetscFunctionReturn(0); 544 } 545 546 /*@ 547 PetscDualSpaceSetNumComponents - Set the number of components for this space 548 549 Input Parameters: 550 + sp - The PetscDualSpace 551 - order - The number of components 552 553 Level: intermediate 554 555 .seealso: PetscDualSpaceGetNumComponents(), PetscDualSpaceCreate(), PetscDualSpace 556 @*/ 557 PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc) 558 { 559 PetscFunctionBegin; 560 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 561 sp->Nc = Nc; 562 PetscFunctionReturn(0); 563 } 564 565 /*@ 566 PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space 567 568 Not collective 569 570 Input Parameters: 571 + sp - The PetscDualSpace 572 - i - The basis number 573 574 Output Parameter: 575 . functional - The basis functional 576 577 Level: intermediate 578 579 .seealso: PetscDualSpaceGetDimension(), PetscDualSpaceCreate() 580 @*/ 581 PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional) 582 { 583 PetscInt dim; 584 PetscErrorCode ierr; 585 586 PetscFunctionBegin; 587 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 588 PetscValidPointer(functional, 3); 589 ierr = PetscDualSpaceGetDimension(sp, &dim);CHKERRQ(ierr); 590 if ((i < 0) || (i >= dim)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %d must be in [0, %d)", i, dim); 591 *functional = sp->functional[i]; 592 PetscFunctionReturn(0); 593 } 594 595 /*@ 596 PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals 597 598 Not collective 599 600 Input Parameter: 601 . sp - The PetscDualSpace 602 603 Output Parameter: 604 . dim - The dimension 605 606 Level: intermediate 607 608 .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate() 609 @*/ 610 PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim) 611 { 612 PetscErrorCode ierr; 613 614 PetscFunctionBegin; 615 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 616 PetscValidPointer(dim, 2); 617 *dim = 0; 618 if (sp->ops->getdimension) {ierr = (*sp->ops->getdimension)(sp, dim);CHKERRQ(ierr);} 619 PetscFunctionReturn(0); 620 } 621 622 /*@C 623 PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension 624 625 Not collective 626 627 Input Parameter: 628 . sp - The PetscDualSpace 629 630 Output Parameter: 631 . numDof - An array of length dim+1 which holds the number of dofs for each dimension 632 633 Level: intermediate 634 635 .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate() 636 @*/ 637 PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt **numDof) 638 { 639 PetscErrorCode ierr; 640 641 PetscFunctionBegin; 642 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 643 PetscValidPointer(numDof, 2); 644 ierr = (*sp->ops->getnumdof)(sp, numDof);CHKERRQ(ierr); 645 if (!*numDof) SETERRQ(PetscObjectComm((PetscObject) sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation"); 646 PetscFunctionReturn(0); 647 } 648 649 PetscErrorCode PetscDualSpaceCreateSection(PetscDualSpace sp, PetscSection *section) 650 { 651 DM dm; 652 PetscInt pStart, pEnd, depth, h, offset; 653 const PetscInt *numDof; 654 PetscErrorCode ierr; 655 656 PetscFunctionBegin; 657 ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr); 658 ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 659 ierr = PetscSectionCreate(PetscObjectComm((PetscObject)sp),section);CHKERRQ(ierr); 660 ierr = PetscSectionSetChart(*section,pStart,pEnd);CHKERRQ(ierr); 661 ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 662 ierr = PetscDualSpaceGetNumDof(sp,&numDof);CHKERRQ(ierr); 663 for (h = 0; h <= depth; h++) { 664 PetscInt hStart, hEnd, p, dof; 665 666 ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr); 667 dof = numDof[depth - h]; 668 for (p = hStart; p < hEnd; p++) { 669 ierr = PetscSectionSetDof(*section,p,dof);CHKERRQ(ierr); 670 } 671 } 672 ierr = PetscSectionSetUp(*section);CHKERRQ(ierr); 673 for (h = 0, offset = 0; h <= depth; h++) { 674 PetscInt hStart, hEnd, p, dof; 675 676 ierr = DMPlexGetHeightStratum(dm,h,&hStart,&hEnd);CHKERRQ(ierr); 677 dof = numDof[depth - h]; 678 for (p = hStart; p < hEnd; p++) { 679 ierr = PetscSectionGetDof(*section,p,&dof);CHKERRQ(ierr); 680 ierr = PetscSectionSetOffset(*section,p,offset);CHKERRQ(ierr); 681 offset += dof; 682 } 683 } 684 PetscFunctionReturn(0); 685 } 686 687 /*@ 688 PetscDualSpaceCreateReferenceCell - Create a DMPLEX with the appropriate FEM reference cell 689 690 Collective on sp 691 692 Input Parameters: 693 + sp - The PetscDualSpace 694 . dim - The spatial dimension 695 - simplex - Flag for simplex, otherwise use a tensor-product cell 696 697 Output Parameter: 698 . refdm - The reference cell 699 700 Level: advanced 701 702 .seealso: PetscDualSpaceCreate(), DMPLEX 703 @*/ 704 PetscErrorCode PetscDualSpaceCreateReferenceCell(PetscDualSpace sp, PetscInt dim, PetscBool simplex, DM *refdm) 705 { 706 PetscErrorCode ierr; 707 708 PetscFunctionBeginUser; 709 ierr = DMPlexCreateReferenceCell(PetscObjectComm((PetscObject) sp), dim, simplex, refdm);CHKERRQ(ierr); 710 PetscFunctionReturn(0); 711 } 712 713 /*@C 714 PetscDualSpaceApply - Apply a functional from the dual space basis to an input function 715 716 Input Parameters: 717 + sp - The PetscDualSpace object 718 . f - The basis functional index 719 . time - The time 720 . 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) 721 . numComp - The number of components for the function 722 . func - The input function 723 - ctx - A context for the function 724 725 Output Parameter: 726 . value - numComp output values 727 728 Note: The calling sequence for the callback func is given by: 729 730 $ func(PetscInt dim, PetscReal time, const PetscReal x[], 731 $ PetscInt numComponents, PetscScalar values[], void *ctx) 732 733 Level: developer 734 735 .seealso: PetscDualSpaceCreate() 736 @*/ 737 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) 738 { 739 PetscErrorCode ierr; 740 741 PetscFunctionBegin; 742 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 743 PetscValidPointer(cgeom, 4); 744 PetscValidPointer(value, 8); 745 ierr = (*sp->ops->apply)(sp, f, time, cgeom, numComp, func, ctx, value);CHKERRQ(ierr); 746 PetscFunctionReturn(0); 747 } 748 749 /*@C 750 PetscDualSpaceApplyAll - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints() 751 752 Input Parameters: 753 + sp - The PetscDualSpace object 754 - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints() 755 756 Output Parameter: 757 . spValue - The values of all dual space functionals 758 759 Level: developer 760 761 .seealso: PetscDualSpaceCreate() 762 @*/ 763 PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue) 764 { 765 PetscErrorCode ierr; 766 767 PetscFunctionBegin; 768 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 769 ierr = (*sp->ops->applyall)(sp, pointEval, spValue);CHKERRQ(ierr); 770 PetscFunctionReturn(0); 771 } 772 773 /*@C 774 PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional. 775 776 Input Parameters: 777 + sp - The PetscDualSpace object 778 . f - The basis functional index 779 . time - The time 780 . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian) 781 . Nc - The number of components for the function 782 . func - The input function 783 - ctx - A context for the function 784 785 Output Parameter: 786 . value - The output value 787 788 Note: The calling sequence for the callback func is given by: 789 790 $ func(PetscInt dim, PetscReal time, const PetscReal x[], 791 $ PetscInt numComponents, PetscScalar values[], void *ctx) 792 793 and the idea is to evaluate the functional as an integral 794 795 $ n(f) = int dx n(x) . f(x) 796 797 where both n and f have Nc components. 798 799 Level: developer 800 801 .seealso: PetscDualSpaceCreate() 802 @*/ 803 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) 804 { 805 DM dm; 806 PetscQuadrature n; 807 const PetscReal *points, *weights; 808 PetscReal x[3]; 809 PetscScalar *val; 810 PetscInt dim, dE, qNc, c, Nq, q; 811 PetscBool isAffine; 812 PetscErrorCode ierr; 813 814 PetscFunctionBegin; 815 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 816 PetscValidPointer(value, 5); 817 ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 818 ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr); 819 ierr = PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights);CHKERRQ(ierr); 820 if (dim != cgeom->dim) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature spatial dimension %D != cell geometry dimension %D", dim, cgeom->dim); 821 if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc); 822 ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr); 823 *value = 0.0; 824 isAffine = cgeom->isAffine; 825 dE = cgeom->dimEmbed; 826 for (q = 0; q < Nq; ++q) { 827 if (isAffine) { 828 CoordinatesRefToReal(dE, cgeom->dim, cgeom->xi, cgeom->v, cgeom->J, &points[q*dim], x); 829 ierr = (*func)(dE, time, x, Nc, val, ctx);CHKERRQ(ierr); 830 } else { 831 ierr = (*func)(dE, time, &cgeom->v[dE*q], Nc, val, ctx);CHKERRQ(ierr); 832 } 833 for (c = 0; c < Nc; ++c) { 834 *value += val[c]*weights[q*Nc+c]; 835 } 836 } 837 ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr); 838 PetscFunctionReturn(0); 839 } 840 841 /*@C 842 PetscDualSpaceApplyAllDefault - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllPoints() 843 844 Input Parameters: 845 + sp - The PetscDualSpace object 846 - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllPoints() 847 848 Output Parameter: 849 . spValue - The values of all dual space functionals 850 851 Level: developer 852 853 .seealso: PetscDualSpaceCreate() 854 @*/ 855 PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue) 856 { 857 PetscQuadrature n; 858 const PetscReal *points, *weights; 859 PetscInt qNc, c, Nq, q, f, spdim, Nc; 860 PetscInt offset; 861 PetscErrorCode ierr; 862 863 PetscFunctionBegin; 864 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 865 PetscValidScalarPointer(pointEval, 2); 866 PetscValidScalarPointer(spValue, 5); 867 ierr = PetscDualSpaceGetDimension(sp, &spdim);CHKERRQ(ierr); 868 ierr = PetscDualSpaceGetNumComponents(sp, &Nc);CHKERRQ(ierr); 869 for (f = 0, offset = 0; f < spdim; f++) { 870 ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr); 871 ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr); 872 if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc); 873 spValue[f] = 0.0; 874 for (q = 0; q < Nq; ++q) { 875 for (c = 0; c < Nc; ++c) { 876 spValue[f] += pointEval[offset++]*weights[q*Nc+c]; 877 } 878 } 879 } 880 PetscFunctionReturn(0); 881 } 882 883 PetscErrorCode PetscDualSpaceGetAllPoints(PetscDualSpace sp, PetscQuadrature *allPoints) 884 { 885 PetscErrorCode ierr; 886 887 PetscFunctionBegin; 888 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 889 PetscValidPointer(allPoints,2); 890 if (!sp->allPoints && sp->ops->createallpoints) { 891 ierr = (*sp->ops->createallpoints)(sp,&sp->allPoints);CHKERRQ(ierr); 892 } 893 *allPoints = sp->allPoints; 894 PetscFunctionReturn(0); 895 } 896 897 PetscErrorCode PetscDualSpaceCreateAllPointsDefault(PetscDualSpace sp, PetscQuadrature *allPoints) 898 { 899 PetscInt spdim; 900 PetscInt numPoints, offset; 901 PetscReal *points; 902 PetscInt f, dim; 903 PetscQuadrature q; 904 PetscErrorCode ierr; 905 906 PetscFunctionBegin; 907 ierr = PetscDualSpaceGetDimension(sp,&spdim);CHKERRQ(ierr); 908 if (!spdim) { 909 ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr); 910 ierr = PetscQuadratureSetData(*allPoints,0,0,0,NULL,NULL);CHKERRQ(ierr); 911 } 912 ierr = PetscDualSpaceGetFunctional(sp,0,&q);CHKERRQ(ierr); 913 ierr = PetscQuadratureGetData(q,&dim,NULL,&numPoints,NULL,NULL);CHKERRQ(ierr); 914 for (f = 1; f < spdim; f++) { 915 PetscInt Np; 916 917 ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr); 918 ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL);CHKERRQ(ierr); 919 numPoints += Np; 920 } 921 ierr = PetscMalloc1(dim*numPoints,&points);CHKERRQ(ierr); 922 for (f = 0, offset = 0; f < spdim; f++) { 923 const PetscReal *p; 924 PetscInt Np, i; 925 926 ierr = PetscDualSpaceGetFunctional(sp,f,&q);CHKERRQ(ierr); 927 ierr = PetscQuadratureGetData(q,NULL,NULL,&Np,&p,NULL);CHKERRQ(ierr); 928 for (i = 0; i < Np * dim; i++) { 929 points[offset + i] = p[i]; 930 } 931 offset += Np * dim; 932 } 933 ierr = PetscQuadratureCreate(PETSC_COMM_SELF,allPoints);CHKERRQ(ierr); 934 ierr = PetscQuadratureSetData(*allPoints,dim,0,numPoints,points,NULL);CHKERRQ(ierr); 935 PetscFunctionReturn(0); 936 } 937 938 /*@C 939 PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid. 940 941 Input Parameters: 942 + sp - The PetscDualSpace object 943 . f - The basis functional index 944 . time - The time 945 . cgeom - A context with geometric information for this cell, we currently just use the centroid 946 . Nc - The number of components for the function 947 . func - The input function 948 - ctx - A context for the function 949 950 Output Parameter: 951 . value - The output value (scalar) 952 953 Note: The calling sequence for the callback func is given by: 954 955 $ func(PetscInt dim, PetscReal time, const PetscReal x[], 956 $ PetscInt numComponents, PetscScalar values[], void *ctx) 957 958 and the idea is to evaluate the functional as an integral 959 960 $ n(f) = int dx n(x) . f(x) 961 962 where both n and f have Nc components. 963 964 Level: developer 965 966 .seealso: PetscDualSpaceCreate() 967 @*/ 968 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) 969 { 970 DM dm; 971 PetscQuadrature n; 972 const PetscReal *points, *weights; 973 PetscScalar *val; 974 PetscInt dimEmbed, qNc, c, Nq, q; 975 PetscErrorCode ierr; 976 977 PetscFunctionBegin; 978 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 979 PetscValidPointer(value, 5); 980 ierr = PetscDualSpaceGetDM(sp, &dm);CHKERRQ(ierr); 981 ierr = DMGetCoordinateDim(dm, &dimEmbed);CHKERRQ(ierr); 982 ierr = PetscDualSpaceGetFunctional(sp, f, &n);CHKERRQ(ierr); 983 ierr = PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights);CHKERRQ(ierr); 984 if (qNc != Nc) SETERRQ2(PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc); 985 ierr = DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr); 986 *value = 0.; 987 for (q = 0; q < Nq; ++q) { 988 ierr = (*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx);CHKERRQ(ierr); 989 for (c = 0; c < Nc; ++c) { 990 *value += val[c]*weights[q*Nc+c]; 991 } 992 } 993 ierr = DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val);CHKERRQ(ierr); 994 PetscFunctionReturn(0); 995 } 996 997 /*@ 998 PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a 999 given height. This assumes that the reference cell is symmetric over points of this height. 1000 1001 If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and 1002 pointwise values are not defined on the element boundaries), or if the implementation of PetscDualSpace does not 1003 support extracting subspaces, then NULL is returned. 1004 1005 This does not increment the reference count on the returned dual space, and the user should not destroy it. 1006 1007 Not collective 1008 1009 Input Parameters: 1010 + sp - the PetscDualSpace object 1011 - height - the height of the mesh point for which the subspace is desired 1012 1013 Output Parameter: 1014 . subsp - the subspace. Note that the functionals in the subspace are with respect to the intrinsic geometry of the 1015 point, which will be of lesser dimension if height > 0. 1016 1017 Level: advanced 1018 1019 .seealso: PetscSpaceGetHeightSubspace(), PetscDualSpace 1020 @*/ 1021 PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp) 1022 { 1023 PetscErrorCode ierr; 1024 1025 PetscFunctionBegin; 1026 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1027 PetscValidPointer(subsp, 3); 1028 *subsp = NULL; 1029 if (sp->ops->getheightsubspace) {ierr = (*sp->ops->getheightsubspace)(sp, height, subsp);CHKERRQ(ierr);} 1030 PetscFunctionReturn(0); 1031 } 1032 1033 /*@ 1034 PetscDualSpaceGetPointSubspace - Get the subset of the dual space basis that is supported on a particular mesh point. 1035 1036 If the dual space is not defined on the mesh point (e.g. if the space is discontinuous and pointwise values are not 1037 defined on the element boundaries), or if the implementation of PetscDualSpace does not support extracting 1038 subspaces, then NULL is returned. 1039 1040 This does not increment the reference count on the returned dual space, and the user should not destroy it. 1041 1042 Not collective 1043 1044 Input Parameters: 1045 + sp - the PetscDualSpace object 1046 - point - the point (in the dual space's DM) for which the subspace is desired 1047 1048 Output Parameters: 1049 bdsp - the subspace. Note that the functionals in the subspace are with respect to the intrinsic geometry of the 1050 point, which will be of lesser dimension if height > 0. 1051 1052 Level: advanced 1053 1054 .seealso: PetscDualSpace 1055 @*/ 1056 PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace sp, PetscInt point, PetscDualSpace *bdsp) 1057 { 1058 PetscErrorCode ierr; 1059 1060 PetscFunctionBegin; 1061 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1062 PetscValidPointer(bdsp,2); 1063 *bdsp = NULL; 1064 if (sp->ops->getpointsubspace) { 1065 ierr = (*sp->ops->getpointsubspace)(sp,point,bdsp);CHKERRQ(ierr); 1066 } else if (sp->ops->getheightsubspace) { 1067 DM dm; 1068 DMLabel label; 1069 PetscInt dim, depth, height; 1070 1071 ierr = PetscDualSpaceGetDM(sp,&dm);CHKERRQ(ierr); 1072 ierr = DMPlexGetDepth(dm,&dim);CHKERRQ(ierr); 1073 ierr = DMPlexGetDepthLabel(dm,&label);CHKERRQ(ierr); 1074 ierr = DMLabelGetValue(label,point,&depth);CHKERRQ(ierr); 1075 height = dim - depth; 1076 ierr = (*sp->ops->getheightsubspace)(sp,height,bdsp);CHKERRQ(ierr); 1077 } 1078 PetscFunctionReturn(0); 1079 } 1080 1081 /*@C 1082 PetscDualSpaceGetSymmetries - Returns a description of the symmetries of this basis 1083 1084 Not collective 1085 1086 Input Parameter: 1087 . sp - the PetscDualSpace object 1088 1089 Output Parameters: 1090 + perms - Permutations of the local degrees of freedom, parameterized by the point orientation 1091 - flips - Sign reversal of the local degrees of freedom, parameterized by the point orientation 1092 1093 Note: The permutation and flip arrays are organized in the following way 1094 $ perms[p][ornt][dof # on point] = new local dof # 1095 $ flips[p][ornt][dof # on point] = reversal or not 1096 1097 Level: developer 1098 1099 .seealso: PetscDualSpaceSetSymmetries() 1100 @*/ 1101 PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips) 1102 { 1103 PetscErrorCode ierr; 1104 1105 PetscFunctionBegin; 1106 PetscValidHeaderSpecific(sp,PETSCDUALSPACE_CLASSID,1); 1107 if (perms) {PetscValidPointer(perms,2); *perms = NULL;} 1108 if (flips) {PetscValidPointer(flips,3); *flips = NULL;} 1109 if (sp->ops->getsymmetries) {ierr = (sp->ops->getsymmetries)(sp,perms,flips);CHKERRQ(ierr);} 1110 PetscFunctionReturn(0); 1111 } 1112 1113 /*@ 1114 PetscDualSpaceGetDeRahm - Get the k-simplex associated with the functionals in this dual space 1115 1116 Input Parameter: 1117 . dsp - The PetscDualSpace 1118 1119 Output Parameter: 1120 . k - The simplex dimension 1121 1122 Level: advanced 1123 1124 Note: Currently supported values are 1125 $ 0: These are H_1 methods that only transform coordinates 1126 $ 1: These are Hcurl methods that transform functions using the covariant Piola transform (COVARIANT_PIOLA_TRANSFORM) 1127 $ 2: These are the same as 1 1128 $ 3: These are Hdiv methods that transform functions using the contravariant Piola transform (CONTRAVARIANT_PIOLA_TRANSFORM) 1129 1130 .seealso: PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceTransformType 1131 @*/ 1132 PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace dsp, PetscInt *k) 1133 { 1134 PetscFunctionBeginHot; 1135 PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 1136 PetscValidPointer(k, 2); 1137 *k = dsp->k; 1138 PetscFunctionReturn(0); 1139 } 1140 1141 /*@C 1142 PetscDualSpaceTransform - Transform the function values 1143 1144 Input Parameters: 1145 + dsp - The PetscDualSpace 1146 . trans - The type of transform 1147 . isInverse - Flag to invert the transform 1148 . fegeom - The cell geometry 1149 . Nv - The number of function samples 1150 . Nc - The number of function components 1151 - vals - The function values 1152 1153 Output Parameter: 1154 . vals - The transformed function values 1155 1156 Level: developer 1157 1158 .seealso: PetscDualSpaceTransformGradient(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransformType 1159 @*/ 1160 PetscErrorCode PetscDualSpaceTransform(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[]) 1161 { 1162 PetscInt dim, v, c; 1163 1164 PetscFunctionBeginHot; 1165 PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 1166 PetscValidPointer(fegeom, 4); 1167 PetscValidPointer(vals, 7); 1168 dim = dsp->dm->dim; 1169 /* Assume its a vector, otherwise assume its a bunch of scalars */ 1170 if (Nc == 1 || Nc != dim) PetscFunctionReturn(0); 1171 switch (trans) { 1172 case IDENTITY_TRANSFORM: break; 1173 case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */ 1174 if (isInverse) { 1175 for (v = 0; v < Nv; ++v) { 1176 switch (dim) 1177 { 1178 case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break; 1179 case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break; 1180 default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 1181 } 1182 } 1183 } else { 1184 for (v = 0; v < Nv; ++v) { 1185 switch (dim) 1186 { 1187 case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break; 1188 case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break; 1189 default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 1190 } 1191 } 1192 } 1193 break; 1194 case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */ 1195 if (isInverse) { 1196 for (v = 0; v < Nv; ++v) { 1197 switch (dim) 1198 { 1199 case 2: DMPlex_Mult2DReal_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break; 1200 case 3: DMPlex_Mult3DReal_Internal(fegeom->invJ, 1, &vals[v*Nc], &vals[v*Nc]);break; 1201 default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 1202 } 1203 for (c = 0; c < Nc; ++c) vals[v*Nc+c] *= fegeom->detJ[0]; 1204 } 1205 } else { 1206 for (v = 0; v < Nv; ++v) { 1207 switch (dim) 1208 { 1209 case 2: DMPlex_Mult2DReal_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break; 1210 case 3: DMPlex_Mult3DReal_Internal(fegeom->J, 1, &vals[v*Nc], &vals[v*Nc]);break; 1211 default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 1212 } 1213 for (c = 0; c < Nc; ++c) vals[v*Nc+c] /= fegeom->detJ[0]; 1214 } 1215 } 1216 break; 1217 } 1218 PetscFunctionReturn(0); 1219 } 1220 /*@C 1221 PetscDualSpaceTransformGradient - Transform the function gradient values 1222 1223 Input Parameters: 1224 + dsp - The PetscDualSpace 1225 . trans - The type of transform 1226 . isInverse - Flag to invert the transform 1227 . fegeom - The cell geometry 1228 . Nv - The number of function gradient samples 1229 . Nc - The number of function components 1230 - vals - The function gradient values 1231 1232 Output Parameter: 1233 . vals - The transformed function values 1234 1235 Level: developer 1236 1237 .seealso: PetscDualSpaceTransform(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransformType 1238 @*/ 1239 PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[]) 1240 { 1241 PetscInt dim, v, c, d; 1242 1243 PetscFunctionBeginHot; 1244 PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 1245 PetscValidPointer(fegeom, 4); 1246 PetscValidPointer(vals, 7); 1247 dim = dsp->dm->dim; 1248 /* Transform gradient */ 1249 for (v = 0; v < Nv; ++v) { 1250 for (c = 0; c < Nc; ++c) { 1251 switch (dim) 1252 { 1253 case 1: vals[(v*Nc+c)*dim] *= fegeom->invJ[0]; 1254 case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break; 1255 case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break; 1256 default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 1257 } 1258 } 1259 } 1260 /* Assume its a vector, otherwise assume its a bunch of scalars */ 1261 if (Nc == 1 || Nc != dim) PetscFunctionReturn(0); 1262 switch (trans) { 1263 case IDENTITY_TRANSFORM: break; 1264 case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */ 1265 if (isInverse) { 1266 for (v = 0; v < Nv; ++v) { 1267 for (d = 0; d < dim; ++d) { 1268 switch (dim) 1269 { 1270 case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 1271 case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 1272 default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 1273 } 1274 } 1275 } 1276 } else { 1277 for (v = 0; v < Nv; ++v) { 1278 for (d = 0; d < dim; ++d) { 1279 switch (dim) 1280 { 1281 case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 1282 case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 1283 default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 1284 } 1285 } 1286 } 1287 } 1288 break; 1289 case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */ 1290 if (isInverse) { 1291 for (v = 0; v < Nv; ++v) { 1292 for (d = 0; d < dim; ++d) { 1293 switch (dim) 1294 { 1295 case 2: DMPlex_Mult2DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 1296 case 3: DMPlex_Mult3DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 1297 default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 1298 } 1299 for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] *= fegeom->detJ[0]; 1300 } 1301 } 1302 } else { 1303 for (v = 0; v < Nv; ++v) { 1304 for (d = 0; d < dim; ++d) { 1305 switch (dim) 1306 { 1307 case 2: DMPlex_Mult2DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 1308 case 3: DMPlex_Mult3DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 1309 default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 1310 } 1311 for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] /= fegeom->detJ[0]; 1312 } 1313 } 1314 } 1315 break; 1316 } 1317 PetscFunctionReturn(0); 1318 } 1319 1320 /*@C 1321 PetscDualSpacePullback - Transform the given functional so that it operates on real space, rather than the reference element. Operationally, this means that we map the function evaluations depending on continuity requirements of our finite element method. 1322 1323 Input Parameters: 1324 + dsp - The PetscDualSpace 1325 . fegeom - The geometry for this cell 1326 . Nq - The number of function samples 1327 . Nc - The number of function components 1328 - pointEval - The function values 1329 1330 Output Parameter: 1331 . pointEval - The transformed function values 1332 1333 Level: advanced 1334 1335 Note: Functions transform in a complementary way (pushforward) to functionals, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex. 1336 1337 .seealso: PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm() 1338 @*/ 1339 PetscErrorCode PetscDualSpacePullback(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[]) 1340 { 1341 PetscDualSpaceTransformType trans; 1342 PetscErrorCode ierr; 1343 1344 PetscFunctionBeginHot; 1345 PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 1346 PetscValidPointer(fegeom, 2); 1347 PetscValidPointer(pointEval, 5); 1348 /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k. 1349 This determines their transformation properties. */ 1350 switch (dsp->k) 1351 { 1352 case 0: /* H^1 point evaluations */ 1353 trans = IDENTITY_TRANSFORM;break; 1354 case 1: /* Hcurl preserves tangential edge traces */ 1355 case 2: 1356 trans = COVARIANT_PIOLA_TRANSFORM;break; 1357 case 3: /* Hdiv preserve normal traces */ 1358 trans = CONTRAVARIANT_PIOLA_TRANSFORM;break; 1359 default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", dsp->k); 1360 } 1361 ierr = PetscDualSpaceTransform(dsp, trans, PETSC_TRUE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr); 1362 PetscFunctionReturn(0); 1363 } 1364 1365 /*@C 1366 PetscDualSpacePushforward - Transform the given function so that it operates on real space, rather than the reference element. Operationally, this means that we map the function evaluations depending on continuity requirements of our finite element method. 1367 1368 Input Parameters: 1369 + dsp - The PetscDualSpace 1370 . fegeom - The geometry for this cell 1371 . Nq - The number of function samples 1372 . Nc - The number of function components 1373 - pointEval - The function values 1374 1375 Output Parameter: 1376 . pointEval - The transformed function values 1377 1378 Level: advanced 1379 1380 Note: Functionals transform in a complementary way (pullback) to functions, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex. 1381 1382 .seealso: PetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm() 1383 @*/ 1384 PetscErrorCode PetscDualSpacePushforward(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[]) 1385 { 1386 PetscDualSpaceTransformType trans; 1387 PetscErrorCode ierr; 1388 1389 PetscFunctionBeginHot; 1390 PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 1391 PetscValidPointer(fegeom, 2); 1392 PetscValidPointer(pointEval, 5); 1393 /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k. 1394 This determines their transformation properties. */ 1395 switch (dsp->k) 1396 { 1397 case 0: /* H^1 point evaluations */ 1398 trans = IDENTITY_TRANSFORM;break; 1399 case 1: /* Hcurl preserves tangential edge traces */ 1400 case 2: 1401 trans = COVARIANT_PIOLA_TRANSFORM;break; 1402 case 3: /* Hdiv preserve normal traces */ 1403 trans = CONTRAVARIANT_PIOLA_TRANSFORM;break; 1404 default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", dsp->k); 1405 } 1406 ierr = PetscDualSpaceTransform(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr); 1407 PetscFunctionReturn(0); 1408 } 1409 1410 /*@C 1411 PetscDualSpacePushforwardGradient - Transform the given function gradient so that it operates on real space, rather than the reference element. Operationally, this means that we map the function evaluations depending on continuity requirements of our finite element method. 1412 1413 Input Parameters: 1414 + dsp - The PetscDualSpace 1415 . fegeom - The geometry for this cell 1416 . Nq - The number of function gradient samples 1417 . Nc - The number of function components 1418 - pointEval - The function gradient values 1419 1420 Output Parameter: 1421 . pointEval - The transformed function gradient values 1422 1423 Level: advanced 1424 1425 Note: Functionals transform in a complementary way (pullback) to functions, so that the scalar product is invariant. The type of transform is dependent on the associated k-simplex from the DeRahm complex. 1426 1427 .seealso: PetscDualSpacePushforward(), PPetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm() 1428 @*/ 1429 PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[]) 1430 { 1431 PetscDualSpaceTransformType trans; 1432 PetscErrorCode ierr; 1433 1434 PetscFunctionBeginHot; 1435 PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 1436 PetscValidPointer(fegeom, 2); 1437 PetscValidPointer(pointEval, 5); 1438 /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k. 1439 This determines their transformation properties. */ 1440 switch (dsp->k) 1441 { 1442 case 0: /* H^1 point evaluations */ 1443 trans = IDENTITY_TRANSFORM;break; 1444 case 1: /* Hcurl preserves tangential edge traces */ 1445 case 2: 1446 trans = COVARIANT_PIOLA_TRANSFORM;break; 1447 case 3: /* Hdiv preserve normal traces */ 1448 trans = CONTRAVARIANT_PIOLA_TRANSFORM;break; 1449 default: SETERRQ1(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", dsp->k); 1450 } 1451 ierr = PetscDualSpaceTransformGradient(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval);CHKERRQ(ierr); 1452 PetscFunctionReturn(0); 1453 } 1454