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