1 /* 2 This file contains routines for basic section object implementation. 3 */ 4 5 #include <petsc/private/sectionimpl.h> /*I "petscsection.h" I*/ 6 #include <petscsf.h> 7 8 PetscClassId PETSC_SECTION_CLASSID; 9 10 /*@ 11 PetscSectionCreate - Allocates PetscSection space and sets the map contents to the default. 12 13 Collective 14 15 Input Parameters: 16 + comm - the MPI communicator 17 - s - pointer to the section 18 19 Level: beginner 20 21 Notes: 22 Typical calling sequence 23 $ PetscSectionCreate(MPI_Comm,PetscSection *); 24 $ PetscSectionSetNumFields(PetscSection, numFields); 25 $ PetscSectionSetChart(PetscSection,low,high); 26 $ PetscSectionSetDof(PetscSection,point,numdof); 27 $ PetscSectionSetUp(PetscSection); 28 $ PetscSectionGetOffset(PetscSection,point,PetscInt *); 29 $ PetscSectionDestroy(PetscSection); 30 31 The PetscSection object and methods are intended to be used in the PETSc Vec and Mat implementions; it is 32 recommended they not be used in user codes unless you really gain something in their use. 33 34 .seealso: PetscSection, PetscSectionDestroy() 35 @*/ 36 PetscErrorCode PetscSectionCreate(MPI_Comm comm, PetscSection *s) 37 { 38 PetscFunctionBegin; 39 PetscValidPointer(s,2); 40 PetscCall(ISInitializePackage()); 41 42 PetscCall(PetscHeaderCreate(*s,PETSC_SECTION_CLASSID,"PetscSection","Section","IS",comm,PetscSectionDestroy,PetscSectionView)); 43 44 (*s)->pStart = -1; 45 (*s)->pEnd = -1; 46 (*s)->perm = NULL; 47 (*s)->pointMajor = PETSC_TRUE; 48 (*s)->includesConstraints = PETSC_TRUE; 49 (*s)->atlasDof = NULL; 50 (*s)->atlasOff = NULL; 51 (*s)->bc = NULL; 52 (*s)->bcIndices = NULL; 53 (*s)->setup = PETSC_FALSE; 54 (*s)->numFields = 0; 55 (*s)->fieldNames = NULL; 56 (*s)->field = NULL; 57 (*s)->useFieldOff = PETSC_FALSE; 58 (*s)->compNames = NULL; 59 (*s)->clObj = NULL; 60 (*s)->clHash = NULL; 61 (*s)->clSection = NULL; 62 (*s)->clPoints = NULL; 63 PetscCall(PetscSectionInvalidateMaxDof_Internal(*s)); 64 PetscFunctionReturn(0); 65 } 66 67 /*@ 68 PetscSectionCopy - Creates a shallow (if possible) copy of the PetscSection 69 70 Collective 71 72 Input Parameter: 73 . section - the PetscSection 74 75 Output Parameter: 76 . newSection - the copy 77 78 Level: intermediate 79 80 .seealso: PetscSection, PetscSectionCreate(), PetscSectionDestroy() 81 @*/ 82 PetscErrorCode PetscSectionCopy(PetscSection section, PetscSection newSection) 83 { 84 PetscSectionSym sym; 85 IS perm; 86 PetscInt numFields, f, c, pStart, pEnd, p; 87 88 PetscFunctionBegin; 89 PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1); 90 PetscValidHeaderSpecific(newSection, PETSC_SECTION_CLASSID, 2); 91 PetscCall(PetscSectionReset(newSection)); 92 PetscCall(PetscSectionGetNumFields(section, &numFields)); 93 if (numFields) PetscCall(PetscSectionSetNumFields(newSection, numFields)); 94 for (f = 0; f < numFields; ++f) { 95 const char *fieldName = NULL, *compName = NULL; 96 PetscInt numComp = 0; 97 98 PetscCall(PetscSectionGetFieldName(section, f, &fieldName)); 99 PetscCall(PetscSectionSetFieldName(newSection, f, fieldName)); 100 PetscCall(PetscSectionGetFieldComponents(section, f, &numComp)); 101 PetscCall(PetscSectionSetFieldComponents(newSection, f, numComp)); 102 for (c = 0; c < numComp; ++c) { 103 PetscCall(PetscSectionGetComponentName(section, f, c, &compName)); 104 PetscCall(PetscSectionSetComponentName(newSection, f, c, compName)); 105 } 106 PetscCall(PetscSectionGetFieldSym(section, f, &sym)); 107 PetscCall(PetscSectionSetFieldSym(newSection, f, sym)); 108 } 109 PetscCall(PetscSectionGetChart(section, &pStart, &pEnd)); 110 PetscCall(PetscSectionSetChart(newSection, pStart, pEnd)); 111 PetscCall(PetscSectionGetPermutation(section, &perm)); 112 PetscCall(PetscSectionSetPermutation(newSection, perm)); 113 PetscCall(PetscSectionGetSym(section, &sym)); 114 PetscCall(PetscSectionSetSym(newSection, sym)); 115 for (p = pStart; p < pEnd; ++p) { 116 PetscInt dof, cdof, fcdof = 0; 117 118 PetscCall(PetscSectionGetDof(section, p, &dof)); 119 PetscCall(PetscSectionSetDof(newSection, p, dof)); 120 PetscCall(PetscSectionGetConstraintDof(section, p, &cdof)); 121 if (cdof) PetscCall(PetscSectionSetConstraintDof(newSection, p, cdof)); 122 for (f = 0; f < numFields; ++f) { 123 PetscCall(PetscSectionGetFieldDof(section, p, f, &dof)); 124 PetscCall(PetscSectionSetFieldDof(newSection, p, f, dof)); 125 if (cdof) { 126 PetscCall(PetscSectionGetFieldConstraintDof(section, p, f, &fcdof)); 127 if (fcdof) PetscCall(PetscSectionSetFieldConstraintDof(newSection, p, f, fcdof)); 128 } 129 } 130 } 131 PetscCall(PetscSectionSetUp(newSection)); 132 for (p = pStart; p < pEnd; ++p) { 133 PetscInt off, cdof, fcdof = 0; 134 const PetscInt *cInd; 135 136 /* Must set offsets in case they do not agree with the prefix sums */ 137 PetscCall(PetscSectionGetOffset(section, p, &off)); 138 PetscCall(PetscSectionSetOffset(newSection, p, off)); 139 PetscCall(PetscSectionGetConstraintDof(section, p, &cdof)); 140 if (cdof) { 141 PetscCall(PetscSectionGetConstraintIndices(section, p, &cInd)); 142 PetscCall(PetscSectionSetConstraintIndices(newSection, p, cInd)); 143 for (f = 0; f < numFields; ++f) { 144 PetscCall(PetscSectionGetFieldConstraintDof(section, p, f, &fcdof)); 145 if (fcdof) { 146 PetscCall(PetscSectionGetFieldConstraintIndices(section, p, f, &cInd)); 147 PetscCall(PetscSectionSetFieldConstraintIndices(newSection, p, f, cInd)); 148 } 149 } 150 } 151 } 152 PetscFunctionReturn(0); 153 } 154 155 /*@ 156 PetscSectionClone - Creates a shallow (if possible) copy of the PetscSection 157 158 Collective 159 160 Input Parameter: 161 . section - the PetscSection 162 163 Output Parameter: 164 . newSection - the copy 165 166 Level: beginner 167 168 .seealso: PetscSection, PetscSectionCreate(), PetscSectionDestroy() 169 @*/ 170 PetscErrorCode PetscSectionClone(PetscSection section, PetscSection *newSection) 171 { 172 PetscFunctionBegin; 173 PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1); 174 PetscValidPointer(newSection, 2); 175 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) section), newSection)); 176 PetscCall(PetscSectionCopy(section, *newSection)); 177 PetscFunctionReturn(0); 178 } 179 180 /*@ 181 PetscSectionSetFromOptions - sets parameters in a PetscSection from the options database 182 183 Collective 184 185 Input Parameter: 186 . section - the PetscSection 187 188 Options Database: 189 . -petscsection_point_major - PETSC_TRUE for point-major order 190 191 Level: intermediate 192 193 .seealso: PetscSection, PetscSectionCreate(), PetscSectionDestroy() 194 @*/ 195 PetscErrorCode PetscSectionSetFromOptions(PetscSection s) 196 { 197 PetscErrorCode ierr; 198 199 PetscFunctionBegin; 200 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 201 ierr = PetscObjectOptionsBegin((PetscObject) s);PetscCall(ierr); 202 PetscCall(PetscOptionsBool("-petscsection_point_major", "The for ordering, either point major or field major", "PetscSectionSetPointMajor", s->pointMajor, &s->pointMajor, NULL)); 203 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 204 PetscCall(PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) s)); 205 ierr = PetscOptionsEnd();PetscCall(ierr); 206 PetscCall(PetscObjectViewFromOptions((PetscObject) s, NULL, "-petscsection_view")); 207 PetscFunctionReturn(0); 208 } 209 210 /*@ 211 PetscSectionCompare - Compares two sections 212 213 Collective 214 215 Input Parameters: 216 + s1 - the first PetscSection 217 - s2 - the second PetscSection 218 219 Output Parameter: 220 . congruent - PETSC_TRUE if the two sections are congruent, PETSC_FALSE otherwise 221 222 Level: intermediate 223 224 Notes: 225 Field names are disregarded. 226 227 .seealso: PetscSection, PetscSectionCreate(), PetscSectionCopy(), PetscSectionClone() 228 @*/ 229 PetscErrorCode PetscSectionCompare(PetscSection s1, PetscSection s2, PetscBool *congruent) 230 { 231 PetscInt pStart, pEnd, nfields, ncdof, nfcdof, p, f, n1, n2; 232 const PetscInt *idx1, *idx2; 233 IS perm1, perm2; 234 PetscBool flg; 235 PetscMPIInt mflg; 236 237 PetscFunctionBegin; 238 PetscValidHeaderSpecific(s1, PETSC_SECTION_CLASSID, 1); 239 PetscValidHeaderSpecific(s2, PETSC_SECTION_CLASSID, 2); 240 PetscValidBoolPointer(congruent,3); 241 flg = PETSC_FALSE; 242 243 PetscCallMPI(MPI_Comm_compare(PetscObjectComm((PetscObject)s1),PetscObjectComm((PetscObject)s2),&mflg)); 244 if (mflg != MPI_CONGRUENT && mflg != MPI_IDENT) { 245 *congruent = PETSC_FALSE; 246 PetscFunctionReturn(0); 247 } 248 249 PetscCall(PetscSectionGetChart(s1, &pStart, &pEnd)); 250 PetscCall(PetscSectionGetChart(s2, &n1, &n2)); 251 if (pStart != n1 || pEnd != n2) goto not_congruent; 252 253 PetscCall(PetscSectionGetPermutation(s1, &perm1)); 254 PetscCall(PetscSectionGetPermutation(s2, &perm2)); 255 if (perm1 && perm2) { 256 PetscCall(ISEqual(perm1, perm2, congruent)); 257 if (!(*congruent)) goto not_congruent; 258 } else if (perm1 != perm2) goto not_congruent; 259 260 for (p = pStart; p < pEnd; ++p) { 261 PetscCall(PetscSectionGetOffset(s1, p, &n1)); 262 PetscCall(PetscSectionGetOffset(s2, p, &n2)); 263 if (n1 != n2) goto not_congruent; 264 265 PetscCall(PetscSectionGetDof(s1, p, &n1)); 266 PetscCall(PetscSectionGetDof(s2, p, &n2)); 267 if (n1 != n2) goto not_congruent; 268 269 PetscCall(PetscSectionGetConstraintDof(s1, p, &ncdof)); 270 PetscCall(PetscSectionGetConstraintDof(s2, p, &n2)); 271 if (ncdof != n2) goto not_congruent; 272 273 PetscCall(PetscSectionGetConstraintIndices(s1, p, &idx1)); 274 PetscCall(PetscSectionGetConstraintIndices(s2, p, &idx2)); 275 PetscCall(PetscArraycmp(idx1, idx2, ncdof, congruent)); 276 if (!(*congruent)) goto not_congruent; 277 } 278 279 PetscCall(PetscSectionGetNumFields(s1, &nfields)); 280 PetscCall(PetscSectionGetNumFields(s2, &n2)); 281 if (nfields != n2) goto not_congruent; 282 283 for (f = 0; f < nfields; ++f) { 284 PetscCall(PetscSectionGetFieldComponents(s1, f, &n1)); 285 PetscCall(PetscSectionGetFieldComponents(s2, f, &n2)); 286 if (n1 != n2) goto not_congruent; 287 288 for (p = pStart; p < pEnd; ++p) { 289 PetscCall(PetscSectionGetFieldOffset(s1, p, f, &n1)); 290 PetscCall(PetscSectionGetFieldOffset(s2, p, f, &n2)); 291 if (n1 != n2) goto not_congruent; 292 293 PetscCall(PetscSectionGetFieldDof(s1, p, f, &n1)); 294 PetscCall(PetscSectionGetFieldDof(s2, p, f, &n2)); 295 if (n1 != n2) goto not_congruent; 296 297 PetscCall(PetscSectionGetFieldConstraintDof(s1, p, f, &nfcdof)); 298 PetscCall(PetscSectionGetFieldConstraintDof(s2, p, f, &n2)); 299 if (nfcdof != n2) goto not_congruent; 300 301 PetscCall(PetscSectionGetFieldConstraintIndices(s1, p, f, &idx1)); 302 PetscCall(PetscSectionGetFieldConstraintIndices(s2, p, f, &idx2)); 303 PetscCall(PetscArraycmp(idx1, idx2, nfcdof, congruent)); 304 if (!(*congruent)) goto not_congruent; 305 } 306 } 307 308 flg = PETSC_TRUE; 309 not_congruent: 310 PetscCall(MPIU_Allreduce(&flg,congruent,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)s1))); 311 PetscFunctionReturn(0); 312 } 313 314 /*@ 315 PetscSectionGetNumFields - Returns the number of fields, or 0 if no fields were defined. 316 317 Not Collective 318 319 Input Parameter: 320 . s - the PetscSection 321 322 Output Parameter: 323 . numFields - the number of fields defined, or 0 if none were defined 324 325 Level: intermediate 326 327 .seealso: PetscSectionSetNumFields() 328 @*/ 329 PetscErrorCode PetscSectionGetNumFields(PetscSection s, PetscInt *numFields) 330 { 331 PetscFunctionBegin; 332 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 333 PetscValidIntPointer(numFields,2); 334 *numFields = s->numFields; 335 PetscFunctionReturn(0); 336 } 337 338 /*@ 339 PetscSectionSetNumFields - Sets the number of fields. 340 341 Not Collective 342 343 Input Parameters: 344 + s - the PetscSection 345 - numFields - the number of fields 346 347 Level: intermediate 348 349 .seealso: PetscSectionGetNumFields() 350 @*/ 351 PetscErrorCode PetscSectionSetNumFields(PetscSection s, PetscInt numFields) 352 { 353 PetscInt f; 354 355 PetscFunctionBegin; 356 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 357 PetscCheck(numFields > 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "The number of fields %" PetscInt_FMT " must be positive", numFields); 358 PetscCall(PetscSectionReset(s)); 359 360 s->numFields = numFields; 361 PetscCall(PetscMalloc1(s->numFields, &s->numFieldComponents)); 362 PetscCall(PetscMalloc1(s->numFields, &s->fieldNames)); 363 PetscCall(PetscMalloc1(s->numFields, &s->compNames)); 364 PetscCall(PetscMalloc1(s->numFields, &s->field)); 365 for (f = 0; f < s->numFields; ++f) { 366 char name[64]; 367 368 s->numFieldComponents[f] = 1; 369 370 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) s), &s->field[f])); 371 PetscCall(PetscSNPrintf(name, 64, "Field_%" PetscInt_FMT, f)); 372 PetscCall(PetscStrallocpy(name, (char **) &s->fieldNames[f])); 373 PetscCall(PetscSNPrintf(name, 64, "Component_0")); 374 PetscCall(PetscMalloc1(s->numFieldComponents[f], &s->compNames[f])); 375 PetscCall(PetscStrallocpy(name, (char **) &s->compNames[f][0])); 376 } 377 PetscFunctionReturn(0); 378 } 379 380 /*@C 381 PetscSectionGetFieldName - Returns the name of a field in the PetscSection 382 383 Not Collective 384 385 Input Parameters: 386 + s - the PetscSection 387 - field - the field number 388 389 Output Parameter: 390 . fieldName - the field name 391 392 Level: intermediate 393 394 .seealso: PetscSectionSetFieldName() 395 @*/ 396 PetscErrorCode PetscSectionGetFieldName(PetscSection s, PetscInt field, const char *fieldName[]) 397 { 398 PetscFunctionBegin; 399 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 400 PetscValidPointer(fieldName, 3); 401 PetscSectionCheckValidField(field,s->numFields); 402 *fieldName = s->fieldNames[field]; 403 PetscFunctionReturn(0); 404 } 405 406 /*@C 407 PetscSectionSetFieldName - Sets the name of a field in the PetscSection 408 409 Not Collective 410 411 Input Parameters: 412 + s - the PetscSection 413 . field - the field number 414 - fieldName - the field name 415 416 Level: intermediate 417 418 .seealso: PetscSectionGetFieldName() 419 @*/ 420 PetscErrorCode PetscSectionSetFieldName(PetscSection s, PetscInt field, const char fieldName[]) 421 { 422 PetscFunctionBegin; 423 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 424 if (fieldName) PetscValidCharPointer(fieldName, 3); 425 PetscSectionCheckValidField(field,s->numFields); 426 PetscCall(PetscFree(s->fieldNames[field])); 427 PetscCall(PetscStrallocpy(fieldName, (char**) &s->fieldNames[field])); 428 PetscFunctionReturn(0); 429 } 430 431 /*@C 432 PetscSectionGetComponentName - Gets the name of a field component in the PetscSection 433 434 Not Collective 435 436 Input Parameters: 437 + s - the PetscSection 438 . field - the field number 439 . comp - the component number 440 - compName - the component name 441 442 Level: intermediate 443 444 .seealso: PetscSectionSetComponentName() 445 @*/ 446 PetscErrorCode PetscSectionGetComponentName(PetscSection s, PetscInt field, PetscInt comp, const char *compName[]) 447 { 448 PetscFunctionBegin; 449 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 450 PetscValidPointer(compName, 4); 451 PetscSectionCheckValidField(field,s->numFields); 452 PetscSectionCheckValidFieldComponent(comp,s->numFieldComponents[field]); 453 *compName = s->compNames[field][comp]; 454 PetscFunctionReturn(0); 455 } 456 457 /*@C 458 PetscSectionSetComponentName - Sets the name of a field component in the PetscSection 459 460 Not Collective 461 462 Input Parameters: 463 + s - the PetscSection 464 . field - the field number 465 . comp - the component number 466 - compName - the component name 467 468 Level: intermediate 469 470 .seealso: PetscSectionGetComponentName() 471 @*/ 472 PetscErrorCode PetscSectionSetComponentName(PetscSection s, PetscInt field, PetscInt comp, const char compName[]) 473 { 474 PetscFunctionBegin; 475 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 476 if (compName) PetscValidCharPointer(compName, 4); 477 PetscSectionCheckValidField(field,s->numFields); 478 PetscSectionCheckValidFieldComponent(comp,s->numFieldComponents[field]); 479 PetscCall(PetscFree(s->compNames[field][comp])); 480 PetscCall(PetscStrallocpy(compName, (char**) &s->compNames[field][comp])); 481 PetscFunctionReturn(0); 482 } 483 484 /*@ 485 PetscSectionGetFieldComponents - Returns the number of field components for the given field. 486 487 Not Collective 488 489 Input Parameters: 490 + s - the PetscSection 491 - field - the field number 492 493 Output Parameter: 494 . numComp - the number of field components 495 496 Level: intermediate 497 498 .seealso: PetscSectionSetFieldComponents(), PetscSectionGetNumFields() 499 @*/ 500 PetscErrorCode PetscSectionGetFieldComponents(PetscSection s, PetscInt field, PetscInt *numComp) 501 { 502 PetscFunctionBegin; 503 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 504 PetscValidIntPointer(numComp, 3); 505 PetscSectionCheckValidField(field,s->numFields); 506 *numComp = s->numFieldComponents[field]; 507 PetscFunctionReturn(0); 508 } 509 510 /*@ 511 PetscSectionSetFieldComponents - Sets the number of field components for the given field. 512 513 Not Collective 514 515 Input Parameters: 516 + s - the PetscSection 517 . field - the field number 518 - numComp - the number of field components 519 520 Level: intermediate 521 522 .seealso: PetscSectionGetFieldComponents(), PetscSectionGetNumFields() 523 @*/ 524 PetscErrorCode PetscSectionSetFieldComponents(PetscSection s, PetscInt field, PetscInt numComp) 525 { 526 PetscInt c; 527 528 PetscFunctionBegin; 529 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 530 PetscSectionCheckValidField(field,s->numFields); 531 if (s->compNames) { 532 for (c = 0; c < s->numFieldComponents[field]; ++c) { 533 PetscCall(PetscFree(s->compNames[field][c])); 534 } 535 PetscCall(PetscFree(s->compNames[field])); 536 } 537 538 s->numFieldComponents[field] = numComp; 539 if (numComp) { 540 PetscCall(PetscMalloc1(numComp, (char ***) &s->compNames[field])); 541 for (c = 0; c < numComp; ++c) { 542 char name[64]; 543 544 PetscCall(PetscSNPrintf(name, 64, "%" PetscInt_FMT, c)); 545 PetscCall(PetscStrallocpy(name, (char **) &s->compNames[field][c])); 546 } 547 } 548 PetscFunctionReturn(0); 549 } 550 551 /*@ 552 PetscSectionGetChart - Returns the range [pStart, pEnd) in which points lie. 553 554 Not Collective 555 556 Input Parameter: 557 . s - the PetscSection 558 559 Output Parameters: 560 + pStart - the first point 561 - pEnd - one past the last point 562 563 Level: intermediate 564 565 .seealso: PetscSectionSetChart(), PetscSectionCreate() 566 @*/ 567 PetscErrorCode PetscSectionGetChart(PetscSection s, PetscInt *pStart, PetscInt *pEnd) 568 { 569 PetscFunctionBegin; 570 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 571 if (pStart) *pStart = s->pStart; 572 if (pEnd) *pEnd = s->pEnd; 573 PetscFunctionReturn(0); 574 } 575 576 /*@ 577 PetscSectionSetChart - Sets the range [pStart, pEnd) in which points lie. 578 579 Not Collective 580 581 Input Parameters: 582 + s - the PetscSection 583 . pStart - the first point 584 - pEnd - one past the last point 585 586 Level: intermediate 587 588 .seealso: PetscSectionGetChart(), PetscSectionCreate() 589 @*/ 590 PetscErrorCode PetscSectionSetChart(PetscSection s, PetscInt pStart, PetscInt pEnd) 591 { 592 PetscInt f; 593 594 PetscFunctionBegin; 595 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 596 if (pStart == s->pStart && pEnd == s->pEnd) PetscFunctionReturn(0); 597 /* Cannot Reset() because it destroys field information */ 598 s->setup = PETSC_FALSE; 599 PetscCall(PetscSectionDestroy(&s->bc)); 600 PetscCall(PetscFree(s->bcIndices)); 601 PetscCall(PetscFree2(s->atlasDof, s->atlasOff)); 602 603 s->pStart = pStart; 604 s->pEnd = pEnd; 605 PetscCall(PetscMalloc2((pEnd - pStart), &s->atlasDof, (pEnd - pStart), &s->atlasOff)); 606 PetscCall(PetscArrayzero(s->atlasDof, pEnd - pStart)); 607 for (f = 0; f < s->numFields; ++f) { 608 PetscCall(PetscSectionSetChart(s->field[f], pStart, pEnd)); 609 } 610 PetscFunctionReturn(0); 611 } 612 613 /*@ 614 PetscSectionGetPermutation - Returns the permutation of [0, pEnd-pStart) or NULL 615 616 Not Collective 617 618 Input Parameter: 619 . s - the PetscSection 620 621 Output Parameters: 622 . perm - The permutation as an IS 623 624 Level: intermediate 625 626 .seealso: PetscSectionSetPermutation(), PetscSectionCreate() 627 @*/ 628 PetscErrorCode PetscSectionGetPermutation(PetscSection s, IS *perm) 629 { 630 PetscFunctionBegin; 631 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 632 if (perm) {PetscValidPointer(perm, 2); *perm = s->perm;} 633 PetscFunctionReturn(0); 634 } 635 636 /*@ 637 PetscSectionSetPermutation - Sets the permutation for [0, pEnd-pStart) 638 639 Not Collective 640 641 Input Parameters: 642 + s - the PetscSection 643 - perm - the permutation of points 644 645 Level: intermediate 646 647 .seealso: PetscSectionGetPermutation(), PetscSectionCreate() 648 @*/ 649 PetscErrorCode PetscSectionSetPermutation(PetscSection s, IS perm) 650 { 651 PetscFunctionBegin; 652 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 653 if (perm) PetscValidHeaderSpecific(perm, IS_CLASSID, 2); 654 PetscCheck(!s->setup,PetscObjectComm((PetscObject) s), PETSC_ERR_ARG_WRONGSTATE, "Cannot set a permutation after the section is setup"); 655 if (s->perm != perm) { 656 PetscCall(ISDestroy(&s->perm)); 657 if (perm) { 658 s->perm = perm; 659 PetscCall(PetscObjectReference((PetscObject) s->perm)); 660 } 661 } 662 PetscFunctionReturn(0); 663 } 664 665 /*@ 666 PetscSectionGetPointMajor - Returns the flag for dof ordering, true if it is point major, otherwise field major 667 668 Not Collective 669 670 Input Parameter: 671 . s - the PetscSection 672 673 Output Parameter: 674 . pm - the flag for point major ordering 675 676 Level: intermediate 677 678 .seealso: PetscSectionSetPointMajor() 679 @*/ 680 PetscErrorCode PetscSectionGetPointMajor(PetscSection s, PetscBool *pm) 681 { 682 PetscFunctionBegin; 683 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 684 PetscValidBoolPointer(pm,2); 685 *pm = s->pointMajor; 686 PetscFunctionReturn(0); 687 } 688 689 /*@ 690 PetscSectionSetPointMajor - Sets the flag for dof ordering, true if it is point major, otherwise field major 691 692 Not Collective 693 694 Input Parameters: 695 + s - the PetscSection 696 - pm - the flag for point major ordering 697 698 Level: intermediate 699 700 .seealso: PetscSectionGetPointMajor() 701 @*/ 702 PetscErrorCode PetscSectionSetPointMajor(PetscSection s, PetscBool pm) 703 { 704 PetscFunctionBegin; 705 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 706 PetscCheck(!s->setup,PetscObjectComm((PetscObject) s), PETSC_ERR_ARG_WRONGSTATE, "Cannot set the dof ordering after the section is setup"); 707 s->pointMajor = pm; 708 PetscFunctionReturn(0); 709 } 710 711 /*@ 712 PetscSectionGetIncludesConstraints - Returns the flag indicating if constrained dofs were included when computing offsets 713 714 Not Collective 715 716 Input Parameter: 717 . s - the PetscSection 718 719 Output Parameter: 720 . includesConstraints - the flag indicating if constrained dofs were included when computing offsets 721 722 Level: intermediate 723 724 .seealso: PetscSectionSetIncludesConstraints() 725 @*/ 726 PetscErrorCode PetscSectionGetIncludesConstraints(PetscSection s, PetscBool *includesConstraints) 727 { 728 PetscFunctionBegin; 729 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 730 PetscValidBoolPointer(includesConstraints,2); 731 *includesConstraints = s->includesConstraints; 732 PetscFunctionReturn(0); 733 } 734 735 /*@ 736 PetscSectionSetIncludesConstraints - Sets the flag indicating if constrained dofs are to be included when computing offsets 737 738 Not Collective 739 740 Input Parameters: 741 + s - the PetscSection 742 - includesConstraints - the flag indicating if constrained dofs are to be included when computing offsets 743 744 Level: intermediate 745 746 .seealso: PetscSectionGetIncludesConstraints() 747 @*/ 748 PetscErrorCode PetscSectionSetIncludesConstraints(PetscSection s, PetscBool includesConstraints) 749 { 750 PetscFunctionBegin; 751 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 752 PetscCheck(!s->setup,PetscObjectComm((PetscObject) s), PETSC_ERR_ARG_WRONGSTATE, "Cannot set includesConstraints after the section is set up"); 753 s->includesConstraints = includesConstraints; 754 PetscFunctionReturn(0); 755 } 756 757 /*@ 758 PetscSectionGetDof - Return the number of degrees of freedom associated with a given point. 759 760 Not Collective 761 762 Input Parameters: 763 + s - the PetscSection 764 - point - the point 765 766 Output Parameter: 767 . numDof - the number of dof 768 769 Level: intermediate 770 771 .seealso: PetscSectionSetDof(), PetscSectionCreate() 772 @*/ 773 PetscErrorCode PetscSectionGetDof(PetscSection s, PetscInt point, PetscInt *numDof) 774 { 775 PetscFunctionBeginHot; 776 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 777 PetscValidIntPointer(numDof, 3); 778 PetscAssert(point >= s->pStart && point < s->pEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %" PetscInt_FMT " should be in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, s->pStart, s->pEnd); 779 *numDof = s->atlasDof[point - s->pStart]; 780 PetscFunctionReturn(0); 781 } 782 783 /*@ 784 PetscSectionSetDof - Sets the number of degrees of freedom associated with a given point. 785 786 Not Collective 787 788 Input Parameters: 789 + s - the PetscSection 790 . point - the point 791 - numDof - the number of dof 792 793 Level: intermediate 794 795 .seealso: PetscSectionGetDof(), PetscSectionAddDof(), PetscSectionCreate() 796 @*/ 797 PetscErrorCode PetscSectionSetDof(PetscSection s, PetscInt point, PetscInt numDof) 798 { 799 PetscFunctionBegin; 800 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 801 PetscAssert(point >= s->pStart && point < s->pEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %" PetscInt_FMT " should be in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, s->pStart, s->pEnd); 802 s->atlasDof[point - s->pStart] = numDof; 803 PetscCall(PetscSectionInvalidateMaxDof_Internal(s)); 804 PetscFunctionReturn(0); 805 } 806 807 /*@ 808 PetscSectionAddDof - Adds to the number of degrees of freedom associated with a given point. 809 810 Not Collective 811 812 Input Parameters: 813 + s - the PetscSection 814 . point - the point 815 - numDof - the number of additional dof 816 817 Level: intermediate 818 819 .seealso: PetscSectionGetDof(), PetscSectionSetDof(), PetscSectionCreate() 820 @*/ 821 PetscErrorCode PetscSectionAddDof(PetscSection s, PetscInt point, PetscInt numDof) 822 { 823 PetscFunctionBeginHot; 824 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 825 PetscAssert(point >= s->pStart && point < s->pEnd, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %" PetscInt_FMT " should be in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, s->pStart, s->pEnd); 826 s->atlasDof[point - s->pStart] += numDof; 827 PetscCall(PetscSectionInvalidateMaxDof_Internal(s)); 828 PetscFunctionReturn(0); 829 } 830 831 /*@ 832 PetscSectionGetFieldDof - Return the number of degrees of freedom associated with a field on a given point. 833 834 Not Collective 835 836 Input Parameters: 837 + s - the PetscSection 838 . point - the point 839 - field - the field 840 841 Output Parameter: 842 . numDof - the number of dof 843 844 Level: intermediate 845 846 .seealso: PetscSectionSetFieldDof(), PetscSectionCreate() 847 @*/ 848 PetscErrorCode PetscSectionGetFieldDof(PetscSection s, PetscInt point, PetscInt field, PetscInt *numDof) 849 { 850 PetscFunctionBegin; 851 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 852 PetscValidIntPointer(numDof,4); 853 PetscSectionCheckValidField(field,s->numFields); 854 PetscCall(PetscSectionGetDof(s->field[field], point, numDof)); 855 PetscFunctionReturn(0); 856 } 857 858 /*@ 859 PetscSectionSetFieldDof - Sets the number of degrees of freedom associated with a field on a given point. 860 861 Not Collective 862 863 Input Parameters: 864 + s - the PetscSection 865 . point - the point 866 . field - the field 867 - numDof - the number of dof 868 869 Level: intermediate 870 871 .seealso: PetscSectionGetFieldDof(), PetscSectionCreate() 872 @*/ 873 PetscErrorCode PetscSectionSetFieldDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof) 874 { 875 PetscFunctionBegin; 876 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 877 PetscSectionCheckValidField(field,s->numFields); 878 PetscCall(PetscSectionSetDof(s->field[field], point, numDof)); 879 PetscFunctionReturn(0); 880 } 881 882 /*@ 883 PetscSectionAddFieldDof - Adds a number of degrees of freedom associated with a field on a given point. 884 885 Not Collective 886 887 Input Parameters: 888 + s - the PetscSection 889 . point - the point 890 . field - the field 891 - numDof - the number of dof 892 893 Level: intermediate 894 895 .seealso: PetscSectionSetFieldDof(), PetscSectionGetFieldDof(), PetscSectionCreate() 896 @*/ 897 PetscErrorCode PetscSectionAddFieldDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof) 898 { 899 PetscFunctionBegin; 900 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 901 PetscSectionCheckValidField(field,s->numFields); 902 PetscCall(PetscSectionAddDof(s->field[field], point, numDof)); 903 PetscFunctionReturn(0); 904 } 905 906 /*@ 907 PetscSectionGetConstraintDof - Return the number of constrained degrees of freedom associated with a given point. 908 909 Not Collective 910 911 Input Parameters: 912 + s - the PetscSection 913 - point - the point 914 915 Output Parameter: 916 . numDof - the number of dof which are fixed by constraints 917 918 Level: intermediate 919 920 .seealso: PetscSectionGetDof(), PetscSectionSetConstraintDof(), PetscSectionCreate() 921 @*/ 922 PetscErrorCode PetscSectionGetConstraintDof(PetscSection s, PetscInt point, PetscInt *numDof) 923 { 924 PetscFunctionBegin; 925 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 926 PetscValidIntPointer(numDof, 3); 927 if (s->bc) { 928 PetscCall(PetscSectionGetDof(s->bc, point, numDof)); 929 } else *numDof = 0; 930 PetscFunctionReturn(0); 931 } 932 933 /*@ 934 PetscSectionSetConstraintDof - Set the number of constrained degrees of freedom associated with a given point. 935 936 Not Collective 937 938 Input Parameters: 939 + s - the PetscSection 940 . point - the point 941 - numDof - the number of dof which are fixed by constraints 942 943 Level: intermediate 944 945 .seealso: PetscSectionSetDof(), PetscSectionGetConstraintDof(), PetscSectionCreate() 946 @*/ 947 PetscErrorCode PetscSectionSetConstraintDof(PetscSection s, PetscInt point, PetscInt numDof) 948 { 949 PetscFunctionBegin; 950 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 951 if (numDof) { 952 PetscCall(PetscSectionCheckConstraints_Private(s)); 953 PetscCall(PetscSectionSetDof(s->bc, point, numDof)); 954 } 955 PetscFunctionReturn(0); 956 } 957 958 /*@ 959 PetscSectionAddConstraintDof - Increment the number of constrained degrees of freedom associated with a given point. 960 961 Not Collective 962 963 Input Parameters: 964 + s - the PetscSection 965 . point - the point 966 - numDof - the number of additional dof which are fixed by constraints 967 968 Level: intermediate 969 970 .seealso: PetscSectionAddDof(), PetscSectionGetConstraintDof(), PetscSectionCreate() 971 @*/ 972 PetscErrorCode PetscSectionAddConstraintDof(PetscSection s, PetscInt point, PetscInt numDof) 973 { 974 PetscFunctionBegin; 975 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 976 if (numDof) { 977 PetscCall(PetscSectionCheckConstraints_Private(s)); 978 PetscCall(PetscSectionAddDof(s->bc, point, numDof)); 979 } 980 PetscFunctionReturn(0); 981 } 982 983 /*@ 984 PetscSectionGetFieldConstraintDof - Return the number of constrained degrees of freedom associated with a given field on a point. 985 986 Not Collective 987 988 Input Parameters: 989 + s - the PetscSection 990 . point - the point 991 - field - the field 992 993 Output Parameter: 994 . numDof - the number of dof which are fixed by constraints 995 996 Level: intermediate 997 998 .seealso: PetscSectionGetDof(), PetscSectionSetFieldConstraintDof(), PetscSectionCreate() 999 @*/ 1000 PetscErrorCode PetscSectionGetFieldConstraintDof(PetscSection s, PetscInt point, PetscInt field, PetscInt *numDof) 1001 { 1002 PetscFunctionBegin; 1003 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1004 PetscValidIntPointer(numDof, 4); 1005 PetscSectionCheckValidField(field,s->numFields); 1006 PetscCall(PetscSectionGetConstraintDof(s->field[field], point, numDof)); 1007 PetscFunctionReturn(0); 1008 } 1009 1010 /*@ 1011 PetscSectionSetFieldConstraintDof - Set the number of constrained degrees of freedom associated with a given field on a point. 1012 1013 Not Collective 1014 1015 Input Parameters: 1016 + s - the PetscSection 1017 . point - the point 1018 . field - the field 1019 - numDof - the number of dof which are fixed by constraints 1020 1021 Level: intermediate 1022 1023 .seealso: PetscSectionSetDof(), PetscSectionGetFieldConstraintDof(), PetscSectionCreate() 1024 @*/ 1025 PetscErrorCode PetscSectionSetFieldConstraintDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof) 1026 { 1027 PetscFunctionBegin; 1028 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1029 PetscSectionCheckValidField(field,s->numFields); 1030 PetscCall(PetscSectionSetConstraintDof(s->field[field], point, numDof)); 1031 PetscFunctionReturn(0); 1032 } 1033 1034 /*@ 1035 PetscSectionAddFieldConstraintDof - Increment the number of constrained degrees of freedom associated with a given field on a point. 1036 1037 Not Collective 1038 1039 Input Parameters: 1040 + s - the PetscSection 1041 . point - the point 1042 . field - the field 1043 - numDof - the number of additional dof which are fixed by constraints 1044 1045 Level: intermediate 1046 1047 .seealso: PetscSectionAddDof(), PetscSectionGetFieldConstraintDof(), PetscSectionCreate() 1048 @*/ 1049 PetscErrorCode PetscSectionAddFieldConstraintDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof) 1050 { 1051 PetscFunctionBegin; 1052 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1053 PetscSectionCheckValidField(field,s->numFields); 1054 PetscCall(PetscSectionAddConstraintDof(s->field[field], point, numDof)); 1055 PetscFunctionReturn(0); 1056 } 1057 1058 /*@ 1059 PetscSectionSetUpBC - Setup the subsections describing boundary conditions. 1060 1061 Not Collective 1062 1063 Input Parameter: 1064 . s - the PetscSection 1065 1066 Level: advanced 1067 1068 .seealso: PetscSectionSetUp(), PetscSectionCreate() 1069 @*/ 1070 PetscErrorCode PetscSectionSetUpBC(PetscSection s) 1071 { 1072 PetscFunctionBegin; 1073 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1074 if (s->bc) { 1075 const PetscInt last = (s->bc->pEnd-s->bc->pStart) - 1; 1076 1077 PetscCall(PetscSectionSetUp(s->bc)); 1078 PetscCall(PetscMalloc1(last >= 0 ? s->bc->atlasOff[last] + s->bc->atlasDof[last] : 0, &s->bcIndices)); 1079 } 1080 PetscFunctionReturn(0); 1081 } 1082 1083 /*@ 1084 PetscSectionSetUp - Calculate offsets based upon the number of degrees of freedom for each point. 1085 1086 Not Collective 1087 1088 Input Parameter: 1089 . s - the PetscSection 1090 1091 Level: intermediate 1092 1093 .seealso: PetscSectionCreate() 1094 @*/ 1095 PetscErrorCode PetscSectionSetUp(PetscSection s) 1096 { 1097 const PetscInt *pind = NULL; 1098 PetscInt offset = 0, foff, p, f; 1099 1100 PetscFunctionBegin; 1101 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1102 if (s->setup) PetscFunctionReturn(0); 1103 s->setup = PETSC_TRUE; 1104 /* Set offsets and field offsets for all points */ 1105 /* Assume that all fields have the same chart */ 1106 PetscCheck(s->includesConstraints,PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSectionSetUp is currently unsupported for includesConstraints = PETSC_TRUE"); 1107 if (s->perm) PetscCall(ISGetIndices(s->perm, &pind)); 1108 if (s->pointMajor) { 1109 for (p = 0; p < s->pEnd - s->pStart; ++p) { 1110 const PetscInt q = pind ? pind[p] : p; 1111 1112 /* Set point offset */ 1113 s->atlasOff[q] = offset; 1114 offset += s->atlasDof[q]; 1115 /* Set field offset */ 1116 for (f = 0, foff = s->atlasOff[q]; f < s->numFields; ++f) { 1117 PetscSection sf = s->field[f]; 1118 1119 sf->atlasOff[q] = foff; 1120 foff += sf->atlasDof[q]; 1121 } 1122 } 1123 } else { 1124 /* Set field offsets for all points */ 1125 for (f = 0; f < s->numFields; ++f) { 1126 PetscSection sf = s->field[f]; 1127 1128 for (p = 0; p < s->pEnd - s->pStart; ++p) { 1129 const PetscInt q = pind ? pind[p] : p; 1130 1131 sf->atlasOff[q] = offset; 1132 offset += sf->atlasDof[q]; 1133 } 1134 } 1135 /* Disable point offsets since these are unused */ 1136 for (p = 0; p < s->pEnd - s->pStart; ++p) { 1137 s->atlasOff[p] = -1; 1138 } 1139 } 1140 if (s->perm) PetscCall(ISRestoreIndices(s->perm, &pind)); 1141 /* Setup BC sections */ 1142 PetscCall(PetscSectionSetUpBC(s)); 1143 for (f = 0; f < s->numFields; ++f) PetscCall(PetscSectionSetUpBC(s->field[f])); 1144 PetscFunctionReturn(0); 1145 } 1146 1147 /*@ 1148 PetscSectionGetMaxDof - Return the maximum number of degrees of freedom on any point in the chart 1149 1150 Not Collective 1151 1152 Input Parameters: 1153 . s - the PetscSection 1154 1155 Output Parameter: 1156 . maxDof - the maximum dof 1157 1158 Level: intermediate 1159 1160 Notes: 1161 The returned number is up-to-date without need for PetscSectionSetUp(). 1162 1163 Developer Notes: 1164 The returned number is calculated lazily and stashed. 1165 A call to PetscSectionInvalidateMaxDof_Internal() invalidates the stashed value. 1166 PetscSectionInvalidateMaxDof_Internal() is called in PetscSectionSetDof(), PetscSectionAddDof() and PetscSectionReset(). 1167 It should also be called every time atlasDof is modified directly. 1168 1169 .seealso: PetscSectionGetDof(), PetscSectionSetDof(), PetscSectionAddDof(), PetscSectionCreate() 1170 @*/ 1171 PetscErrorCode PetscSectionGetMaxDof(PetscSection s, PetscInt *maxDof) 1172 { 1173 PetscInt p; 1174 1175 PetscFunctionBegin; 1176 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1177 PetscValidIntPointer(maxDof, 2); 1178 if (s->maxDof == PETSC_MIN_INT) { 1179 s->maxDof = 0; 1180 for (p = 0; p < s->pEnd - s->pStart; ++p) { 1181 s->maxDof = PetscMax(s->maxDof, s->atlasDof[p]); 1182 } 1183 } 1184 *maxDof = s->maxDof; 1185 PetscFunctionReturn(0); 1186 } 1187 1188 /*@ 1189 PetscSectionGetStorageSize - Return the size of an array or local Vec capable of holding all the degrees of freedom. 1190 1191 Not Collective 1192 1193 Input Parameter: 1194 . s - the PetscSection 1195 1196 Output Parameter: 1197 . size - the size of an array which can hold all the dofs 1198 1199 Level: intermediate 1200 1201 .seealso: PetscSectionGetOffset(), PetscSectionGetConstrainedStorageSize(), PetscSectionCreate() 1202 @*/ 1203 PetscErrorCode PetscSectionGetStorageSize(PetscSection s, PetscInt *size) 1204 { 1205 PetscInt p, n = 0; 1206 1207 PetscFunctionBegin; 1208 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1209 PetscValidIntPointer(size, 2); 1210 for (p = 0; p < s->pEnd - s->pStart; ++p) n += s->atlasDof[p] > 0 ? s->atlasDof[p] : 0; 1211 *size = n; 1212 PetscFunctionReturn(0); 1213 } 1214 1215 /*@ 1216 PetscSectionGetConstrainedStorageSize - Return the size of an array or local Vec capable of holding all unconstrained degrees of freedom. 1217 1218 Not Collective 1219 1220 Input Parameter: 1221 . s - the PetscSection 1222 1223 Output Parameter: 1224 . size - the size of an array which can hold all unconstrained dofs 1225 1226 Level: intermediate 1227 1228 .seealso: PetscSectionGetStorageSize(), PetscSectionGetOffset(), PetscSectionCreate() 1229 @*/ 1230 PetscErrorCode PetscSectionGetConstrainedStorageSize(PetscSection s, PetscInt *size) 1231 { 1232 PetscInt p, n = 0; 1233 1234 PetscFunctionBegin; 1235 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1236 PetscValidIntPointer(size, 2); 1237 for (p = 0; p < s->pEnd - s->pStart; ++p) { 1238 const PetscInt cdof = s->bc ? s->bc->atlasDof[p] : 0; 1239 n += s->atlasDof[p] > 0 ? s->atlasDof[p] - cdof : 0; 1240 } 1241 *size = n; 1242 PetscFunctionReturn(0); 1243 } 1244 1245 /*@ 1246 PetscSectionCreateGlobalSection - Create a section describing the global field layout using 1247 the local section and an SF describing the section point overlap. 1248 1249 Input Parameters: 1250 + s - The PetscSection for the local field layout 1251 . sf - The SF describing parallel layout of the section points (leaves are unowned local points) 1252 . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs 1253 - localOffsets - If PETSC_TRUE, use local rather than global offsets for the points 1254 1255 Output Parameter: 1256 . gsection - The PetscSection for the global field layout 1257 1258 Note: This gives negative sizes and offsets to points not owned by this process 1259 1260 Level: intermediate 1261 1262 .seealso: PetscSectionCreate() 1263 @*/ 1264 PetscErrorCode PetscSectionCreateGlobalSection(PetscSection s, PetscSF sf, PetscBool includeConstraints, PetscBool localOffsets, PetscSection *gsection) 1265 { 1266 PetscSection gs; 1267 const PetscInt *pind = NULL; 1268 PetscInt *recv = NULL, *neg = NULL; 1269 PetscInt pStart, pEnd, p, dof, cdof, off, foff, globalOff = 0, nroots, nlocal, maxleaf; 1270 PetscInt numFields, f, numComponents; 1271 1272 PetscFunctionBegin; 1273 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1274 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 1275 PetscValidLogicalCollectiveBool(s, includeConstraints, 3); 1276 PetscValidLogicalCollectiveBool(s, localOffsets, 4); 1277 PetscValidPointer(gsection, 5); 1278 PetscCheck(s->pointMajor,PETSC_COMM_SELF,PETSC_ERR_SUP, "No support for field major ordering"); 1279 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) s), &gs)); 1280 PetscCall(PetscSectionGetNumFields(s, &numFields)); 1281 if (numFields > 0) PetscCall(PetscSectionSetNumFields(gs, numFields)); 1282 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 1283 PetscCall(PetscSectionSetChart(gs, pStart, pEnd)); 1284 gs->includesConstraints = includeConstraints; 1285 PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 1286 nlocal = nroots; /* The local/leaf space matches global/root space */ 1287 /* Must allocate for all points visible to SF, which may be more than this section */ 1288 if (nroots >= 0) { /* nroots < 0 means that the graph has not been set, only happens in serial */ 1289 PetscCall(PetscSFGetLeafRange(sf, NULL, &maxleaf)); 1290 PetscCheck(nroots >= pEnd,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "SF roots %" PetscInt_FMT " < pEnd %" PetscInt_FMT, nroots, pEnd); 1291 PetscCheck(maxleaf < nroots,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Max local leaf %" PetscInt_FMT " >= nroots %" PetscInt_FMT, maxleaf, nroots); 1292 PetscCall(PetscMalloc2(nroots,&neg,nlocal,&recv)); 1293 PetscCall(PetscArrayzero(neg,nroots)); 1294 } 1295 /* Mark all local points with negative dof */ 1296 for (p = pStart; p < pEnd; ++p) { 1297 PetscCall(PetscSectionGetDof(s, p, &dof)); 1298 PetscCall(PetscSectionSetDof(gs, p, dof)); 1299 PetscCall(PetscSectionGetConstraintDof(s, p, &cdof)); 1300 if (!includeConstraints && cdof > 0) PetscCall(PetscSectionSetConstraintDof(gs, p, cdof)); 1301 if (neg) neg[p] = -(dof+1); 1302 } 1303 PetscCall(PetscSectionSetUpBC(gs)); 1304 if (gs->bcIndices) PetscCall(PetscArraycpy(gs->bcIndices, s->bcIndices, gs->bc->atlasOff[gs->bc->pEnd-gs->bc->pStart-1] + gs->bc->atlasDof[gs->bc->pEnd-gs->bc->pStart-1])); 1305 if (nroots >= 0) { 1306 PetscCall(PetscArrayzero(recv,nlocal)); 1307 PetscCall(PetscSFBcastBegin(sf, MPIU_INT, neg, recv,MPI_REPLACE)); 1308 PetscCall(PetscSFBcastEnd(sf, MPIU_INT, neg, recv,MPI_REPLACE)); 1309 for (p = pStart; p < pEnd; ++p) { 1310 if (recv[p] < 0) { 1311 gs->atlasDof[p-pStart] = recv[p]; 1312 PetscCall(PetscSectionGetDof(s, p, &dof)); 1313 PetscCheck(-(recv[p]+1) == dof,PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global dof %" PetscInt_FMT " for point %" PetscInt_FMT " is not the unconstrained %" PetscInt_FMT, -(recv[p]+1), p, dof); 1314 } 1315 } 1316 } 1317 /* Calculate new sizes, get process offset, and calculate point offsets */ 1318 if (s->perm) PetscCall(ISGetIndices(s->perm, &pind)); 1319 for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1320 const PetscInt q = pind ? pind[p] : p; 1321 1322 cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[q] : 0; 1323 gs->atlasOff[q] = off; 1324 off += gs->atlasDof[q] > 0 ? gs->atlasDof[q]-cdof : 0; 1325 } 1326 if (!localOffsets) { 1327 PetscCallMPI(MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s))); 1328 globalOff -= off; 1329 } 1330 for (p = pStart, off = 0; p < pEnd; ++p) { 1331 gs->atlasOff[p-pStart] += globalOff; 1332 if (neg) neg[p] = -(gs->atlasOff[p-pStart]+1); 1333 } 1334 if (s->perm) PetscCall(ISRestoreIndices(s->perm, &pind)); 1335 /* Put in negative offsets for ghost points */ 1336 if (nroots >= 0) { 1337 PetscCall(PetscArrayzero(recv,nlocal)); 1338 PetscCall(PetscSFBcastBegin(sf, MPIU_INT, neg, recv,MPI_REPLACE)); 1339 PetscCall(PetscSFBcastEnd(sf, MPIU_INT, neg, recv,MPI_REPLACE)); 1340 for (p = pStart; p < pEnd; ++p) { 1341 if (recv[p] < 0) gs->atlasOff[p-pStart] = recv[p]; 1342 } 1343 } 1344 PetscCall(PetscFree2(neg,recv)); 1345 /* Set field dofs/offsets/constraints */ 1346 for (f = 0; f < numFields; ++f) { 1347 gs->field[f]->includesConstraints = includeConstraints; 1348 PetscCall(PetscSectionGetFieldComponents(s, f, &numComponents)); 1349 PetscCall(PetscSectionSetFieldComponents(gs, f, numComponents)); 1350 } 1351 for (p = pStart; p < pEnd; ++p) { 1352 PetscCall(PetscSectionGetOffset(gs, p, &off)); 1353 for (f = 0, foff = off; f < numFields; ++f) { 1354 PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cdof)); 1355 if (!includeConstraints && cdof > 0) PetscCall(PetscSectionSetFieldConstraintDof(gs, p, f, cdof)); 1356 PetscCall(PetscSectionGetFieldDof(s, p, f, &dof)); 1357 PetscCall(PetscSectionSetFieldDof(gs, p, f, off < 0 ? -(dof + 1) : dof)); 1358 PetscCall(PetscSectionSetFieldOffset(gs, p, f, foff)); 1359 PetscCall(PetscSectionGetFieldConstraintDof(gs, p, f, &cdof)); 1360 foff = off < 0 ? foff - (dof - cdof) : foff + (dof - cdof); 1361 } 1362 } 1363 for (f = 0; f < numFields; ++f) { 1364 PetscSection gfs = gs->field[f]; 1365 1366 PetscCall(PetscSectionSetUpBC(gfs)); 1367 if (gfs->bcIndices) PetscCall(PetscArraycpy(gfs->bcIndices, s->field[f]->bcIndices, gfs->bc->atlasOff[gfs->bc->pEnd-gfs->bc->pStart-1] + gfs->bc->atlasDof[gfs->bc->pEnd-gfs->bc->pStart-1])); 1368 } 1369 gs->setup = PETSC_TRUE; 1370 PetscCall(PetscSectionViewFromOptions(gs, NULL, "-global_section_view")); 1371 *gsection = gs; 1372 PetscFunctionReturn(0); 1373 } 1374 1375 /*@ 1376 PetscSectionCreateGlobalSectionCensored - Create a section describing the global field layout using 1377 the local section and an SF describing the section point overlap. 1378 1379 Input Parameters: 1380 + s - The PetscSection for the local field layout 1381 . sf - The SF describing parallel layout of the section points 1382 . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs 1383 . numExcludes - The number of exclusion ranges 1384 - excludes - An array [start_0, end_0, start_1, end_1, ...] where there are numExcludes pairs 1385 1386 Output Parameter: 1387 . gsection - The PetscSection for the global field layout 1388 1389 Note: This gives negative sizes and offsets to points not owned by this process 1390 1391 Level: advanced 1392 1393 .seealso: PetscSectionCreate() 1394 @*/ 1395 PetscErrorCode PetscSectionCreateGlobalSectionCensored(PetscSection s, PetscSF sf, PetscBool includeConstraints, PetscInt numExcludes, const PetscInt excludes[], PetscSection *gsection) 1396 { 1397 const PetscInt *pind = NULL; 1398 PetscInt *neg = NULL, *tmpOff = NULL; 1399 PetscInt pStart, pEnd, p, e, dof, cdof, off, globalOff = 0, nroots; 1400 1401 PetscFunctionBegin; 1402 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1403 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 1404 PetscValidPointer(gsection, 6); 1405 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection)); 1406 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 1407 PetscCall(PetscSectionSetChart(*gsection, pStart, pEnd)); 1408 PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 1409 if (nroots >= 0) { 1410 PetscCheck(nroots >= pEnd-pStart,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %" PetscInt_FMT " < %" PetscInt_FMT " section size", nroots, pEnd-pStart); 1411 PetscCall(PetscCalloc1(nroots, &neg)); 1412 if (nroots > pEnd-pStart) { 1413 PetscCall(PetscCalloc1(nroots, &tmpOff)); 1414 } else { 1415 tmpOff = &(*gsection)->atlasDof[-pStart]; 1416 } 1417 } 1418 /* Mark ghost points with negative dof */ 1419 for (p = pStart; p < pEnd; ++p) { 1420 for (e = 0; e < numExcludes; ++e) { 1421 if ((p >= excludes[e*2+0]) && (p < excludes[e*2+1])) { 1422 PetscCall(PetscSectionSetDof(*gsection, p, 0)); 1423 break; 1424 } 1425 } 1426 if (e < numExcludes) continue; 1427 PetscCall(PetscSectionGetDof(s, p, &dof)); 1428 PetscCall(PetscSectionSetDof(*gsection, p, dof)); 1429 PetscCall(PetscSectionGetConstraintDof(s, p, &cdof)); 1430 if (!includeConstraints && cdof > 0) PetscCall(PetscSectionSetConstraintDof(*gsection, p, cdof)); 1431 if (neg) neg[p] = -(dof+1); 1432 } 1433 PetscCall(PetscSectionSetUpBC(*gsection)); 1434 if (nroots >= 0) { 1435 PetscCall(PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff,MPI_REPLACE)); 1436 PetscCall(PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff,MPI_REPLACE)); 1437 if (nroots > pEnd - pStart) { 1438 for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];} 1439 } 1440 } 1441 /* Calculate new sizes, get proccess offset, and calculate point offsets */ 1442 if (s->perm) PetscCall(ISGetIndices(s->perm, &pind)); 1443 for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1444 const PetscInt q = pind ? pind[p] : p; 1445 1446 cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[q] : 0; 1447 (*gsection)->atlasOff[q] = off; 1448 off += (*gsection)->atlasDof[q] > 0 ? (*gsection)->atlasDof[q]-cdof : 0; 1449 } 1450 PetscCallMPI(MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s))); 1451 globalOff -= off; 1452 for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1453 (*gsection)->atlasOff[p] += globalOff; 1454 if (neg) neg[p+pStart] = -((*gsection)->atlasOff[p]+1); 1455 } 1456 if (s->perm) PetscCall(ISRestoreIndices(s->perm, &pind)); 1457 /* Put in negative offsets for ghost points */ 1458 if (nroots >= 0) { 1459 if (nroots == pEnd-pStart) tmpOff = &(*gsection)->atlasOff[-pStart]; 1460 PetscCall(PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff,MPI_REPLACE)); 1461 PetscCall(PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff,MPI_REPLACE)); 1462 if (nroots > pEnd - pStart) { 1463 for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];} 1464 } 1465 } 1466 if (nroots >= 0 && nroots > pEnd-pStart) PetscCall(PetscFree(tmpOff)); 1467 PetscCall(PetscFree(neg)); 1468 PetscFunctionReturn(0); 1469 } 1470 1471 /*@ 1472 PetscSectionGetPointLayout - Get the PetscLayout associated with the section points 1473 1474 Collective on comm 1475 1476 Input Parameters: 1477 + comm - The MPI_Comm 1478 - s - The PetscSection 1479 1480 Output Parameter: 1481 . layout - The point layout for the section 1482 1483 Note: This is usually called for the default global section. 1484 1485 Level: advanced 1486 1487 .seealso: PetscSectionGetValueLayout(), PetscSectionCreate() 1488 @*/ 1489 PetscErrorCode PetscSectionGetPointLayout(MPI_Comm comm, PetscSection s, PetscLayout *layout) 1490 { 1491 PetscInt pStart, pEnd, p, localSize = 0; 1492 1493 PetscFunctionBegin; 1494 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 1495 for (p = pStart; p < pEnd; ++p) { 1496 PetscInt dof; 1497 1498 PetscCall(PetscSectionGetDof(s, p, &dof)); 1499 if (dof >= 0) ++localSize; 1500 } 1501 PetscCall(PetscLayoutCreate(comm, layout)); 1502 PetscCall(PetscLayoutSetLocalSize(*layout, localSize)); 1503 PetscCall(PetscLayoutSetBlockSize(*layout, 1)); 1504 PetscCall(PetscLayoutSetUp(*layout)); 1505 PetscFunctionReturn(0); 1506 } 1507 1508 /*@ 1509 PetscSectionGetValueLayout - Get the PetscLayout associated with the section dofs. 1510 1511 Collective on comm 1512 1513 Input Parameters: 1514 + comm - The MPI_Comm 1515 - s - The PetscSection 1516 1517 Output Parameter: 1518 . layout - The dof layout for the section 1519 1520 Note: This is usually called for the default global section. 1521 1522 Level: advanced 1523 1524 .seealso: PetscSectionGetPointLayout(), PetscSectionCreate() 1525 @*/ 1526 PetscErrorCode PetscSectionGetValueLayout(MPI_Comm comm, PetscSection s, PetscLayout *layout) 1527 { 1528 PetscInt pStart, pEnd, p, localSize = 0; 1529 1530 PetscFunctionBegin; 1531 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 2); 1532 PetscValidPointer(layout, 3); 1533 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 1534 for (p = pStart; p < pEnd; ++p) { 1535 PetscInt dof,cdof; 1536 1537 PetscCall(PetscSectionGetDof(s, p, &dof)); 1538 PetscCall(PetscSectionGetConstraintDof(s, p, &cdof)); 1539 if (dof-cdof > 0) localSize += dof-cdof; 1540 } 1541 PetscCall(PetscLayoutCreate(comm, layout)); 1542 PetscCall(PetscLayoutSetLocalSize(*layout, localSize)); 1543 PetscCall(PetscLayoutSetBlockSize(*layout, 1)); 1544 PetscCall(PetscLayoutSetUp(*layout)); 1545 PetscFunctionReturn(0); 1546 } 1547 1548 /*@ 1549 PetscSectionGetOffset - Return the offset into an array or local Vec for the dof associated with the given point. 1550 1551 Not Collective 1552 1553 Input Parameters: 1554 + s - the PetscSection 1555 - point - the point 1556 1557 Output Parameter: 1558 . offset - the offset 1559 1560 Level: intermediate 1561 1562 .seealso: PetscSectionGetFieldOffset(), PetscSectionCreate() 1563 @*/ 1564 PetscErrorCode PetscSectionGetOffset(PetscSection s, PetscInt point, PetscInt *offset) 1565 { 1566 PetscFunctionBegin; 1567 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1568 PetscValidIntPointer(offset, 3); 1569 if (PetscDefined(USE_DEBUG)) { 1570 PetscCheck(!(point < s->pStart) && !(point >= s->pEnd),PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %" PetscInt_FMT " should be in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, s->pStart, s->pEnd); 1571 } 1572 *offset = s->atlasOff[point - s->pStart]; 1573 PetscFunctionReturn(0); 1574 } 1575 1576 /*@ 1577 PetscSectionSetOffset - Set the offset into an array or local Vec for the dof associated with the given point. 1578 1579 Not Collective 1580 1581 Input Parameters: 1582 + s - the PetscSection 1583 . point - the point 1584 - offset - the offset 1585 1586 Note: The user usually does not call this function, but uses PetscSectionSetUp() 1587 1588 Level: intermediate 1589 1590 .seealso: PetscSectionGetFieldOffset(), PetscSectionCreate(), PetscSectionSetUp() 1591 @*/ 1592 PetscErrorCode PetscSectionSetOffset(PetscSection s, PetscInt point, PetscInt offset) 1593 { 1594 PetscFunctionBegin; 1595 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1596 PetscCheck(!(point < s->pStart) && !(point >= s->pEnd),PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %" PetscInt_FMT " should be in [%" PetscInt_FMT ", %" PetscInt_FMT ")", point, s->pStart, s->pEnd); 1597 s->atlasOff[point - s->pStart] = offset; 1598 PetscFunctionReturn(0); 1599 } 1600 1601 /*@ 1602 PetscSectionGetFieldOffset - Return the offset into an array or local Vec for the dof associated with the given point. 1603 1604 Not Collective 1605 1606 Input Parameters: 1607 + s - the PetscSection 1608 . point - the point 1609 - field - the field 1610 1611 Output Parameter: 1612 . offset - the offset 1613 1614 Level: intermediate 1615 1616 .seealso: PetscSectionGetOffset(), PetscSectionCreate() 1617 @*/ 1618 PetscErrorCode PetscSectionGetFieldOffset(PetscSection s, PetscInt point, PetscInt field, PetscInt *offset) 1619 { 1620 PetscFunctionBegin; 1621 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1622 PetscValidIntPointer(offset, 4); 1623 PetscSectionCheckValidField(field,s->numFields); 1624 PetscCall(PetscSectionGetOffset(s->field[field], point, offset)); 1625 PetscFunctionReturn(0); 1626 } 1627 1628 /*@ 1629 PetscSectionSetFieldOffset - Set the offset into an array or local Vec for the dof associated with the given field at a point. 1630 1631 Not Collective 1632 1633 Input Parameters: 1634 + s - the PetscSection 1635 . point - the point 1636 . field - the field 1637 - offset - the offset 1638 1639 Note: The user usually does not call this function, but uses PetscSectionSetUp() 1640 1641 Level: intermediate 1642 1643 .seealso: PetscSectionGetFieldOffset(), PetscSectionSetOffset(), PetscSectionCreate(), PetscSectionSetUp() 1644 @*/ 1645 PetscErrorCode PetscSectionSetFieldOffset(PetscSection s, PetscInt point, PetscInt field, PetscInt offset) 1646 { 1647 PetscFunctionBegin; 1648 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1649 PetscSectionCheckValidField(field,s->numFields); 1650 PetscCall(PetscSectionSetOffset(s->field[field], point, offset)); 1651 PetscFunctionReturn(0); 1652 } 1653 1654 /*@ 1655 PetscSectionGetFieldPointOffset - Return the offset on the given point for the dof associated with the given point. 1656 1657 Not Collective 1658 1659 Input Parameters: 1660 + s - the PetscSection 1661 . point - the point 1662 - field - the field 1663 1664 Output Parameter: 1665 . offset - the offset 1666 1667 Note: This gives the offset on a point of the field, ignoring constraints, meaning starting at the first dof for 1668 this point, what number is the first dof with this field. 1669 1670 Level: advanced 1671 1672 .seealso: PetscSectionGetOffset(), PetscSectionCreate() 1673 @*/ 1674 PetscErrorCode PetscSectionGetFieldPointOffset(PetscSection s, PetscInt point, PetscInt field, PetscInt *offset) 1675 { 1676 PetscInt off, foff; 1677 1678 PetscFunctionBegin; 1679 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1680 PetscValidIntPointer(offset, 4); 1681 PetscSectionCheckValidField(field,s->numFields); 1682 PetscCall(PetscSectionGetOffset(s, point, &off)); 1683 PetscCall(PetscSectionGetOffset(s->field[field], point, &foff)); 1684 *offset = foff - off; 1685 PetscFunctionReturn(0); 1686 } 1687 1688 /*@ 1689 PetscSectionGetOffsetRange - Return the full range of offsets [start, end) 1690 1691 Not Collective 1692 1693 Input Parameter: 1694 . s - the PetscSection 1695 1696 Output Parameters: 1697 + start - the minimum offset 1698 - end - one more than the maximum offset 1699 1700 Level: intermediate 1701 1702 .seealso: PetscSectionGetOffset(), PetscSectionCreate() 1703 @*/ 1704 PetscErrorCode PetscSectionGetOffsetRange(PetscSection s, PetscInt *start, PetscInt *end) 1705 { 1706 PetscInt os = 0, oe = 0, pStart, pEnd, p; 1707 1708 PetscFunctionBegin; 1709 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1710 if (s->atlasOff) {os = s->atlasOff[0]; oe = s->atlasOff[0];} 1711 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 1712 for (p = 0; p < pEnd-pStart; ++p) { 1713 PetscInt dof = s->atlasDof[p], off = s->atlasOff[p]; 1714 1715 if (off >= 0) { 1716 os = PetscMin(os, off); 1717 oe = PetscMax(oe, off+dof); 1718 } 1719 } 1720 if (start) *start = os; 1721 if (end) *end = oe; 1722 PetscFunctionReturn(0); 1723 } 1724 1725 /*@ 1726 PetscSectionCreateSubsection - Create a new, smaller section composed of only the selected fields 1727 1728 Collective 1729 1730 Input Parameters: 1731 + s - the PetscSection 1732 . len - the number of subfields 1733 - fields - the subfield numbers 1734 1735 Output Parameter: 1736 . subs - the subsection 1737 1738 Note: The section offsets now refer to a new, smaller vector. 1739 1740 Level: advanced 1741 1742 .seealso: PetscSectionCreateSupersection(), PetscSectionCreate() 1743 @*/ 1744 PetscErrorCode PetscSectionCreateSubsection(PetscSection s, PetscInt len, const PetscInt fields[], PetscSection *subs) 1745 { 1746 PetscInt nF, f, c, pStart, pEnd, p, maxCdof = 0; 1747 1748 PetscFunctionBegin; 1749 if (!len) PetscFunctionReturn(0); 1750 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1751 PetscValidIntPointer(fields, 3); 1752 PetscValidPointer(subs, 4); 1753 PetscCall(PetscSectionGetNumFields(s, &nF)); 1754 PetscCheck(len <= nF,PetscObjectComm((PetscObject) s), PETSC_ERR_ARG_WRONG, "Number of requested fields %" PetscInt_FMT " greater than number of fields %" PetscInt_FMT, len, nF); 1755 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) s), subs)); 1756 PetscCall(PetscSectionSetNumFields(*subs, len)); 1757 for (f = 0; f < len; ++f) { 1758 const char *name = NULL; 1759 PetscInt numComp = 0; 1760 1761 PetscCall(PetscSectionGetFieldName(s, fields[f], &name)); 1762 PetscCall(PetscSectionSetFieldName(*subs, f, name)); 1763 PetscCall(PetscSectionGetFieldComponents(s, fields[f], &numComp)); 1764 PetscCall(PetscSectionSetFieldComponents(*subs, f, numComp)); 1765 for (c = 0; c < s->numFieldComponents[fields[f]]; ++c) { 1766 PetscCall(PetscSectionGetComponentName(s, fields[f], c, &name)); 1767 PetscCall(PetscSectionSetComponentName(*subs, f, c, name)); 1768 } 1769 } 1770 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 1771 PetscCall(PetscSectionSetChart(*subs, pStart, pEnd)); 1772 for (p = pStart; p < pEnd; ++p) { 1773 PetscInt dof = 0, cdof = 0, fdof = 0, cfdof = 0; 1774 1775 for (f = 0; f < len; ++f) { 1776 PetscCall(PetscSectionGetFieldDof(s, p, fields[f], &fdof)); 1777 PetscCall(PetscSectionSetFieldDof(*subs, p, f, fdof)); 1778 PetscCall(PetscSectionGetFieldConstraintDof(s, p, fields[f], &cfdof)); 1779 if (cfdof) PetscCall(PetscSectionSetFieldConstraintDof(*subs, p, f, cfdof)); 1780 dof += fdof; 1781 cdof += cfdof; 1782 } 1783 PetscCall(PetscSectionSetDof(*subs, p, dof)); 1784 if (cdof) PetscCall(PetscSectionSetConstraintDof(*subs, p, cdof)); 1785 maxCdof = PetscMax(cdof, maxCdof); 1786 } 1787 PetscCall(PetscSectionSetUp(*subs)); 1788 if (maxCdof) { 1789 PetscInt *indices; 1790 1791 PetscCall(PetscMalloc1(maxCdof, &indices)); 1792 for (p = pStart; p < pEnd; ++p) { 1793 PetscInt cdof; 1794 1795 PetscCall(PetscSectionGetConstraintDof(*subs, p, &cdof)); 1796 if (cdof) { 1797 const PetscInt *oldIndices = NULL; 1798 PetscInt fdof = 0, cfdof = 0, fc, numConst = 0, fOff = 0; 1799 1800 for (f = 0; f < len; ++f) { 1801 PetscCall(PetscSectionGetFieldDof(s, p, fields[f], &fdof)); 1802 PetscCall(PetscSectionGetFieldConstraintDof(s, p, fields[f], &cfdof)); 1803 PetscCall(PetscSectionGetFieldConstraintIndices(s, p, fields[f], &oldIndices)); 1804 PetscCall(PetscSectionSetFieldConstraintIndices(*subs, p, f, oldIndices)); 1805 for (fc = 0; fc < cfdof; ++fc) indices[numConst+fc] = oldIndices[fc] + fOff; 1806 numConst += cfdof; 1807 fOff += fdof; 1808 } 1809 PetscCall(PetscSectionSetConstraintIndices(*subs, p, indices)); 1810 } 1811 } 1812 PetscCall(PetscFree(indices)); 1813 } 1814 PetscFunctionReturn(0); 1815 } 1816 1817 /*@ 1818 PetscSectionCreateSupersection - Create a new, larger section composed of the input sections 1819 1820 Collective 1821 1822 Input Parameters: 1823 + s - the input sections 1824 - len - the number of input sections 1825 1826 Output Parameter: 1827 . supers - the supersection 1828 1829 Note: The section offsets now refer to a new, larger vector. 1830 1831 Level: advanced 1832 1833 .seealso: PetscSectionCreateSubsection(), PetscSectionCreate() 1834 @*/ 1835 PetscErrorCode PetscSectionCreateSupersection(PetscSection s[], PetscInt len, PetscSection *supers) 1836 { 1837 PetscInt Nf = 0, f, pStart = PETSC_MAX_INT, pEnd = 0, p, maxCdof = 0, i; 1838 1839 PetscFunctionBegin; 1840 if (!len) PetscFunctionReturn(0); 1841 for (i = 0; i < len; ++i) { 1842 PetscInt nf, pStarti, pEndi; 1843 1844 PetscCall(PetscSectionGetNumFields(s[i], &nf)); 1845 PetscCall(PetscSectionGetChart(s[i], &pStarti, &pEndi)); 1846 pStart = PetscMin(pStart, pStarti); 1847 pEnd = PetscMax(pEnd, pEndi); 1848 Nf += nf; 1849 } 1850 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) s[0]), supers)); 1851 PetscCall(PetscSectionSetNumFields(*supers, Nf)); 1852 for (i = 0, f = 0; i < len; ++i) { 1853 PetscInt nf, fi, ci; 1854 1855 PetscCall(PetscSectionGetNumFields(s[i], &nf)); 1856 for (fi = 0; fi < nf; ++fi, ++f) { 1857 const char *name = NULL; 1858 PetscInt numComp = 0; 1859 1860 PetscCall(PetscSectionGetFieldName(s[i], fi, &name)); 1861 PetscCall(PetscSectionSetFieldName(*supers, f, name)); 1862 PetscCall(PetscSectionGetFieldComponents(s[i], fi, &numComp)); 1863 PetscCall(PetscSectionSetFieldComponents(*supers, f, numComp)); 1864 for (ci = 0; ci < s[i]->numFieldComponents[fi]; ++ci) { 1865 PetscCall(PetscSectionGetComponentName(s[i], fi, ci, &name)); 1866 PetscCall(PetscSectionSetComponentName(*supers, f, ci, name)); 1867 } 1868 } 1869 } 1870 PetscCall(PetscSectionSetChart(*supers, pStart, pEnd)); 1871 for (p = pStart; p < pEnd; ++p) { 1872 PetscInt dof = 0, cdof = 0; 1873 1874 for (i = 0, f = 0; i < len; ++i) { 1875 PetscInt nf, fi, pStarti, pEndi; 1876 PetscInt fdof = 0, cfdof = 0; 1877 1878 PetscCall(PetscSectionGetNumFields(s[i], &nf)); 1879 PetscCall(PetscSectionGetChart(s[i], &pStarti, &pEndi)); 1880 if ((p < pStarti) || (p >= pEndi)) continue; 1881 for (fi = 0; fi < nf; ++fi, ++f) { 1882 PetscCall(PetscSectionGetFieldDof(s[i], p, fi, &fdof)); 1883 PetscCall(PetscSectionAddFieldDof(*supers, p, f, fdof)); 1884 PetscCall(PetscSectionGetFieldConstraintDof(s[i], p, fi, &cfdof)); 1885 if (cfdof) PetscCall(PetscSectionAddFieldConstraintDof(*supers, p, f, cfdof)); 1886 dof += fdof; 1887 cdof += cfdof; 1888 } 1889 } 1890 PetscCall(PetscSectionSetDof(*supers, p, dof)); 1891 if (cdof) PetscCall(PetscSectionSetConstraintDof(*supers, p, cdof)); 1892 maxCdof = PetscMax(cdof, maxCdof); 1893 } 1894 PetscCall(PetscSectionSetUp(*supers)); 1895 if (maxCdof) { 1896 PetscInt *indices; 1897 1898 PetscCall(PetscMalloc1(maxCdof, &indices)); 1899 for (p = pStart; p < pEnd; ++p) { 1900 PetscInt cdof; 1901 1902 PetscCall(PetscSectionGetConstraintDof(*supers, p, &cdof)); 1903 if (cdof) { 1904 PetscInt dof, numConst = 0, fOff = 0; 1905 1906 for (i = 0, f = 0; i < len; ++i) { 1907 const PetscInt *oldIndices = NULL; 1908 PetscInt nf, fi, pStarti, pEndi, fdof, cfdof, fc; 1909 1910 PetscCall(PetscSectionGetNumFields(s[i], &nf)); 1911 PetscCall(PetscSectionGetChart(s[i], &pStarti, &pEndi)); 1912 if ((p < pStarti) || (p >= pEndi)) continue; 1913 for (fi = 0; fi < nf; ++fi, ++f) { 1914 PetscCall(PetscSectionGetFieldDof(s[i], p, fi, &fdof)); 1915 PetscCall(PetscSectionGetFieldConstraintDof(s[i], p, fi, &cfdof)); 1916 PetscCall(PetscSectionGetFieldConstraintIndices(s[i], p, fi, &oldIndices)); 1917 for (fc = 0; fc < cfdof; ++fc) indices[numConst+fc] = oldIndices[fc]; 1918 PetscCall(PetscSectionSetFieldConstraintIndices(*supers, p, f, &indices[numConst])); 1919 for (fc = 0; fc < cfdof; ++fc) indices[numConst+fc] += fOff; 1920 numConst += cfdof; 1921 } 1922 PetscCall(PetscSectionGetDof(s[i], p, &dof)); 1923 fOff += dof; 1924 } 1925 PetscCall(PetscSectionSetConstraintIndices(*supers, p, indices)); 1926 } 1927 } 1928 PetscCall(PetscFree(indices)); 1929 } 1930 PetscFunctionReturn(0); 1931 } 1932 1933 PetscErrorCode PetscSectionCreateSubplexSection_Internal(PetscSection s, IS subpointMap, PetscBool renumberPoints, PetscSection *subs) 1934 { 1935 const PetscInt *points = NULL, *indices = NULL; 1936 PetscInt numFields, f, c, numSubpoints = 0, pStart, pEnd, p, spStart, spEnd, subp; 1937 1938 PetscFunctionBegin; 1939 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1940 PetscValidHeaderSpecific(subpointMap, IS_CLASSID, 2); 1941 PetscValidPointer(subs, 4); 1942 PetscCall(PetscSectionGetNumFields(s, &numFields)); 1943 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) s), subs)); 1944 if (numFields) PetscCall(PetscSectionSetNumFields(*subs, numFields)); 1945 for (f = 0; f < numFields; ++f) { 1946 const char *name = NULL; 1947 PetscInt numComp = 0; 1948 1949 PetscCall(PetscSectionGetFieldName(s, f, &name)); 1950 PetscCall(PetscSectionSetFieldName(*subs, f, name)); 1951 PetscCall(PetscSectionGetFieldComponents(s, f, &numComp)); 1952 PetscCall(PetscSectionSetFieldComponents(*subs, f, numComp)); 1953 for (c = 0; c < s->numFieldComponents[f]; ++c) { 1954 PetscCall(PetscSectionGetComponentName(s, f, c, &name)); 1955 PetscCall(PetscSectionSetComponentName(*subs, f, c, name)); 1956 } 1957 } 1958 /* For right now, we do not try to squeeze the subchart */ 1959 if (subpointMap) { 1960 PetscCall(ISGetSize(subpointMap, &numSubpoints)); 1961 PetscCall(ISGetIndices(subpointMap, &points)); 1962 } 1963 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 1964 if (renumberPoints) { 1965 spStart = 0; 1966 spEnd = numSubpoints; 1967 } else { 1968 PetscCall(ISGetMinMax(subpointMap, &spStart, &spEnd)); 1969 ++spEnd; 1970 } 1971 PetscCall(PetscSectionSetChart(*subs, spStart, spEnd)); 1972 for (p = pStart; p < pEnd; ++p) { 1973 PetscInt dof, cdof, fdof = 0, cfdof = 0; 1974 1975 PetscCall(PetscFindInt(p, numSubpoints, points, &subp)); 1976 if (subp < 0) continue; 1977 if (!renumberPoints) subp = p; 1978 for (f = 0; f < numFields; ++f) { 1979 PetscCall(PetscSectionGetFieldDof(s, p, f, &fdof)); 1980 PetscCall(PetscSectionSetFieldDof(*subs, subp, f, fdof)); 1981 PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cfdof)); 1982 if (cfdof) PetscCall(PetscSectionSetFieldConstraintDof(*subs, subp, f, cfdof)); 1983 } 1984 PetscCall(PetscSectionGetDof(s, p, &dof)); 1985 PetscCall(PetscSectionSetDof(*subs, subp, dof)); 1986 PetscCall(PetscSectionGetConstraintDof(s, p, &cdof)); 1987 if (cdof) PetscCall(PetscSectionSetConstraintDof(*subs, subp, cdof)); 1988 } 1989 PetscCall(PetscSectionSetUp(*subs)); 1990 /* Change offsets to original offsets */ 1991 for (p = pStart; p < pEnd; ++p) { 1992 PetscInt off, foff = 0; 1993 1994 PetscCall(PetscFindInt(p, numSubpoints, points, &subp)); 1995 if (subp < 0) continue; 1996 if (!renumberPoints) subp = p; 1997 for (f = 0; f < numFields; ++f) { 1998 PetscCall(PetscSectionGetFieldOffset(s, p, f, &foff)); 1999 PetscCall(PetscSectionSetFieldOffset(*subs, subp, f, foff)); 2000 } 2001 PetscCall(PetscSectionGetOffset(s, p, &off)); 2002 PetscCall(PetscSectionSetOffset(*subs, subp, off)); 2003 } 2004 /* Copy constraint indices */ 2005 for (subp = spStart; subp < spEnd; ++subp) { 2006 PetscInt cdof; 2007 2008 PetscCall(PetscSectionGetConstraintDof(*subs, subp, &cdof)); 2009 if (cdof) { 2010 for (f = 0; f < numFields; ++f) { 2011 PetscCall(PetscSectionGetFieldConstraintIndices(s, points[subp-spStart], f, &indices)); 2012 PetscCall(PetscSectionSetFieldConstraintIndices(*subs, subp, f, indices)); 2013 } 2014 PetscCall(PetscSectionGetConstraintIndices(s, points[subp-spStart], &indices)); 2015 PetscCall(PetscSectionSetConstraintIndices(*subs, subp, indices)); 2016 } 2017 } 2018 if (subpointMap) PetscCall(ISRestoreIndices(subpointMap, &points)); 2019 PetscFunctionReturn(0); 2020 } 2021 2022 /*@ 2023 PetscSectionCreateSubmeshSection - Create a new, smaller section with support on the submesh 2024 2025 Collective on s 2026 2027 Input Parameters: 2028 + s - the PetscSection 2029 - subpointMap - a sorted list of points in the original mesh which are in the submesh 2030 2031 Output Parameter: 2032 . subs - the subsection 2033 2034 Note: 2035 The points are renumbered from 0, and the section offsets now refer to a new, smaller vector. 2036 2037 Level: advanced 2038 2039 .seealso: PetscSectionCreateSubdomainSection(), PetscSectionCreateSubsection(), DMPlexGetSubpointMap(), PetscSectionCreate() 2040 @*/ 2041 PetscErrorCode PetscSectionCreateSubmeshSection(PetscSection s, IS subpointMap, PetscSection *subs) 2042 { 2043 PetscFunctionBegin; 2044 PetscCall(PetscSectionCreateSubplexSection_Internal(s, subpointMap, PETSC_TRUE, subs)); 2045 PetscFunctionReturn(0); 2046 } 2047 2048 /*@ 2049 PetscSectionCreateSubdomainSection - Create a new, smaller section with support on a subdomain of the mesh 2050 2051 Collective on s 2052 2053 Input Parameters: 2054 + s - the PetscSection 2055 - subpointMap - a sorted list of points in the original mesh which are in the subdomain 2056 2057 Output Parameter: 2058 . subs - the subsection 2059 2060 Note: 2061 The point numbers remain the same, but the section offsets now refer to a new, smaller vector. 2062 2063 Level: advanced 2064 2065 .seealso: PetscSectionCreateSubmeshSection(), PetscSectionCreateSubsection(), DMPlexGetSubpointMap(), PetscSectionCreate() 2066 @*/ 2067 PetscErrorCode PetscSectionCreateSubdomainSection(PetscSection s, IS subpointMap, PetscSection *subs) 2068 { 2069 PetscFunctionBegin; 2070 PetscCall(PetscSectionCreateSubplexSection_Internal(s, subpointMap, PETSC_FALSE, subs)); 2071 PetscFunctionReturn(0); 2072 } 2073 2074 static PetscErrorCode PetscSectionView_ASCII(PetscSection s, PetscViewer viewer) 2075 { 2076 PetscInt p; 2077 PetscMPIInt rank; 2078 2079 PetscFunctionBegin; 2080 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank)); 2081 PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 2082 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Process %d:\n", rank)); 2083 for (p = 0; p < s->pEnd - s->pStart; ++p) { 2084 if ((s->bc) && (s->bc->atlasDof[p] > 0)) { 2085 PetscInt b; 2086 2087 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " (%4" PetscInt_FMT ") dim %2" PetscInt_FMT " offset %3" PetscInt_FMT " constrained", p+s->pStart, s->atlasDof[p], s->atlasOff[p])); 2088 if (s->bcIndices) { 2089 for (b = 0; b < s->bc->atlasDof[p]; ++b) { 2090 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %" PetscInt_FMT, s->bcIndices[s->bc->atlasOff[p]+b])); 2091 } 2092 } 2093 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n")); 2094 } else { 2095 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " (%4" PetscInt_FMT ") dim %2" PetscInt_FMT " offset %3" PetscInt_FMT "\n", p+s->pStart, s->atlasDof[p], s->atlasOff[p])); 2096 } 2097 } 2098 PetscCall(PetscViewerFlush(viewer)); 2099 PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 2100 if (s->sym) { 2101 PetscCall(PetscViewerASCIIPushTab(viewer)); 2102 PetscCall(PetscSectionSymView(s->sym,viewer)); 2103 PetscCall(PetscViewerASCIIPopTab(viewer)); 2104 } 2105 PetscFunctionReturn(0); 2106 } 2107 2108 /*@C 2109 PetscSectionViewFromOptions - View from Options 2110 2111 Collective 2112 2113 Input Parameters: 2114 + A - the PetscSection object to view 2115 . obj - Optional object 2116 - name - command line option 2117 2118 Level: intermediate 2119 .seealso: PetscSection, PetscSectionView, PetscObjectViewFromOptions(), PetscSectionCreate() 2120 @*/ 2121 PetscErrorCode PetscSectionViewFromOptions(PetscSection A,PetscObject obj,const char name[]) 2122 { 2123 PetscFunctionBegin; 2124 PetscValidHeaderSpecific(A,PETSC_SECTION_CLASSID,1); 2125 PetscCall(PetscObjectViewFromOptions((PetscObject)A,obj,name)); 2126 PetscFunctionReturn(0); 2127 } 2128 2129 /*@C 2130 PetscSectionView - Views a PetscSection 2131 2132 Collective 2133 2134 Input Parameters: 2135 + s - the PetscSection object to view 2136 - v - the viewer 2137 2138 Note: 2139 PetscSectionView(), when viewer is of type PETSCVIEWERHDF5, only saves 2140 distribution independent data, such as dofs, offsets, constraint dofs, 2141 and constraint indices. Points that have negative dofs, for instance, 2142 are not saved as they represent points owned by other processes. 2143 Point numbering and rank assignment is currently not stored. 2144 The saved section can be loaded with PetscSectionLoad(). 2145 2146 Level: beginner 2147 2148 .seealso PetscSectionCreate(), PetscSectionDestroy(), PetscSectionLoad() 2149 @*/ 2150 PetscErrorCode PetscSectionView(PetscSection s, PetscViewer viewer) 2151 { 2152 PetscBool isascii, ishdf5; 2153 PetscInt f; 2154 2155 PetscFunctionBegin; 2156 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2157 if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)s), &viewer)); 2158 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 2159 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii)); 2160 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5)); 2161 if (isascii) { 2162 PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)s,viewer)); 2163 if (s->numFields) { 2164 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " fields\n", s->numFields)); 2165 for (f = 0; f < s->numFields; ++f) { 2166 PetscCall(PetscViewerASCIIPrintf(viewer, " field %" PetscInt_FMT " with %" PetscInt_FMT " components\n", f, s->numFieldComponents[f])); 2167 PetscCall(PetscSectionView_ASCII(s->field[f], viewer)); 2168 } 2169 } else { 2170 PetscCall(PetscSectionView_ASCII(s, viewer)); 2171 } 2172 } else if (ishdf5) { 2173 #if PetscDefined(HAVE_HDF5) 2174 PetscCall(PetscSectionView_HDF5_Internal(s, viewer)); 2175 #else 2176 SETERRQ(PetscObjectComm((PetscObject) s), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5"); 2177 #endif 2178 } 2179 PetscFunctionReturn(0); 2180 } 2181 2182 /*@C 2183 PetscSectionLoad - Loads a PetscSection 2184 2185 Collective 2186 2187 Input Parameters: 2188 + s - the PetscSection object to load 2189 - v - the viewer 2190 2191 Note: 2192 PetscSectionLoad(), when viewer is of type PETSCVIEWERHDF5, loads 2193 a section saved with PetscSectionView(). The number of processes 2194 used here (N) does not need to be the same as that used when saving. 2195 After calling this function, the chart of s on rank i will be set 2196 to [0, E_i), where \sum_{i=0}^{N-1}E_i equals to the total number of 2197 saved section points. 2198 2199 Level: beginner 2200 2201 .seealso PetscSectionCreate(), PetscSectionDestroy(), PetscSectionView() 2202 @*/ 2203 PetscErrorCode PetscSectionLoad(PetscSection s, PetscViewer viewer) 2204 { 2205 PetscBool ishdf5; 2206 2207 PetscFunctionBegin; 2208 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2209 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 2210 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5)); 2211 if (ishdf5) { 2212 #if PetscDefined(HAVE_HDF5) 2213 PetscCall(PetscSectionLoad_HDF5_Internal(s, viewer)); 2214 PetscFunctionReturn(0); 2215 #else 2216 SETERRQ(PetscObjectComm((PetscObject) s), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5"); 2217 #endif 2218 } else SETERRQ(PetscObjectComm((PetscObject) s), PETSC_ERR_SUP, "Viewer type %s not yet supported for PetscSection loading", ((PetscObject)viewer)->type_name); 2219 } 2220 2221 static PetscErrorCode PetscSectionResetClosurePermutation(PetscSection section) 2222 { 2223 PetscSectionClosurePermVal clVal; 2224 2225 PetscFunctionBegin; 2226 if (!section->clHash) PetscFunctionReturn(0); 2227 kh_foreach_value(section->clHash, clVal, { 2228 PetscCall(PetscFree(clVal.perm)); 2229 PetscCall(PetscFree(clVal.invPerm)); 2230 }); 2231 kh_destroy(ClPerm, section->clHash); 2232 section->clHash = NULL; 2233 PetscFunctionReturn(0); 2234 } 2235 2236 /*@ 2237 PetscSectionReset - Frees all section data. 2238 2239 Not Collective 2240 2241 Input Parameters: 2242 . s - the PetscSection 2243 2244 Level: beginner 2245 2246 .seealso: PetscSection, PetscSectionCreate() 2247 @*/ 2248 PetscErrorCode PetscSectionReset(PetscSection s) 2249 { 2250 PetscInt f, c; 2251 2252 PetscFunctionBegin; 2253 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2254 for (f = 0; f < s->numFields; ++f) { 2255 PetscCall(PetscSectionDestroy(&s->field[f])); 2256 PetscCall(PetscFree(s->fieldNames[f])); 2257 for (c = 0; c < s->numFieldComponents[f]; ++c) { 2258 PetscCall(PetscFree(s->compNames[f][c])); 2259 } 2260 PetscCall(PetscFree(s->compNames[f])); 2261 } 2262 PetscCall(PetscFree(s->numFieldComponents)); 2263 PetscCall(PetscFree(s->fieldNames)); 2264 PetscCall(PetscFree(s->compNames)); 2265 PetscCall(PetscFree(s->field)); 2266 PetscCall(PetscSectionDestroy(&s->bc)); 2267 PetscCall(PetscFree(s->bcIndices)); 2268 PetscCall(PetscFree2(s->atlasDof, s->atlasOff)); 2269 PetscCall(PetscSectionDestroy(&s->clSection)); 2270 PetscCall(ISDestroy(&s->clPoints)); 2271 PetscCall(ISDestroy(&s->perm)); 2272 PetscCall(PetscSectionResetClosurePermutation(s)); 2273 PetscCall(PetscSectionSymDestroy(&s->sym)); 2274 PetscCall(PetscSectionDestroy(&s->clSection)); 2275 PetscCall(ISDestroy(&s->clPoints)); 2276 PetscCall(PetscSectionInvalidateMaxDof_Internal(s)); 2277 s->pStart = -1; 2278 s->pEnd = -1; 2279 s->maxDof = 0; 2280 s->setup = PETSC_FALSE; 2281 s->numFields = 0; 2282 s->clObj = NULL; 2283 PetscFunctionReturn(0); 2284 } 2285 2286 /*@ 2287 PetscSectionDestroy - Frees a section object and frees its range if that exists. 2288 2289 Not Collective 2290 2291 Input Parameters: 2292 . s - the PetscSection 2293 2294 Level: beginner 2295 2296 .seealso: PetscSection, PetscSectionCreate() 2297 @*/ 2298 PetscErrorCode PetscSectionDestroy(PetscSection *s) 2299 { 2300 PetscFunctionBegin; 2301 if (!*s) PetscFunctionReturn(0); 2302 PetscValidHeaderSpecific(*s,PETSC_SECTION_CLASSID,1); 2303 if (--((PetscObject)(*s))->refct > 0) { 2304 *s = NULL; 2305 PetscFunctionReturn(0); 2306 } 2307 PetscCall(PetscSectionReset(*s)); 2308 PetscCall(PetscHeaderDestroy(s)); 2309 PetscFunctionReturn(0); 2310 } 2311 2312 PetscErrorCode VecIntGetValuesSection(PetscInt *baseArray, PetscSection s, PetscInt point, const PetscInt **values) 2313 { 2314 const PetscInt p = point - s->pStart; 2315 2316 PetscFunctionBegin; 2317 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 2); 2318 *values = &baseArray[s->atlasOff[p]]; 2319 PetscFunctionReturn(0); 2320 } 2321 2322 PetscErrorCode VecIntSetValuesSection(PetscInt *baseArray, PetscSection s, PetscInt point, const PetscInt values[], InsertMode mode) 2323 { 2324 PetscInt *array; 2325 const PetscInt p = point - s->pStart; 2326 const PetscInt orientation = 0; /* Needs to be included for use in closure operations */ 2327 PetscInt cDim = 0; 2328 2329 PetscFunctionBegin; 2330 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 2); 2331 PetscCall(PetscSectionGetConstraintDof(s, p, &cDim)); 2332 array = &baseArray[s->atlasOff[p]]; 2333 if (!cDim) { 2334 if (orientation >= 0) { 2335 const PetscInt dim = s->atlasDof[p]; 2336 PetscInt i; 2337 2338 if (mode == INSERT_VALUES) { 2339 for (i = 0; i < dim; ++i) array[i] = values[i]; 2340 } else { 2341 for (i = 0; i < dim; ++i) array[i] += values[i]; 2342 } 2343 } else { 2344 PetscInt offset = 0; 2345 PetscInt j = -1, field, i; 2346 2347 for (field = 0; field < s->numFields; ++field) { 2348 const PetscInt dim = s->field[field]->atlasDof[p]; 2349 2350 for (i = dim-1; i >= 0; --i) array[++j] = values[i+offset]; 2351 offset += dim; 2352 } 2353 } 2354 } else { 2355 if (orientation >= 0) { 2356 const PetscInt dim = s->atlasDof[p]; 2357 PetscInt cInd = 0, i; 2358 const PetscInt *cDof; 2359 2360 PetscCall(PetscSectionGetConstraintIndices(s, point, &cDof)); 2361 if (mode == INSERT_VALUES) { 2362 for (i = 0; i < dim; ++i) { 2363 if ((cInd < cDim) && (i == cDof[cInd])) {++cInd; continue;} 2364 array[i] = values[i]; 2365 } 2366 } else { 2367 for (i = 0; i < dim; ++i) { 2368 if ((cInd < cDim) && (i == cDof[cInd])) {++cInd; continue;} 2369 array[i] += values[i]; 2370 } 2371 } 2372 } else { 2373 const PetscInt *cDof; 2374 PetscInt offset = 0; 2375 PetscInt cOffset = 0; 2376 PetscInt j = 0, field; 2377 2378 PetscCall(PetscSectionGetConstraintIndices(s, point, &cDof)); 2379 for (field = 0; field < s->numFields; ++field) { 2380 const PetscInt dim = s->field[field]->atlasDof[p]; /* PetscSectionGetFieldDof() */ 2381 const PetscInt tDim = s->field[field]->bc->atlasDof[p]; /* PetscSectionGetFieldConstraintDof() */ 2382 const PetscInt sDim = dim - tDim; 2383 PetscInt cInd = 0, i,k; 2384 2385 for (i = 0, k = dim+offset-1; i < dim; ++i, ++j, --k) { 2386 if ((cInd < sDim) && (j == cDof[cInd+cOffset])) {++cInd; continue;} 2387 array[j] = values[k]; 2388 } 2389 offset += dim; 2390 cOffset += dim - tDim; 2391 } 2392 } 2393 } 2394 PetscFunctionReturn(0); 2395 } 2396 2397 /*@C 2398 PetscSectionHasConstraints - Determine whether a section has constrained dofs 2399 2400 Not Collective 2401 2402 Input Parameter: 2403 . s - The PetscSection 2404 2405 Output Parameter: 2406 . hasConstraints - flag indicating that the section has constrained dofs 2407 2408 Level: intermediate 2409 2410 .seealso: PetscSectionSetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection 2411 @*/ 2412 PetscErrorCode PetscSectionHasConstraints(PetscSection s, PetscBool *hasConstraints) 2413 { 2414 PetscFunctionBegin; 2415 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2416 PetscValidBoolPointer(hasConstraints, 2); 2417 *hasConstraints = s->bc ? PETSC_TRUE : PETSC_FALSE; 2418 PetscFunctionReturn(0); 2419 } 2420 2421 /*@C 2422 PetscSectionGetConstraintIndices - Get the point dof numbers, in [0, dof), which are constrained 2423 2424 Not Collective 2425 2426 Input Parameters: 2427 + s - The PetscSection 2428 - point - The point 2429 2430 Output Parameter: 2431 . indices - The constrained dofs 2432 2433 Note: In Fortran, you call PetscSectionGetConstraintIndicesF90() and PetscSectionRestoreConstraintIndicesF90() 2434 2435 Level: intermediate 2436 2437 .seealso: PetscSectionSetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection 2438 @*/ 2439 PetscErrorCode PetscSectionGetConstraintIndices(PetscSection s, PetscInt point, const PetscInt **indices) 2440 { 2441 PetscFunctionBegin; 2442 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2443 if (s->bc) { 2444 PetscCall(VecIntGetValuesSection(s->bcIndices, s->bc, point, indices)); 2445 } else *indices = NULL; 2446 PetscFunctionReturn(0); 2447 } 2448 2449 /*@C 2450 PetscSectionSetConstraintIndices - Set the point dof numbers, in [0, dof), which are constrained 2451 2452 Not Collective 2453 2454 Input Parameters: 2455 + s - The PetscSection 2456 . point - The point 2457 - indices - The constrained dofs 2458 2459 Note: The Fortran is PetscSectionSetConstraintIndicesF90() 2460 2461 Level: intermediate 2462 2463 .seealso: PetscSectionGetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection 2464 @*/ 2465 PetscErrorCode PetscSectionSetConstraintIndices(PetscSection s, PetscInt point, const PetscInt indices[]) 2466 { 2467 PetscFunctionBegin; 2468 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2469 if (s->bc) { 2470 const PetscInt dof = s->atlasDof[point]; 2471 const PetscInt cdof = s->bc->atlasDof[point]; 2472 PetscInt d; 2473 2474 for (d = 0; d < cdof; ++d) { 2475 PetscCheck(indices[d] < dof,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %" PetscInt_FMT " dof %" PetscInt_FMT ", invalid constraint index[%" PetscInt_FMT "]: %" PetscInt_FMT, point, dof, d, indices[d]); 2476 } 2477 PetscCall(VecIntSetValuesSection(s->bcIndices, s->bc, point, indices, INSERT_VALUES)); 2478 } 2479 PetscFunctionReturn(0); 2480 } 2481 2482 /*@C 2483 PetscSectionGetFieldConstraintIndices - Get the field dof numbers, in [0, fdof), which are constrained 2484 2485 Not Collective 2486 2487 Input Parameters: 2488 + s - The PetscSection 2489 . field - The field number 2490 - point - The point 2491 2492 Output Parameter: 2493 . indices - The constrained dofs sorted in ascending order 2494 2495 Notes: 2496 The indices array, which is provided by the caller, must have capacity to hold the number of constrained dofs, e.g., as returned by PetscSectionGetConstraintDof(). 2497 2498 Fortran Note: 2499 In Fortran, you call PetscSectionGetFieldConstraintIndicesF90() and PetscSectionRestoreFieldConstraintIndicesF90() 2500 2501 Level: intermediate 2502 2503 .seealso: PetscSectionSetFieldConstraintIndices(), PetscSectionGetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection 2504 @*/ 2505 PetscErrorCode PetscSectionGetFieldConstraintIndices(PetscSection s, PetscInt point, PetscInt field, const PetscInt **indices) 2506 { 2507 PetscFunctionBegin; 2508 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2509 PetscValidPointer(indices,4); 2510 PetscSectionCheckValidField(field,s->numFields); 2511 PetscCall(PetscSectionGetConstraintIndices(s->field[field], point, indices)); 2512 PetscFunctionReturn(0); 2513 } 2514 2515 /*@C 2516 PetscSectionSetFieldConstraintIndices - Set the field dof numbers, in [0, fdof), which are constrained 2517 2518 Not Collective 2519 2520 Input Parameters: 2521 + s - The PetscSection 2522 . point - The point 2523 . field - The field number 2524 - indices - The constrained dofs 2525 2526 Note: The Fortran is PetscSectionSetFieldConstraintIndicesF90() 2527 2528 Level: intermediate 2529 2530 .seealso: PetscSectionSetConstraintIndices(), PetscSectionGetFieldConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection 2531 @*/ 2532 PetscErrorCode PetscSectionSetFieldConstraintIndices(PetscSection s, PetscInt point, PetscInt field, const PetscInt indices[]) 2533 { 2534 PetscFunctionBegin; 2535 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2536 if (PetscDefined(USE_DEBUG)) { 2537 PetscInt nfdof; 2538 2539 PetscCall(PetscSectionGetFieldConstraintDof(s, point, field, &nfdof)); 2540 if (nfdof) PetscValidIntPointer(indices, 4); 2541 } 2542 PetscSectionCheckValidField(field,s->numFields); 2543 PetscCall(PetscSectionSetConstraintIndices(s->field[field], point, indices)); 2544 PetscFunctionReturn(0); 2545 } 2546 2547 /*@ 2548 PetscSectionPermute - Reorder the section according to the input point permutation 2549 2550 Collective 2551 2552 Input Parameters: 2553 + section - The PetscSection object 2554 - perm - The point permutation, old point p becomes new point perm[p] 2555 2556 Output Parameter: 2557 . sectionNew - The permuted PetscSection 2558 2559 Level: intermediate 2560 2561 .seealso: MatPermute() 2562 @*/ 2563 PetscErrorCode PetscSectionPermute(PetscSection section, IS permutation, PetscSection *sectionNew) 2564 { 2565 PetscSection s = section, sNew; 2566 const PetscInt *perm; 2567 PetscInt numFields, f, c, numPoints, pStart, pEnd, p; 2568 2569 PetscFunctionBegin; 2570 PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1); 2571 PetscValidHeaderSpecific(permutation, IS_CLASSID, 2); 2572 PetscValidPointer(sectionNew, 3); 2573 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) s), &sNew)); 2574 PetscCall(PetscSectionGetNumFields(s, &numFields)); 2575 if (numFields) PetscCall(PetscSectionSetNumFields(sNew, numFields)); 2576 for (f = 0; f < numFields; ++f) { 2577 const char *name; 2578 PetscInt numComp; 2579 2580 PetscCall(PetscSectionGetFieldName(s, f, &name)); 2581 PetscCall(PetscSectionSetFieldName(sNew, f, name)); 2582 PetscCall(PetscSectionGetFieldComponents(s, f, &numComp)); 2583 PetscCall(PetscSectionSetFieldComponents(sNew, f, numComp)); 2584 for (c = 0; c < s->numFieldComponents[f]; ++c) { 2585 PetscCall(PetscSectionGetComponentName(s, f, c, &name)); 2586 PetscCall(PetscSectionSetComponentName(sNew, f, c, name)); 2587 } 2588 } 2589 PetscCall(ISGetLocalSize(permutation, &numPoints)); 2590 PetscCall(ISGetIndices(permutation, &perm)); 2591 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 2592 PetscCall(PetscSectionSetChart(sNew, pStart, pEnd)); 2593 PetscCheck(numPoints >= pEnd,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Permutation size %" PetscInt_FMT " is less than largest Section point %" PetscInt_FMT, numPoints, pEnd); 2594 for (p = pStart; p < pEnd; ++p) { 2595 PetscInt dof, cdof; 2596 2597 PetscCall(PetscSectionGetDof(s, p, &dof)); 2598 PetscCall(PetscSectionSetDof(sNew, perm[p], dof)); 2599 PetscCall(PetscSectionGetConstraintDof(s, p, &cdof)); 2600 if (cdof) PetscCall(PetscSectionSetConstraintDof(sNew, perm[p], cdof)); 2601 for (f = 0; f < numFields; ++f) { 2602 PetscCall(PetscSectionGetFieldDof(s, p, f, &dof)); 2603 PetscCall(PetscSectionSetFieldDof(sNew, perm[p], f, dof)); 2604 PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cdof)); 2605 if (cdof) PetscCall(PetscSectionSetFieldConstraintDof(sNew, perm[p], f, cdof)); 2606 } 2607 } 2608 PetscCall(PetscSectionSetUp(sNew)); 2609 for (p = pStart; p < pEnd; ++p) { 2610 const PetscInt *cind; 2611 PetscInt cdof; 2612 2613 PetscCall(PetscSectionGetConstraintDof(s, p, &cdof)); 2614 if (cdof) { 2615 PetscCall(PetscSectionGetConstraintIndices(s, p, &cind)); 2616 PetscCall(PetscSectionSetConstraintIndices(sNew, perm[p], cind)); 2617 } 2618 for (f = 0; f < numFields; ++f) { 2619 PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cdof)); 2620 if (cdof) { 2621 PetscCall(PetscSectionGetFieldConstraintIndices(s, p, f, &cind)); 2622 PetscCall(PetscSectionSetFieldConstraintIndices(sNew, perm[p], f, cind)); 2623 } 2624 } 2625 } 2626 PetscCall(ISRestoreIndices(permutation, &perm)); 2627 *sectionNew = sNew; 2628 PetscFunctionReturn(0); 2629 } 2630 2631 /*@ 2632 PetscSectionSetClosureIndex - Set a cache of points in the closure of each point in the section 2633 2634 Collective 2635 2636 Input Parameters: 2637 + section - The PetscSection 2638 . obj - A PetscObject which serves as the key for this index 2639 . clSection - Section giving the size of the closure of each point 2640 - clPoints - IS giving the points in each closure 2641 2642 Note: We compress out closure points with no dofs in this section 2643 2644 Level: advanced 2645 2646 .seealso: PetscSectionGetClosureIndex(), DMPlexCreateClosureIndex() 2647 @*/ 2648 PetscErrorCode PetscSectionSetClosureIndex(PetscSection section, PetscObject obj, PetscSection clSection, IS clPoints) 2649 { 2650 PetscFunctionBegin; 2651 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 2652 PetscValidHeaderSpecific(clSection,PETSC_SECTION_CLASSID,3); 2653 PetscValidHeaderSpecific(clPoints,IS_CLASSID,4); 2654 if (section->clObj != obj) PetscCall(PetscSectionResetClosurePermutation(section)); 2655 section->clObj = obj; 2656 PetscCall(PetscObjectReference((PetscObject)clSection)); 2657 PetscCall(PetscObjectReference((PetscObject)clPoints)); 2658 PetscCall(PetscSectionDestroy(§ion->clSection)); 2659 PetscCall(ISDestroy(§ion->clPoints)); 2660 section->clSection = clSection; 2661 section->clPoints = clPoints; 2662 PetscFunctionReturn(0); 2663 } 2664 2665 /*@ 2666 PetscSectionGetClosureIndex - Get the cache of points in the closure of each point in the section 2667 2668 Collective 2669 2670 Input Parameters: 2671 + section - The PetscSection 2672 - obj - A PetscObject which serves as the key for this index 2673 2674 Output Parameters: 2675 + clSection - Section giving the size of the closure of each point 2676 - clPoints - IS giving the points in each closure 2677 2678 Note: We compress out closure points with no dofs in this section 2679 2680 Level: advanced 2681 2682 .seealso: PetscSectionSetClosureIndex(), DMPlexCreateClosureIndex() 2683 @*/ 2684 PetscErrorCode PetscSectionGetClosureIndex(PetscSection section, PetscObject obj, PetscSection *clSection, IS *clPoints) 2685 { 2686 PetscFunctionBegin; 2687 if (section->clObj == obj) { 2688 if (clSection) *clSection = section->clSection; 2689 if (clPoints) *clPoints = section->clPoints; 2690 } else { 2691 if (clSection) *clSection = NULL; 2692 if (clPoints) *clPoints = NULL; 2693 } 2694 PetscFunctionReturn(0); 2695 } 2696 2697 PetscErrorCode PetscSectionSetClosurePermutation_Internal(PetscSection section, PetscObject obj, PetscInt depth, PetscInt clSize, PetscCopyMode mode, PetscInt *clPerm) 2698 { 2699 PetscInt i; 2700 khiter_t iter; 2701 int new_entry; 2702 PetscSectionClosurePermKey key = {depth, clSize}; 2703 PetscSectionClosurePermVal *val; 2704 2705 PetscFunctionBegin; 2706 if (section->clObj != obj) { 2707 PetscCall(PetscSectionDestroy(§ion->clSection)); 2708 PetscCall(ISDestroy(§ion->clPoints)); 2709 } 2710 section->clObj = obj; 2711 if (!section->clHash) PetscCall(PetscClPermCreate(§ion->clHash)); 2712 iter = kh_put(ClPerm, section->clHash, key, &new_entry); 2713 val = &kh_val(section->clHash, iter); 2714 if (!new_entry) { 2715 PetscCall(PetscFree(val->perm)); 2716 PetscCall(PetscFree(val->invPerm)); 2717 } 2718 if (mode == PETSC_COPY_VALUES) { 2719 PetscCall(PetscMalloc1(clSize, &val->perm)); 2720 PetscCall(PetscLogObjectMemory((PetscObject) obj, clSize*sizeof(PetscInt))); 2721 PetscCall(PetscArraycpy(val->perm, clPerm, clSize)); 2722 } else if (mode == PETSC_OWN_POINTER) { 2723 val->perm = clPerm; 2724 } else SETERRQ(PetscObjectComm(obj), PETSC_ERR_SUP, "Do not support borrowed arrays"); 2725 PetscCall(PetscMalloc1(clSize, &val->invPerm)); 2726 for (i = 0; i < clSize; ++i) val->invPerm[clPerm[i]] = i; 2727 PetscFunctionReturn(0); 2728 } 2729 2730 /*@ 2731 PetscSectionSetClosurePermutation - Set the dof permutation for the closure of each cell in the section, meaning clPerm[newIndex] = oldIndex. 2732 2733 Not Collective 2734 2735 Input Parameters: 2736 + section - The PetscSection 2737 . obj - A PetscObject which serves as the key for this index (usually a DM) 2738 . depth - Depth of points on which to apply the given permutation 2739 - perm - Permutation of the cell dof closure 2740 2741 Note: 2742 The specified permutation will only be applied to points at depth whose closure size matches the length of perm. In a 2743 mixed-topology or variable-degree finite element space, this function can be called multiple times at each depth for 2744 each topology and degree. 2745 2746 This approach assumes that (depth, len(perm)) uniquely identifies the desired permutation; this might not be true for 2747 exotic/enriched spaces on mixed topology meshes. 2748 2749 Level: intermediate 2750 2751 .seealso: PetscSectionGetClosurePermutation(), PetscSectionGetClosureIndex(), DMPlexCreateClosureIndex(), PetscCopyMode 2752 @*/ 2753 PetscErrorCode PetscSectionSetClosurePermutation(PetscSection section, PetscObject obj, PetscInt depth, IS perm) 2754 { 2755 const PetscInt *clPerm = NULL; 2756 PetscInt clSize = 0; 2757 2758 PetscFunctionBegin; 2759 if (perm) { 2760 PetscCall(ISGetLocalSize(perm, &clSize)); 2761 PetscCall(ISGetIndices(perm, &clPerm)); 2762 } 2763 PetscCall(PetscSectionSetClosurePermutation_Internal(section, obj, depth, clSize, PETSC_COPY_VALUES, (PetscInt *) clPerm)); 2764 if (perm) PetscCall(ISRestoreIndices(perm, &clPerm)); 2765 PetscFunctionReturn(0); 2766 } 2767 2768 PetscErrorCode PetscSectionGetClosurePermutation_Internal(PetscSection section, PetscObject obj, PetscInt depth, PetscInt size, const PetscInt *perm[]) 2769 { 2770 PetscFunctionBegin; 2771 if (section->clObj == obj) { 2772 PetscSectionClosurePermKey k = {depth, size}; 2773 PetscSectionClosurePermVal v; 2774 PetscCall(PetscClPermGet(section->clHash, k, &v)); 2775 if (perm) *perm = v.perm; 2776 } else { 2777 if (perm) *perm = NULL; 2778 } 2779 PetscFunctionReturn(0); 2780 } 2781 2782 /*@ 2783 PetscSectionGetClosurePermutation - Get the dof permutation for the closure of each cell in the section, meaning clPerm[newIndex] = oldIndex. 2784 2785 Not Collective 2786 2787 Input Parameters: 2788 + section - The PetscSection 2789 . obj - A PetscObject which serves as the key for this index (usually a DM) 2790 . depth - Depth stratum on which to obtain closure permutation 2791 - clSize - Closure size to be permuted (e.g., may vary with element topology and degree) 2792 2793 Output Parameter: 2794 . perm - The dof closure permutation 2795 2796 Note: 2797 The user must destroy the IS that is returned. 2798 2799 Level: intermediate 2800 2801 .seealso: PetscSectionSetClosurePermutation(), PetscSectionGetClosureInversePermutation(), PetscSectionGetClosureIndex(), PetscSectionSetClosureIndex(), DMPlexCreateClosureIndex() 2802 @*/ 2803 PetscErrorCode PetscSectionGetClosurePermutation(PetscSection section, PetscObject obj, PetscInt depth, PetscInt clSize, IS *perm) 2804 { 2805 const PetscInt *clPerm; 2806 2807 PetscFunctionBegin; 2808 PetscCall(PetscSectionGetClosurePermutation_Internal(section, obj, depth, clSize, &clPerm)); 2809 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, clSize, clPerm, PETSC_USE_POINTER, perm)); 2810 PetscFunctionReturn(0); 2811 } 2812 2813 PetscErrorCode PetscSectionGetClosureInversePermutation_Internal(PetscSection section, PetscObject obj, PetscInt depth, PetscInt size, const PetscInt *perm[]) 2814 { 2815 PetscFunctionBegin; 2816 if (section->clObj == obj && section->clHash) { 2817 PetscSectionClosurePermKey k = {depth, size}; 2818 PetscSectionClosurePermVal v; 2819 PetscCall(PetscClPermGet(section->clHash, k, &v)); 2820 if (perm) *perm = v.invPerm; 2821 } else { 2822 if (perm) *perm = NULL; 2823 } 2824 PetscFunctionReturn(0); 2825 } 2826 2827 /*@ 2828 PetscSectionGetClosureInversePermutation - Get the inverse dof permutation for the closure of each cell in the section, meaning clPerm[oldIndex] = newIndex. 2829 2830 Not Collective 2831 2832 Input Parameters: 2833 + section - The PetscSection 2834 . obj - A PetscObject which serves as the key for this index (usually a DM) 2835 . depth - Depth stratum on which to obtain closure permutation 2836 - clSize - Closure size to be permuted (e.g., may vary with element topology and degree) 2837 2838 Output Parameters: 2839 . perm - The dof closure permutation 2840 2841 Note: 2842 The user must destroy the IS that is returned. 2843 2844 Level: intermediate 2845 2846 .seealso: PetscSectionSetClosurePermutation(), PetscSectionGetClosureIndex(), PetscSectionSetClosureIndex(), DMPlexCreateClosureIndex() 2847 @*/ 2848 PetscErrorCode PetscSectionGetClosureInversePermutation(PetscSection section, PetscObject obj, PetscInt depth, PetscInt clSize, IS *perm) 2849 { 2850 const PetscInt *clPerm; 2851 2852 PetscFunctionBegin; 2853 PetscCall(PetscSectionGetClosureInversePermutation_Internal(section, obj, depth, clSize, &clPerm)); 2854 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, clSize, clPerm, PETSC_USE_POINTER, perm)); 2855 PetscFunctionReturn(0); 2856 } 2857 2858 /*@ 2859 PetscSectionGetField - Get the subsection associated with a single field 2860 2861 Input Parameters: 2862 + s - The PetscSection 2863 - field - The field number 2864 2865 Output Parameter: 2866 . subs - The subsection for the given field 2867 2868 Level: intermediate 2869 2870 .seealso: PetscSectionSetNumFields() 2871 @*/ 2872 PetscErrorCode PetscSectionGetField(PetscSection s, PetscInt field, PetscSection *subs) 2873 { 2874 PetscFunctionBegin; 2875 PetscValidHeaderSpecific(s,PETSC_SECTION_CLASSID,1); 2876 PetscValidPointer(subs,3); 2877 PetscSectionCheckValidField(field,s->numFields); 2878 *subs = s->field[field]; 2879 PetscFunctionReturn(0); 2880 } 2881 2882 PetscClassId PETSC_SECTION_SYM_CLASSID; 2883 PetscFunctionList PetscSectionSymList = NULL; 2884 2885 /*@ 2886 PetscSectionSymCreate - Creates an empty PetscSectionSym object. 2887 2888 Collective 2889 2890 Input Parameter: 2891 . comm - the MPI communicator 2892 2893 Output Parameter: 2894 . sym - pointer to the new set of symmetries 2895 2896 Level: developer 2897 2898 .seealso: PetscSectionSym, PetscSectionSymDestroy() 2899 @*/ 2900 PetscErrorCode PetscSectionSymCreate(MPI_Comm comm, PetscSectionSym *sym) 2901 { 2902 PetscFunctionBegin; 2903 PetscValidPointer(sym,2); 2904 PetscCall(ISInitializePackage()); 2905 PetscCall(PetscHeaderCreate(*sym,PETSC_SECTION_SYM_CLASSID,"PetscSectionSym","Section Symmetry","IS",comm,PetscSectionSymDestroy,PetscSectionSymView)); 2906 PetscFunctionReturn(0); 2907 } 2908 2909 /*@C 2910 PetscSectionSymSetType - Builds a PetscSection symmetry, for a particular implementation. 2911 2912 Collective 2913 2914 Input Parameters: 2915 + sym - The section symmetry object 2916 - method - The name of the section symmetry type 2917 2918 Level: developer 2919 2920 .seealso: PetscSectionSymGetType(), PetscSectionSymCreate() 2921 @*/ 2922 PetscErrorCode PetscSectionSymSetType(PetscSectionSym sym, PetscSectionSymType method) 2923 { 2924 PetscErrorCode (*r)(PetscSectionSym); 2925 PetscBool match; 2926 2927 PetscFunctionBegin; 2928 PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID,1); 2929 PetscCall(PetscObjectTypeCompare((PetscObject) sym, method, &match)); 2930 if (match) PetscFunctionReturn(0); 2931 2932 PetscCall(PetscFunctionListFind(PetscSectionSymList,method,&r)); 2933 PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscSectionSym type: %s", method); 2934 if (sym->ops->destroy) { 2935 PetscCall((*sym->ops->destroy)(sym)); 2936 sym->ops->destroy = NULL; 2937 } 2938 PetscCall((*r)(sym)); 2939 PetscCall(PetscObjectChangeTypeName((PetscObject)sym,method)); 2940 PetscFunctionReturn(0); 2941 } 2942 2943 /*@C 2944 PetscSectionSymGetType - Gets the section symmetry type name (as a string) from the PetscSectionSym. 2945 2946 Not Collective 2947 2948 Input Parameter: 2949 . sym - The section symmetry 2950 2951 Output Parameter: 2952 . type - The index set type name 2953 2954 Level: developer 2955 2956 .seealso: PetscSectionSymSetType(), PetscSectionSymCreate() 2957 @*/ 2958 PetscErrorCode PetscSectionSymGetType(PetscSectionSym sym, PetscSectionSymType *type) 2959 { 2960 PetscFunctionBegin; 2961 PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID,1); 2962 PetscValidPointer(type,2); 2963 *type = ((PetscObject)sym)->type_name; 2964 PetscFunctionReturn(0); 2965 } 2966 2967 /*@C 2968 PetscSectionSymRegister - Adds a new section symmetry implementation 2969 2970 Not Collective 2971 2972 Input Parameters: 2973 + name - The name of a new user-defined creation routine 2974 - create_func - The creation routine itself 2975 2976 Notes: 2977 PetscSectionSymRegister() may be called multiple times to add several user-defined vectors 2978 2979 Level: developer 2980 2981 .seealso: PetscSectionSymCreate(), PetscSectionSymSetType() 2982 @*/ 2983 PetscErrorCode PetscSectionSymRegister(const char sname[], PetscErrorCode (*function)(PetscSectionSym)) 2984 { 2985 PetscFunctionBegin; 2986 PetscCall(ISInitializePackage()); 2987 PetscCall(PetscFunctionListAdd(&PetscSectionSymList,sname,function)); 2988 PetscFunctionReturn(0); 2989 } 2990 2991 /*@ 2992 PetscSectionSymDestroy - Destroys a section symmetry. 2993 2994 Collective 2995 2996 Input Parameters: 2997 . sym - the section symmetry 2998 2999 Level: developer 3000 3001 .seealso: PetscSectionSymCreate(), PetscSectionSymDestroy() 3002 @*/ 3003 PetscErrorCode PetscSectionSymDestroy(PetscSectionSym *sym) 3004 { 3005 SymWorkLink link,next; 3006 3007 PetscFunctionBegin; 3008 if (!*sym) PetscFunctionReturn(0); 3009 PetscValidHeaderSpecific((*sym),PETSC_SECTION_SYM_CLASSID,1); 3010 if (--((PetscObject)(*sym))->refct > 0) {*sym = NULL; PetscFunctionReturn(0);} 3011 if ((*sym)->ops->destroy) { 3012 PetscCall((*(*sym)->ops->destroy)(*sym)); 3013 } 3014 PetscCheckFalse((*sym)->workout,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Work array still checked out"); 3015 for (link=(*sym)->workin; link; link=next) { 3016 next = link->next; 3017 PetscCall(PetscFree2(link->perms,link->rots)); 3018 PetscCall(PetscFree(link)); 3019 } 3020 (*sym)->workin = NULL; 3021 PetscCall(PetscHeaderDestroy(sym)); 3022 PetscFunctionReturn(0); 3023 } 3024 3025 /*@C 3026 PetscSectionSymView - Displays a section symmetry 3027 3028 Collective 3029 3030 Input Parameters: 3031 + sym - the index set 3032 - viewer - viewer used to display the set, for example PETSC_VIEWER_STDOUT_SELF. 3033 3034 Level: developer 3035 3036 .seealso: PetscViewerASCIIOpen() 3037 @*/ 3038 PetscErrorCode PetscSectionSymView(PetscSectionSym sym,PetscViewer viewer) 3039 { 3040 PetscFunctionBegin; 3041 PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1); 3042 if (!viewer) { 3043 PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sym),&viewer)); 3044 } 3045 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 3046 PetscCheckSameComm(sym,1,viewer,2); 3047 PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)sym,viewer)); 3048 if (sym->ops->view) { 3049 PetscCall((*sym->ops->view)(sym,viewer)); 3050 } 3051 PetscFunctionReturn(0); 3052 } 3053 3054 /*@ 3055 PetscSectionSetSym - Set the symmetries for the data referred to by the section 3056 3057 Collective 3058 3059 Input Parameters: 3060 + section - the section describing data layout 3061 - sym - the symmetry describing the affect of orientation on the access of the data 3062 3063 Level: developer 3064 3065 .seealso: PetscSectionGetSym(), PetscSectionSymCreate() 3066 @*/ 3067 PetscErrorCode PetscSectionSetSym(PetscSection section, PetscSectionSym sym) 3068 { 3069 PetscFunctionBegin; 3070 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3071 PetscCall(PetscSectionSymDestroy(&(section->sym))); 3072 if (sym) { 3073 PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,2); 3074 PetscCheckSameComm(section,1,sym,2); 3075 PetscCall(PetscObjectReference((PetscObject) sym)); 3076 } 3077 section->sym = sym; 3078 PetscFunctionReturn(0); 3079 } 3080 3081 /*@ 3082 PetscSectionGetSym - Get the symmetries for the data referred to by the section 3083 3084 Not Collective 3085 3086 Input Parameters: 3087 . section - the section describing data layout 3088 3089 Output Parameters: 3090 . sym - the symmetry describing the affect of orientation on the access of the data 3091 3092 Level: developer 3093 3094 .seealso: PetscSectionSetSym(), PetscSectionSymCreate() 3095 @*/ 3096 PetscErrorCode PetscSectionGetSym(PetscSection section, PetscSectionSym *sym) 3097 { 3098 PetscFunctionBegin; 3099 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3100 *sym = section->sym; 3101 PetscFunctionReturn(0); 3102 } 3103 3104 /*@ 3105 PetscSectionSetFieldSym - Set the symmetries for the data referred to by a field of the section 3106 3107 Collective 3108 3109 Input Parameters: 3110 + section - the section describing data layout 3111 . field - the field number 3112 - sym - the symmetry describing the affect of orientation on the access of the data 3113 3114 Level: developer 3115 3116 .seealso: PetscSectionGetFieldSym(), PetscSectionSymCreate() 3117 @*/ 3118 PetscErrorCode PetscSectionSetFieldSym(PetscSection section, PetscInt field, PetscSectionSym sym) 3119 { 3120 PetscFunctionBegin; 3121 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3122 PetscSectionCheckValidField(field,section->numFields); 3123 PetscCall(PetscSectionSetSym(section->field[field],sym)); 3124 PetscFunctionReturn(0); 3125 } 3126 3127 /*@ 3128 PetscSectionGetFieldSym - Get the symmetries for the data referred to by a field of the section 3129 3130 Collective 3131 3132 Input Parameters: 3133 + section - the section describing data layout 3134 - field - the field number 3135 3136 Output Parameters: 3137 . sym - the symmetry describing the affect of orientation on the access of the data 3138 3139 Level: developer 3140 3141 .seealso: PetscSectionSetFieldSym(), PetscSectionSymCreate() 3142 @*/ 3143 PetscErrorCode PetscSectionGetFieldSym(PetscSection section, PetscInt field, PetscSectionSym *sym) 3144 { 3145 PetscFunctionBegin; 3146 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3147 PetscSectionCheckValidField(field,section->numFields); 3148 *sym = section->field[field]->sym; 3149 PetscFunctionReturn(0); 3150 } 3151 3152 /*@C 3153 PetscSectionGetPointSyms - Get the symmetries for a set of points in a PetscSection under specific orientations. 3154 3155 Not Collective 3156 3157 Input Parameters: 3158 + section - the section 3159 . numPoints - the number of points 3160 - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an 3161 arbitrary integer: its interpretation is up to sym. Orientations are used by DM: for their interpretation in that 3162 context, see DMPlexGetConeOrientation()). 3163 3164 Output Parameters: 3165 + perms - The permutations for the given orientations (or NULL if there is no symmetry or the permutation is the identity). 3166 - rots - The field rotations symmetries for the given orientations (or NULL if there is no symmetry or the rotations are all 3167 identity). 3168 3169 Example of usage, gathering dofs into a local array (lArray) from a section array (sArray): 3170 .vb 3171 const PetscInt **perms; 3172 const PetscScalar **rots; 3173 PetscInt lOffset; 3174 3175 PetscSectionGetPointSyms(section,numPoints,points,&perms,&rots); 3176 for (i = 0, lOffset = 0; i < numPoints; i++) { 3177 PetscInt point = points[2*i], dof, sOffset; 3178 const PetscInt *perm = perms ? perms[i] : NULL; 3179 const PetscScalar *rot = rots ? rots[i] : NULL; 3180 3181 PetscSectionGetDof(section,point,&dof); 3182 PetscSectionGetOffset(section,point,&sOffset); 3183 3184 if (perm) {for (j = 0; j < dof; j++) {lArray[lOffset + perm[j]] = sArray[sOffset + j];}} 3185 else {for (j = 0; j < dof; j++) {lArray[lOffset + j ] = sArray[sOffset + j];}} 3186 if (rot) {for (j = 0; j < dof; j++) {lArray[lOffset + j ] *= rot[j]; }} 3187 lOffset += dof; 3188 } 3189 PetscSectionRestorePointSyms(section,numPoints,points,&perms,&rots); 3190 .ve 3191 3192 Example of usage, adding dofs into a section array (sArray) from a local array (lArray): 3193 .vb 3194 const PetscInt **perms; 3195 const PetscScalar **rots; 3196 PetscInt lOffset; 3197 3198 PetscSectionGetPointSyms(section,numPoints,points,&perms,&rots); 3199 for (i = 0, lOffset = 0; i < numPoints; i++) { 3200 PetscInt point = points[2*i], dof, sOffset; 3201 const PetscInt *perm = perms ? perms[i] : NULL; 3202 const PetscScalar *rot = rots ? rots[i] : NULL; 3203 3204 PetscSectionGetDof(section,point,&dof); 3205 PetscSectionGetOffset(section,point,&sOff); 3206 3207 if (perm) {for (j = 0; j < dof; j++) {sArray[sOffset + j] += lArray[lOffset + perm[j]] * (rot ? PetscConj(rot[perm[j]]) : 1.);}} 3208 else {for (j = 0; j < dof; j++) {sArray[sOffset + j] += lArray[lOffset + j ] * (rot ? PetscConj(rot[ j ]) : 1.);}} 3209 offset += dof; 3210 } 3211 PetscSectionRestorePointSyms(section,numPoints,points,&perms,&rots); 3212 .ve 3213 3214 Level: developer 3215 3216 .seealso: PetscSectionRestorePointSyms(), PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym() 3217 @*/ 3218 PetscErrorCode PetscSectionGetPointSyms(PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots) 3219 { 3220 PetscSectionSym sym; 3221 3222 PetscFunctionBegin; 3223 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3224 if (numPoints) PetscValidIntPointer(points,3); 3225 if (perms) *perms = NULL; 3226 if (rots) *rots = NULL; 3227 sym = section->sym; 3228 if (sym && (perms || rots)) { 3229 SymWorkLink link; 3230 3231 if (sym->workin) { 3232 link = sym->workin; 3233 sym->workin = sym->workin->next; 3234 } else { 3235 PetscCall(PetscNewLog(sym,&link)); 3236 } 3237 if (numPoints > link->numPoints) { 3238 PetscCall(PetscFree2(link->perms,link->rots)); 3239 PetscCall(PetscMalloc2(numPoints,&link->perms,numPoints,&link->rots)); 3240 link->numPoints = numPoints; 3241 } 3242 link->next = sym->workout; 3243 sym->workout = link; 3244 PetscCall(PetscArrayzero((PetscInt**)link->perms,numPoints)); 3245 PetscCall(PetscArrayzero((PetscInt**)link->rots,numPoints)); 3246 PetscCall((*sym->ops->getpoints) (sym, section, numPoints, points, link->perms, link->rots)); 3247 if (perms) *perms = link->perms; 3248 if (rots) *rots = link->rots; 3249 } 3250 PetscFunctionReturn(0); 3251 } 3252 3253 /*@C 3254 PetscSectionRestorePointSyms - Restore the symmetries returned by PetscSectionGetPointSyms() 3255 3256 Not Collective 3257 3258 Input Parameters: 3259 + section - the section 3260 . numPoints - the number of points 3261 - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an 3262 arbitrary integer: its interpretation is up to sym. Orientations are used by DM: for their interpretation in that 3263 context, see DMPlexGetConeOrientation()). 3264 3265 Output Parameters: 3266 + perms - The permutations for the given orientations: set to NULL at conclusion 3267 - rots - The field rotations symmetries for the given orientations: set to NULL at conclusion 3268 3269 Level: developer 3270 3271 .seealso: PetscSectionGetPointSyms(), PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym() 3272 @*/ 3273 PetscErrorCode PetscSectionRestorePointSyms(PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots) 3274 { 3275 PetscSectionSym sym; 3276 3277 PetscFunctionBegin; 3278 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3279 sym = section->sym; 3280 if (sym && (perms || rots)) { 3281 SymWorkLink *p,link; 3282 3283 for (p=&sym->workout; (link=*p); p=&link->next) { 3284 if ((perms && link->perms == *perms) || (rots && link->rots == *rots)) { 3285 *p = link->next; 3286 link->next = sym->workin; 3287 sym->workin = link; 3288 if (perms) *perms = NULL; 3289 if (rots) *rots = NULL; 3290 PetscFunctionReturn(0); 3291 } 3292 } 3293 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Array was not checked out"); 3294 } 3295 PetscFunctionReturn(0); 3296 } 3297 3298 /*@C 3299 PetscSectionGetFieldPointSyms - Get the symmetries for a set of points in a field of a PetscSection under specific orientations. 3300 3301 Not Collective 3302 3303 Input Parameters: 3304 + section - the section 3305 . field - the field of the section 3306 . numPoints - the number of points 3307 - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an 3308 arbitrary integer: its interpretation is up to sym. Orientations are used by DM: for their interpretation in that 3309 context, see DMPlexGetConeOrientation()). 3310 3311 Output Parameters: 3312 + perms - The permutations for the given orientations (or NULL if there is no symmetry or the permutation is the identity). 3313 - rots - The field rotations symmetries for the given orientations (or NULL if there is no symmetry or the rotations are all 3314 identity). 3315 3316 Level: developer 3317 3318 .seealso: PetscSectionGetPointSyms(), PetscSectionRestoreFieldPointSyms() 3319 @*/ 3320 PetscErrorCode PetscSectionGetFieldPointSyms(PetscSection section, PetscInt field, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots) 3321 { 3322 PetscFunctionBegin; 3323 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3324 PetscCheck(field <= section->numFields,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"field %" PetscInt_FMT " greater than number of fields (%" PetscInt_FMT ") in section",field,section->numFields); 3325 PetscCall(PetscSectionGetPointSyms(section->field[field],numPoints,points,perms,rots)); 3326 PetscFunctionReturn(0); 3327 } 3328 3329 /*@C 3330 PetscSectionRestoreFieldPointSyms - Restore the symmetries returned by PetscSectionGetFieldPointSyms() 3331 3332 Not Collective 3333 3334 Input Parameters: 3335 + section - the section 3336 . field - the field number 3337 . numPoints - the number of points 3338 - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an 3339 arbitrary integer: its interpretation is up to sym. Orientations are used by DM: for their interpretation in that 3340 context, see DMPlexGetConeOrientation()). 3341 3342 Output Parameters: 3343 + perms - The permutations for the given orientations: set to NULL at conclusion 3344 - rots - The field rotations symmetries for the given orientations: set to NULL at conclusion 3345 3346 Level: developer 3347 3348 .seealso: PetscSectionRestorePointSyms(), petscSectionGetFieldPointSyms(), PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym() 3349 @*/ 3350 PetscErrorCode PetscSectionRestoreFieldPointSyms(PetscSection section, PetscInt field, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots) 3351 { 3352 PetscFunctionBegin; 3353 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3354 PetscCheck(field <= section->numFields,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"field %" PetscInt_FMT " greater than number of fields (%" PetscInt_FMT ") in section",field,section->numFields); 3355 PetscCall(PetscSectionRestorePointSyms(section->field[field],numPoints,points,perms,rots)); 3356 PetscFunctionReturn(0); 3357 } 3358 3359 /*@ 3360 PetscSectionSymCopy - Copy the symmetries, assuming that the point structure is compatible 3361 3362 Not Collective 3363 3364 Input Parameter: 3365 . sym - the PetscSectionSym 3366 3367 Output Parameter: 3368 . nsym - the equivalent symmetries 3369 3370 Level: developer 3371 3372 .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms() 3373 @*/ 3374 PetscErrorCode PetscSectionSymCopy(PetscSectionSym sym, PetscSectionSym nsym) 3375 { 3376 PetscFunctionBegin; 3377 PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 1); 3378 PetscValidHeaderSpecific(nsym, PETSC_SECTION_SYM_CLASSID, 2); 3379 if (sym->ops->copy) PetscCall((*sym->ops->copy)(sym, nsym)); 3380 PetscFunctionReturn(0); 3381 } 3382 3383 /*@ 3384 PetscSectionSymDistribute - Distribute the symmetries in accordance with the input SF 3385 3386 Collective 3387 3388 Input Parameters: 3389 + sym - the PetscSectionSym 3390 - migrationSF - the distribution map from roots to leaves 3391 3392 Output Parameters: 3393 . dsym - the redistributed symmetries 3394 3395 Level: developer 3396 3397 .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms() 3398 @*/ 3399 PetscErrorCode PetscSectionSymDistribute(PetscSectionSym sym, PetscSF migrationSF, PetscSectionSym *dsym) 3400 { 3401 PetscFunctionBegin; 3402 PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 1); 3403 PetscValidHeaderSpecific(migrationSF, PETSCSF_CLASSID, 2); 3404 PetscValidPointer(dsym, 3); 3405 if (sym->ops->distribute) PetscCall((*sym->ops->distribute)(sym, migrationSF, dsym)); 3406 PetscFunctionReturn(0); 3407 } 3408 3409 /*@ 3410 PetscSectionGetUseFieldOffsets - Get the flag to use field offsets directly in a global section, rather than just the point offset 3411 3412 Not Collective 3413 3414 Input Parameter: 3415 . s - the global PetscSection 3416 3417 Output Parameters: 3418 . flg - the flag 3419 3420 Level: developer 3421 3422 .seealso: PetscSectionSetChart(), PetscSectionCreate() 3423 @*/ 3424 PetscErrorCode PetscSectionGetUseFieldOffsets(PetscSection s, PetscBool *flg) 3425 { 3426 PetscFunctionBegin; 3427 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 3428 *flg = s->useFieldOff; 3429 PetscFunctionReturn(0); 3430 } 3431 3432 /*@ 3433 PetscSectionSetUseFieldOffsets - Set the flag to use field offsets directly in a global section, rather than just the point offset 3434 3435 Not Collective 3436 3437 Input Parameters: 3438 + s - the global PetscSection 3439 - flg - the flag 3440 3441 Level: developer 3442 3443 .seealso: PetscSectionGetUseFieldOffsets(), PetscSectionSetChart(), PetscSectionCreate() 3444 @*/ 3445 PetscErrorCode PetscSectionSetUseFieldOffsets(PetscSection s, PetscBool flg) 3446 { 3447 PetscFunctionBegin; 3448 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 3449 s->useFieldOff = flg; 3450 PetscFunctionReturn(0); 3451 } 3452 3453 #define PetscSectionExpandPoints_Loop(TYPE) \ 3454 { \ 3455 PetscInt i, n, o0, o1, size; \ 3456 TYPE *a0 = (TYPE*)origArray, *a1; \ 3457 PetscCall(PetscSectionGetStorageSize(s, &size)); \ 3458 PetscCall(PetscMalloc1(size, &a1)); \ 3459 for (i=0; i<npoints; i++) { \ 3460 PetscCall(PetscSectionGetOffset(origSection, points_[i], &o0)); \ 3461 PetscCall(PetscSectionGetOffset(s, i, &o1)); \ 3462 PetscCall(PetscSectionGetDof(s, i, &n)); \ 3463 PetscCall(PetscMemcpy(&a1[o1], &a0[o0], n*unitsize)); \ 3464 } \ 3465 *newArray = (void*)a1; \ 3466 } 3467 3468 /*@ 3469 PetscSectionExtractDofsFromArray - Extracts elements of an array corresponding to DOFs of specified points. 3470 3471 Not Collective 3472 3473 Input Parameters: 3474 + origSection - the PetscSection describing the layout of the array 3475 . dataType - MPI_Datatype describing the data type of the array (currently only MPIU_INT, MPIU_SCALAR, MPIU_REAL) 3476 . origArray - the array; its size must be equal to the storage size of origSection 3477 - points - IS with points to extract; its indices must lie in the chart of origSection 3478 3479 Output Parameters: 3480 + newSection - the new PetscSection desribing the layout of the new array (with points renumbered 0,1,... but preserving numbers of DOFs) 3481 - newArray - the array of the extracted DOFs; its size is the storage size of newSection 3482 3483 Level: developer 3484 3485 .seealso: PetscSectionGetChart(), PetscSectionGetDof(), PetscSectionGetStorageSize(), PetscSectionCreate() 3486 @*/ 3487 PetscErrorCode PetscSectionExtractDofsFromArray(PetscSection origSection, MPI_Datatype dataType, const void *origArray, IS points, PetscSection *newSection, void *newArray[]) 3488 { 3489 PetscSection s; 3490 const PetscInt *points_; 3491 PetscInt i, n, npoints, pStart, pEnd; 3492 PetscMPIInt unitsize; 3493 3494 PetscFunctionBegin; 3495 PetscValidHeaderSpecific(origSection, PETSC_SECTION_CLASSID, 1); 3496 PetscValidPointer(origArray, 3); 3497 PetscValidHeaderSpecific(points, IS_CLASSID, 4); 3498 if (newSection) PetscValidPointer(newSection, 5); 3499 if (newArray) PetscValidPointer(newArray, 6); 3500 PetscCallMPI(MPI_Type_size(dataType, &unitsize)); 3501 PetscCall(ISGetLocalSize(points, &npoints)); 3502 PetscCall(ISGetIndices(points, &points_)); 3503 PetscCall(PetscSectionGetChart(origSection, &pStart, &pEnd)); 3504 PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &s)); 3505 PetscCall(PetscSectionSetChart(s, 0, npoints)); 3506 for (i=0; i<npoints; i++) { 3507 PetscCheckFalse(points_[i] < pStart || points_[i] >= pEnd,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "point %" PetscInt_FMT " (index %" PetscInt_FMT ") in input IS out of input section's chart", points_[i], i); 3508 PetscCall(PetscSectionGetDof(origSection, points_[i], &n)); 3509 PetscCall(PetscSectionSetDof(s, i, n)); 3510 } 3511 PetscCall(PetscSectionSetUp(s)); 3512 if (newArray) { 3513 if (dataType == MPIU_INT) {PetscSectionExpandPoints_Loop(PetscInt);} 3514 else if (dataType == MPIU_SCALAR) {PetscSectionExpandPoints_Loop(PetscScalar);} 3515 else if (dataType == MPIU_REAL) {PetscSectionExpandPoints_Loop(PetscReal);} 3516 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "not implemented for this MPI_Datatype"); 3517 } 3518 if (newSection) { 3519 *newSection = s; 3520 } else { 3521 PetscCall(PetscSectionDestroy(&s)); 3522 } 3523 PetscCall(ISRestoreIndices(points, &points_)); 3524 PetscFunctionReturn(0); 3525 } 3526