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