1 #include <petsc/private/petscfeimpl.h> /*I "petscfe.h" I*/ 2 #include <petscdmplex.h> 3 4 PetscClassId PETSCDUALSPACE_CLASSID = 0; 5 6 PetscLogEvent PETSCDUALSPACE_SetUp; 7 8 PetscFunctionList PetscDualSpaceList = NULL; 9 PetscBool PetscDualSpaceRegisterAllCalled = PETSC_FALSE; 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}, {0,2}. 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 PetscFunctionBegin; 110 PetscCall(PetscFunctionListAdd(&PetscDualSpaceList, sname, function)); 111 PetscFunctionReturn(0); 112 } 113 114 /*@C 115 PetscDualSpaceSetType - Builds a particular PetscDualSpace 116 117 Collective on sp 118 119 Input Parameters: 120 + sp - The PetscDualSpace object 121 - name - The kind of space 122 123 Options Database Key: 124 . -petscdualspace_type <type> - Sets the PetscDualSpace type; use -help for a list of available types 125 126 Level: intermediate 127 128 .seealso: PetscDualSpaceGetType(), PetscDualSpaceCreate() 129 @*/ 130 PetscErrorCode PetscDualSpaceSetType(PetscDualSpace sp, PetscDualSpaceType name) 131 { 132 PetscErrorCode (*r)(PetscDualSpace); 133 PetscBool match; 134 135 PetscFunctionBegin; 136 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 137 PetscCall(PetscObjectTypeCompare((PetscObject) sp, name, &match)); 138 if (match) PetscFunctionReturn(0); 139 140 if (!PetscDualSpaceRegisterAllCalled) PetscCall(PetscDualSpaceRegisterAll()); 141 PetscCall(PetscFunctionListFind(PetscDualSpaceList, name, &r)); 142 PetscCheck(r,PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscDualSpace type: %s", name); 143 144 if (sp->ops->destroy) { 145 PetscCall((*sp->ops->destroy)(sp)); 146 sp->ops->destroy = NULL; 147 } 148 PetscCall((*r)(sp)); 149 PetscCall(PetscObjectChangeTypeName((PetscObject) sp, name)); 150 PetscFunctionReturn(0); 151 } 152 153 /*@C 154 PetscDualSpaceGetType - Gets the PetscDualSpace type name (as a string) from the object. 155 156 Not Collective 157 158 Input Parameter: 159 . sp - The PetscDualSpace 160 161 Output Parameter: 162 . name - The PetscDualSpace type name 163 164 Level: intermediate 165 166 .seealso: PetscDualSpaceSetType(), PetscDualSpaceCreate() 167 @*/ 168 PetscErrorCode PetscDualSpaceGetType(PetscDualSpace sp, PetscDualSpaceType *name) 169 { 170 PetscFunctionBegin; 171 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 172 PetscValidPointer(name, 2); 173 if (!PetscDualSpaceRegisterAllCalled) { 174 PetscCall(PetscDualSpaceRegisterAll()); 175 } 176 *name = ((PetscObject) sp)->type_name; 177 PetscFunctionReturn(0); 178 } 179 180 static PetscErrorCode PetscDualSpaceView_ASCII(PetscDualSpace sp, PetscViewer v) 181 { 182 PetscViewerFormat format; 183 PetscInt pdim, f; 184 185 PetscFunctionBegin; 186 PetscCall(PetscDualSpaceGetDimension(sp, &pdim)); 187 PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject) sp, v)); 188 PetscCall(PetscViewerASCIIPushTab(v)); 189 if (sp->k) { 190 PetscCall(PetscViewerASCIIPrintf(v, "Dual space for %D-forms %swith %D components, size %D\n", PetscAbsInt(sp->k), sp->k < 0 ? "(stored in dual form) ": "", sp->Nc, pdim)); 191 } else { 192 PetscCall(PetscViewerASCIIPrintf(v, "Dual space with %D components, size %D\n", sp->Nc, pdim)); 193 } 194 if (sp->ops->view) PetscCall((*sp->ops->view)(sp, v)); 195 PetscCall(PetscViewerGetFormat(v, &format)); 196 if (format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 197 PetscCall(PetscViewerASCIIPushTab(v)); 198 for (f = 0; f < pdim; ++f) { 199 PetscCall(PetscViewerASCIIPrintf(v, "Dual basis vector %D\n", f)); 200 PetscCall(PetscViewerASCIIPushTab(v)); 201 PetscCall(PetscQuadratureView(sp->functional[f], v)); 202 PetscCall(PetscViewerASCIIPopTab(v)); 203 } 204 PetscCall(PetscViewerASCIIPopTab(v)); 205 } 206 PetscCall(PetscViewerASCIIPopTab(v)); 207 PetscFunctionReturn(0); 208 } 209 210 /*@C 211 PetscDualSpaceViewFromOptions - View from Options 212 213 Collective on PetscDualSpace 214 215 Input Parameters: 216 + A - the PetscDualSpace object 217 . obj - Optional object, proivides prefix 218 - name - command line option 219 220 Level: intermediate 221 .seealso: PetscDualSpace, PetscDualSpaceView(), PetscObjectViewFromOptions(), PetscDualSpaceCreate() 222 @*/ 223 PetscErrorCode PetscDualSpaceViewFromOptions(PetscDualSpace A,PetscObject obj,const char name[]) 224 { 225 PetscFunctionBegin; 226 PetscValidHeaderSpecific(A,PETSCDUALSPACE_CLASSID,1); 227 PetscCall(PetscObjectViewFromOptions((PetscObject)A,obj,name)); 228 PetscFunctionReturn(0); 229 } 230 231 /*@ 232 PetscDualSpaceView - Views a PetscDualSpace 233 234 Collective on sp 235 236 Input Parameters: 237 + sp - the PetscDualSpace object to view 238 - v - the viewer 239 240 Level: beginner 241 242 .seealso PetscDualSpaceDestroy(), PetscDualSpace 243 @*/ 244 PetscErrorCode PetscDualSpaceView(PetscDualSpace sp, PetscViewer v) 245 { 246 PetscBool iascii; 247 248 PetscFunctionBegin; 249 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 250 if (v) PetscValidHeaderSpecific(v, PETSC_VIEWER_CLASSID, 2); 251 if (!v) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject) sp), &v)); 252 PetscCall(PetscObjectTypeCompare((PetscObject) v, PETSCVIEWERASCII, &iascii)); 253 if (iascii) PetscCall(PetscDualSpaceView_ASCII(sp, v)); 254 PetscFunctionReturn(0); 255 } 256 257 /*@ 258 PetscDualSpaceSetFromOptions - sets parameters in a PetscDualSpace from the options database 259 260 Collective on sp 261 262 Input Parameter: 263 . sp - the PetscDualSpace object to set options for 264 265 Options Database: 266 + -petscdualspace_order <order> - the approximation order of the space 267 . -petscdualspace_form_degree <deg> - the form degree, say 0 for point evaluations, or 2 for area integrals 268 . -petscdualspace_components <c> - the number of components, say d for a vector field 269 - -petscdualspace_refcell <celltype> - Reference cell type name 270 271 Level: intermediate 272 273 .seealso PetscDualSpaceView(), PetscDualSpace, PetscObjectSetFromOptions() 274 @*/ 275 PetscErrorCode PetscDualSpaceSetFromOptions(PetscDualSpace sp) 276 { 277 DMPolytopeType refCell = DM_POLYTOPE_TRIANGLE; 278 const char *defaultType; 279 char name[256]; 280 PetscBool flg; 281 PetscErrorCode ierr; 282 283 PetscFunctionBegin; 284 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 285 if (!((PetscObject) sp)->type_name) { 286 defaultType = PETSCDUALSPACELAGRANGE; 287 } else { 288 defaultType = ((PetscObject) sp)->type_name; 289 } 290 if (!PetscSpaceRegisterAllCalled) PetscCall(PetscSpaceRegisterAll()); 291 292 ierr = PetscObjectOptionsBegin((PetscObject) sp);PetscCall(ierr); 293 PetscCall(PetscOptionsFList("-petscdualspace_type", "Dual space", "PetscDualSpaceSetType", PetscDualSpaceList, defaultType, name, 256, &flg)); 294 if (flg) { 295 PetscCall(PetscDualSpaceSetType(sp, name)); 296 } else if (!((PetscObject) sp)->type_name) { 297 PetscCall(PetscDualSpaceSetType(sp, defaultType)); 298 } 299 PetscCall(PetscOptionsBoundedInt("-petscdualspace_order", "The approximation order", "PetscDualSpaceSetOrder", sp->order, &sp->order, NULL,0)); 300 PetscCall(PetscOptionsInt("-petscdualspace_form_degree", "The form degree of the dofs", "PetscDualSpaceSetFormDegree", sp->k, &sp->k, NULL)); 301 PetscCall(PetscOptionsBoundedInt("-petscdualspace_components", "The number of components", "PetscDualSpaceSetNumComponents", sp->Nc, &sp->Nc, NULL,1)); 302 if (sp->ops->setfromoptions) { 303 PetscCall((*sp->ops->setfromoptions)(PetscOptionsObject,sp)); 304 } 305 PetscCall(PetscOptionsEnum("-petscdualspace_refcell", "Reference cell shape", "PetscDualSpaceSetReferenceCell", DMPolytopeTypes, (PetscEnum) refCell, (PetscEnum *) &refCell, &flg)); 306 if (flg) { 307 DM K; 308 309 PetscCall(DMPlexCreateReferenceCell(PETSC_COMM_SELF, refCell, &K)); 310 PetscCall(PetscDualSpaceSetDM(sp, K)); 311 PetscCall(DMDestroy(&K)); 312 } 313 314 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 315 PetscCall(PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) sp)); 316 ierr = PetscOptionsEnd();PetscCall(ierr); 317 sp->setfromoptionscalled = PETSC_TRUE; 318 PetscFunctionReturn(0); 319 } 320 321 /*@ 322 PetscDualSpaceSetUp - Construct a basis for the PetscDualSpace 323 324 Collective on sp 325 326 Input Parameter: 327 . sp - the PetscDualSpace object to setup 328 329 Level: intermediate 330 331 .seealso PetscDualSpaceView(), PetscDualSpaceDestroy(), PetscDualSpace 332 @*/ 333 PetscErrorCode PetscDualSpaceSetUp(PetscDualSpace sp) 334 { 335 PetscFunctionBegin; 336 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 337 if (sp->setupcalled) PetscFunctionReturn(0); 338 PetscCall(PetscLogEventBegin(PETSCDUALSPACE_SetUp, sp, 0, 0, 0)); 339 sp->setupcalled = PETSC_TRUE; 340 if (sp->ops->setup) PetscCall((*sp->ops->setup)(sp)); 341 PetscCall(PetscLogEventEnd(PETSCDUALSPACE_SetUp, sp, 0, 0, 0)); 342 if (sp->setfromoptionscalled) PetscCall(PetscDualSpaceViewFromOptions(sp, NULL, "-petscdualspace_view")); 343 PetscFunctionReturn(0); 344 } 345 346 static PetscErrorCode PetscDualSpaceClearDMData_Internal(PetscDualSpace sp, DM dm) 347 { 348 PetscInt pStart = -1, pEnd = -1, depth = -1; 349 350 PetscFunctionBegin; 351 if (!dm) PetscFunctionReturn(0); 352 PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 353 PetscCall(DMPlexGetDepth(dm, &depth)); 354 355 if (sp->pointSpaces) { 356 PetscInt i; 357 358 for (i = 0; i < pEnd - pStart; i++) { 359 PetscCall(PetscDualSpaceDestroy(&(sp->pointSpaces[i]))); 360 } 361 } 362 PetscCall(PetscFree(sp->pointSpaces)); 363 364 if (sp->heightSpaces) { 365 PetscInt i; 366 367 for (i = 0; i <= depth; i++) { 368 PetscCall(PetscDualSpaceDestroy(&(sp->heightSpaces[i]))); 369 } 370 } 371 PetscCall(PetscFree(sp->heightSpaces)); 372 373 PetscCall(PetscSectionDestroy(&(sp->pointSection))); 374 PetscCall(PetscQuadratureDestroy(&(sp->intNodes))); 375 PetscCall(VecDestroy(&(sp->intDofValues))); 376 PetscCall(VecDestroy(&(sp->intNodeValues))); 377 PetscCall(MatDestroy(&(sp->intMat))); 378 PetscCall(PetscQuadratureDestroy(&(sp->allNodes))); 379 PetscCall(VecDestroy(&(sp->allDofValues))); 380 PetscCall(VecDestroy(&(sp->allNodeValues))); 381 PetscCall(MatDestroy(&(sp->allMat))); 382 PetscCall(PetscFree(sp->numDof)); 383 PetscFunctionReturn(0); 384 } 385 386 /*@ 387 PetscDualSpaceDestroy - Destroys a PetscDualSpace object 388 389 Collective on sp 390 391 Input Parameter: 392 . sp - the PetscDualSpace object to destroy 393 394 Level: beginner 395 396 .seealso PetscDualSpaceView(), PetscDualSpace(), PetscDualSpaceCreate() 397 @*/ 398 PetscErrorCode PetscDualSpaceDestroy(PetscDualSpace *sp) 399 { 400 PetscInt dim, f; 401 DM dm; 402 403 PetscFunctionBegin; 404 if (!*sp) PetscFunctionReturn(0); 405 PetscValidHeaderSpecific((*sp), PETSCDUALSPACE_CLASSID, 1); 406 407 if (--((PetscObject)(*sp))->refct > 0) {*sp = NULL; PetscFunctionReturn(0);} 408 ((PetscObject) (*sp))->refct = 0; 409 410 PetscCall(PetscDualSpaceGetDimension(*sp, &dim)); 411 dm = (*sp)->dm; 412 413 if ((*sp)->ops->destroy) PetscCall((*(*sp)->ops->destroy)(*sp)); 414 PetscCall(PetscDualSpaceClearDMData_Internal(*sp, dm)); 415 416 for (f = 0; f < dim; ++f) { 417 PetscCall(PetscQuadratureDestroy(&(*sp)->functional[f])); 418 } 419 PetscCall(PetscFree((*sp)->functional)); 420 PetscCall(DMDestroy(&(*sp)->dm)); 421 PetscCall(PetscHeaderDestroy(sp)); 422 PetscFunctionReturn(0); 423 } 424 425 /*@ 426 PetscDualSpaceCreate - Creates an empty PetscDualSpace object. The type can then be set with PetscDualSpaceSetType(). 427 428 Collective 429 430 Input Parameter: 431 . comm - The communicator for the PetscDualSpace object 432 433 Output Parameter: 434 . sp - The PetscDualSpace object 435 436 Level: beginner 437 438 .seealso: PetscDualSpaceSetType(), PETSCDUALSPACELAGRANGE 439 @*/ 440 PetscErrorCode PetscDualSpaceCreate(MPI_Comm comm, PetscDualSpace *sp) 441 { 442 PetscDualSpace s; 443 444 PetscFunctionBegin; 445 PetscValidPointer(sp, 2); 446 PetscCall(PetscCitationsRegister(FECitation,&FEcite)); 447 *sp = NULL; 448 PetscCall(PetscFEInitializePackage()); 449 450 PetscCall(PetscHeaderCreate(s, PETSCDUALSPACE_CLASSID, "PetscDualSpace", "Dual Space", "PetscDualSpace", comm, PetscDualSpaceDestroy, PetscDualSpaceView)); 451 452 s->order = 0; 453 s->Nc = 1; 454 s->k = 0; 455 s->spdim = -1; 456 s->spintdim = -1; 457 s->uniform = PETSC_TRUE; 458 s->setupcalled = PETSC_FALSE; 459 460 *sp = s; 461 PetscFunctionReturn(0); 462 } 463 464 /*@ 465 PetscDualSpaceDuplicate - Creates a duplicate PetscDualSpace object, however it is not setup. 466 467 Collective on sp 468 469 Input Parameter: 470 . sp - The original PetscDualSpace 471 472 Output Parameter: 473 . spNew - The duplicate PetscDualSpace 474 475 Level: beginner 476 477 .seealso: PetscDualSpaceCreate(), PetscDualSpaceSetType() 478 @*/ 479 PetscErrorCode PetscDualSpaceDuplicate(PetscDualSpace sp, PetscDualSpace *spNew) 480 { 481 DM dm; 482 PetscDualSpaceType type; 483 const char *name; 484 485 PetscFunctionBegin; 486 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 487 PetscValidPointer(spNew, 2); 488 PetscCall(PetscDualSpaceCreate(PetscObjectComm((PetscObject)sp), spNew)); 489 PetscCall(PetscObjectGetName((PetscObject) sp, &name)); 490 PetscCall(PetscObjectSetName((PetscObject) *spNew, name)); 491 PetscCall(PetscDualSpaceGetType(sp, &type)); 492 PetscCall(PetscDualSpaceSetType(*spNew, type)); 493 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 494 PetscCall(PetscDualSpaceSetDM(*spNew, dm)); 495 496 (*spNew)->order = sp->order; 497 (*spNew)->k = sp->k; 498 (*spNew)->Nc = sp->Nc; 499 (*spNew)->uniform = sp->uniform; 500 if (sp->ops->duplicate) PetscCall((*sp->ops->duplicate)(sp, *spNew)); 501 PetscFunctionReturn(0); 502 } 503 504 /*@ 505 PetscDualSpaceGetDM - Get the DM representing the reference cell 506 507 Not collective 508 509 Input Parameter: 510 . sp - The PetscDualSpace 511 512 Output Parameter: 513 . dm - The reference cell 514 515 Level: intermediate 516 517 .seealso: PetscDualSpaceSetDM(), PetscDualSpaceCreate() 518 @*/ 519 PetscErrorCode PetscDualSpaceGetDM(PetscDualSpace sp, DM *dm) 520 { 521 PetscFunctionBegin; 522 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 523 PetscValidPointer(dm, 2); 524 *dm = sp->dm; 525 PetscFunctionReturn(0); 526 } 527 528 /*@ 529 PetscDualSpaceSetDM - Get the DM representing the reference cell 530 531 Not collective 532 533 Input Parameters: 534 + sp - The PetscDualSpace 535 - dm - The reference cell 536 537 Level: intermediate 538 539 .seealso: PetscDualSpaceGetDM(), PetscDualSpaceCreate() 540 @*/ 541 PetscErrorCode PetscDualSpaceSetDM(PetscDualSpace sp, DM dm) 542 { 543 PetscFunctionBegin; 544 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 545 PetscValidHeaderSpecific(dm, DM_CLASSID, 2); 546 PetscCheck(!sp->setupcalled,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change DM after dualspace is set up"); 547 PetscCall(PetscObjectReference((PetscObject) dm)); 548 if (sp->dm && sp->dm != dm) { 549 PetscCall(PetscDualSpaceClearDMData_Internal(sp, sp->dm)); 550 } 551 PetscCall(DMDestroy(&sp->dm)); 552 sp->dm = dm; 553 PetscFunctionReturn(0); 554 } 555 556 /*@ 557 PetscDualSpaceGetOrder - Get the order of the dual space 558 559 Not collective 560 561 Input Parameter: 562 . sp - The PetscDualSpace 563 564 Output Parameter: 565 . order - The order 566 567 Level: intermediate 568 569 .seealso: PetscDualSpaceSetOrder(), PetscDualSpaceCreate() 570 @*/ 571 PetscErrorCode PetscDualSpaceGetOrder(PetscDualSpace sp, PetscInt *order) 572 { 573 PetscFunctionBegin; 574 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 575 PetscValidIntPointer(order, 2); 576 *order = sp->order; 577 PetscFunctionReturn(0); 578 } 579 580 /*@ 581 PetscDualSpaceSetOrder - Set the order of the dual space 582 583 Not collective 584 585 Input Parameters: 586 + sp - The PetscDualSpace 587 - order - The order 588 589 Level: intermediate 590 591 .seealso: PetscDualSpaceGetOrder(), PetscDualSpaceCreate() 592 @*/ 593 PetscErrorCode PetscDualSpaceSetOrder(PetscDualSpace sp, PetscInt order) 594 { 595 PetscFunctionBegin; 596 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 597 PetscCheck(!sp->setupcalled,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change order after dualspace is set up"); 598 sp->order = order; 599 PetscFunctionReturn(0); 600 } 601 602 /*@ 603 PetscDualSpaceGetNumComponents - Return the number of components for this space 604 605 Input Parameter: 606 . sp - The PetscDualSpace 607 608 Output Parameter: 609 . Nc - The number of components 610 611 Note: A vector space, for example, will have d components, where d is the spatial dimension 612 613 Level: intermediate 614 615 .seealso: PetscDualSpaceSetNumComponents(), PetscDualSpaceGetDimension(), PetscDualSpaceCreate(), PetscDualSpace 616 @*/ 617 PetscErrorCode PetscDualSpaceGetNumComponents(PetscDualSpace sp, PetscInt *Nc) 618 { 619 PetscFunctionBegin; 620 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 621 PetscValidIntPointer(Nc, 2); 622 *Nc = sp->Nc; 623 PetscFunctionReturn(0); 624 } 625 626 /*@ 627 PetscDualSpaceSetNumComponents - Set the number of components for this space 628 629 Input Parameters: 630 + sp - The PetscDualSpace 631 - order - The number of components 632 633 Level: intermediate 634 635 .seealso: PetscDualSpaceGetNumComponents(), PetscDualSpaceCreate(), PetscDualSpace 636 @*/ 637 PetscErrorCode PetscDualSpaceSetNumComponents(PetscDualSpace sp, PetscInt Nc) 638 { 639 PetscFunctionBegin; 640 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 641 PetscCheck(!sp->setupcalled,PetscObjectComm((PetscObject)sp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change number of components after dualspace is set up"); 642 sp->Nc = Nc; 643 PetscFunctionReturn(0); 644 } 645 646 /*@ 647 PetscDualSpaceGetFunctional - Get the i-th basis functional in the dual space 648 649 Not collective 650 651 Input Parameters: 652 + sp - The PetscDualSpace 653 - i - The basis number 654 655 Output Parameter: 656 . functional - The basis functional 657 658 Level: intermediate 659 660 .seealso: PetscDualSpaceGetDimension(), PetscDualSpaceCreate() 661 @*/ 662 PetscErrorCode PetscDualSpaceGetFunctional(PetscDualSpace sp, PetscInt i, PetscQuadrature *functional) 663 { 664 PetscInt dim; 665 666 PetscFunctionBegin; 667 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 668 PetscValidPointer(functional, 3); 669 PetscCall(PetscDualSpaceGetDimension(sp, &dim)); 670 PetscCheckFalse((i < 0) || (i >= dim),PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Functional index %d must be in [0, %d)", i, dim); 671 *functional = sp->functional[i]; 672 PetscFunctionReturn(0); 673 } 674 675 /*@ 676 PetscDualSpaceGetDimension - Get the dimension of the dual space, i.e. the number of basis functionals 677 678 Not collective 679 680 Input Parameter: 681 . sp - The PetscDualSpace 682 683 Output Parameter: 684 . dim - The dimension 685 686 Level: intermediate 687 688 .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate() 689 @*/ 690 PetscErrorCode PetscDualSpaceGetDimension(PetscDualSpace sp, PetscInt *dim) 691 { 692 PetscFunctionBegin; 693 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 694 PetscValidIntPointer(dim, 2); 695 if (sp->spdim < 0) { 696 PetscSection section; 697 698 PetscCall(PetscDualSpaceGetSection(sp, §ion)); 699 if (section) { 700 PetscCall(PetscSectionGetStorageSize(section, &(sp->spdim))); 701 } else sp->spdim = 0; 702 } 703 *dim = sp->spdim; 704 PetscFunctionReturn(0); 705 } 706 707 /*@ 708 PetscDualSpaceGetInteriorDimension - Get the interior dimension of the dual space, i.e. the number of basis functionals assigned to the interior of the reference domain 709 710 Not collective 711 712 Input Parameter: 713 . sp - The PetscDualSpace 714 715 Output Parameter: 716 . dim - The dimension 717 718 Level: intermediate 719 720 .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate() 721 @*/ 722 PetscErrorCode PetscDualSpaceGetInteriorDimension(PetscDualSpace sp, PetscInt *intdim) 723 { 724 PetscFunctionBegin; 725 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 726 PetscValidIntPointer(intdim, 2); 727 if (sp->spintdim < 0) { 728 PetscSection section; 729 730 PetscCall(PetscDualSpaceGetSection(sp, §ion)); 731 if (section) { 732 PetscCall(PetscSectionGetConstrainedStorageSize(section, &(sp->spintdim))); 733 } else sp->spintdim = 0; 734 } 735 *intdim = sp->spintdim; 736 PetscFunctionReturn(0); 737 } 738 739 /*@ 740 PetscDualSpaceGetUniform - Whether this dual space is uniform 741 742 Not collective 743 744 Input Parameters: 745 . sp - A dual space 746 747 Output Parameters: 748 . uniform - PETSC_TRUE if (a) the dual space is the same for each point in a stratum of the reference DMPlex, and 749 (b) every symmetry of each point in the reference DMPlex is also a symmetry of the point's dual space. 750 751 Level: advanced 752 753 Note: all of the usual spaces on simplex or tensor-product elements will be uniform, only reference cells 754 with non-uniform strata (like trianguar-prisms) or anisotropic hp dual spaces will not be uniform. 755 756 .seealso: PetscDualSpaceGetPointSubspace(), PetscDualSpaceGetSymmetries() 757 @*/ 758 PetscErrorCode PetscDualSpaceGetUniform(PetscDualSpace sp, PetscBool *uniform) 759 { 760 PetscFunctionBegin; 761 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 762 PetscValidBoolPointer(uniform, 2); 763 *uniform = sp->uniform; 764 PetscFunctionReturn(0); 765 } 766 767 /*@C 768 PetscDualSpaceGetNumDof - Get the number of degrees of freedom for each spatial (topological) dimension 769 770 Not collective 771 772 Input Parameter: 773 . sp - The PetscDualSpace 774 775 Output Parameter: 776 . numDof - An array of length dim+1 which holds the number of dofs for each dimension 777 778 Level: intermediate 779 780 .seealso: PetscDualSpaceGetFunctional(), PetscDualSpaceCreate() 781 @*/ 782 PetscErrorCode PetscDualSpaceGetNumDof(PetscDualSpace sp, const PetscInt **numDof) 783 { 784 PetscFunctionBegin; 785 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 786 PetscValidPointer(numDof, 2); 787 PetscCheck(sp->uniform,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "A non-uniform space does not have a fixed number of dofs for each height"); 788 if (!sp->numDof) { 789 DM dm; 790 PetscInt depth, d; 791 PetscSection section; 792 793 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 794 PetscCall(DMPlexGetDepth(dm, &depth)); 795 PetscCall(PetscCalloc1(depth+1,&(sp->numDof))); 796 PetscCall(PetscDualSpaceGetSection(sp, §ion)); 797 for (d = 0; d <= depth; d++) { 798 PetscInt dStart, dEnd; 799 800 PetscCall(DMPlexGetDepthStratum(dm, d, &dStart, &dEnd)); 801 if (dEnd <= dStart) continue; 802 PetscCall(PetscSectionGetDof(section, dStart, &(sp->numDof[d]))); 803 804 } 805 } 806 *numDof = sp->numDof; 807 PetscCheckFalse(!*numDof,PetscObjectComm((PetscObject) sp), PETSC_ERR_LIB, "Empty numDof[] returned from dual space implementation"); 808 PetscFunctionReturn(0); 809 } 810 811 /* create the section of the right size and set a permutation for topological ordering */ 812 PetscErrorCode PetscDualSpaceSectionCreate_Internal(PetscDualSpace sp, PetscSection *topSection) 813 { 814 DM dm; 815 PetscInt pStart, pEnd, cStart, cEnd, c, depth, count, i; 816 PetscInt *seen, *perm; 817 PetscSection section; 818 819 PetscFunctionBegin; 820 dm = sp->dm; 821 PetscCall(PetscSectionCreate(PETSC_COMM_SELF, §ion)); 822 PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 823 PetscCall(PetscSectionSetChart(section, pStart, pEnd)); 824 PetscCall(PetscCalloc1(pEnd - pStart, &seen)); 825 PetscCall(PetscMalloc1(pEnd - pStart, &perm)); 826 PetscCall(DMPlexGetDepth(dm, &depth)); 827 PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd)); 828 for (c = cStart, count = 0; c < cEnd; c++) { 829 PetscInt closureSize = -1, e; 830 PetscInt *closure = NULL; 831 832 perm[count++] = c; 833 seen[c-pStart] = 1; 834 PetscCall(DMPlexGetTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 835 for (e = 0; e < closureSize; e++) { 836 PetscInt point = closure[2*e]; 837 838 if (seen[point-pStart]) continue; 839 perm[count++] = point; 840 seen[point-pStart] = 1; 841 } 842 PetscCall(DMPlexRestoreTransitiveClosure(dm, c, PETSC_TRUE, &closureSize, &closure)); 843 } 844 PetscCheckFalse(count != pEnd - pStart,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Bad topological ordering"); 845 for (i = 0; i < pEnd - pStart; i++) if (perm[i] != i) break; 846 if (i < pEnd - pStart) { 847 IS permIS; 848 849 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, pEnd - pStart, perm, PETSC_OWN_POINTER, &permIS)); 850 PetscCall(ISSetPermutation(permIS)); 851 PetscCall(PetscSectionSetPermutation(section, permIS)); 852 PetscCall(ISDestroy(&permIS)); 853 } else { 854 PetscCall(PetscFree(perm)); 855 } 856 PetscCall(PetscFree(seen)); 857 *topSection = section; 858 PetscFunctionReturn(0); 859 } 860 861 /* mark boundary points and set up */ 862 PetscErrorCode PetscDualSpaceSectionSetUp_Internal(PetscDualSpace sp, PetscSection section) 863 { 864 DM dm; 865 DMLabel boundary; 866 PetscInt pStart, pEnd, p; 867 868 PetscFunctionBegin; 869 dm = sp->dm; 870 PetscCall(DMLabelCreate(PETSC_COMM_SELF,"boundary",&boundary)); 871 PetscCall(PetscDualSpaceGetDM(sp,&dm)); 872 PetscCall(DMPlexMarkBoundaryFaces(dm,1,boundary)); 873 PetscCall(DMPlexLabelComplete(dm,boundary)); 874 PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 875 for (p = pStart; p < pEnd; p++) { 876 PetscInt bval; 877 878 PetscCall(DMLabelGetValue(boundary, p, &bval)); 879 if (bval == 1) { 880 PetscInt dof; 881 882 PetscCall(PetscSectionGetDof(section, p, &dof)); 883 PetscCall(PetscSectionSetConstraintDof(section, p, dof)); 884 } 885 } 886 PetscCall(DMLabelDestroy(&boundary)); 887 PetscCall(PetscSectionSetUp(section)); 888 PetscFunctionReturn(0); 889 } 890 891 /*@ 892 PetscDualSpaceGetSection - Create a PetscSection over the reference cell with the layout from this space 893 894 Collective on sp 895 896 Input Parameters: 897 . sp - The PetscDualSpace 898 899 Output Parameter: 900 . section - The section 901 902 Level: advanced 903 904 .seealso: PetscDualSpaceCreate(), DMPLEX 905 @*/ 906 PetscErrorCode PetscDualSpaceGetSection(PetscDualSpace sp, PetscSection *section) 907 { 908 PetscInt pStart, pEnd, p; 909 910 PetscFunctionBegin; 911 if (!sp->pointSection) { 912 /* mark the boundary */ 913 PetscCall(PetscDualSpaceSectionCreate_Internal(sp, &(sp->pointSection))); 914 PetscCall(DMPlexGetChart(sp->dm,&pStart,&pEnd)); 915 for (p = pStart; p < pEnd; p++) { 916 PetscDualSpace psp; 917 918 PetscCall(PetscDualSpaceGetPointSubspace(sp, p, &psp)); 919 if (psp) { 920 PetscInt dof; 921 922 PetscCall(PetscDualSpaceGetInteriorDimension(psp, &dof)); 923 PetscCall(PetscSectionSetDof(sp->pointSection,p,dof)); 924 } 925 } 926 PetscCall(PetscDualSpaceSectionSetUp_Internal(sp,sp->pointSection)); 927 } 928 *section = sp->pointSection; 929 PetscFunctionReturn(0); 930 } 931 932 /* this assumes that all of the point dual spaces store their interior dofs first, which is true when the point DMs 933 * have one cell */ 934 PetscErrorCode PetscDualSpacePushForwardSubspaces_Internal(PetscDualSpace sp, PetscInt sStart, PetscInt sEnd) 935 { 936 PetscReal *sv0, *v0, *J; 937 PetscSection section; 938 PetscInt dim, s, k; 939 DM dm; 940 941 PetscFunctionBegin; 942 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 943 PetscCall(DMGetDimension(dm, &dim)); 944 PetscCall(PetscDualSpaceGetSection(sp, §ion)); 945 PetscCall(PetscMalloc3(dim, &v0, dim, &sv0, dim*dim, &J)); 946 PetscCall(PetscDualSpaceGetFormDegree(sp, &k)); 947 for (s = sStart; s < sEnd; s++) { 948 PetscReal detJ, hdetJ; 949 PetscDualSpace ssp; 950 PetscInt dof, off, f, sdim; 951 PetscInt i, j; 952 DM sdm; 953 954 PetscCall(PetscDualSpaceGetPointSubspace(sp, s, &ssp)); 955 if (!ssp) continue; 956 PetscCall(PetscSectionGetDof(section, s, &dof)); 957 PetscCall(PetscSectionGetOffset(section, s, &off)); 958 /* get the first vertex of the reference cell */ 959 PetscCall(PetscDualSpaceGetDM(ssp, &sdm)); 960 PetscCall(DMGetDimension(sdm, &sdim)); 961 PetscCall(DMPlexComputeCellGeometryAffineFEM(sdm, 0, sv0, NULL, NULL, &hdetJ)); 962 PetscCall(DMPlexComputeCellGeometryAffineFEM(dm, s, v0, J, NULL, &detJ)); 963 /* compactify Jacobian */ 964 for (i = 0; i < dim; i++) for (j = 0; j < sdim; j++) J[i* sdim + j] = J[i * dim + j]; 965 for (f = 0; f < dof; f++) { 966 PetscQuadrature fn; 967 968 PetscCall(PetscDualSpaceGetFunctional(ssp, f, &fn)); 969 PetscCall(PetscQuadraturePushForward(fn, dim, sv0, v0, J, k, &(sp->functional[off+f]))); 970 } 971 } 972 PetscCall(PetscFree3(v0, sv0, J)); 973 PetscFunctionReturn(0); 974 } 975 976 /*@C 977 PetscDualSpaceApply - Apply a functional from the dual space basis to an input function 978 979 Input Parameters: 980 + sp - The PetscDualSpace object 981 . f - The basis functional index 982 . time - The time 983 . 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) 984 . numComp - The number of components for the function 985 . func - The input function 986 - ctx - A context for the function 987 988 Output Parameter: 989 . value - numComp output values 990 991 Note: The calling sequence for the callback func is given by: 992 993 $ func(PetscInt dim, PetscReal time, const PetscReal x[], 994 $ PetscInt numComponents, PetscScalar values[], void *ctx) 995 996 Level: beginner 997 998 .seealso: PetscDualSpaceCreate() 999 @*/ 1000 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) 1001 { 1002 PetscFunctionBegin; 1003 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1004 PetscValidPointer(cgeom, 4); 1005 PetscValidScalarPointer(value, 8); 1006 PetscCall((*sp->ops->apply)(sp, f, time, cgeom, numComp, func, ctx, value)); 1007 PetscFunctionReturn(0); 1008 } 1009 1010 /*@C 1011 PetscDualSpaceApplyAll - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllData() 1012 1013 Input Parameters: 1014 + sp - The PetscDualSpace object 1015 - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllData() 1016 1017 Output Parameter: 1018 . spValue - The values of all dual space functionals 1019 1020 Level: beginner 1021 1022 .seealso: PetscDualSpaceCreate() 1023 @*/ 1024 PetscErrorCode PetscDualSpaceApplyAll(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue) 1025 { 1026 PetscFunctionBegin; 1027 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1028 PetscCall((*sp->ops->applyall)(sp, pointEval, spValue)); 1029 PetscFunctionReturn(0); 1030 } 1031 1032 /*@C 1033 PetscDualSpaceApplyInterior - Apply interior functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetInteriorData() 1034 1035 Input Parameters: 1036 + sp - The PetscDualSpace object 1037 - pointEval - Evaluation at the points returned by PetscDualSpaceGetInteriorData() 1038 1039 Output Parameter: 1040 . spValue - The values of interior dual space functionals 1041 1042 Level: beginner 1043 1044 .seealso: PetscDualSpaceCreate() 1045 @*/ 1046 PetscErrorCode PetscDualSpaceApplyInterior(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue) 1047 { 1048 PetscFunctionBegin; 1049 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1050 PetscCall((*sp->ops->applyint)(sp, pointEval, spValue)); 1051 PetscFunctionReturn(0); 1052 } 1053 1054 /*@C 1055 PetscDualSpaceApplyDefault - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional. 1056 1057 Input Parameters: 1058 + sp - The PetscDualSpace object 1059 . f - The basis functional index 1060 . time - The time 1061 . cgeom - A context with geometric information for this cell, we use v0 (the initial vertex) and J (the Jacobian) 1062 . Nc - The number of components for the function 1063 . func - The input function 1064 - ctx - A context for the function 1065 1066 Output Parameter: 1067 . value - The output value 1068 1069 Note: The calling sequence for the callback func is given by: 1070 1071 $ func(PetscInt dim, PetscReal time, const PetscReal x[], 1072 $ PetscInt numComponents, PetscScalar values[], void *ctx) 1073 1074 and the idea is to evaluate the functional as an integral 1075 1076 $ n(f) = int dx n(x) . f(x) 1077 1078 where both n and f have Nc components. 1079 1080 Level: beginner 1081 1082 .seealso: PetscDualSpaceCreate() 1083 @*/ 1084 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) 1085 { 1086 DM dm; 1087 PetscQuadrature n; 1088 const PetscReal *points, *weights; 1089 PetscReal x[3]; 1090 PetscScalar *val; 1091 PetscInt dim, dE, qNc, c, Nq, q; 1092 PetscBool isAffine; 1093 1094 PetscFunctionBegin; 1095 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1096 PetscValidScalarPointer(value, 8); 1097 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 1098 PetscCall(PetscDualSpaceGetFunctional(sp, f, &n)); 1099 PetscCall(PetscQuadratureGetData(n, &dim, &qNc, &Nq, &points, &weights)); 1100 PetscCheckFalse(dim != cgeom->dim,PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature spatial dimension %D != cell geometry dimension %D", dim, cgeom->dim); 1101 PetscCheckFalse(qNc != Nc,PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc); 1102 PetscCall(DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val)); 1103 *value = 0.0; 1104 isAffine = cgeom->isAffine; 1105 dE = cgeom->dimEmbed; 1106 for (q = 0; q < Nq; ++q) { 1107 if (isAffine) { 1108 CoordinatesRefToReal(dE, cgeom->dim, cgeom->xi, cgeom->v, cgeom->J, &points[q*dim], x); 1109 PetscCall((*func)(dE, time, x, Nc, val, ctx)); 1110 } else { 1111 PetscCall((*func)(dE, time, &cgeom->v[dE*q], Nc, val, ctx)); 1112 } 1113 for (c = 0; c < Nc; ++c) { 1114 *value += val[c]*weights[q*Nc+c]; 1115 } 1116 } 1117 PetscCall(DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val)); 1118 PetscFunctionReturn(0); 1119 } 1120 1121 /*@C 1122 PetscDualSpaceApplyAllDefault - Apply all functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetAllData() 1123 1124 Input Parameters: 1125 + sp - The PetscDualSpace object 1126 - pointEval - Evaluation at the points returned by PetscDualSpaceGetAllData() 1127 1128 Output Parameter: 1129 . spValue - The values of all dual space functionals 1130 1131 Level: beginner 1132 1133 .seealso: PetscDualSpaceCreate() 1134 @*/ 1135 PetscErrorCode PetscDualSpaceApplyAllDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue) 1136 { 1137 Vec pointValues, dofValues; 1138 Mat allMat; 1139 1140 PetscFunctionBegin; 1141 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1142 PetscValidScalarPointer(pointEval, 2); 1143 PetscValidScalarPointer(spValue, 3); 1144 PetscCall(PetscDualSpaceGetAllData(sp, NULL, &allMat)); 1145 if (!(sp->allNodeValues)) { 1146 PetscCall(MatCreateVecs(allMat, &(sp->allNodeValues), NULL)); 1147 } 1148 pointValues = sp->allNodeValues; 1149 if (!(sp->allDofValues)) { 1150 PetscCall(MatCreateVecs(allMat, NULL, &(sp->allDofValues))); 1151 } 1152 dofValues = sp->allDofValues; 1153 PetscCall(VecPlaceArray(pointValues, pointEval)); 1154 PetscCall(VecPlaceArray(dofValues, spValue)); 1155 PetscCall(MatMult(allMat, pointValues, dofValues)); 1156 PetscCall(VecResetArray(dofValues)); 1157 PetscCall(VecResetArray(pointValues)); 1158 PetscFunctionReturn(0); 1159 } 1160 1161 /*@C 1162 PetscDualSpaceApplyInteriorDefault - Apply interior functionals from the dual space basis to the result of an evaluation at the points returned by PetscDualSpaceGetInteriorData() 1163 1164 Input Parameters: 1165 + sp - The PetscDualSpace object 1166 - pointEval - Evaluation at the points returned by PetscDualSpaceGetInteriorData() 1167 1168 Output Parameter: 1169 . spValue - The values of interior dual space functionals 1170 1171 Level: beginner 1172 1173 .seealso: PetscDualSpaceCreate() 1174 @*/ 1175 PetscErrorCode PetscDualSpaceApplyInteriorDefault(PetscDualSpace sp, const PetscScalar *pointEval, PetscScalar *spValue) 1176 { 1177 Vec pointValues, dofValues; 1178 Mat intMat; 1179 1180 PetscFunctionBegin; 1181 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1182 PetscValidScalarPointer(pointEval, 2); 1183 PetscValidScalarPointer(spValue, 3); 1184 PetscCall(PetscDualSpaceGetInteriorData(sp, NULL, &intMat)); 1185 if (!(sp->intNodeValues)) { 1186 PetscCall(MatCreateVecs(intMat, &(sp->intNodeValues), NULL)); 1187 } 1188 pointValues = sp->intNodeValues; 1189 if (!(sp->intDofValues)) { 1190 PetscCall(MatCreateVecs(intMat, NULL, &(sp->intDofValues))); 1191 } 1192 dofValues = sp->intDofValues; 1193 PetscCall(VecPlaceArray(pointValues, pointEval)); 1194 PetscCall(VecPlaceArray(dofValues, spValue)); 1195 PetscCall(MatMult(intMat, pointValues, dofValues)); 1196 PetscCall(VecResetArray(dofValues)); 1197 PetscCall(VecResetArray(pointValues)); 1198 PetscFunctionReturn(0); 1199 } 1200 1201 /*@ 1202 PetscDualSpaceGetAllData - Get all quadrature nodes from this space, and the matrix that sends quadrature node values to degree-of-freedom values 1203 1204 Input Parameter: 1205 . sp - The dualspace 1206 1207 Output Parameters: 1208 + allNodes - A PetscQuadrature object containing all evaluation nodes 1209 - allMat - A Mat for the node-to-dof transformation 1210 1211 Level: advanced 1212 1213 .seealso: PetscDualSpaceCreate() 1214 @*/ 1215 PetscErrorCode PetscDualSpaceGetAllData(PetscDualSpace sp, PetscQuadrature *allNodes, Mat *allMat) 1216 { 1217 PetscFunctionBegin; 1218 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1219 if (allNodes) PetscValidPointer(allNodes,2); 1220 if (allMat) PetscValidPointer(allMat,3); 1221 if ((!sp->allNodes || !sp->allMat) && sp->ops->createalldata) { 1222 PetscQuadrature qpoints; 1223 Mat amat; 1224 1225 PetscCall((*sp->ops->createalldata)(sp,&qpoints,&amat)); 1226 PetscCall(PetscQuadratureDestroy(&(sp->allNodes))); 1227 PetscCall(MatDestroy(&(sp->allMat))); 1228 sp->allNodes = qpoints; 1229 sp->allMat = amat; 1230 } 1231 if (allNodes) *allNodes = sp->allNodes; 1232 if (allMat) *allMat = sp->allMat; 1233 PetscFunctionReturn(0); 1234 } 1235 1236 /*@ 1237 PetscDualSpaceCreateAllDataDefault - Create all evaluation nodes and the node-to-dof matrix by examining functionals 1238 1239 Input Parameter: 1240 . sp - The dualspace 1241 1242 Output Parameters: 1243 + allNodes - A PetscQuadrature object containing all evaluation nodes 1244 - allMat - A Mat for the node-to-dof transformation 1245 1246 Level: advanced 1247 1248 .seealso: PetscDualSpaceCreate() 1249 @*/ 1250 PetscErrorCode PetscDualSpaceCreateAllDataDefault(PetscDualSpace sp, PetscQuadrature *allNodes, Mat *allMat) 1251 { 1252 PetscInt spdim; 1253 PetscInt numPoints, offset; 1254 PetscReal *points; 1255 PetscInt f, dim; 1256 PetscInt Nc, nrows, ncols; 1257 PetscInt maxNumPoints; 1258 PetscQuadrature q; 1259 Mat A; 1260 1261 PetscFunctionBegin; 1262 PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc)); 1263 PetscCall(PetscDualSpaceGetDimension(sp,&spdim)); 1264 if (!spdim) { 1265 PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF,allNodes)); 1266 PetscCall(PetscQuadratureSetData(*allNodes,0,0,0,NULL,NULL)); 1267 } 1268 nrows = spdim; 1269 PetscCall(PetscDualSpaceGetFunctional(sp,0,&q)); 1270 PetscCall(PetscQuadratureGetData(q,&dim,NULL,&numPoints,NULL,NULL)); 1271 maxNumPoints = numPoints; 1272 for (f = 1; f < spdim; f++) { 1273 PetscInt Np; 1274 1275 PetscCall(PetscDualSpaceGetFunctional(sp,f,&q)); 1276 PetscCall(PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL)); 1277 numPoints += Np; 1278 maxNumPoints = PetscMax(maxNumPoints,Np); 1279 } 1280 ncols = numPoints * Nc; 1281 PetscCall(PetscMalloc1(dim*numPoints,&points)); 1282 PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, nrows, ncols, maxNumPoints * Nc, NULL, &A)); 1283 for (f = 0, offset = 0; f < spdim; f++) { 1284 const PetscReal *p, *w; 1285 PetscInt Np, i; 1286 PetscInt fnc; 1287 1288 PetscCall(PetscDualSpaceGetFunctional(sp,f,&q)); 1289 PetscCall(PetscQuadratureGetData(q,NULL,&fnc,&Np,&p,&w)); 1290 PetscCheckFalse(fnc != Nc,PETSC_COMM_SELF, PETSC_ERR_PLIB, "functional component mismatch"); 1291 for (i = 0; i < Np * dim; i++) { 1292 points[offset* dim + i] = p[i]; 1293 } 1294 for (i = 0; i < Np * Nc; i++) { 1295 PetscCall(MatSetValue(A, f, offset * Nc, w[i], INSERT_VALUES)); 1296 } 1297 offset += Np; 1298 } 1299 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1300 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1301 PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF,allNodes)); 1302 PetscCall(PetscQuadratureSetData(*allNodes,dim,0,numPoints,points,NULL)); 1303 *allMat = A; 1304 PetscFunctionReturn(0); 1305 } 1306 1307 /*@ 1308 PetscDualSpaceGetInteriorData - Get all quadrature points necessary to compute the interior degrees of freedom from 1309 this space, as well as the matrix that computes the degrees of freedom from the quadrature values. Degrees of 1310 freedom are interior degrees of freedom if they belong (by PetscDualSpaceGetSection()) to interior points in the 1311 reference DMPlex: complementary boundary degrees of freedom are marked as constrained in the section returned by 1312 PetscDualSpaceGetSection()). 1313 1314 Input Parameter: 1315 . sp - The dualspace 1316 1317 Output Parameters: 1318 + intNodes - A PetscQuadrature object containing all evaluation points needed to evaluate interior degrees of freedom 1319 - intMat - A matrix that computes dual space values from point values: size [spdim0 x (npoints * nc)], where spdim0 is 1320 the size of the constrained layout (PetscSectionGetConstrainStorageSize()) of the dual space section, 1321 npoints is the number of points in intNodes and nc is PetscDualSpaceGetNumComponents(). 1322 1323 Level: advanced 1324 1325 .seealso: PetscDualSpaceCreate(), PetscDualSpaceGetDimension(), PetscDualSpaceGetNumComponents(), PetscQuadratureGetData() 1326 @*/ 1327 PetscErrorCode PetscDualSpaceGetInteriorData(PetscDualSpace sp, PetscQuadrature *intNodes, Mat *intMat) 1328 { 1329 PetscFunctionBegin; 1330 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1331 if (intNodes) PetscValidPointer(intNodes,2); 1332 if (intMat) PetscValidPointer(intMat,3); 1333 if ((!sp->intNodes || !sp->intMat) && sp->ops->createintdata) { 1334 PetscQuadrature qpoints; 1335 Mat imat; 1336 1337 PetscCall((*sp->ops->createintdata)(sp,&qpoints,&imat)); 1338 PetscCall(PetscQuadratureDestroy(&(sp->intNodes))); 1339 PetscCall(MatDestroy(&(sp->intMat))); 1340 sp->intNodes = qpoints; 1341 sp->intMat = imat; 1342 } 1343 if (intNodes) *intNodes = sp->intNodes; 1344 if (intMat) *intMat = sp->intMat; 1345 PetscFunctionReturn(0); 1346 } 1347 1348 /*@ 1349 PetscDualSpaceCreateInteriorDataDefault - Create quadrature points by examining interior functionals and create the matrix mapping quadrature point values to interior dual space values 1350 1351 Input Parameter: 1352 . sp - The dualspace 1353 1354 Output Parameters: 1355 + intNodes - A PetscQuadrature object containing all evaluation points needed to evaluate interior degrees of freedom 1356 - intMat - A matrix that computes dual space values from point values: size [spdim0 x (npoints * nc)], where spdim0 is 1357 the size of the constrained layout (PetscSectionGetConstrainStorageSize()) of the dual space section, 1358 npoints is the number of points in allNodes and nc is PetscDualSpaceGetNumComponents(). 1359 1360 Level: advanced 1361 1362 .seealso: PetscDualSpaceCreate(), PetscDualSpaceGetInteriorData() 1363 @*/ 1364 PetscErrorCode PetscDualSpaceCreateInteriorDataDefault(PetscDualSpace sp, PetscQuadrature *intNodes, Mat *intMat) 1365 { 1366 DM dm; 1367 PetscInt spdim0; 1368 PetscInt Nc; 1369 PetscInt pStart, pEnd, p, f; 1370 PetscSection section; 1371 PetscInt numPoints, offset, matoffset; 1372 PetscReal *points; 1373 PetscInt dim; 1374 PetscInt *nnz; 1375 PetscQuadrature q; 1376 Mat imat; 1377 1378 PetscFunctionBegin; 1379 PetscValidHeaderSpecific(sp,PETSCDUALSPACE_CLASSID,1); 1380 PetscCall(PetscDualSpaceGetSection(sp, §ion)); 1381 PetscCall(PetscSectionGetConstrainedStorageSize(section, &spdim0)); 1382 if (!spdim0) { 1383 *intNodes = NULL; 1384 *intMat = NULL; 1385 PetscFunctionReturn(0); 1386 } 1387 PetscCall(PetscDualSpaceGetNumComponents(sp, &Nc)); 1388 PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 1389 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 1390 PetscCall(DMGetDimension(dm, &dim)); 1391 PetscCall(PetscMalloc1(spdim0, &nnz)); 1392 for (p = pStart, f = 0, numPoints = 0; p < pEnd; p++) { 1393 PetscInt dof, cdof, off, d; 1394 1395 PetscCall(PetscSectionGetDof(section, p, &dof)); 1396 PetscCall(PetscSectionGetConstraintDof(section, p, &cdof)); 1397 if (!(dof - cdof)) continue; 1398 PetscCall(PetscSectionGetOffset(section, p, &off)); 1399 for (d = 0; d < dof; d++, off++, f++) { 1400 PetscInt Np; 1401 1402 PetscCall(PetscDualSpaceGetFunctional(sp,off,&q)); 1403 PetscCall(PetscQuadratureGetData(q,NULL,NULL,&Np,NULL,NULL)); 1404 nnz[f] = Np * Nc; 1405 numPoints += Np; 1406 } 1407 } 1408 PetscCall(MatCreateSeqAIJ(PETSC_COMM_SELF, spdim0, numPoints * Nc, 0, nnz, &imat)); 1409 PetscCall(PetscFree(nnz)); 1410 PetscCall(PetscMalloc1(dim*numPoints,&points)); 1411 for (p = pStart, f = 0, offset = 0, matoffset = 0; p < pEnd; p++) { 1412 PetscInt dof, cdof, off, d; 1413 1414 PetscCall(PetscSectionGetDof(section, p, &dof)); 1415 PetscCall(PetscSectionGetConstraintDof(section, p, &cdof)); 1416 if (!(dof - cdof)) continue; 1417 PetscCall(PetscSectionGetOffset(section, p, &off)); 1418 for (d = 0; d < dof; d++, off++, f++) { 1419 const PetscReal *p; 1420 const PetscReal *w; 1421 PetscInt Np, i; 1422 1423 PetscCall(PetscDualSpaceGetFunctional(sp,off,&q)); 1424 PetscCall(PetscQuadratureGetData(q,NULL,NULL,&Np,&p,&w)); 1425 for (i = 0; i < Np * dim; i++) { 1426 points[offset + i] = p[i]; 1427 } 1428 for (i = 0; i < Np * Nc; i++) { 1429 PetscCall(MatSetValue(imat, f, matoffset + i, w[i],INSERT_VALUES)); 1430 } 1431 offset += Np * dim; 1432 matoffset += Np * Nc; 1433 } 1434 } 1435 PetscCall(PetscQuadratureCreate(PETSC_COMM_SELF,intNodes)); 1436 PetscCall(PetscQuadratureSetData(*intNodes,dim,0,numPoints,points,NULL)); 1437 PetscCall(MatAssemblyBegin(imat, MAT_FINAL_ASSEMBLY)); 1438 PetscCall(MatAssemblyEnd(imat, MAT_FINAL_ASSEMBLY)); 1439 *intMat = imat; 1440 PetscFunctionReturn(0); 1441 } 1442 1443 /*@ 1444 PetscDualSpaceEqual - Determine if a dual space is equivalent 1445 1446 Input Parameters: 1447 + A - A PetscDualSpace object 1448 - B - Another PetscDualSpace object 1449 1450 Output Parameter: 1451 . equal - PETSC_TRUE if the dual spaces are equivalent 1452 1453 Level: advanced 1454 1455 .seealso: PetscDualSpaceCreate() 1456 @*/ 1457 PetscErrorCode PetscDualSpaceEqual(PetscDualSpace A, PetscDualSpace B, PetscBool *equal) 1458 { 1459 PetscInt sizeA, sizeB, dimA, dimB; 1460 const PetscInt *dofA, *dofB; 1461 PetscQuadrature quadA, quadB; 1462 Mat matA, matB; 1463 1464 PetscFunctionBegin; 1465 PetscValidHeaderSpecific(A,PETSCDUALSPACE_CLASSID,1); 1466 PetscValidHeaderSpecific(B,PETSCDUALSPACE_CLASSID,2); 1467 PetscValidBoolPointer(equal,3); 1468 *equal = PETSC_FALSE; 1469 PetscCall(PetscDualSpaceGetDimension(A, &sizeA)); 1470 PetscCall(PetscDualSpaceGetDimension(B, &sizeB)); 1471 if (sizeB != sizeA) { 1472 PetscFunctionReturn(0); 1473 } 1474 PetscCall(DMGetDimension(A->dm, &dimA)); 1475 PetscCall(DMGetDimension(B->dm, &dimB)); 1476 if (dimA != dimB) { 1477 PetscFunctionReturn(0); 1478 } 1479 1480 PetscCall(PetscDualSpaceGetNumDof(A, &dofA)); 1481 PetscCall(PetscDualSpaceGetNumDof(B, &dofB)); 1482 for (PetscInt d=0; d<dimA; d++) { 1483 if (dofA[d] != dofB[d]) { 1484 PetscFunctionReturn(0); 1485 } 1486 } 1487 1488 PetscCall(PetscDualSpaceGetInteriorData(A, &quadA, &matA)); 1489 PetscCall(PetscDualSpaceGetInteriorData(B, &quadB, &matB)); 1490 if (!quadA && !quadB) { 1491 *equal = PETSC_TRUE; 1492 } else if (quadA && quadB) { 1493 PetscCall(PetscQuadratureEqual(quadA, quadB, equal)); 1494 if (*equal == PETSC_FALSE) PetscFunctionReturn(0); 1495 if (!matA && !matB) PetscFunctionReturn(0); 1496 if (matA && matB) PetscCall(MatEqual(matA, matB, equal)); 1497 else *equal = PETSC_FALSE; 1498 } 1499 PetscFunctionReturn(0); 1500 } 1501 1502 /*@C 1503 PetscDualSpaceApplyFVM - Apply a functional from the dual space basis to an input function by assuming a point evaluation functional at the cell centroid. 1504 1505 Input Parameters: 1506 + sp - The PetscDualSpace object 1507 . f - The basis functional index 1508 . time - The time 1509 . cgeom - A context with geometric information for this cell, we currently just use the centroid 1510 . Nc - The number of components for the function 1511 . func - The input function 1512 - ctx - A context for the function 1513 1514 Output Parameter: 1515 . value - The output value (scalar) 1516 1517 Note: The calling sequence for the callback func is given by: 1518 1519 $ func(PetscInt dim, PetscReal time, const PetscReal x[], 1520 $ PetscInt numComponents, PetscScalar values[], void *ctx) 1521 1522 and the idea is to evaluate the functional as an integral 1523 1524 $ n(f) = int dx n(x) . f(x) 1525 1526 where both n and f have Nc components. 1527 1528 Level: beginner 1529 1530 .seealso: PetscDualSpaceCreate() 1531 @*/ 1532 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) 1533 { 1534 DM dm; 1535 PetscQuadrature n; 1536 const PetscReal *points, *weights; 1537 PetscScalar *val; 1538 PetscInt dimEmbed, qNc, c, Nq, q; 1539 1540 PetscFunctionBegin; 1541 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1542 PetscValidScalarPointer(value, 8); 1543 PetscCall(PetscDualSpaceGetDM(sp, &dm)); 1544 PetscCall(DMGetCoordinateDim(dm, &dimEmbed)); 1545 PetscCall(PetscDualSpaceGetFunctional(sp, f, &n)); 1546 PetscCall(PetscQuadratureGetData(n, NULL, &qNc, &Nq, &points, &weights)); 1547 PetscCheckFalse(qNc != Nc,PetscObjectComm((PetscObject) sp), PETSC_ERR_ARG_SIZ, "The quadrature components %D != function components %D", qNc, Nc); 1548 PetscCall(DMGetWorkArray(dm, Nc, MPIU_SCALAR, &val)); 1549 *value = 0.; 1550 for (q = 0; q < Nq; ++q) { 1551 PetscCall((*func)(dimEmbed, time, cgeom->centroid, Nc, val, ctx)); 1552 for (c = 0; c < Nc; ++c) { 1553 *value += val[c]*weights[q*Nc+c]; 1554 } 1555 } 1556 PetscCall(DMRestoreWorkArray(dm, Nc, MPIU_SCALAR, &val)); 1557 PetscFunctionReturn(0); 1558 } 1559 1560 /*@ 1561 PetscDualSpaceGetHeightSubspace - Get the subset of the dual space basis that is supported on a mesh point of a 1562 given height. This assumes that the reference cell is symmetric over points of this height. 1563 1564 If the dual space is not defined on mesh points of the given height (e.g. if the space is discontinuous and 1565 pointwise values are not defined on the element boundaries), or if the implementation of PetscDualSpace does not 1566 support extracting subspaces, then NULL is returned. 1567 1568 This does not increment the reference count on the returned dual space, and the user should not destroy it. 1569 1570 Not collective 1571 1572 Input Parameters: 1573 + sp - the PetscDualSpace object 1574 - height - the height of the mesh point for which the subspace is desired 1575 1576 Output Parameter: 1577 . subsp - the subspace. Note that the functionals in the subspace are with respect to the intrinsic geometry of the 1578 point, which will be of lesser dimension if height > 0. 1579 1580 Level: advanced 1581 1582 .seealso: PetscSpaceGetHeightSubspace(), PetscDualSpace 1583 @*/ 1584 PetscErrorCode PetscDualSpaceGetHeightSubspace(PetscDualSpace sp, PetscInt height, PetscDualSpace *subsp) 1585 { 1586 PetscInt depth = -1, cStart, cEnd; 1587 DM dm; 1588 1589 PetscFunctionBegin; 1590 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1591 PetscValidPointer(subsp,3); 1592 PetscCheckFalse(!(sp->uniform),PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "A non-uniform dual space does not have a single dual space at each height"); 1593 *subsp = NULL; 1594 dm = sp->dm; 1595 PetscCall(DMPlexGetDepth(dm, &depth)); 1596 PetscCheckFalse(height < 0 || height > depth,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid height"); 1597 PetscCall(DMPlexGetHeightStratum(dm,0,&cStart,&cEnd)); 1598 if (height == 0 && cEnd == cStart + 1) { 1599 *subsp = sp; 1600 PetscFunctionReturn(0); 1601 } 1602 if (!sp->heightSpaces) { 1603 PetscInt h; 1604 PetscCall(PetscCalloc1(depth+1, &(sp->heightSpaces))); 1605 1606 for (h = 0; h <= depth; h++) { 1607 if (h == 0 && cEnd == cStart + 1) continue; 1608 if (sp->ops->createheightsubspace) PetscCall((*sp->ops->createheightsubspace)(sp,height,&(sp->heightSpaces[h]))); 1609 else if (sp->pointSpaces) { 1610 PetscInt hStart, hEnd; 1611 1612 PetscCall(DMPlexGetHeightStratum(dm,h,&hStart,&hEnd)); 1613 if (hEnd > hStart) { 1614 const char *name; 1615 1616 PetscCall(PetscObjectReference((PetscObject)(sp->pointSpaces[hStart]))); 1617 if (sp->pointSpaces[hStart]) { 1618 PetscCall(PetscObjectGetName((PetscObject) sp, &name)); 1619 PetscCall(PetscObjectSetName((PetscObject) sp->pointSpaces[hStart], name)); 1620 } 1621 sp->heightSpaces[h] = sp->pointSpaces[hStart]; 1622 } 1623 } 1624 } 1625 } 1626 *subsp = sp->heightSpaces[height]; 1627 PetscFunctionReturn(0); 1628 } 1629 1630 /*@ 1631 PetscDualSpaceGetPointSubspace - Get the subset of the dual space basis that is supported on a particular mesh point. 1632 1633 If the dual space is not defined on the mesh point (e.g. if the space is discontinuous and pointwise values are not 1634 defined on the element boundaries), or if the implementation of PetscDualSpace does not support extracting 1635 subspaces, then NULL is returned. 1636 1637 This does not increment the reference count on the returned dual space, and the user should not destroy it. 1638 1639 Not collective 1640 1641 Input Parameters: 1642 + sp - the PetscDualSpace object 1643 - point - the point (in the dual space's DM) for which the subspace is desired 1644 1645 Output Parameters: 1646 bdsp - the subspace. Note that the functionals in the subspace are with respect to the intrinsic geometry of the 1647 point, which will be of lesser dimension if height > 0. 1648 1649 Level: advanced 1650 1651 .seealso: PetscDualSpace 1652 @*/ 1653 PetscErrorCode PetscDualSpaceGetPointSubspace(PetscDualSpace sp, PetscInt point, PetscDualSpace *bdsp) 1654 { 1655 PetscInt pStart = 0, pEnd = 0, cStart, cEnd; 1656 DM dm; 1657 1658 PetscFunctionBegin; 1659 PetscValidHeaderSpecific(sp, PETSCDUALSPACE_CLASSID, 1); 1660 PetscValidPointer(bdsp,3); 1661 *bdsp = NULL; 1662 dm = sp->dm; 1663 PetscCall(DMPlexGetChart(dm, &pStart, &pEnd)); 1664 PetscCheckFalse(point < pStart || point > pEnd,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid point"); 1665 PetscCall(DMPlexGetHeightStratum(dm,0,&cStart,&cEnd)); 1666 if (point == cStart && cEnd == cStart + 1) { /* the dual space is only equivalent to the dual space on a cell if the reference mesh has just one cell */ 1667 *bdsp = sp; 1668 PetscFunctionReturn(0); 1669 } 1670 if (!sp->pointSpaces) { 1671 PetscInt p; 1672 PetscCall(PetscCalloc1(pEnd - pStart, &(sp->pointSpaces))); 1673 1674 for (p = 0; p < pEnd - pStart; p++) { 1675 if (p + pStart == cStart && cEnd == cStart + 1) continue; 1676 if (sp->ops->createpointsubspace) PetscCall((*sp->ops->createpointsubspace)(sp,p+pStart,&(sp->pointSpaces[p]))); 1677 else if (sp->heightSpaces || sp->ops->createheightsubspace) { 1678 PetscInt dim, depth, height; 1679 DMLabel label; 1680 1681 PetscCall(DMPlexGetDepth(dm,&dim)); 1682 PetscCall(DMPlexGetDepthLabel(dm,&label)); 1683 PetscCall(DMLabelGetValue(label,p+pStart,&depth)); 1684 height = dim - depth; 1685 PetscCall(PetscDualSpaceGetHeightSubspace(sp, height, &(sp->pointSpaces[p]))); 1686 PetscCall(PetscObjectReference((PetscObject)sp->pointSpaces[p])); 1687 } 1688 } 1689 } 1690 *bdsp = sp->pointSpaces[point - pStart]; 1691 PetscFunctionReturn(0); 1692 } 1693 1694 /*@C 1695 PetscDualSpaceGetSymmetries - Returns a description of the symmetries of this basis 1696 1697 Not collective 1698 1699 Input Parameter: 1700 . sp - the PetscDualSpace object 1701 1702 Output Parameters: 1703 + perms - Permutations of the interior degrees of freedom, parameterized by the point orientation 1704 - flips - Sign reversal of the interior degrees of freedom, parameterized by the point orientation 1705 1706 Note: The permutation and flip arrays are organized in the following way 1707 $ perms[p][ornt][dof # on point] = new local dof # 1708 $ flips[p][ornt][dof # on point] = reversal or not 1709 1710 Level: developer 1711 1712 @*/ 1713 PetscErrorCode PetscDualSpaceGetSymmetries(PetscDualSpace sp, const PetscInt ****perms, const PetscScalar ****flips) 1714 { 1715 PetscFunctionBegin; 1716 PetscValidHeaderSpecific(sp,PETSCDUALSPACE_CLASSID,1); 1717 if (perms) {PetscValidPointer(perms,2); *perms = NULL;} 1718 if (flips) {PetscValidPointer(flips,3); *flips = NULL;} 1719 if (sp->ops->getsymmetries) PetscCall((sp->ops->getsymmetries)(sp,perms,flips)); 1720 PetscFunctionReturn(0); 1721 } 1722 1723 /*@ 1724 PetscDualSpaceGetFormDegree - Get the form degree k for the k-form the describes the pushforwards/pullbacks of this 1725 dual space's functionals. 1726 1727 Input Parameter: 1728 . dsp - The PetscDualSpace 1729 1730 Output Parameter: 1731 . k - The *signed* degree k of the k. If k >= 0, this means that the degrees of freedom are k-forms, and are stored 1732 in lexicographic order according to the basis of k-forms constructed from the wedge product of 1-forms. So for example, 1733 the 1-form basis in 3-D is (dx, dy, dz), and the 2-form basis in 3-D is (dx wedge dy, dx wedge dz, dy wedge dz). 1734 If k < 0, this means that the degrees transform as k-forms, but are stored as (N-k) forms according to the 1735 Hodge star map. So for example if k = -2 and N = 3, this means that the degrees of freedom transform as 2-forms 1736 but are stored as 1-forms. 1737 1738 Level: developer 1739 1740 .seealso: PetscDTAltV, PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceTransformType 1741 @*/ 1742 PetscErrorCode PetscDualSpaceGetFormDegree(PetscDualSpace dsp, PetscInt *k) 1743 { 1744 PetscFunctionBeginHot; 1745 PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 1746 PetscValidIntPointer(k, 2); 1747 *k = dsp->k; 1748 PetscFunctionReturn(0); 1749 } 1750 1751 /*@ 1752 PetscDualSpaceSetFormDegree - Set the form degree k for the k-form the describes the pushforwards/pullbacks of this 1753 dual space's functionals. 1754 1755 Input Parameters: 1756 + dsp - The PetscDualSpace 1757 - k - The *signed* degree k of the k. If k >= 0, this means that the degrees of freedom are k-forms, and are stored 1758 in lexicographic order according to the basis of k-forms constructed from the wedge product of 1-forms. So for example, 1759 the 1-form basis in 3-D is (dx, dy, dz), and the 2-form basis in 3-D is (dx wedge dy, dx wedge dz, dy wedge dz). 1760 If k < 0, this means that the degrees transform as k-forms, but are stored as (N-k) forms according to the 1761 Hodge star map. So for example if k = -2 and N = 3, this means that the degrees of freedom transform as 2-forms 1762 but are stored as 1-forms. 1763 1764 Level: developer 1765 1766 .seealso: PetscDTAltV, PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceTransformType 1767 @*/ 1768 PetscErrorCode PetscDualSpaceSetFormDegree(PetscDualSpace dsp, PetscInt k) 1769 { 1770 PetscInt dim; 1771 1772 PetscFunctionBeginHot; 1773 PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 1774 PetscCheck(!dsp->setupcalled,PetscObjectComm((PetscObject)dsp), PETSC_ERR_ARG_WRONGSTATE, "Cannot change number of components after dualspace is set up"); 1775 dim = dsp->dm->dim; 1776 PetscCheckFalse(k < -dim || k > dim,PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported %D-form on %D-dimensional reference cell", PetscAbsInt(k), dim); 1777 dsp->k = k; 1778 PetscFunctionReturn(0); 1779 } 1780 1781 /*@ 1782 PetscDualSpaceGetDeRahm - Get the k-simplex associated with the functionals in this dual space 1783 1784 Input Parameter: 1785 . dsp - The PetscDualSpace 1786 1787 Output Parameter: 1788 . k - The simplex dimension 1789 1790 Level: developer 1791 1792 Note: Currently supported values are 1793 $ 0: These are H_1 methods that only transform coordinates 1794 $ 1: These are Hcurl methods that transform functions using the covariant Piola transform (COVARIANT_PIOLA_TRANSFORM) 1795 $ 2: These are the same as 1 1796 $ 3: These are Hdiv methods that transform functions using the contravariant Piola transform (CONTRAVARIANT_PIOLA_TRANSFORM) 1797 1798 .seealso: PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceTransformType 1799 @*/ 1800 PetscErrorCode PetscDualSpaceGetDeRahm(PetscDualSpace dsp, PetscInt *k) 1801 { 1802 PetscInt dim; 1803 1804 PetscFunctionBeginHot; 1805 PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 1806 PetscValidIntPointer(k, 2); 1807 dim = dsp->dm->dim; 1808 if (!dsp->k) *k = IDENTITY_TRANSFORM; 1809 else if (dsp->k == 1) *k = COVARIANT_PIOLA_TRANSFORM; 1810 else if (dsp->k == -(dim - 1)) *k = CONTRAVARIANT_PIOLA_TRANSFORM; 1811 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_PLIB, "Unsupported transformation"); 1812 PetscFunctionReturn(0); 1813 } 1814 1815 /*@C 1816 PetscDualSpaceTransform - Transform the function values 1817 1818 Input Parameters: 1819 + dsp - The PetscDualSpace 1820 . trans - The type of transform 1821 . isInverse - Flag to invert the transform 1822 . fegeom - The cell geometry 1823 . Nv - The number of function samples 1824 . Nc - The number of function components 1825 - vals - The function values 1826 1827 Output Parameter: 1828 . vals - The transformed function values 1829 1830 Level: intermediate 1831 1832 Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 1833 1834 .seealso: PetscDualSpaceTransformGradient(), PetscDualSpaceTransformHessian(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransformType 1835 @*/ 1836 PetscErrorCode PetscDualSpaceTransform(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[]) 1837 { 1838 PetscReal Jstar[9] = {0}; 1839 PetscInt dim, v, c, Nk; 1840 1841 PetscFunctionBeginHot; 1842 PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 1843 PetscValidPointer(fegeom, 4); 1844 PetscValidScalarPointer(vals, 7); 1845 /* TODO: not handling dimEmbed != dim right now */ 1846 dim = dsp->dm->dim; 1847 /* No change needed for 0-forms */ 1848 if (!dsp->k) PetscFunctionReturn(0); 1849 PetscCall(PetscDTBinomialInt(dim, PetscAbsInt(dsp->k), &Nk)); 1850 /* TODO: use fegeom->isAffine */ 1851 PetscCall(PetscDTAltVPullbackMatrix(dim, dim, isInverse ? fegeom->J : fegeom->invJ, dsp->k, Jstar)); 1852 for (v = 0; v < Nv; ++v) { 1853 switch (Nk) { 1854 case 1: 1855 for (c = 0; c < Nc; c++) vals[v*Nc + c] *= Jstar[0]; 1856 break; 1857 case 2: 1858 for (c = 0; c < Nc; c += 2) DMPlex_Mult2DReal_Internal(Jstar, 1, &vals[v*Nc + c], &vals[v*Nc + c]); 1859 break; 1860 case 3: 1861 for (c = 0; c < Nc; c += 3) DMPlex_Mult3DReal_Internal(Jstar, 1, &vals[v*Nc + c], &vals[v*Nc + c]); 1862 break; 1863 default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported form size %D for transformation", Nk); 1864 } 1865 } 1866 PetscFunctionReturn(0); 1867 } 1868 1869 /*@C 1870 PetscDualSpaceTransformGradient - Transform the function gradient values 1871 1872 Input Parameters: 1873 + dsp - The PetscDualSpace 1874 . trans - The type of transform 1875 . isInverse - Flag to invert the transform 1876 . fegeom - The cell geometry 1877 . Nv - The number of function gradient samples 1878 . Nc - The number of function components 1879 - vals - The function gradient values 1880 1881 Output Parameter: 1882 . vals - The transformed function gradient values 1883 1884 Level: intermediate 1885 1886 Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 1887 1888 .seealso: PetscDualSpaceTransform(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransformType 1889 @*/ 1890 PetscErrorCode PetscDualSpaceTransformGradient(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[]) 1891 { 1892 const PetscInt dim = dsp->dm->dim, dE = fegeom->dimEmbed; 1893 PetscInt v, c, d; 1894 1895 PetscFunctionBeginHot; 1896 PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 1897 PetscValidPointer(fegeom, 4); 1898 PetscValidScalarPointer(vals, 7); 1899 #ifdef PETSC_USE_DEBUG 1900 PetscCheckFalse(dE <= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid embedding dimension %D", dE); 1901 #endif 1902 /* Transform gradient */ 1903 if (dim == dE) { 1904 for (v = 0; v < Nv; ++v) { 1905 for (c = 0; c < Nc; ++c) { 1906 switch (dim) 1907 { 1908 case 1: vals[(v*Nc+c)*dim] *= fegeom->invJ[0];break; 1909 case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break; 1910 case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, 1, &vals[(v*Nc+c)*dim], &vals[(v*Nc+c)*dim]);break; 1911 default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 1912 } 1913 } 1914 } 1915 } else { 1916 for (v = 0; v < Nv; ++v) { 1917 for (c = 0; c < Nc; ++c) { 1918 DMPlex_MultTransposeReal_Internal(fegeom->invJ, dim, dE, 1, &vals[(v*Nc+c)*dE], &vals[(v*Nc+c)*dE]); 1919 } 1920 } 1921 } 1922 /* Assume its a vector, otherwise assume its a bunch of scalars */ 1923 if (Nc == 1 || Nc != dim) PetscFunctionReturn(0); 1924 switch (trans) { 1925 case IDENTITY_TRANSFORM: break; 1926 case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */ 1927 if (isInverse) { 1928 for (v = 0; v < Nv; ++v) { 1929 for (d = 0; d < dim; ++d) { 1930 switch (dim) 1931 { 1932 case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 1933 case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 1934 default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 1935 } 1936 } 1937 } 1938 } else { 1939 for (v = 0; v < Nv; ++v) { 1940 for (d = 0; d < dim; ++d) { 1941 switch (dim) 1942 { 1943 case 2: DMPlex_MultTranspose2DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 1944 case 3: DMPlex_MultTranspose3DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 1945 default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 1946 } 1947 } 1948 } 1949 } 1950 break; 1951 case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */ 1952 if (isInverse) { 1953 for (v = 0; v < Nv; ++v) { 1954 for (d = 0; d < dim; ++d) { 1955 switch (dim) 1956 { 1957 case 2: DMPlex_Mult2DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 1958 case 3: DMPlex_Mult3DReal_Internal(fegeom->invJ, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 1959 default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 1960 } 1961 for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] *= fegeom->detJ[0]; 1962 } 1963 } 1964 } else { 1965 for (v = 0; v < Nv; ++v) { 1966 for (d = 0; d < dim; ++d) { 1967 switch (dim) 1968 { 1969 case 2: DMPlex_Mult2DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 1970 case 3: DMPlex_Mult3DReal_Internal(fegeom->J, dim, &vals[v*Nc*dim+d], &vals[v*Nc*dim+d]);break; 1971 default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 1972 } 1973 for (c = 0; c < Nc; ++c) vals[(v*Nc+c)*dim+d] /= fegeom->detJ[0]; 1974 } 1975 } 1976 } 1977 break; 1978 } 1979 PetscFunctionReturn(0); 1980 } 1981 1982 /*@C 1983 PetscDualSpaceTransformHessian - Transform the function Hessian values 1984 1985 Input Parameters: 1986 + dsp - The PetscDualSpace 1987 . trans - The type of transform 1988 . isInverse - Flag to invert the transform 1989 . fegeom - The cell geometry 1990 . Nv - The number of function Hessian samples 1991 . Nc - The number of function components 1992 - vals - The function gradient values 1993 1994 Output Parameter: 1995 . vals - The transformed function Hessian values 1996 1997 Level: intermediate 1998 1999 Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 2000 2001 .seealso: PetscDualSpaceTransform(), PetscDualSpacePullback(), PetscDualSpacePushforward(), PetscDualSpaceTransformType 2002 @*/ 2003 PetscErrorCode PetscDualSpaceTransformHessian(PetscDualSpace dsp, PetscDualSpaceTransformType trans, PetscBool isInverse, PetscFEGeom *fegeom, PetscInt Nv, PetscInt Nc, PetscScalar vals[]) 2004 { 2005 const PetscInt dim = dsp->dm->dim, dE = fegeom->dimEmbed; 2006 PetscInt v, c; 2007 2008 PetscFunctionBeginHot; 2009 PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 2010 PetscValidPointer(fegeom, 4); 2011 PetscValidScalarPointer(vals, 7); 2012 #ifdef PETSC_USE_DEBUG 2013 PetscCheckFalse(dE <= 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Invalid embedding dimension %D", dE); 2014 #endif 2015 /* Transform Hessian: J^{-T}_{ik} J^{-T}_{jl} H(f)_{kl} = J^{-T}_{ik} H(f)_{kl} J^{-1}_{lj} */ 2016 if (dim == dE) { 2017 for (v = 0; v < Nv; ++v) { 2018 for (c = 0; c < Nc; ++c) { 2019 switch (dim) 2020 { 2021 case 1: vals[(v*Nc+c)*dim*dim] *= PetscSqr(fegeom->invJ[0]);break; 2022 case 2: DMPlex_PTAP2DReal_Internal(fegeom->invJ, &vals[(v*Nc+c)*dim*dim], &vals[(v*Nc+c)*dim*dim]);break; 2023 case 3: DMPlex_PTAP3DReal_Internal(fegeom->invJ, &vals[(v*Nc+c)*dim*dim], &vals[(v*Nc+c)*dim*dim]);break; 2024 default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported dim %D for transformation", dim); 2025 } 2026 } 2027 } 2028 } else { 2029 for (v = 0; v < Nv; ++v) { 2030 for (c = 0; c < Nc; ++c) { 2031 DMPlex_PTAPReal_Internal(fegeom->invJ, dim, dE, &vals[(v*Nc+c)*dE*dE], &vals[(v*Nc+c)*dE*dE]); 2032 } 2033 } 2034 } 2035 /* Assume its a vector, otherwise assume its a bunch of scalars */ 2036 if (Nc == 1 || Nc != dim) PetscFunctionReturn(0); 2037 switch (trans) { 2038 case IDENTITY_TRANSFORM: break; 2039 case COVARIANT_PIOLA_TRANSFORM: /* Covariant Piola mapping $\sigma^*(F) = J^{-T} F \circ \phi^{-1)$ */ 2040 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Piola mapping for Hessians not yet supported"); 2041 case CONTRAVARIANT_PIOLA_TRANSFORM: /* Contravariant Piola mapping $\sigma^*(F) = \frac{1}{|\det J|} J F \circ \phi^{-1}$ */ 2042 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Piola mapping for Hessians not yet supported"); 2043 } 2044 PetscFunctionReturn(0); 2045 } 2046 2047 /*@C 2048 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. 2049 2050 Input Parameters: 2051 + dsp - The PetscDualSpace 2052 . fegeom - The geometry for this cell 2053 . Nq - The number of function samples 2054 . Nc - The number of function components 2055 - pointEval - The function values 2056 2057 Output Parameter: 2058 . pointEval - The transformed function values 2059 2060 Level: advanced 2061 2062 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. 2063 2064 Note: This only handles tranformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 2065 2066 .seealso: PetscDualSpacePushforward(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm() 2067 @*/ 2068 PetscErrorCode PetscDualSpacePullback(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[]) 2069 { 2070 PetscDualSpaceTransformType trans; 2071 PetscInt k; 2072 2073 PetscFunctionBeginHot; 2074 PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 2075 PetscValidPointer(fegeom, 2); 2076 PetscValidScalarPointer(pointEval, 5); 2077 /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k. 2078 This determines their transformation properties. */ 2079 PetscCall(PetscDualSpaceGetDeRahm(dsp, &k)); 2080 switch (k) 2081 { 2082 case 0: /* H^1 point evaluations */ 2083 trans = IDENTITY_TRANSFORM;break; 2084 case 1: /* Hcurl preserves tangential edge traces */ 2085 trans = COVARIANT_PIOLA_TRANSFORM;break; 2086 case 2: 2087 case 3: /* Hdiv preserve normal traces */ 2088 trans = CONTRAVARIANT_PIOLA_TRANSFORM;break; 2089 default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", k); 2090 } 2091 PetscCall(PetscDualSpaceTransform(dsp, trans, PETSC_TRUE, fegeom, Nq, Nc, pointEval)); 2092 PetscFunctionReturn(0); 2093 } 2094 2095 /*@C 2096 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. 2097 2098 Input Parameters: 2099 + dsp - The PetscDualSpace 2100 . fegeom - The geometry for this cell 2101 . Nq - The number of function samples 2102 . Nc - The number of function components 2103 - pointEval - The function values 2104 2105 Output Parameter: 2106 . pointEval - The transformed function values 2107 2108 Level: advanced 2109 2110 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. 2111 2112 Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 2113 2114 .seealso: PetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm() 2115 @*/ 2116 PetscErrorCode PetscDualSpacePushforward(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[]) 2117 { 2118 PetscDualSpaceTransformType trans; 2119 PetscInt k; 2120 2121 PetscFunctionBeginHot; 2122 PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 2123 PetscValidPointer(fegeom, 2); 2124 PetscValidScalarPointer(pointEval, 5); 2125 /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k. 2126 This determines their transformation properties. */ 2127 PetscCall(PetscDualSpaceGetDeRahm(dsp, &k)); 2128 switch (k) 2129 { 2130 case 0: /* H^1 point evaluations */ 2131 trans = IDENTITY_TRANSFORM;break; 2132 case 1: /* Hcurl preserves tangential edge traces */ 2133 trans = COVARIANT_PIOLA_TRANSFORM;break; 2134 case 2: 2135 case 3: /* Hdiv preserve normal traces */ 2136 trans = CONTRAVARIANT_PIOLA_TRANSFORM;break; 2137 default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", k); 2138 } 2139 PetscCall(PetscDualSpaceTransform(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval)); 2140 PetscFunctionReturn(0); 2141 } 2142 2143 /*@C 2144 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. 2145 2146 Input Parameters: 2147 + dsp - The PetscDualSpace 2148 . fegeom - The geometry for this cell 2149 . Nq - The number of function gradient samples 2150 . Nc - The number of function components 2151 - pointEval - The function gradient values 2152 2153 Output Parameter: 2154 . pointEval - The transformed function gradient values 2155 2156 Level: advanced 2157 2158 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. 2159 2160 Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 2161 2162 .seealso: PetscDualSpacePushforward(), PPetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm() 2163 @*/ 2164 PetscErrorCode PetscDualSpacePushforwardGradient(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[]) 2165 { 2166 PetscDualSpaceTransformType trans; 2167 PetscInt k; 2168 2169 PetscFunctionBeginHot; 2170 PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 2171 PetscValidPointer(fegeom, 2); 2172 PetscValidScalarPointer(pointEval, 5); 2173 /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k. 2174 This determines their transformation properties. */ 2175 PetscCall(PetscDualSpaceGetDeRahm(dsp, &k)); 2176 switch (k) 2177 { 2178 case 0: /* H^1 point evaluations */ 2179 trans = IDENTITY_TRANSFORM;break; 2180 case 1: /* Hcurl preserves tangential edge traces */ 2181 trans = COVARIANT_PIOLA_TRANSFORM;break; 2182 case 2: 2183 case 3: /* Hdiv preserve normal traces */ 2184 trans = CONTRAVARIANT_PIOLA_TRANSFORM;break; 2185 default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", k); 2186 } 2187 PetscCall(PetscDualSpaceTransformGradient(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval)); 2188 PetscFunctionReturn(0); 2189 } 2190 2191 /*@C 2192 PetscDualSpacePushforwardHessian - Transform the given function Hessian 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. 2193 2194 Input Parameters: 2195 + dsp - The PetscDualSpace 2196 . fegeom - The geometry for this cell 2197 . Nq - The number of function Hessian samples 2198 . Nc - The number of function components 2199 - pointEval - The function gradient values 2200 2201 Output Parameter: 2202 . pointEval - The transformed function Hessian values 2203 2204 Level: advanced 2205 2206 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. 2207 2208 Note: This only handles transformations when the embedding dimension of the geometry in fegeom is the same as the reference dimension. 2209 2210 .seealso: PetscDualSpacePushforward(), PPetscDualSpacePullback(), PetscDualSpaceTransform(), PetscDualSpaceGetDeRahm() 2211 @*/ 2212 PetscErrorCode PetscDualSpacePushforwardHessian(PetscDualSpace dsp, PetscFEGeom *fegeom, PetscInt Nq, PetscInt Nc, PetscScalar pointEval[]) 2213 { 2214 PetscDualSpaceTransformType trans; 2215 PetscInt k; 2216 2217 PetscFunctionBeginHot; 2218 PetscValidHeaderSpecific(dsp, PETSCDUALSPACE_CLASSID, 1); 2219 PetscValidPointer(fegeom, 2); 2220 PetscValidScalarPointer(pointEval, 5); 2221 /* The dualspace dofs correspond to some simplex in the DeRahm complex, which we label by k. 2222 This determines their transformation properties. */ 2223 PetscCall(PetscDualSpaceGetDeRahm(dsp, &k)); 2224 switch (k) 2225 { 2226 case 0: /* H^1 point evaluations */ 2227 trans = IDENTITY_TRANSFORM;break; 2228 case 1: /* Hcurl preserves tangential edge traces */ 2229 trans = COVARIANT_PIOLA_TRANSFORM;break; 2230 case 2: 2231 case 3: /* Hdiv preserve normal traces */ 2232 trans = CONTRAVARIANT_PIOLA_TRANSFORM;break; 2233 default: SETERRQ(PetscObjectComm((PetscObject) dsp), PETSC_ERR_ARG_OUTOFRANGE, "Unsupported simplex dim %D for transformation", k); 2234 } 2235 PetscCall(PetscDualSpaceTransformHessian(dsp, trans, PETSC_FALSE, fegeom, Nq, Nc, pointEval)); 2236 PetscFunctionReturn(0); 2237 } 2238