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