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