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)->maxDof = 0; 50 (*s)->atlasDof = NULL; 51 (*s)->atlasOff = NULL; 52 (*s)->bc = NULL; 53 (*s)->bcIndices = NULL; 54 (*s)->setup = PETSC_FALSE; 55 (*s)->numFields = 0; 56 (*s)->fieldNames = NULL; 57 (*s)->field = NULL; 58 (*s)->useFieldOff = PETSC_FALSE; 59 (*s)->compNames = NULL; 60 (*s)->clObj = NULL; 61 (*s)->clHash = NULL; 62 (*s)->clSection = NULL; 63 (*s)->clPoints = NULL; 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 PetscCallMPI(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 PetscCheckFalse(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 PetscFunctionReturn(0); 804 } 805 806 /*@ 807 PetscSectionAddDof - Adds to the number of degrees of freedom associated with a given point. 808 809 Not Collective 810 811 Input Parameters: 812 + s - the PetscSection 813 . point - the point 814 - numDof - the number of additional dof 815 816 Level: intermediate 817 818 .seealso: PetscSectionGetDof(), PetscSectionSetDof(), PetscSectionCreate() 819 @*/ 820 PetscErrorCode PetscSectionAddDof(PetscSection s, PetscInt point, PetscInt numDof) 821 { 822 PetscFunctionBeginHot; 823 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 824 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); 825 s->atlasDof[point - s->pStart] += numDof; 826 PetscFunctionReturn(0); 827 } 828 829 /*@ 830 PetscSectionGetFieldDof - Return the number of degrees of freedom associated with a field on a given point. 831 832 Not Collective 833 834 Input Parameters: 835 + s - the PetscSection 836 . point - the point 837 - field - the field 838 839 Output Parameter: 840 . numDof - the number of dof 841 842 Level: intermediate 843 844 .seealso: PetscSectionSetFieldDof(), PetscSectionCreate() 845 @*/ 846 PetscErrorCode PetscSectionGetFieldDof(PetscSection s, PetscInt point, PetscInt field, PetscInt *numDof) 847 { 848 PetscFunctionBegin; 849 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 850 PetscValidIntPointer(numDof,4); 851 PetscSectionCheckValidField(field,s->numFields); 852 PetscCall(PetscSectionGetDof(s->field[field], point, numDof)); 853 PetscFunctionReturn(0); 854 } 855 856 /*@ 857 PetscSectionSetFieldDof - Sets the number of degrees of freedom associated with a field on a given point. 858 859 Not Collective 860 861 Input Parameters: 862 + s - the PetscSection 863 . point - the point 864 . field - the field 865 - numDof - the number of dof 866 867 Level: intermediate 868 869 .seealso: PetscSectionGetFieldDof(), PetscSectionCreate() 870 @*/ 871 PetscErrorCode PetscSectionSetFieldDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof) 872 { 873 PetscFunctionBegin; 874 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 875 PetscSectionCheckValidField(field,s->numFields); 876 PetscCall(PetscSectionSetDof(s->field[field], point, numDof)); 877 PetscFunctionReturn(0); 878 } 879 880 /*@ 881 PetscSectionAddFieldDof - Adds a number of degrees of freedom associated with a field on a given point. 882 883 Not Collective 884 885 Input Parameters: 886 + s - the PetscSection 887 . point - the point 888 . field - the field 889 - numDof - the number of dof 890 891 Level: intermediate 892 893 .seealso: PetscSectionSetFieldDof(), PetscSectionGetFieldDof(), PetscSectionCreate() 894 @*/ 895 PetscErrorCode PetscSectionAddFieldDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof) 896 { 897 PetscFunctionBegin; 898 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 899 PetscSectionCheckValidField(field,s->numFields); 900 PetscCall(PetscSectionAddDof(s->field[field], point, numDof)); 901 PetscFunctionReturn(0); 902 } 903 904 /*@ 905 PetscSectionGetConstraintDof - Return the number of constrained degrees of freedom associated with a given point. 906 907 Not Collective 908 909 Input Parameters: 910 + s - the PetscSection 911 - point - the point 912 913 Output Parameter: 914 . numDof - the number of dof which are fixed by constraints 915 916 Level: intermediate 917 918 .seealso: PetscSectionGetDof(), PetscSectionSetConstraintDof(), PetscSectionCreate() 919 @*/ 920 PetscErrorCode PetscSectionGetConstraintDof(PetscSection s, PetscInt point, PetscInt *numDof) 921 { 922 PetscFunctionBegin; 923 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 924 PetscValidIntPointer(numDof, 3); 925 if (s->bc) { 926 PetscCall(PetscSectionGetDof(s->bc, point, numDof)); 927 } else *numDof = 0; 928 PetscFunctionReturn(0); 929 } 930 931 /*@ 932 PetscSectionSetConstraintDof - Set the number of constrained degrees of freedom associated with a given point. 933 934 Not Collective 935 936 Input Parameters: 937 + s - the PetscSection 938 . point - the point 939 - numDof - the number of dof which are fixed by constraints 940 941 Level: intermediate 942 943 .seealso: PetscSectionSetDof(), PetscSectionGetConstraintDof(), PetscSectionCreate() 944 @*/ 945 PetscErrorCode PetscSectionSetConstraintDof(PetscSection s, PetscInt point, PetscInt numDof) 946 { 947 PetscFunctionBegin; 948 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 949 if (numDof) { 950 PetscCall(PetscSectionCheckConstraints_Private(s)); 951 PetscCall(PetscSectionSetDof(s->bc, point, numDof)); 952 } 953 PetscFunctionReturn(0); 954 } 955 956 /*@ 957 PetscSectionAddConstraintDof - Increment the number of constrained degrees of freedom associated with a given point. 958 959 Not Collective 960 961 Input Parameters: 962 + s - the PetscSection 963 . point - the point 964 - numDof - the number of additional dof which are fixed by constraints 965 966 Level: intermediate 967 968 .seealso: PetscSectionAddDof(), PetscSectionGetConstraintDof(), PetscSectionCreate() 969 @*/ 970 PetscErrorCode PetscSectionAddConstraintDof(PetscSection s, PetscInt point, PetscInt numDof) 971 { 972 PetscFunctionBegin; 973 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 974 if (numDof) { 975 PetscCall(PetscSectionCheckConstraints_Private(s)); 976 PetscCall(PetscSectionAddDof(s->bc, point, numDof)); 977 } 978 PetscFunctionReturn(0); 979 } 980 981 /*@ 982 PetscSectionGetFieldConstraintDof - Return the number of constrained degrees of freedom associated with a given field on a point. 983 984 Not Collective 985 986 Input Parameters: 987 + s - the PetscSection 988 . point - the point 989 - field - the field 990 991 Output Parameter: 992 . numDof - the number of dof which are fixed by constraints 993 994 Level: intermediate 995 996 .seealso: PetscSectionGetDof(), PetscSectionSetFieldConstraintDof(), PetscSectionCreate() 997 @*/ 998 PetscErrorCode PetscSectionGetFieldConstraintDof(PetscSection s, PetscInt point, PetscInt field, PetscInt *numDof) 999 { 1000 PetscFunctionBegin; 1001 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1002 PetscValidIntPointer(numDof, 4); 1003 PetscSectionCheckValidField(field,s->numFields); 1004 PetscCall(PetscSectionGetConstraintDof(s->field[field], point, numDof)); 1005 PetscFunctionReturn(0); 1006 } 1007 1008 /*@ 1009 PetscSectionSetFieldConstraintDof - Set the number of constrained degrees of freedom associated with a given field on a point. 1010 1011 Not Collective 1012 1013 Input Parameters: 1014 + s - the PetscSection 1015 . point - the point 1016 . field - the field 1017 - numDof - the number of dof which are fixed by constraints 1018 1019 Level: intermediate 1020 1021 .seealso: PetscSectionSetDof(), PetscSectionGetFieldConstraintDof(), PetscSectionCreate() 1022 @*/ 1023 PetscErrorCode PetscSectionSetFieldConstraintDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof) 1024 { 1025 PetscFunctionBegin; 1026 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1027 PetscSectionCheckValidField(field,s->numFields); 1028 PetscCall(PetscSectionSetConstraintDof(s->field[field], point, numDof)); 1029 PetscFunctionReturn(0); 1030 } 1031 1032 /*@ 1033 PetscSectionAddFieldConstraintDof - Increment the number of constrained degrees of freedom associated with a given field on a point. 1034 1035 Not Collective 1036 1037 Input Parameters: 1038 + s - the PetscSection 1039 . point - the point 1040 . field - the field 1041 - numDof - the number of additional dof which are fixed by constraints 1042 1043 Level: intermediate 1044 1045 .seealso: PetscSectionAddDof(), PetscSectionGetFieldConstraintDof(), PetscSectionCreate() 1046 @*/ 1047 PetscErrorCode PetscSectionAddFieldConstraintDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof) 1048 { 1049 PetscFunctionBegin; 1050 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1051 PetscSectionCheckValidField(field,s->numFields); 1052 PetscCall(PetscSectionAddConstraintDof(s->field[field], point, numDof)); 1053 PetscFunctionReturn(0); 1054 } 1055 1056 /*@ 1057 PetscSectionSetUpBC - Setup the subsections describing boundary conditions. 1058 1059 Not Collective 1060 1061 Input Parameter: 1062 . s - the PetscSection 1063 1064 Level: advanced 1065 1066 .seealso: PetscSectionSetUp(), PetscSectionCreate() 1067 @*/ 1068 PetscErrorCode PetscSectionSetUpBC(PetscSection s) 1069 { 1070 PetscFunctionBegin; 1071 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1072 if (s->bc) { 1073 const PetscInt last = (s->bc->pEnd-s->bc->pStart) - 1; 1074 1075 PetscCall(PetscSectionSetUp(s->bc)); 1076 PetscCall(PetscMalloc1(last >= 0 ? s->bc->atlasOff[last] + s->bc->atlasDof[last] : 0, &s->bcIndices)); 1077 } 1078 PetscFunctionReturn(0); 1079 } 1080 1081 /*@ 1082 PetscSectionSetUp - Calculate offsets based upon the number of degrees of freedom for each point. 1083 1084 Not Collective 1085 1086 Input Parameter: 1087 . s - the PetscSection 1088 1089 Level: intermediate 1090 1091 .seealso: PetscSectionCreate() 1092 @*/ 1093 PetscErrorCode PetscSectionSetUp(PetscSection s) 1094 { 1095 const PetscInt *pind = NULL; 1096 PetscInt offset = 0, foff, p, f; 1097 1098 PetscFunctionBegin; 1099 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1100 if (s->setup) PetscFunctionReturn(0); 1101 s->setup = PETSC_TRUE; 1102 /* Set offsets and field offsets for all points */ 1103 /* Assume that all fields have the same chart */ 1104 PetscCheck(s->includesConstraints,PETSC_COMM_SELF,PETSC_ERR_SUP,"PetscSectionSetUp is currently unsupported for includesConstraints = PETSC_TRUE"); 1105 if (s->perm) PetscCall(ISGetIndices(s->perm, &pind)); 1106 if (s->pointMajor) { 1107 for (p = 0; p < s->pEnd - s->pStart; ++p) { 1108 const PetscInt q = pind ? pind[p] : p; 1109 1110 /* Set point offset */ 1111 s->atlasOff[q] = offset; 1112 offset += s->atlasDof[q]; 1113 s->maxDof = PetscMax(s->maxDof, s->atlasDof[q]); 1114 /* Set field offset */ 1115 for (f = 0, foff = s->atlasOff[q]; f < s->numFields; ++f) { 1116 PetscSection sf = s->field[f]; 1117 1118 sf->atlasOff[q] = foff; 1119 foff += sf->atlasDof[q]; 1120 } 1121 } 1122 } else { 1123 /* Set field offsets for all points */ 1124 for (f = 0; f < s->numFields; ++f) { 1125 PetscSection sf = s->field[f]; 1126 1127 for (p = 0; p < s->pEnd - s->pStart; ++p) { 1128 const PetscInt q = pind ? pind[p] : p; 1129 1130 sf->atlasOff[q] = offset; 1131 offset += sf->atlasDof[q]; 1132 } 1133 } 1134 /* Disable point offsets since these are unused */ 1135 for (p = 0; p < s->pEnd - s->pStart; ++p) { 1136 s->atlasOff[p] = -1; 1137 s->maxDof = PetscMax(s->maxDof, s->atlasDof[p]); 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 .seealso: PetscSectionGetDof(), PetscSectionSetDof(), PetscSectionCreate() 1161 @*/ 1162 PetscErrorCode PetscSectionGetMaxDof(PetscSection s, PetscInt *maxDof) 1163 { 1164 PetscFunctionBegin; 1165 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1166 PetscValidIntPointer(maxDof, 2); 1167 *maxDof = s->maxDof; 1168 PetscFunctionReturn(0); 1169 } 1170 1171 /*@ 1172 PetscSectionGetStorageSize - Return the size of an array or local Vec capable of holding all the degrees of freedom. 1173 1174 Not Collective 1175 1176 Input Parameter: 1177 . s - the PetscSection 1178 1179 Output Parameter: 1180 . size - the size of an array which can hold all the dofs 1181 1182 Level: intermediate 1183 1184 .seealso: PetscSectionGetOffset(), PetscSectionGetConstrainedStorageSize(), PetscSectionCreate() 1185 @*/ 1186 PetscErrorCode PetscSectionGetStorageSize(PetscSection s, PetscInt *size) 1187 { 1188 PetscInt p, n = 0; 1189 1190 PetscFunctionBegin; 1191 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1192 PetscValidIntPointer(size, 2); 1193 for (p = 0; p < s->pEnd - s->pStart; ++p) n += s->atlasDof[p] > 0 ? s->atlasDof[p] : 0; 1194 *size = n; 1195 PetscFunctionReturn(0); 1196 } 1197 1198 /*@ 1199 PetscSectionGetConstrainedStorageSize - Return the size of an array or local Vec capable of holding all unconstrained degrees of freedom. 1200 1201 Not Collective 1202 1203 Input Parameter: 1204 . s - the PetscSection 1205 1206 Output Parameter: 1207 . size - the size of an array which can hold all unconstrained dofs 1208 1209 Level: intermediate 1210 1211 .seealso: PetscSectionGetStorageSize(), PetscSectionGetOffset(), PetscSectionCreate() 1212 @*/ 1213 PetscErrorCode PetscSectionGetConstrainedStorageSize(PetscSection s, PetscInt *size) 1214 { 1215 PetscInt p, n = 0; 1216 1217 PetscFunctionBegin; 1218 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1219 PetscValidIntPointer(size, 2); 1220 for (p = 0; p < s->pEnd - s->pStart; ++p) { 1221 const PetscInt cdof = s->bc ? s->bc->atlasDof[p] : 0; 1222 n += s->atlasDof[p] > 0 ? s->atlasDof[p] - cdof : 0; 1223 } 1224 *size = n; 1225 PetscFunctionReturn(0); 1226 } 1227 1228 /*@ 1229 PetscSectionCreateGlobalSection - Create a section describing the global field layout using 1230 the local section and an SF describing the section point overlap. 1231 1232 Input Parameters: 1233 + s - The PetscSection for the local field layout 1234 . sf - The SF describing parallel layout of the section points (leaves are unowned local points) 1235 . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs 1236 - localOffsets - If PETSC_TRUE, use local rather than global offsets for the points 1237 1238 Output Parameter: 1239 . gsection - The PetscSection for the global field layout 1240 1241 Note: This gives negative sizes and offsets to points not owned by this process 1242 1243 Level: intermediate 1244 1245 .seealso: PetscSectionCreate() 1246 @*/ 1247 PetscErrorCode PetscSectionCreateGlobalSection(PetscSection s, PetscSF sf, PetscBool includeConstraints, PetscBool localOffsets, PetscSection *gsection) 1248 { 1249 PetscSection gs; 1250 const PetscInt *pind = NULL; 1251 PetscInt *recv = NULL, *neg = NULL; 1252 PetscInt pStart, pEnd, p, dof, cdof, off, foff, globalOff = 0, nroots, nlocal, maxleaf; 1253 PetscInt numFields, f, numComponents; 1254 1255 PetscFunctionBegin; 1256 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1257 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 1258 PetscValidLogicalCollectiveBool(s, includeConstraints, 3); 1259 PetscValidLogicalCollectiveBool(s, localOffsets, 4); 1260 PetscValidPointer(gsection, 5); 1261 PetscCheck(s->pointMajor,PETSC_COMM_SELF,PETSC_ERR_SUP, "No support for field major ordering"); 1262 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) s), &gs)); 1263 PetscCall(PetscSectionGetNumFields(s, &numFields)); 1264 if (numFields > 0) PetscCall(PetscSectionSetNumFields(gs, numFields)); 1265 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 1266 PetscCall(PetscSectionSetChart(gs, pStart, pEnd)); 1267 gs->includesConstraints = includeConstraints; 1268 PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 1269 nlocal = nroots; /* The local/leaf space matches global/root space */ 1270 /* Must allocate for all points visible to SF, which may be more than this section */ 1271 if (nroots >= 0) { /* nroots < 0 means that the graph has not been set, only happens in serial */ 1272 PetscCall(PetscSFGetLeafRange(sf, NULL, &maxleaf)); 1273 PetscCheckFalse(nroots < pEnd,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "SF roots %" PetscInt_FMT " < pEnd %" PetscInt_FMT, nroots, pEnd); 1274 PetscCheckFalse(maxleaf >= nroots,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Max local leaf %" PetscInt_FMT " >= nroots %" PetscInt_FMT, maxleaf, nroots); 1275 PetscCall(PetscMalloc2(nroots,&neg,nlocal,&recv)); 1276 PetscCall(PetscArrayzero(neg,nroots)); 1277 } 1278 /* Mark all local points with negative dof */ 1279 for (p = pStart; p < pEnd; ++p) { 1280 PetscCall(PetscSectionGetDof(s, p, &dof)); 1281 PetscCall(PetscSectionSetDof(gs, p, dof)); 1282 PetscCall(PetscSectionGetConstraintDof(s, p, &cdof)); 1283 if (!includeConstraints && cdof > 0) PetscCall(PetscSectionSetConstraintDof(gs, p, cdof)); 1284 if (neg) neg[p] = -(dof+1); 1285 } 1286 PetscCall(PetscSectionSetUpBC(gs)); 1287 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])); 1288 if (nroots >= 0) { 1289 PetscCall(PetscArrayzero(recv,nlocal)); 1290 PetscCall(PetscSFBcastBegin(sf, MPIU_INT, neg, recv,MPI_REPLACE)); 1291 PetscCall(PetscSFBcastEnd(sf, MPIU_INT, neg, recv,MPI_REPLACE)); 1292 for (p = pStart; p < pEnd; ++p) { 1293 if (recv[p] < 0) { 1294 gs->atlasDof[p-pStart] = recv[p]; 1295 PetscCall(PetscSectionGetDof(s, p, &dof)); 1296 PetscCheckFalse(-(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); 1297 } 1298 } 1299 } 1300 /* Calculate new sizes, get process offset, and calculate point offsets */ 1301 if (s->perm) PetscCall(ISGetIndices(s->perm, &pind)); 1302 for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1303 const PetscInt q = pind ? pind[p] : p; 1304 1305 cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[q] : 0; 1306 gs->atlasOff[q] = off; 1307 off += gs->atlasDof[q] > 0 ? gs->atlasDof[q]-cdof : 0; 1308 } 1309 if (!localOffsets) { 1310 PetscCallMPI(MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s))); 1311 globalOff -= off; 1312 } 1313 for (p = pStart, off = 0; p < pEnd; ++p) { 1314 gs->atlasOff[p-pStart] += globalOff; 1315 if (neg) neg[p] = -(gs->atlasOff[p-pStart]+1); 1316 } 1317 if (s->perm) PetscCall(ISRestoreIndices(s->perm, &pind)); 1318 /* Put in negative offsets for ghost points */ 1319 if (nroots >= 0) { 1320 PetscCall(PetscArrayzero(recv,nlocal)); 1321 PetscCall(PetscSFBcastBegin(sf, MPIU_INT, neg, recv,MPI_REPLACE)); 1322 PetscCall(PetscSFBcastEnd(sf, MPIU_INT, neg, recv,MPI_REPLACE)); 1323 for (p = pStart; p < pEnd; ++p) { 1324 if (recv[p] < 0) gs->atlasOff[p-pStart] = recv[p]; 1325 } 1326 } 1327 PetscCall(PetscFree2(neg,recv)); 1328 /* Set field dofs/offsets/constraints */ 1329 for (f = 0; f < numFields; ++f) { 1330 gs->field[f]->includesConstraints = includeConstraints; 1331 PetscCall(PetscSectionGetFieldComponents(s, f, &numComponents)); 1332 PetscCall(PetscSectionSetFieldComponents(gs, f, numComponents)); 1333 } 1334 for (p = pStart; p < pEnd; ++p) { 1335 PetscCall(PetscSectionGetOffset(gs, p, &off)); 1336 for (f = 0, foff = off; f < numFields; ++f) { 1337 PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cdof)); 1338 if (!includeConstraints && cdof > 0) PetscCall(PetscSectionSetFieldConstraintDof(gs, p, f, cdof)); 1339 PetscCall(PetscSectionGetFieldDof(s, p, f, &dof)); 1340 PetscCall(PetscSectionSetFieldDof(gs, p, f, off < 0 ? -(dof + 1) : dof)); 1341 PetscCall(PetscSectionSetFieldOffset(gs, p, f, foff)); 1342 PetscCall(PetscSectionGetFieldConstraintDof(gs, p, f, &cdof)); 1343 foff = off < 0 ? foff - (dof - cdof) : foff + (dof - cdof); 1344 } 1345 } 1346 for (f = 0; f < numFields; ++f) { 1347 PetscSection gfs = gs->field[f]; 1348 1349 PetscCall(PetscSectionSetUpBC(gfs)); 1350 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])); 1351 } 1352 gs->setup = PETSC_TRUE; 1353 PetscCall(PetscSectionViewFromOptions(gs, NULL, "-global_section_view")); 1354 *gsection = gs; 1355 PetscFunctionReturn(0); 1356 } 1357 1358 /*@ 1359 PetscSectionCreateGlobalSectionCensored - Create a section describing the global field layout using 1360 the local section and an SF describing the section point overlap. 1361 1362 Input Parameters: 1363 + s - The PetscSection for the local field layout 1364 . sf - The SF describing parallel layout of the section points 1365 . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs 1366 . numExcludes - The number of exclusion ranges 1367 - excludes - An array [start_0, end_0, start_1, end_1, ...] where there are numExcludes pairs 1368 1369 Output Parameter: 1370 . gsection - The PetscSection for the global field layout 1371 1372 Note: This gives negative sizes and offsets to points not owned by this process 1373 1374 Level: advanced 1375 1376 .seealso: PetscSectionCreate() 1377 @*/ 1378 PetscErrorCode PetscSectionCreateGlobalSectionCensored(PetscSection s, PetscSF sf, PetscBool includeConstraints, PetscInt numExcludes, const PetscInt excludes[], PetscSection *gsection) 1379 { 1380 const PetscInt *pind = NULL; 1381 PetscInt *neg = NULL, *tmpOff = NULL; 1382 PetscInt pStart, pEnd, p, e, dof, cdof, off, globalOff = 0, nroots; 1383 1384 PetscFunctionBegin; 1385 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1386 PetscValidHeaderSpecific(sf, PETSCSF_CLASSID, 2); 1387 PetscValidPointer(gsection, 6); 1388 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) s), gsection)); 1389 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 1390 PetscCall(PetscSectionSetChart(*gsection, pStart, pEnd)); 1391 PetscCall(PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL)); 1392 if (nroots >= 0) { 1393 PetscCheckFalse(nroots < pEnd-pStart,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %" PetscInt_FMT " < %" PetscInt_FMT " section size", nroots, pEnd-pStart); 1394 PetscCall(PetscCalloc1(nroots, &neg)); 1395 if (nroots > pEnd-pStart) { 1396 PetscCall(PetscCalloc1(nroots, &tmpOff)); 1397 } else { 1398 tmpOff = &(*gsection)->atlasDof[-pStart]; 1399 } 1400 } 1401 /* Mark ghost points with negative dof */ 1402 for (p = pStart; p < pEnd; ++p) { 1403 for (e = 0; e < numExcludes; ++e) { 1404 if ((p >= excludes[e*2+0]) && (p < excludes[e*2+1])) { 1405 PetscCall(PetscSectionSetDof(*gsection, p, 0)); 1406 break; 1407 } 1408 } 1409 if (e < numExcludes) continue; 1410 PetscCall(PetscSectionGetDof(s, p, &dof)); 1411 PetscCall(PetscSectionSetDof(*gsection, p, dof)); 1412 PetscCall(PetscSectionGetConstraintDof(s, p, &cdof)); 1413 if (!includeConstraints && cdof > 0) PetscCall(PetscSectionSetConstraintDof(*gsection, p, cdof)); 1414 if (neg) neg[p] = -(dof+1); 1415 } 1416 PetscCall(PetscSectionSetUpBC(*gsection)); 1417 if (nroots >= 0) { 1418 PetscCall(PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff,MPI_REPLACE)); 1419 PetscCall(PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff,MPI_REPLACE)); 1420 if (nroots > pEnd - pStart) { 1421 for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];} 1422 } 1423 } 1424 /* Calculate new sizes, get proccess offset, and calculate point offsets */ 1425 if (s->perm) PetscCall(ISGetIndices(s->perm, &pind)); 1426 for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1427 const PetscInt q = pind ? pind[p] : p; 1428 1429 cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[q] : 0; 1430 (*gsection)->atlasOff[q] = off; 1431 off += (*gsection)->atlasDof[q] > 0 ? (*gsection)->atlasDof[q]-cdof : 0; 1432 } 1433 PetscCallMPI(MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject) s))); 1434 globalOff -= off; 1435 for (p = 0, off = 0; p < pEnd-pStart; ++p) { 1436 (*gsection)->atlasOff[p] += globalOff; 1437 if (neg) neg[p+pStart] = -((*gsection)->atlasOff[p]+1); 1438 } 1439 if (s->perm) PetscCall(ISRestoreIndices(s->perm, &pind)); 1440 /* Put in negative offsets for ghost points */ 1441 if (nroots >= 0) { 1442 if (nroots == pEnd-pStart) tmpOff = &(*gsection)->atlasOff[-pStart]; 1443 PetscCall(PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff,MPI_REPLACE)); 1444 PetscCall(PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff,MPI_REPLACE)); 1445 if (nroots > pEnd - pStart) { 1446 for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];} 1447 } 1448 } 1449 if (nroots >= 0 && nroots > pEnd-pStart) PetscCall(PetscFree(tmpOff)); 1450 PetscCall(PetscFree(neg)); 1451 PetscFunctionReturn(0); 1452 } 1453 1454 /*@ 1455 PetscSectionGetPointLayout - Get the PetscLayout associated with the section points 1456 1457 Collective on comm 1458 1459 Input Parameters: 1460 + comm - The MPI_Comm 1461 - s - The PetscSection 1462 1463 Output Parameter: 1464 . layout - The point layout for the section 1465 1466 Note: This is usually called for the default global section. 1467 1468 Level: advanced 1469 1470 .seealso: PetscSectionGetValueLayout(), PetscSectionCreate() 1471 @*/ 1472 PetscErrorCode PetscSectionGetPointLayout(MPI_Comm comm, PetscSection s, PetscLayout *layout) 1473 { 1474 PetscInt pStart, pEnd, p, localSize = 0; 1475 1476 PetscFunctionBegin; 1477 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 1478 for (p = pStart; p < pEnd; ++p) { 1479 PetscInt dof; 1480 1481 PetscCall(PetscSectionGetDof(s, p, &dof)); 1482 if (dof >= 0) ++localSize; 1483 } 1484 PetscCall(PetscLayoutCreate(comm, layout)); 1485 PetscCall(PetscLayoutSetLocalSize(*layout, localSize)); 1486 PetscCall(PetscLayoutSetBlockSize(*layout, 1)); 1487 PetscCall(PetscLayoutSetUp(*layout)); 1488 PetscFunctionReturn(0); 1489 } 1490 1491 /*@ 1492 PetscSectionGetValueLayout - Get the PetscLayout associated with the section dofs. 1493 1494 Collective on comm 1495 1496 Input Parameters: 1497 + comm - The MPI_Comm 1498 - s - The PetscSection 1499 1500 Output Parameter: 1501 . layout - The dof layout for the section 1502 1503 Note: This is usually called for the default global section. 1504 1505 Level: advanced 1506 1507 .seealso: PetscSectionGetPointLayout(), PetscSectionCreate() 1508 @*/ 1509 PetscErrorCode PetscSectionGetValueLayout(MPI_Comm comm, PetscSection s, PetscLayout *layout) 1510 { 1511 PetscInt pStart, pEnd, p, localSize = 0; 1512 1513 PetscFunctionBegin; 1514 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 2); 1515 PetscValidPointer(layout, 3); 1516 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 1517 for (p = pStart; p < pEnd; ++p) { 1518 PetscInt dof,cdof; 1519 1520 PetscCall(PetscSectionGetDof(s, p, &dof)); 1521 PetscCall(PetscSectionGetConstraintDof(s, p, &cdof)); 1522 if (dof-cdof > 0) localSize += dof-cdof; 1523 } 1524 PetscCall(PetscLayoutCreate(comm, layout)); 1525 PetscCall(PetscLayoutSetLocalSize(*layout, localSize)); 1526 PetscCall(PetscLayoutSetBlockSize(*layout, 1)); 1527 PetscCall(PetscLayoutSetUp(*layout)); 1528 PetscFunctionReturn(0); 1529 } 1530 1531 /*@ 1532 PetscSectionGetOffset - Return the offset into an array or local Vec for the dof associated with the given point. 1533 1534 Not Collective 1535 1536 Input Parameters: 1537 + s - the PetscSection 1538 - point - the point 1539 1540 Output Parameter: 1541 . offset - the offset 1542 1543 Level: intermediate 1544 1545 .seealso: PetscSectionGetFieldOffset(), PetscSectionCreate() 1546 @*/ 1547 PetscErrorCode PetscSectionGetOffset(PetscSection s, PetscInt point, PetscInt *offset) 1548 { 1549 PetscFunctionBegin; 1550 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1551 PetscValidIntPointer(offset, 3); 1552 if (PetscDefined(USE_DEBUG)) { 1553 PetscCheckFalse((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); 1554 } 1555 *offset = s->atlasOff[point - s->pStart]; 1556 PetscFunctionReturn(0); 1557 } 1558 1559 /*@ 1560 PetscSectionSetOffset - Set the offset into an array or local Vec for the dof associated with the given point. 1561 1562 Not Collective 1563 1564 Input Parameters: 1565 + s - the PetscSection 1566 . point - the point 1567 - offset - the offset 1568 1569 Note: The user usually does not call this function, but uses PetscSectionSetUp() 1570 1571 Level: intermediate 1572 1573 .seealso: PetscSectionGetFieldOffset(), PetscSectionCreate(), PetscSectionSetUp() 1574 @*/ 1575 PetscErrorCode PetscSectionSetOffset(PetscSection s, PetscInt point, PetscInt offset) 1576 { 1577 PetscFunctionBegin; 1578 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1579 PetscCheckFalse((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); 1580 s->atlasOff[point - s->pStart] = offset; 1581 PetscFunctionReturn(0); 1582 } 1583 1584 /*@ 1585 PetscSectionGetFieldOffset - Return the offset into an array or local Vec for the dof associated with the given point. 1586 1587 Not Collective 1588 1589 Input Parameters: 1590 + s - the PetscSection 1591 . point - the point 1592 - field - the field 1593 1594 Output Parameter: 1595 . offset - the offset 1596 1597 Level: intermediate 1598 1599 .seealso: PetscSectionGetOffset(), PetscSectionCreate() 1600 @*/ 1601 PetscErrorCode PetscSectionGetFieldOffset(PetscSection s, PetscInt point, PetscInt field, PetscInt *offset) 1602 { 1603 PetscFunctionBegin; 1604 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1605 PetscValidIntPointer(offset, 4); 1606 PetscSectionCheckValidField(field,s->numFields); 1607 PetscCall(PetscSectionGetOffset(s->field[field], point, offset)); 1608 PetscFunctionReturn(0); 1609 } 1610 1611 /*@ 1612 PetscSectionSetFieldOffset - Set the offset into an array or local Vec for the dof associated with the given field at a point. 1613 1614 Not Collective 1615 1616 Input Parameters: 1617 + s - the PetscSection 1618 . point - the point 1619 . field - the field 1620 - offset - the offset 1621 1622 Note: The user usually does not call this function, but uses PetscSectionSetUp() 1623 1624 Level: intermediate 1625 1626 .seealso: PetscSectionGetFieldOffset(), PetscSectionSetOffset(), PetscSectionCreate(), PetscSectionSetUp() 1627 @*/ 1628 PetscErrorCode PetscSectionSetFieldOffset(PetscSection s, PetscInt point, PetscInt field, PetscInt offset) 1629 { 1630 PetscFunctionBegin; 1631 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1632 PetscSectionCheckValidField(field,s->numFields); 1633 PetscCall(PetscSectionSetOffset(s->field[field], point, offset)); 1634 PetscFunctionReturn(0); 1635 } 1636 1637 /*@ 1638 PetscSectionGetFieldPointOffset - Return the offset on the given point for the dof associated with the given point. 1639 1640 Not Collective 1641 1642 Input Parameters: 1643 + s - the PetscSection 1644 . point - the point 1645 - field - the field 1646 1647 Output Parameter: 1648 . offset - the offset 1649 1650 Note: This gives the offset on a point of the field, ignoring constraints, meaning starting at the first dof for 1651 this point, what number is the first dof with this field. 1652 1653 Level: advanced 1654 1655 .seealso: PetscSectionGetOffset(), PetscSectionCreate() 1656 @*/ 1657 PetscErrorCode PetscSectionGetFieldPointOffset(PetscSection s, PetscInt point, PetscInt field, PetscInt *offset) 1658 { 1659 PetscInt off, foff; 1660 1661 PetscFunctionBegin; 1662 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1663 PetscValidIntPointer(offset, 4); 1664 PetscSectionCheckValidField(field,s->numFields); 1665 PetscCall(PetscSectionGetOffset(s, point, &off)); 1666 PetscCall(PetscSectionGetOffset(s->field[field], point, &foff)); 1667 *offset = foff - off; 1668 PetscFunctionReturn(0); 1669 } 1670 1671 /*@ 1672 PetscSectionGetOffsetRange - Return the full range of offsets [start, end) 1673 1674 Not Collective 1675 1676 Input Parameter: 1677 . s - the PetscSection 1678 1679 Output Parameters: 1680 + start - the minimum offset 1681 - end - one more than the maximum offset 1682 1683 Level: intermediate 1684 1685 .seealso: PetscSectionGetOffset(), PetscSectionCreate() 1686 @*/ 1687 PetscErrorCode PetscSectionGetOffsetRange(PetscSection s, PetscInt *start, PetscInt *end) 1688 { 1689 PetscInt os = 0, oe = 0, pStart, pEnd, p; 1690 1691 PetscFunctionBegin; 1692 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1693 if (s->atlasOff) {os = s->atlasOff[0]; oe = s->atlasOff[0];} 1694 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 1695 for (p = 0; p < pEnd-pStart; ++p) { 1696 PetscInt dof = s->atlasDof[p], off = s->atlasOff[p]; 1697 1698 if (off >= 0) { 1699 os = PetscMin(os, off); 1700 oe = PetscMax(oe, off+dof); 1701 } 1702 } 1703 if (start) *start = os; 1704 if (end) *end = oe; 1705 PetscFunctionReturn(0); 1706 } 1707 1708 /*@ 1709 PetscSectionCreateSubsection - Create a new, smaller section composed of only the selected fields 1710 1711 Collective 1712 1713 Input Parameters: 1714 + s - the PetscSection 1715 . len - the number of subfields 1716 - fields - the subfield numbers 1717 1718 Output Parameter: 1719 . subs - the subsection 1720 1721 Note: The section offsets now refer to a new, smaller vector. 1722 1723 Level: advanced 1724 1725 .seealso: PetscSectionCreateSupersection(), PetscSectionCreate() 1726 @*/ 1727 PetscErrorCode PetscSectionCreateSubsection(PetscSection s, PetscInt len, const PetscInt fields[], PetscSection *subs) 1728 { 1729 PetscInt nF, f, c, pStart, pEnd, p, maxCdof = 0; 1730 1731 PetscFunctionBegin; 1732 if (!len) PetscFunctionReturn(0); 1733 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1734 PetscValidIntPointer(fields, 3); 1735 PetscValidPointer(subs, 4); 1736 PetscCall(PetscSectionGetNumFields(s, &nF)); 1737 PetscCheckFalse(len > nF,PetscObjectComm((PetscObject) s), PETSC_ERR_ARG_WRONG, "Number of requested fields %" PetscInt_FMT " greater than number of fields %" PetscInt_FMT, len, nF); 1738 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) s), subs)); 1739 PetscCall(PetscSectionSetNumFields(*subs, len)); 1740 for (f = 0; f < len; ++f) { 1741 const char *name = NULL; 1742 PetscInt numComp = 0; 1743 1744 PetscCall(PetscSectionGetFieldName(s, fields[f], &name)); 1745 PetscCall(PetscSectionSetFieldName(*subs, f, name)); 1746 PetscCall(PetscSectionGetFieldComponents(s, fields[f], &numComp)); 1747 PetscCall(PetscSectionSetFieldComponents(*subs, f, numComp)); 1748 for (c = 0; c < s->numFieldComponents[fields[f]]; ++c) { 1749 PetscCall(PetscSectionGetComponentName(s, fields[f], c, &name)); 1750 PetscCall(PetscSectionSetComponentName(*subs, f, c, name)); 1751 } 1752 } 1753 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 1754 PetscCall(PetscSectionSetChart(*subs, pStart, pEnd)); 1755 for (p = pStart; p < pEnd; ++p) { 1756 PetscInt dof = 0, cdof = 0, fdof = 0, cfdof = 0; 1757 1758 for (f = 0; f < len; ++f) { 1759 PetscCall(PetscSectionGetFieldDof(s, p, fields[f], &fdof)); 1760 PetscCall(PetscSectionSetFieldDof(*subs, p, f, fdof)); 1761 PetscCall(PetscSectionGetFieldConstraintDof(s, p, fields[f], &cfdof)); 1762 if (cfdof) PetscCall(PetscSectionSetFieldConstraintDof(*subs, p, f, cfdof)); 1763 dof += fdof; 1764 cdof += cfdof; 1765 } 1766 PetscCall(PetscSectionSetDof(*subs, p, dof)); 1767 if (cdof) PetscCall(PetscSectionSetConstraintDof(*subs, p, cdof)); 1768 maxCdof = PetscMax(cdof, maxCdof); 1769 } 1770 PetscCall(PetscSectionSetUp(*subs)); 1771 if (maxCdof) { 1772 PetscInt *indices; 1773 1774 PetscCall(PetscMalloc1(maxCdof, &indices)); 1775 for (p = pStart; p < pEnd; ++p) { 1776 PetscInt cdof; 1777 1778 PetscCall(PetscSectionGetConstraintDof(*subs, p, &cdof)); 1779 if (cdof) { 1780 const PetscInt *oldIndices = NULL; 1781 PetscInt fdof = 0, cfdof = 0, fc, numConst = 0, fOff = 0; 1782 1783 for (f = 0; f < len; ++f) { 1784 PetscCall(PetscSectionGetFieldDof(s, p, fields[f], &fdof)); 1785 PetscCall(PetscSectionGetFieldConstraintDof(s, p, fields[f], &cfdof)); 1786 PetscCall(PetscSectionGetFieldConstraintIndices(s, p, fields[f], &oldIndices)); 1787 PetscCall(PetscSectionSetFieldConstraintIndices(*subs, p, f, oldIndices)); 1788 for (fc = 0; fc < cfdof; ++fc) indices[numConst+fc] = oldIndices[fc] + fOff; 1789 numConst += cfdof; 1790 fOff += fdof; 1791 } 1792 PetscCall(PetscSectionSetConstraintIndices(*subs, p, indices)); 1793 } 1794 } 1795 PetscCall(PetscFree(indices)); 1796 } 1797 PetscFunctionReturn(0); 1798 } 1799 1800 /*@ 1801 PetscSectionCreateSupersection - Create a new, larger section composed of the input sections 1802 1803 Collective 1804 1805 Input Parameters: 1806 + s - the input sections 1807 - len - the number of input sections 1808 1809 Output Parameter: 1810 . supers - the supersection 1811 1812 Note: The section offsets now refer to a new, larger vector. 1813 1814 Level: advanced 1815 1816 .seealso: PetscSectionCreateSubsection(), PetscSectionCreate() 1817 @*/ 1818 PetscErrorCode PetscSectionCreateSupersection(PetscSection s[], PetscInt len, PetscSection *supers) 1819 { 1820 PetscInt Nf = 0, f, pStart = PETSC_MAX_INT, pEnd = 0, p, maxCdof = 0, i; 1821 1822 PetscFunctionBegin; 1823 if (!len) PetscFunctionReturn(0); 1824 for (i = 0; i < len; ++i) { 1825 PetscInt nf, pStarti, pEndi; 1826 1827 PetscCall(PetscSectionGetNumFields(s[i], &nf)); 1828 PetscCall(PetscSectionGetChart(s[i], &pStarti, &pEndi)); 1829 pStart = PetscMin(pStart, pStarti); 1830 pEnd = PetscMax(pEnd, pEndi); 1831 Nf += nf; 1832 } 1833 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) s[0]), supers)); 1834 PetscCall(PetscSectionSetNumFields(*supers, Nf)); 1835 for (i = 0, f = 0; i < len; ++i) { 1836 PetscInt nf, fi, ci; 1837 1838 PetscCall(PetscSectionGetNumFields(s[i], &nf)); 1839 for (fi = 0; fi < nf; ++fi, ++f) { 1840 const char *name = NULL; 1841 PetscInt numComp = 0; 1842 1843 PetscCall(PetscSectionGetFieldName(s[i], fi, &name)); 1844 PetscCall(PetscSectionSetFieldName(*supers, f, name)); 1845 PetscCall(PetscSectionGetFieldComponents(s[i], fi, &numComp)); 1846 PetscCall(PetscSectionSetFieldComponents(*supers, f, numComp)); 1847 for (ci = 0; ci < s[i]->numFieldComponents[fi]; ++ci) { 1848 PetscCall(PetscSectionGetComponentName(s[i], fi, ci, &name)); 1849 PetscCall(PetscSectionSetComponentName(*supers, f, ci, name)); 1850 } 1851 } 1852 } 1853 PetscCall(PetscSectionSetChart(*supers, pStart, pEnd)); 1854 for (p = pStart; p < pEnd; ++p) { 1855 PetscInt dof = 0, cdof = 0; 1856 1857 for (i = 0, f = 0; i < len; ++i) { 1858 PetscInt nf, fi, pStarti, pEndi; 1859 PetscInt fdof = 0, cfdof = 0; 1860 1861 PetscCall(PetscSectionGetNumFields(s[i], &nf)); 1862 PetscCall(PetscSectionGetChart(s[i], &pStarti, &pEndi)); 1863 if ((p < pStarti) || (p >= pEndi)) continue; 1864 for (fi = 0; fi < nf; ++fi, ++f) { 1865 PetscCall(PetscSectionGetFieldDof(s[i], p, fi, &fdof)); 1866 PetscCall(PetscSectionAddFieldDof(*supers, p, f, fdof)); 1867 PetscCall(PetscSectionGetFieldConstraintDof(s[i], p, fi, &cfdof)); 1868 if (cfdof) PetscCall(PetscSectionAddFieldConstraintDof(*supers, p, f, cfdof)); 1869 dof += fdof; 1870 cdof += cfdof; 1871 } 1872 } 1873 PetscCall(PetscSectionSetDof(*supers, p, dof)); 1874 if (cdof) PetscCall(PetscSectionSetConstraintDof(*supers, p, cdof)); 1875 maxCdof = PetscMax(cdof, maxCdof); 1876 } 1877 PetscCall(PetscSectionSetUp(*supers)); 1878 if (maxCdof) { 1879 PetscInt *indices; 1880 1881 PetscCall(PetscMalloc1(maxCdof, &indices)); 1882 for (p = pStart; p < pEnd; ++p) { 1883 PetscInt cdof; 1884 1885 PetscCall(PetscSectionGetConstraintDof(*supers, p, &cdof)); 1886 if (cdof) { 1887 PetscInt dof, numConst = 0, fOff = 0; 1888 1889 for (i = 0, f = 0; i < len; ++i) { 1890 const PetscInt *oldIndices = NULL; 1891 PetscInt nf, fi, pStarti, pEndi, fdof, cfdof, fc; 1892 1893 PetscCall(PetscSectionGetNumFields(s[i], &nf)); 1894 PetscCall(PetscSectionGetChart(s[i], &pStarti, &pEndi)); 1895 if ((p < pStarti) || (p >= pEndi)) continue; 1896 for (fi = 0; fi < nf; ++fi, ++f) { 1897 PetscCall(PetscSectionGetFieldDof(s[i], p, fi, &fdof)); 1898 PetscCall(PetscSectionGetFieldConstraintDof(s[i], p, fi, &cfdof)); 1899 PetscCall(PetscSectionGetFieldConstraintIndices(s[i], p, fi, &oldIndices)); 1900 for (fc = 0; fc < cfdof; ++fc) indices[numConst+fc] = oldIndices[fc]; 1901 PetscCall(PetscSectionSetFieldConstraintIndices(*supers, p, f, &indices[numConst])); 1902 for (fc = 0; fc < cfdof; ++fc) indices[numConst+fc] += fOff; 1903 numConst += cfdof; 1904 } 1905 PetscCall(PetscSectionGetDof(s[i], p, &dof)); 1906 fOff += dof; 1907 } 1908 PetscCall(PetscSectionSetConstraintIndices(*supers, p, indices)); 1909 } 1910 } 1911 PetscCall(PetscFree(indices)); 1912 } 1913 PetscFunctionReturn(0); 1914 } 1915 1916 /*@ 1917 PetscSectionCreateSubmeshSection - Create a new, smaller section with support on the submesh 1918 1919 Collective 1920 1921 Input Parameters: 1922 + s - the PetscSection 1923 - subpointMap - a sorted list of points in the original mesh which are in the submesh 1924 1925 Output Parameter: 1926 . subs - the subsection 1927 1928 Note: The section offsets now refer to a new, smaller vector. 1929 1930 Level: advanced 1931 1932 .seealso: PetscSectionCreateSubsection(), DMPlexGetSubpointMap(), PetscSectionCreate() 1933 @*/ 1934 PetscErrorCode PetscSectionCreateSubmeshSection(PetscSection s, IS subpointMap, PetscSection *subs) 1935 { 1936 const PetscInt *points = NULL, *indices = NULL; 1937 PetscInt numFields, f, c, numSubpoints = 0, pStart, pEnd, p, subp; 1938 1939 PetscFunctionBegin; 1940 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 1941 PetscValidHeaderSpecific(subpointMap, IS_CLASSID, 2); 1942 PetscValidPointer(subs, 3); 1943 PetscCall(PetscSectionGetNumFields(s, &numFields)); 1944 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) s), subs)); 1945 if (numFields) PetscCall(PetscSectionSetNumFields(*subs, numFields)); 1946 for (f = 0; f < numFields; ++f) { 1947 const char *name = NULL; 1948 PetscInt numComp = 0; 1949 1950 PetscCall(PetscSectionGetFieldName(s, f, &name)); 1951 PetscCall(PetscSectionSetFieldName(*subs, f, name)); 1952 PetscCall(PetscSectionGetFieldComponents(s, f, &numComp)); 1953 PetscCall(PetscSectionSetFieldComponents(*subs, f, numComp)); 1954 for (c = 0; c < s->numFieldComponents[f]; ++c) { 1955 PetscCall(PetscSectionGetComponentName(s, f, c, &name)); 1956 PetscCall(PetscSectionSetComponentName(*subs, f, c, name)); 1957 } 1958 } 1959 /* For right now, we do not try to squeeze the subchart */ 1960 if (subpointMap) { 1961 PetscCall(ISGetSize(subpointMap, &numSubpoints)); 1962 PetscCall(ISGetIndices(subpointMap, &points)); 1963 } 1964 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 1965 PetscCall(PetscSectionSetChart(*subs, 0, numSubpoints)); 1966 for (p = pStart; p < pEnd; ++p) { 1967 PetscInt dof, cdof, fdof = 0, cfdof = 0; 1968 1969 PetscCall(PetscFindInt(p, numSubpoints, points, &subp)); 1970 if (subp < 0) continue; 1971 for (f = 0; f < numFields; ++f) { 1972 PetscCall(PetscSectionGetFieldDof(s, p, f, &fdof)); 1973 PetscCall(PetscSectionSetFieldDof(*subs, subp, f, fdof)); 1974 PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cfdof)); 1975 if (cfdof) PetscCall(PetscSectionSetFieldConstraintDof(*subs, subp, f, cfdof)); 1976 } 1977 PetscCall(PetscSectionGetDof(s, p, &dof)); 1978 PetscCall(PetscSectionSetDof(*subs, subp, dof)); 1979 PetscCall(PetscSectionGetConstraintDof(s, p, &cdof)); 1980 if (cdof) PetscCall(PetscSectionSetConstraintDof(*subs, subp, cdof)); 1981 } 1982 PetscCall(PetscSectionSetUp(*subs)); 1983 /* Change offsets to original offsets */ 1984 for (p = pStart; p < pEnd; ++p) { 1985 PetscInt off, foff = 0; 1986 1987 PetscCall(PetscFindInt(p, numSubpoints, points, &subp)); 1988 if (subp < 0) continue; 1989 for (f = 0; f < numFields; ++f) { 1990 PetscCall(PetscSectionGetFieldOffset(s, p, f, &foff)); 1991 PetscCall(PetscSectionSetFieldOffset(*subs, subp, f, foff)); 1992 } 1993 PetscCall(PetscSectionGetOffset(s, p, &off)); 1994 PetscCall(PetscSectionSetOffset(*subs, subp, off)); 1995 } 1996 /* Copy constraint indices */ 1997 for (subp = 0; subp < numSubpoints; ++subp) { 1998 PetscInt cdof; 1999 2000 PetscCall(PetscSectionGetConstraintDof(*subs, subp, &cdof)); 2001 if (cdof) { 2002 for (f = 0; f < numFields; ++f) { 2003 PetscCall(PetscSectionGetFieldConstraintIndices(s, points[subp], f, &indices)); 2004 PetscCall(PetscSectionSetFieldConstraintIndices(*subs, subp, f, indices)); 2005 } 2006 PetscCall(PetscSectionGetConstraintIndices(s, points[subp], &indices)); 2007 PetscCall(PetscSectionSetConstraintIndices(*subs, subp, indices)); 2008 } 2009 } 2010 if (subpointMap) PetscCall(ISRestoreIndices(subpointMap, &points)); 2011 PetscFunctionReturn(0); 2012 } 2013 2014 static PetscErrorCode PetscSectionView_ASCII(PetscSection s, PetscViewer viewer) 2015 { 2016 PetscInt p; 2017 PetscMPIInt rank; 2018 2019 PetscFunctionBegin; 2020 PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank)); 2021 PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 2022 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "Process %d:\n", rank)); 2023 for (p = 0; p < s->pEnd - s->pStart; ++p) { 2024 if ((s->bc) && (s->bc->atlasDof[p] > 0)) { 2025 PetscInt b; 2026 2027 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " (%4" PetscInt_FMT ") dim %2" PetscInt_FMT " offset %3" PetscInt_FMT " constrained", p+s->pStart, s->atlasDof[p], s->atlasOff[p])); 2028 if (s->bcIndices) { 2029 for (b = 0; b < s->bc->atlasDof[p]; ++b) { 2030 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " %" PetscInt_FMT, s->bcIndices[s->bc->atlasOff[p]+b])); 2031 } 2032 } 2033 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "\n")); 2034 } else { 2035 PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " (%4" PetscInt_FMT ") dim %2" PetscInt_FMT " offset %3" PetscInt_FMT "\n", p+s->pStart, s->atlasDof[p], s->atlasOff[p])); 2036 } 2037 } 2038 PetscCall(PetscViewerFlush(viewer)); 2039 PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 2040 if (s->sym) { 2041 PetscCall(PetscViewerASCIIPushTab(viewer)); 2042 PetscCall(PetscSectionSymView(s->sym,viewer)); 2043 PetscCall(PetscViewerASCIIPopTab(viewer)); 2044 } 2045 PetscFunctionReturn(0); 2046 } 2047 2048 /*@C 2049 PetscSectionViewFromOptions - View from Options 2050 2051 Collective 2052 2053 Input Parameters: 2054 + A - the PetscSection object to view 2055 . obj - Optional object 2056 - name - command line option 2057 2058 Level: intermediate 2059 .seealso: PetscSection, PetscSectionView, PetscObjectViewFromOptions(), PetscSectionCreate() 2060 @*/ 2061 PetscErrorCode PetscSectionViewFromOptions(PetscSection A,PetscObject obj,const char name[]) 2062 { 2063 PetscFunctionBegin; 2064 PetscValidHeaderSpecific(A,PETSC_SECTION_CLASSID,1); 2065 PetscCall(PetscObjectViewFromOptions((PetscObject)A,obj,name)); 2066 PetscFunctionReturn(0); 2067 } 2068 2069 /*@C 2070 PetscSectionView - Views a PetscSection 2071 2072 Collective 2073 2074 Input Parameters: 2075 + s - the PetscSection object to view 2076 - v - the viewer 2077 2078 Note: 2079 PetscSectionView(), when viewer is of type PETSCVIEWERHDF5, only saves 2080 distribution independent data, such as dofs, offsets, constraint dofs, 2081 and constraint indices. Points that have negative dofs, for instance, 2082 are not saved as they represent points owned by other processes. 2083 Point numbering and rank assignment is currently not stored. 2084 The saved section can be loaded with PetscSectionLoad(). 2085 2086 Level: beginner 2087 2088 .seealso PetscSectionCreate(), PetscSectionDestroy(), PetscSectionLoad() 2089 @*/ 2090 PetscErrorCode PetscSectionView(PetscSection s, PetscViewer viewer) 2091 { 2092 PetscBool isascii, ishdf5; 2093 PetscInt f; 2094 2095 PetscFunctionBegin; 2096 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2097 if (!viewer) PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)s), &viewer)); 2098 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 2099 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii)); 2100 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5)); 2101 if (isascii) { 2102 PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)s,viewer)); 2103 if (s->numFields) { 2104 PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " fields\n", s->numFields)); 2105 for (f = 0; f < s->numFields; ++f) { 2106 PetscCall(PetscViewerASCIIPrintf(viewer, " field %" PetscInt_FMT " with %" PetscInt_FMT " components\n", f, s->numFieldComponents[f])); 2107 PetscCall(PetscSectionView_ASCII(s->field[f], viewer)); 2108 } 2109 } else { 2110 PetscCall(PetscSectionView_ASCII(s, viewer)); 2111 } 2112 } else if (ishdf5) { 2113 #if PetscDefined(HAVE_HDF5) 2114 PetscCall(PetscSectionView_HDF5_Internal(s, viewer)); 2115 #else 2116 SETERRQ(PetscObjectComm((PetscObject) s), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5"); 2117 #endif 2118 } 2119 PetscFunctionReturn(0); 2120 } 2121 2122 /*@C 2123 PetscSectionLoad - Loads a PetscSection 2124 2125 Collective 2126 2127 Input Parameters: 2128 + s - the PetscSection object to load 2129 - v - the viewer 2130 2131 Note: 2132 PetscSectionLoad(), when viewer is of type PETSCVIEWERHDF5, loads 2133 a section saved with PetscSectionView(). The number of processes 2134 used here (N) does not need to be the same as that used when saving. 2135 After calling this function, the chart of s on rank i will be set 2136 to [0, E_i), where \sum_{i=0}^{N-1}E_i equals to the total number of 2137 saved section points. 2138 2139 Level: beginner 2140 2141 .seealso PetscSectionCreate(), PetscSectionDestroy(), PetscSectionView() 2142 @*/ 2143 PetscErrorCode PetscSectionLoad(PetscSection s, PetscViewer viewer) 2144 { 2145 PetscBool ishdf5; 2146 2147 PetscFunctionBegin; 2148 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2149 PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 2150 PetscCall(PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERHDF5, &ishdf5)); 2151 if (ishdf5) { 2152 #if PetscDefined(HAVE_HDF5) 2153 PetscCall(PetscSectionLoad_HDF5_Internal(s, viewer)); 2154 PetscFunctionReturn(0); 2155 #else 2156 SETERRQ(PetscObjectComm((PetscObject) s), PETSC_ERR_SUP, "HDF5 not supported in this build.\nPlease reconfigure using --download-hdf5"); 2157 #endif 2158 } else SETERRQ(PetscObjectComm((PetscObject) s), PETSC_ERR_SUP, "Viewer type %s not yet supported for PetscSection loading", ((PetscObject)viewer)->type_name); 2159 } 2160 2161 static PetscErrorCode PetscSectionResetClosurePermutation(PetscSection section) 2162 { 2163 PetscSectionClosurePermVal clVal; 2164 2165 PetscFunctionBegin; 2166 if (!section->clHash) PetscFunctionReturn(0); 2167 kh_foreach_value(section->clHash, clVal, { 2168 PetscCall(PetscFree(clVal.perm)); 2169 PetscCall(PetscFree(clVal.invPerm)); 2170 }); 2171 kh_destroy(ClPerm, section->clHash); 2172 section->clHash = NULL; 2173 PetscFunctionReturn(0); 2174 } 2175 2176 /*@ 2177 PetscSectionReset - Frees all section data. 2178 2179 Not Collective 2180 2181 Input Parameters: 2182 . s - the PetscSection 2183 2184 Level: beginner 2185 2186 .seealso: PetscSection, PetscSectionCreate() 2187 @*/ 2188 PetscErrorCode PetscSectionReset(PetscSection s) 2189 { 2190 PetscInt f, c; 2191 2192 PetscFunctionBegin; 2193 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2194 for (f = 0; f < s->numFields; ++f) { 2195 PetscCall(PetscSectionDestroy(&s->field[f])); 2196 PetscCall(PetscFree(s->fieldNames[f])); 2197 for (c = 0; c < s->numFieldComponents[f]; ++c) { 2198 PetscCall(PetscFree(s->compNames[f][c])); 2199 } 2200 PetscCall(PetscFree(s->compNames[f])); 2201 } 2202 PetscCall(PetscFree(s->numFieldComponents)); 2203 PetscCall(PetscFree(s->fieldNames)); 2204 PetscCall(PetscFree(s->compNames)); 2205 PetscCall(PetscFree(s->field)); 2206 PetscCall(PetscSectionDestroy(&s->bc)); 2207 PetscCall(PetscFree(s->bcIndices)); 2208 PetscCall(PetscFree2(s->atlasDof, s->atlasOff)); 2209 PetscCall(PetscSectionDestroy(&s->clSection)); 2210 PetscCall(ISDestroy(&s->clPoints)); 2211 PetscCall(ISDestroy(&s->perm)); 2212 PetscCall(PetscSectionResetClosurePermutation(s)); 2213 PetscCall(PetscSectionSymDestroy(&s->sym)); 2214 PetscCall(PetscSectionDestroy(&s->clSection)); 2215 PetscCall(ISDestroy(&s->clPoints)); 2216 2217 s->pStart = -1; 2218 s->pEnd = -1; 2219 s->maxDof = 0; 2220 s->setup = PETSC_FALSE; 2221 s->numFields = 0; 2222 s->clObj = NULL; 2223 PetscFunctionReturn(0); 2224 } 2225 2226 /*@ 2227 PetscSectionDestroy - Frees a section object and frees its range if that exists. 2228 2229 Not Collective 2230 2231 Input Parameters: 2232 . s - the PetscSection 2233 2234 Level: beginner 2235 2236 .seealso: PetscSection, PetscSectionCreate() 2237 @*/ 2238 PetscErrorCode PetscSectionDestroy(PetscSection *s) 2239 { 2240 PetscFunctionBegin; 2241 if (!*s) PetscFunctionReturn(0); 2242 PetscValidHeaderSpecific(*s,PETSC_SECTION_CLASSID,1); 2243 if (--((PetscObject)(*s))->refct > 0) { 2244 *s = NULL; 2245 PetscFunctionReturn(0); 2246 } 2247 PetscCall(PetscSectionReset(*s)); 2248 PetscCall(PetscHeaderDestroy(s)); 2249 PetscFunctionReturn(0); 2250 } 2251 2252 PetscErrorCode VecIntGetValuesSection(PetscInt *baseArray, PetscSection s, PetscInt point, const PetscInt **values) 2253 { 2254 const PetscInt p = point - s->pStart; 2255 2256 PetscFunctionBegin; 2257 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 2); 2258 *values = &baseArray[s->atlasOff[p]]; 2259 PetscFunctionReturn(0); 2260 } 2261 2262 PetscErrorCode VecIntSetValuesSection(PetscInt *baseArray, PetscSection s, PetscInt point, const PetscInt values[], InsertMode mode) 2263 { 2264 PetscInt *array; 2265 const PetscInt p = point - s->pStart; 2266 const PetscInt orientation = 0; /* Needs to be included for use in closure operations */ 2267 PetscInt cDim = 0; 2268 2269 PetscFunctionBegin; 2270 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 2); 2271 PetscCall(PetscSectionGetConstraintDof(s, p, &cDim)); 2272 array = &baseArray[s->atlasOff[p]]; 2273 if (!cDim) { 2274 if (orientation >= 0) { 2275 const PetscInt dim = s->atlasDof[p]; 2276 PetscInt i; 2277 2278 if (mode == INSERT_VALUES) { 2279 for (i = 0; i < dim; ++i) array[i] = values[i]; 2280 } else { 2281 for (i = 0; i < dim; ++i) array[i] += values[i]; 2282 } 2283 } else { 2284 PetscInt offset = 0; 2285 PetscInt j = -1, field, i; 2286 2287 for (field = 0; field < s->numFields; ++field) { 2288 const PetscInt dim = s->field[field]->atlasDof[p]; 2289 2290 for (i = dim-1; i >= 0; --i) array[++j] = values[i+offset]; 2291 offset += dim; 2292 } 2293 } 2294 } else { 2295 if (orientation >= 0) { 2296 const PetscInt dim = s->atlasDof[p]; 2297 PetscInt cInd = 0, i; 2298 const PetscInt *cDof; 2299 2300 PetscCall(PetscSectionGetConstraintIndices(s, point, &cDof)); 2301 if (mode == INSERT_VALUES) { 2302 for (i = 0; i < dim; ++i) { 2303 if ((cInd < cDim) && (i == cDof[cInd])) {++cInd; continue;} 2304 array[i] = values[i]; 2305 } 2306 } else { 2307 for (i = 0; i < dim; ++i) { 2308 if ((cInd < cDim) && (i == cDof[cInd])) {++cInd; continue;} 2309 array[i] += values[i]; 2310 } 2311 } 2312 } else { 2313 const PetscInt *cDof; 2314 PetscInt offset = 0; 2315 PetscInt cOffset = 0; 2316 PetscInt j = 0, field; 2317 2318 PetscCall(PetscSectionGetConstraintIndices(s, point, &cDof)); 2319 for (field = 0; field < s->numFields; ++field) { 2320 const PetscInt dim = s->field[field]->atlasDof[p]; /* PetscSectionGetFieldDof() */ 2321 const PetscInt tDim = s->field[field]->bc->atlasDof[p]; /* PetscSectionGetFieldConstraintDof() */ 2322 const PetscInt sDim = dim - tDim; 2323 PetscInt cInd = 0, i,k; 2324 2325 for (i = 0, k = dim+offset-1; i < dim; ++i, ++j, --k) { 2326 if ((cInd < sDim) && (j == cDof[cInd+cOffset])) {++cInd; continue;} 2327 array[j] = values[k]; 2328 } 2329 offset += dim; 2330 cOffset += dim - tDim; 2331 } 2332 } 2333 } 2334 PetscFunctionReturn(0); 2335 } 2336 2337 /*@C 2338 PetscSectionHasConstraints - Determine whether a section has constrained dofs 2339 2340 Not Collective 2341 2342 Input Parameter: 2343 . s - The PetscSection 2344 2345 Output Parameter: 2346 . hasConstraints - flag indicating that the section has constrained dofs 2347 2348 Level: intermediate 2349 2350 .seealso: PetscSectionSetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection 2351 @*/ 2352 PetscErrorCode PetscSectionHasConstraints(PetscSection s, PetscBool *hasConstraints) 2353 { 2354 PetscFunctionBegin; 2355 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2356 PetscValidBoolPointer(hasConstraints, 2); 2357 *hasConstraints = s->bc ? PETSC_TRUE : PETSC_FALSE; 2358 PetscFunctionReturn(0); 2359 } 2360 2361 /*@C 2362 PetscSectionGetConstraintIndices - Get the point dof numbers, in [0, dof), which are constrained 2363 2364 Not Collective 2365 2366 Input Parameters: 2367 + s - The PetscSection 2368 - point - The point 2369 2370 Output Parameter: 2371 . indices - The constrained dofs 2372 2373 Note: In Fortran, you call PetscSectionGetConstraintIndicesF90() and PetscSectionRestoreConstraintIndicesF90() 2374 2375 Level: intermediate 2376 2377 .seealso: PetscSectionSetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection 2378 @*/ 2379 PetscErrorCode PetscSectionGetConstraintIndices(PetscSection s, PetscInt point, const PetscInt **indices) 2380 { 2381 PetscFunctionBegin; 2382 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2383 if (s->bc) { 2384 PetscCall(VecIntGetValuesSection(s->bcIndices, s->bc, point, indices)); 2385 } else *indices = NULL; 2386 PetscFunctionReturn(0); 2387 } 2388 2389 /*@C 2390 PetscSectionSetConstraintIndices - Set the point dof numbers, in [0, dof), which are constrained 2391 2392 Not Collective 2393 2394 Input Parameters: 2395 + s - The PetscSection 2396 . point - The point 2397 - indices - The constrained dofs 2398 2399 Note: The Fortran is PetscSectionSetConstraintIndicesF90() 2400 2401 Level: intermediate 2402 2403 .seealso: PetscSectionGetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection 2404 @*/ 2405 PetscErrorCode PetscSectionSetConstraintIndices(PetscSection s, PetscInt point, const PetscInt indices[]) 2406 { 2407 PetscFunctionBegin; 2408 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2409 if (s->bc) { 2410 const PetscInt dof = s->atlasDof[point]; 2411 const PetscInt cdof = s->bc->atlasDof[point]; 2412 PetscInt d; 2413 2414 for (d = 0; d < cdof; ++d) { 2415 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]); 2416 } 2417 PetscCall(VecIntSetValuesSection(s->bcIndices, s->bc, point, indices, INSERT_VALUES)); 2418 } 2419 PetscFunctionReturn(0); 2420 } 2421 2422 /*@C 2423 PetscSectionGetFieldConstraintIndices - Get the field dof numbers, in [0, fdof), which are constrained 2424 2425 Not Collective 2426 2427 Input Parameters: 2428 + s - The PetscSection 2429 . field - The field number 2430 - point - The point 2431 2432 Output Parameter: 2433 . indices - The constrained dofs sorted in ascending order 2434 2435 Notes: 2436 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(). 2437 2438 Fortran Note: 2439 In Fortran, you call PetscSectionGetFieldConstraintIndicesF90() and PetscSectionRestoreFieldConstraintIndicesF90() 2440 2441 Level: intermediate 2442 2443 .seealso: PetscSectionSetFieldConstraintIndices(), PetscSectionGetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection 2444 @*/ 2445 PetscErrorCode PetscSectionGetFieldConstraintIndices(PetscSection s, PetscInt point, PetscInt field, const PetscInt **indices) 2446 { 2447 PetscFunctionBegin; 2448 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2449 PetscValidPointer(indices,4); 2450 PetscSectionCheckValidField(field,s->numFields); 2451 PetscCall(PetscSectionGetConstraintIndices(s->field[field], point, indices)); 2452 PetscFunctionReturn(0); 2453 } 2454 2455 /*@C 2456 PetscSectionSetFieldConstraintIndices - Set the field dof numbers, in [0, fdof), which are constrained 2457 2458 Not Collective 2459 2460 Input Parameters: 2461 + s - The PetscSection 2462 . point - The point 2463 . field - The field number 2464 - indices - The constrained dofs 2465 2466 Note: The Fortran is PetscSectionSetFieldConstraintIndicesF90() 2467 2468 Level: intermediate 2469 2470 .seealso: PetscSectionSetConstraintIndices(), PetscSectionGetFieldConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection 2471 @*/ 2472 PetscErrorCode PetscSectionSetFieldConstraintIndices(PetscSection s, PetscInt point, PetscInt field, const PetscInt indices[]) 2473 { 2474 PetscFunctionBegin; 2475 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2476 if (PetscDefined(USE_DEBUG)) { 2477 PetscInt nfdof; 2478 2479 PetscCall(PetscSectionGetFieldConstraintDof(s, point, field, &nfdof)); 2480 if (nfdof) PetscValidIntPointer(indices, 4); 2481 } 2482 PetscSectionCheckValidField(field,s->numFields); 2483 PetscCall(PetscSectionSetConstraintIndices(s->field[field], point, indices)); 2484 PetscFunctionReturn(0); 2485 } 2486 2487 /*@ 2488 PetscSectionPermute - Reorder the section according to the input point permutation 2489 2490 Collective 2491 2492 Input Parameters: 2493 + section - The PetscSection object 2494 - perm - The point permutation, old point p becomes new point perm[p] 2495 2496 Output Parameter: 2497 . sectionNew - The permuted PetscSection 2498 2499 Level: intermediate 2500 2501 .seealso: MatPermute() 2502 @*/ 2503 PetscErrorCode PetscSectionPermute(PetscSection section, IS permutation, PetscSection *sectionNew) 2504 { 2505 PetscSection s = section, sNew; 2506 const PetscInt *perm; 2507 PetscInt numFields, f, c, numPoints, pStart, pEnd, p; 2508 2509 PetscFunctionBegin; 2510 PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1); 2511 PetscValidHeaderSpecific(permutation, IS_CLASSID, 2); 2512 PetscValidPointer(sectionNew, 3); 2513 PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject) s), &sNew)); 2514 PetscCall(PetscSectionGetNumFields(s, &numFields)); 2515 if (numFields) PetscCall(PetscSectionSetNumFields(sNew, numFields)); 2516 for (f = 0; f < numFields; ++f) { 2517 const char *name; 2518 PetscInt numComp; 2519 2520 PetscCall(PetscSectionGetFieldName(s, f, &name)); 2521 PetscCall(PetscSectionSetFieldName(sNew, f, name)); 2522 PetscCall(PetscSectionGetFieldComponents(s, f, &numComp)); 2523 PetscCall(PetscSectionSetFieldComponents(sNew, f, numComp)); 2524 for (c = 0; c < s->numFieldComponents[f]; ++c) { 2525 PetscCall(PetscSectionGetComponentName(s, f, c, &name)); 2526 PetscCall(PetscSectionSetComponentName(sNew, f, c, name)); 2527 } 2528 } 2529 PetscCall(ISGetLocalSize(permutation, &numPoints)); 2530 PetscCall(ISGetIndices(permutation, &perm)); 2531 PetscCall(PetscSectionGetChart(s, &pStart, &pEnd)); 2532 PetscCall(PetscSectionSetChart(sNew, pStart, pEnd)); 2533 PetscCheckFalse(numPoints < pEnd,PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Permutation size %" PetscInt_FMT " is less than largest Section point %" PetscInt_FMT, numPoints, pEnd); 2534 for (p = pStart; p < pEnd; ++p) { 2535 PetscInt dof, cdof; 2536 2537 PetscCall(PetscSectionGetDof(s, p, &dof)); 2538 PetscCall(PetscSectionSetDof(sNew, perm[p], dof)); 2539 PetscCall(PetscSectionGetConstraintDof(s, p, &cdof)); 2540 if (cdof) PetscCall(PetscSectionSetConstraintDof(sNew, perm[p], cdof)); 2541 for (f = 0; f < numFields; ++f) { 2542 PetscCall(PetscSectionGetFieldDof(s, p, f, &dof)); 2543 PetscCall(PetscSectionSetFieldDof(sNew, perm[p], f, dof)); 2544 PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cdof)); 2545 if (cdof) PetscCall(PetscSectionSetFieldConstraintDof(sNew, perm[p], f, cdof)); 2546 } 2547 } 2548 PetscCall(PetscSectionSetUp(sNew)); 2549 for (p = pStart; p < pEnd; ++p) { 2550 const PetscInt *cind; 2551 PetscInt cdof; 2552 2553 PetscCall(PetscSectionGetConstraintDof(s, p, &cdof)); 2554 if (cdof) { 2555 PetscCall(PetscSectionGetConstraintIndices(s, p, &cind)); 2556 PetscCall(PetscSectionSetConstraintIndices(sNew, perm[p], cind)); 2557 } 2558 for (f = 0; f < numFields; ++f) { 2559 PetscCall(PetscSectionGetFieldConstraintDof(s, p, f, &cdof)); 2560 if (cdof) { 2561 PetscCall(PetscSectionGetFieldConstraintIndices(s, p, f, &cind)); 2562 PetscCall(PetscSectionSetFieldConstraintIndices(sNew, perm[p], f, cind)); 2563 } 2564 } 2565 } 2566 PetscCall(ISRestoreIndices(permutation, &perm)); 2567 *sectionNew = sNew; 2568 PetscFunctionReturn(0); 2569 } 2570 2571 /*@ 2572 PetscSectionSetClosureIndex - Set a cache of points in the closure of each point in the section 2573 2574 Collective 2575 2576 Input Parameters: 2577 + section - The PetscSection 2578 . obj - A PetscObject which serves as the key for this index 2579 . clSection - Section giving the size of the closure of each point 2580 - clPoints - IS giving the points in each closure 2581 2582 Note: We compress out closure points with no dofs in this section 2583 2584 Level: advanced 2585 2586 .seealso: PetscSectionGetClosureIndex(), DMPlexCreateClosureIndex() 2587 @*/ 2588 PetscErrorCode PetscSectionSetClosureIndex(PetscSection section, PetscObject obj, PetscSection clSection, IS clPoints) 2589 { 2590 PetscFunctionBegin; 2591 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 2592 PetscValidHeaderSpecific(clSection,PETSC_SECTION_CLASSID,3); 2593 PetscValidHeaderSpecific(clPoints,IS_CLASSID,4); 2594 if (section->clObj != obj) PetscCall(PetscSectionResetClosurePermutation(section)); 2595 section->clObj = obj; 2596 PetscCall(PetscObjectReference((PetscObject)clSection)); 2597 PetscCall(PetscObjectReference((PetscObject)clPoints)); 2598 PetscCall(PetscSectionDestroy(§ion->clSection)); 2599 PetscCall(ISDestroy(§ion->clPoints)); 2600 section->clSection = clSection; 2601 section->clPoints = clPoints; 2602 PetscFunctionReturn(0); 2603 } 2604 2605 /*@ 2606 PetscSectionGetClosureIndex - Get the cache of points in the closure of each point in the section 2607 2608 Collective 2609 2610 Input Parameters: 2611 + section - The PetscSection 2612 - obj - A PetscObject which serves as the key for this index 2613 2614 Output Parameters: 2615 + clSection - Section giving the size of the closure of each point 2616 - clPoints - IS giving the points in each closure 2617 2618 Note: We compress out closure points with no dofs in this section 2619 2620 Level: advanced 2621 2622 .seealso: PetscSectionSetClosureIndex(), DMPlexCreateClosureIndex() 2623 @*/ 2624 PetscErrorCode PetscSectionGetClosureIndex(PetscSection section, PetscObject obj, PetscSection *clSection, IS *clPoints) 2625 { 2626 PetscFunctionBegin; 2627 if (section->clObj == obj) { 2628 if (clSection) *clSection = section->clSection; 2629 if (clPoints) *clPoints = section->clPoints; 2630 } else { 2631 if (clSection) *clSection = NULL; 2632 if (clPoints) *clPoints = NULL; 2633 } 2634 PetscFunctionReturn(0); 2635 } 2636 2637 PetscErrorCode PetscSectionSetClosurePermutation_Internal(PetscSection section, PetscObject obj, PetscInt depth, PetscInt clSize, PetscCopyMode mode, PetscInt *clPerm) 2638 { 2639 PetscInt i; 2640 khiter_t iter; 2641 int new_entry; 2642 PetscSectionClosurePermKey key = {depth, clSize}; 2643 PetscSectionClosurePermVal *val; 2644 2645 PetscFunctionBegin; 2646 if (section->clObj != obj) { 2647 PetscCall(PetscSectionDestroy(§ion->clSection)); 2648 PetscCall(ISDestroy(§ion->clPoints)); 2649 } 2650 section->clObj = obj; 2651 if (!section->clHash) PetscCall(PetscClPermCreate(§ion->clHash)); 2652 iter = kh_put(ClPerm, section->clHash, key, &new_entry); 2653 val = &kh_val(section->clHash, iter); 2654 if (!new_entry) { 2655 PetscCall(PetscFree(val->perm)); 2656 PetscCall(PetscFree(val->invPerm)); 2657 } 2658 if (mode == PETSC_COPY_VALUES) { 2659 PetscCall(PetscMalloc1(clSize, &val->perm)); 2660 PetscCall(PetscLogObjectMemory((PetscObject) obj, clSize*sizeof(PetscInt))); 2661 PetscCall(PetscArraycpy(val->perm, clPerm, clSize)); 2662 } else if (mode == PETSC_OWN_POINTER) { 2663 val->perm = clPerm; 2664 } else SETERRQ(PetscObjectComm(obj), PETSC_ERR_SUP, "Do not support borrowed arrays"); 2665 PetscCall(PetscMalloc1(clSize, &val->invPerm)); 2666 for (i = 0; i < clSize; ++i) val->invPerm[clPerm[i]] = i; 2667 PetscFunctionReturn(0); 2668 } 2669 2670 /*@ 2671 PetscSectionSetClosurePermutation - Set the dof permutation for the closure of each cell in the section, meaning clPerm[newIndex] = oldIndex. 2672 2673 Not Collective 2674 2675 Input Parameters: 2676 + section - The PetscSection 2677 . obj - A PetscObject which serves as the key for this index (usually a DM) 2678 . depth - Depth of points on which to apply the given permutation 2679 - perm - Permutation of the cell dof closure 2680 2681 Note: 2682 The specified permutation will only be applied to points at depth whose closure size matches the length of perm. In a 2683 mixed-topology or variable-degree finite element space, this function can be called multiple times at each depth for 2684 each topology and degree. 2685 2686 This approach assumes that (depth, len(perm)) uniquely identifies the desired permutation; this might not be true for 2687 exotic/enriched spaces on mixed topology meshes. 2688 2689 Level: intermediate 2690 2691 .seealso: PetscSectionGetClosurePermutation(), PetscSectionGetClosureIndex(), DMPlexCreateClosureIndex(), PetscCopyMode 2692 @*/ 2693 PetscErrorCode PetscSectionSetClosurePermutation(PetscSection section, PetscObject obj, PetscInt depth, IS perm) 2694 { 2695 const PetscInt *clPerm = NULL; 2696 PetscInt clSize = 0; 2697 2698 PetscFunctionBegin; 2699 if (perm) { 2700 PetscCall(ISGetLocalSize(perm, &clSize)); 2701 PetscCall(ISGetIndices(perm, &clPerm)); 2702 } 2703 PetscCall(PetscSectionSetClosurePermutation_Internal(section, obj, depth, clSize, PETSC_COPY_VALUES, (PetscInt *) clPerm)); 2704 if (perm) PetscCall(ISRestoreIndices(perm, &clPerm)); 2705 PetscFunctionReturn(0); 2706 } 2707 2708 PetscErrorCode PetscSectionGetClosurePermutation_Internal(PetscSection section, PetscObject obj, PetscInt depth, PetscInt size, const PetscInt *perm[]) 2709 { 2710 PetscFunctionBegin; 2711 if (section->clObj == obj) { 2712 PetscSectionClosurePermKey k = {depth, size}; 2713 PetscSectionClosurePermVal v; 2714 PetscCall(PetscClPermGet(section->clHash, k, &v)); 2715 if (perm) *perm = v.perm; 2716 } else { 2717 if (perm) *perm = NULL; 2718 } 2719 PetscFunctionReturn(0); 2720 } 2721 2722 /*@ 2723 PetscSectionGetClosurePermutation - Get the dof permutation for the closure of each cell in the section, meaning clPerm[newIndex] = oldIndex. 2724 2725 Not Collective 2726 2727 Input Parameters: 2728 + section - The PetscSection 2729 . obj - A PetscObject which serves as the key for this index (usually a DM) 2730 . depth - Depth stratum on which to obtain closure permutation 2731 - clSize - Closure size to be permuted (e.g., may vary with element topology and degree) 2732 2733 Output Parameter: 2734 . perm - The dof closure permutation 2735 2736 Note: 2737 The user must destroy the IS that is returned. 2738 2739 Level: intermediate 2740 2741 .seealso: PetscSectionSetClosurePermutation(), PetscSectionGetClosureInversePermutation(), PetscSectionGetClosureIndex(), PetscSectionSetClosureIndex(), DMPlexCreateClosureIndex() 2742 @*/ 2743 PetscErrorCode PetscSectionGetClosurePermutation(PetscSection section, PetscObject obj, PetscInt depth, PetscInt clSize, IS *perm) 2744 { 2745 const PetscInt *clPerm; 2746 2747 PetscFunctionBegin; 2748 PetscCall(PetscSectionGetClosurePermutation_Internal(section, obj, depth, clSize, &clPerm)); 2749 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, clSize, clPerm, PETSC_USE_POINTER, perm)); 2750 PetscFunctionReturn(0); 2751 } 2752 2753 PetscErrorCode PetscSectionGetClosureInversePermutation_Internal(PetscSection section, PetscObject obj, PetscInt depth, PetscInt size, const PetscInt *perm[]) 2754 { 2755 PetscFunctionBegin; 2756 if (section->clObj == obj && section->clHash) { 2757 PetscSectionClosurePermKey k = {depth, size}; 2758 PetscSectionClosurePermVal v; 2759 PetscCall(PetscClPermGet(section->clHash, k, &v)); 2760 if (perm) *perm = v.invPerm; 2761 } else { 2762 if (perm) *perm = NULL; 2763 } 2764 PetscFunctionReturn(0); 2765 } 2766 2767 /*@ 2768 PetscSectionGetClosureInversePermutation - Get the inverse dof permutation for the closure of each cell in the section, meaning clPerm[oldIndex] = newIndex. 2769 2770 Not Collective 2771 2772 Input Parameters: 2773 + section - The PetscSection 2774 . obj - A PetscObject which serves as the key for this index (usually a DM) 2775 . depth - Depth stratum on which to obtain closure permutation 2776 - clSize - Closure size to be permuted (e.g., may vary with element topology and degree) 2777 2778 Output Parameters: 2779 . perm - The dof closure permutation 2780 2781 Note: 2782 The user must destroy the IS that is returned. 2783 2784 Level: intermediate 2785 2786 .seealso: PetscSectionSetClosurePermutation(), PetscSectionGetClosureIndex(), PetscSectionSetClosureIndex(), DMPlexCreateClosureIndex() 2787 @*/ 2788 PetscErrorCode PetscSectionGetClosureInversePermutation(PetscSection section, PetscObject obj, PetscInt depth, PetscInt clSize, IS *perm) 2789 { 2790 const PetscInt *clPerm; 2791 2792 PetscFunctionBegin; 2793 PetscCall(PetscSectionGetClosureInversePermutation_Internal(section, obj, depth, clSize, &clPerm)); 2794 PetscCall(ISCreateGeneral(PETSC_COMM_SELF, clSize, clPerm, PETSC_USE_POINTER, perm)); 2795 PetscFunctionReturn(0); 2796 } 2797 2798 /*@ 2799 PetscSectionGetField - Get the subsection associated with a single field 2800 2801 Input Parameters: 2802 + s - The PetscSection 2803 - field - The field number 2804 2805 Output Parameter: 2806 . subs - The subsection for the given field 2807 2808 Level: intermediate 2809 2810 .seealso: PetscSectionSetNumFields() 2811 @*/ 2812 PetscErrorCode PetscSectionGetField(PetscSection s, PetscInt field, PetscSection *subs) 2813 { 2814 PetscFunctionBegin; 2815 PetscValidHeaderSpecific(s,PETSC_SECTION_CLASSID,1); 2816 PetscValidPointer(subs,3); 2817 PetscSectionCheckValidField(field,s->numFields); 2818 *subs = s->field[field]; 2819 PetscFunctionReturn(0); 2820 } 2821 2822 PetscClassId PETSC_SECTION_SYM_CLASSID; 2823 PetscFunctionList PetscSectionSymList = NULL; 2824 2825 /*@ 2826 PetscSectionSymCreate - Creates an empty PetscSectionSym object. 2827 2828 Collective 2829 2830 Input Parameter: 2831 . comm - the MPI communicator 2832 2833 Output Parameter: 2834 . sym - pointer to the new set of symmetries 2835 2836 Level: developer 2837 2838 .seealso: PetscSectionSym, PetscSectionSymDestroy() 2839 @*/ 2840 PetscErrorCode PetscSectionSymCreate(MPI_Comm comm, PetscSectionSym *sym) 2841 { 2842 PetscFunctionBegin; 2843 PetscValidPointer(sym,2); 2844 PetscCall(ISInitializePackage()); 2845 PetscCall(PetscHeaderCreate(*sym,PETSC_SECTION_SYM_CLASSID,"PetscSectionSym","Section Symmetry","IS",comm,PetscSectionSymDestroy,PetscSectionSymView)); 2846 PetscFunctionReturn(0); 2847 } 2848 2849 /*@C 2850 PetscSectionSymSetType - Builds a PetscSection symmetry, for a particular implementation. 2851 2852 Collective 2853 2854 Input Parameters: 2855 + sym - The section symmetry object 2856 - method - The name of the section symmetry type 2857 2858 Level: developer 2859 2860 .seealso: PetscSectionSymGetType(), PetscSectionSymCreate() 2861 @*/ 2862 PetscErrorCode PetscSectionSymSetType(PetscSectionSym sym, PetscSectionSymType method) 2863 { 2864 PetscErrorCode (*r)(PetscSectionSym); 2865 PetscBool match; 2866 2867 PetscFunctionBegin; 2868 PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID,1); 2869 PetscCall(PetscObjectTypeCompare((PetscObject) sym, method, &match)); 2870 if (match) PetscFunctionReturn(0); 2871 2872 PetscCall(PetscFunctionListFind(PetscSectionSymList,method,&r)); 2873 PetscCheck(r,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscSectionSym type: %s", method); 2874 if (sym->ops->destroy) { 2875 PetscCall((*sym->ops->destroy)(sym)); 2876 sym->ops->destroy = NULL; 2877 } 2878 PetscCall((*r)(sym)); 2879 PetscCall(PetscObjectChangeTypeName((PetscObject)sym,method)); 2880 PetscFunctionReturn(0); 2881 } 2882 2883 /*@C 2884 PetscSectionSymGetType - Gets the section symmetry type name (as a string) from the PetscSectionSym. 2885 2886 Not Collective 2887 2888 Input Parameter: 2889 . sym - The section symmetry 2890 2891 Output Parameter: 2892 . type - The index set type name 2893 2894 Level: developer 2895 2896 .seealso: PetscSectionSymSetType(), PetscSectionSymCreate() 2897 @*/ 2898 PetscErrorCode PetscSectionSymGetType(PetscSectionSym sym, PetscSectionSymType *type) 2899 { 2900 PetscFunctionBegin; 2901 PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID,1); 2902 PetscValidPointer(type,2); 2903 *type = ((PetscObject)sym)->type_name; 2904 PetscFunctionReturn(0); 2905 } 2906 2907 /*@C 2908 PetscSectionSymRegister - Adds a new section symmetry implementation 2909 2910 Not Collective 2911 2912 Input Parameters: 2913 + name - The name of a new user-defined creation routine 2914 - create_func - The creation routine itself 2915 2916 Notes: 2917 PetscSectionSymRegister() may be called multiple times to add several user-defined vectors 2918 2919 Level: developer 2920 2921 .seealso: PetscSectionSymCreate(), PetscSectionSymSetType() 2922 @*/ 2923 PetscErrorCode PetscSectionSymRegister(const char sname[], PetscErrorCode (*function)(PetscSectionSym)) 2924 { 2925 PetscFunctionBegin; 2926 PetscCall(ISInitializePackage()); 2927 PetscCall(PetscFunctionListAdd(&PetscSectionSymList,sname,function)); 2928 PetscFunctionReturn(0); 2929 } 2930 2931 /*@ 2932 PetscSectionSymDestroy - Destroys a section symmetry. 2933 2934 Collective 2935 2936 Input Parameters: 2937 . sym - the section symmetry 2938 2939 Level: developer 2940 2941 .seealso: PetscSectionSymCreate(), PetscSectionSymDestroy() 2942 @*/ 2943 PetscErrorCode PetscSectionSymDestroy(PetscSectionSym *sym) 2944 { 2945 SymWorkLink link,next; 2946 2947 PetscFunctionBegin; 2948 if (!*sym) PetscFunctionReturn(0); 2949 PetscValidHeaderSpecific((*sym),PETSC_SECTION_SYM_CLASSID,1); 2950 if (--((PetscObject)(*sym))->refct > 0) {*sym = NULL; PetscFunctionReturn(0);} 2951 if ((*sym)->ops->destroy) { 2952 PetscCall((*(*sym)->ops->destroy)(*sym)); 2953 } 2954 PetscCheckFalse((*sym)->workout,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Work array still checked out"); 2955 for (link=(*sym)->workin; link; link=next) { 2956 next = link->next; 2957 PetscCall(PetscFree2(link->perms,link->rots)); 2958 PetscCall(PetscFree(link)); 2959 } 2960 (*sym)->workin = NULL; 2961 PetscCall(PetscHeaderDestroy(sym)); 2962 PetscFunctionReturn(0); 2963 } 2964 2965 /*@C 2966 PetscSectionSymView - Displays a section symmetry 2967 2968 Collective 2969 2970 Input Parameters: 2971 + sym - the index set 2972 - viewer - viewer used to display the set, for example PETSC_VIEWER_STDOUT_SELF. 2973 2974 Level: developer 2975 2976 .seealso: PetscViewerASCIIOpen() 2977 @*/ 2978 PetscErrorCode PetscSectionSymView(PetscSectionSym sym,PetscViewer viewer) 2979 { 2980 PetscFunctionBegin; 2981 PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1); 2982 if (!viewer) { 2983 PetscCall(PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sym),&viewer)); 2984 } 2985 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 2986 PetscCheckSameComm(sym,1,viewer,2); 2987 PetscCall(PetscObjectPrintClassNamePrefixType((PetscObject)sym,viewer)); 2988 if (sym->ops->view) { 2989 PetscCall((*sym->ops->view)(sym,viewer)); 2990 } 2991 PetscFunctionReturn(0); 2992 } 2993 2994 /*@ 2995 PetscSectionSetSym - Set the symmetries for the data referred to by the section 2996 2997 Collective 2998 2999 Input Parameters: 3000 + section - the section describing data layout 3001 - sym - the symmetry describing the affect of orientation on the access of the data 3002 3003 Level: developer 3004 3005 .seealso: PetscSectionGetSym(), PetscSectionSymCreate() 3006 @*/ 3007 PetscErrorCode PetscSectionSetSym(PetscSection section, PetscSectionSym sym) 3008 { 3009 PetscFunctionBegin; 3010 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3011 PetscCall(PetscSectionSymDestroy(&(section->sym))); 3012 if (sym) { 3013 PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,2); 3014 PetscCheckSameComm(section,1,sym,2); 3015 PetscCall(PetscObjectReference((PetscObject) sym)); 3016 } 3017 section->sym = sym; 3018 PetscFunctionReturn(0); 3019 } 3020 3021 /*@ 3022 PetscSectionGetSym - Get the symmetries for the data referred to by the section 3023 3024 Not Collective 3025 3026 Input Parameters: 3027 . section - the section describing data layout 3028 3029 Output Parameters: 3030 . sym - the symmetry describing the affect of orientation on the access of the data 3031 3032 Level: developer 3033 3034 .seealso: PetscSectionSetSym(), PetscSectionSymCreate() 3035 @*/ 3036 PetscErrorCode PetscSectionGetSym(PetscSection section, PetscSectionSym *sym) 3037 { 3038 PetscFunctionBegin; 3039 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3040 *sym = section->sym; 3041 PetscFunctionReturn(0); 3042 } 3043 3044 /*@ 3045 PetscSectionSetFieldSym - Set the symmetries for the data referred to by a field of the section 3046 3047 Collective 3048 3049 Input Parameters: 3050 + section - the section describing data layout 3051 . field - the field number 3052 - sym - the symmetry describing the affect of orientation on the access of the data 3053 3054 Level: developer 3055 3056 .seealso: PetscSectionGetFieldSym(), PetscSectionSymCreate() 3057 @*/ 3058 PetscErrorCode PetscSectionSetFieldSym(PetscSection section, PetscInt field, PetscSectionSym sym) 3059 { 3060 PetscFunctionBegin; 3061 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3062 PetscSectionCheckValidField(field,section->numFields); 3063 PetscCall(PetscSectionSetSym(section->field[field],sym)); 3064 PetscFunctionReturn(0); 3065 } 3066 3067 /*@ 3068 PetscSectionGetFieldSym - Get the symmetries for the data referred to by a field of the section 3069 3070 Collective 3071 3072 Input Parameters: 3073 + section - the section describing data layout 3074 - field - the field number 3075 3076 Output Parameters: 3077 . sym - the symmetry describing the affect of orientation on the access of the data 3078 3079 Level: developer 3080 3081 .seealso: PetscSectionSetFieldSym(), PetscSectionSymCreate() 3082 @*/ 3083 PetscErrorCode PetscSectionGetFieldSym(PetscSection section, PetscInt field, PetscSectionSym *sym) 3084 { 3085 PetscFunctionBegin; 3086 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3087 PetscSectionCheckValidField(field,section->numFields); 3088 *sym = section->field[field]->sym; 3089 PetscFunctionReturn(0); 3090 } 3091 3092 /*@C 3093 PetscSectionGetPointSyms - Get the symmetries for a set of points in a PetscSection under specific orientations. 3094 3095 Not Collective 3096 3097 Input Parameters: 3098 + section - the section 3099 . numPoints - the number of points 3100 - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an 3101 arbitrary integer: its interpretation is up to sym. Orientations are used by DM: for their interpretation in that 3102 context, see DMPlexGetConeOrientation()). 3103 3104 Output Parameters: 3105 + perms - The permutations for the given orientations (or NULL if there is no symmetry or the permutation is the identity). 3106 - rots - The field rotations symmetries for the given orientations (or NULL if there is no symmetry or the rotations are all 3107 identity). 3108 3109 Example of usage, gathering dofs into a local array (lArray) from a section array (sArray): 3110 .vb 3111 const PetscInt **perms; 3112 const PetscScalar **rots; 3113 PetscInt lOffset; 3114 3115 PetscSectionGetPointSyms(section,numPoints,points,&perms,&rots); 3116 for (i = 0, lOffset = 0; i < numPoints; i++) { 3117 PetscInt point = points[2*i], dof, sOffset; 3118 const PetscInt *perm = perms ? perms[i] : NULL; 3119 const PetscScalar *rot = rots ? rots[i] : NULL; 3120 3121 PetscSectionGetDof(section,point,&dof); 3122 PetscSectionGetOffset(section,point,&sOffset); 3123 3124 if (perm) {for (j = 0; j < dof; j++) {lArray[lOffset + perm[j]] = sArray[sOffset + j];}} 3125 else {for (j = 0; j < dof; j++) {lArray[lOffset + j ] = sArray[sOffset + j];}} 3126 if (rot) {for (j = 0; j < dof; j++) {lArray[lOffset + j ] *= rot[j]; }} 3127 lOffset += dof; 3128 } 3129 PetscSectionRestorePointSyms(section,numPoints,points,&perms,&rots); 3130 .ve 3131 3132 Example of usage, adding dofs into a section array (sArray) from a local array (lArray): 3133 .vb 3134 const PetscInt **perms; 3135 const PetscScalar **rots; 3136 PetscInt lOffset; 3137 3138 PetscSectionGetPointSyms(section,numPoints,points,&perms,&rots); 3139 for (i = 0, lOffset = 0; i < numPoints; i++) { 3140 PetscInt point = points[2*i], dof, sOffset; 3141 const PetscInt *perm = perms ? perms[i] : NULL; 3142 const PetscScalar *rot = rots ? rots[i] : NULL; 3143 3144 PetscSectionGetDof(section,point,&dof); 3145 PetscSectionGetOffset(section,point,&sOff); 3146 3147 if (perm) {for (j = 0; j < dof; j++) {sArray[sOffset + j] += lArray[lOffset + perm[j]] * (rot ? PetscConj(rot[perm[j]]) : 1.);}} 3148 else {for (j = 0; j < dof; j++) {sArray[sOffset + j] += lArray[lOffset + j ] * (rot ? PetscConj(rot[ j ]) : 1.);}} 3149 offset += dof; 3150 } 3151 PetscSectionRestorePointSyms(section,numPoints,points,&perms,&rots); 3152 .ve 3153 3154 Level: developer 3155 3156 .seealso: PetscSectionRestorePointSyms(), PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym() 3157 @*/ 3158 PetscErrorCode PetscSectionGetPointSyms(PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots) 3159 { 3160 PetscSectionSym sym; 3161 3162 PetscFunctionBegin; 3163 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3164 if (numPoints) PetscValidIntPointer(points,3); 3165 if (perms) *perms = NULL; 3166 if (rots) *rots = NULL; 3167 sym = section->sym; 3168 if (sym && (perms || rots)) { 3169 SymWorkLink link; 3170 3171 if (sym->workin) { 3172 link = sym->workin; 3173 sym->workin = sym->workin->next; 3174 } else { 3175 PetscCall(PetscNewLog(sym,&link)); 3176 } 3177 if (numPoints > link->numPoints) { 3178 PetscCall(PetscFree2(link->perms,link->rots)); 3179 PetscCall(PetscMalloc2(numPoints,&link->perms,numPoints,&link->rots)); 3180 link->numPoints = numPoints; 3181 } 3182 link->next = sym->workout; 3183 sym->workout = link; 3184 PetscCall(PetscArrayzero((PetscInt**)link->perms,numPoints)); 3185 PetscCall(PetscArrayzero((PetscInt**)link->rots,numPoints)); 3186 PetscCall((*sym->ops->getpoints) (sym, section, numPoints, points, link->perms, link->rots)); 3187 if (perms) *perms = link->perms; 3188 if (rots) *rots = link->rots; 3189 } 3190 PetscFunctionReturn(0); 3191 } 3192 3193 /*@C 3194 PetscSectionRestorePointSyms - Restore the symmetries returned by PetscSectionGetPointSyms() 3195 3196 Not Collective 3197 3198 Input Parameters: 3199 + section - the section 3200 . numPoints - the number of points 3201 - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an 3202 arbitrary integer: its interpretation is up to sym. Orientations are used by DM: for their interpretation in that 3203 context, see DMPlexGetConeOrientation()). 3204 3205 Output Parameters: 3206 + perms - The permutations for the given orientations: set to NULL at conclusion 3207 - rots - The field rotations symmetries for the given orientations: set to NULL at conclusion 3208 3209 Level: developer 3210 3211 .seealso: PetscSectionGetPointSyms(), PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym() 3212 @*/ 3213 PetscErrorCode PetscSectionRestorePointSyms(PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots) 3214 { 3215 PetscSectionSym sym; 3216 3217 PetscFunctionBegin; 3218 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3219 sym = section->sym; 3220 if (sym && (perms || rots)) { 3221 SymWorkLink *p,link; 3222 3223 for (p=&sym->workout; (link=*p); p=&link->next) { 3224 if ((perms && link->perms == *perms) || (rots && link->rots == *rots)) { 3225 *p = link->next; 3226 link->next = sym->workin; 3227 sym->workin = link; 3228 if (perms) *perms = NULL; 3229 if (rots) *rots = NULL; 3230 PetscFunctionReturn(0); 3231 } 3232 } 3233 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Array was not checked out"); 3234 } 3235 PetscFunctionReturn(0); 3236 } 3237 3238 /*@C 3239 PetscSectionGetFieldPointSyms - Get the symmetries for a set of points in a field of a PetscSection under specific orientations. 3240 3241 Not Collective 3242 3243 Input Parameters: 3244 + section - the section 3245 . field - the field of the section 3246 . numPoints - the number of points 3247 - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an 3248 arbitrary integer: its interpretation is up to sym. Orientations are used by DM: for their interpretation in that 3249 context, see DMPlexGetConeOrientation()). 3250 3251 Output Parameters: 3252 + perms - The permutations for the given orientations (or NULL if there is no symmetry or the permutation is the identity). 3253 - rots - The field rotations symmetries for the given orientations (or NULL if there is no symmetry or the rotations are all 3254 identity). 3255 3256 Level: developer 3257 3258 .seealso: PetscSectionGetPointSyms(), PetscSectionRestoreFieldPointSyms() 3259 @*/ 3260 PetscErrorCode PetscSectionGetFieldPointSyms(PetscSection section, PetscInt field, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots) 3261 { 3262 PetscFunctionBegin; 3263 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3264 PetscCheckFalse(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); 3265 PetscCall(PetscSectionGetPointSyms(section->field[field],numPoints,points,perms,rots)); 3266 PetscFunctionReturn(0); 3267 } 3268 3269 /*@C 3270 PetscSectionRestoreFieldPointSyms - Restore the symmetries returned by PetscSectionGetFieldPointSyms() 3271 3272 Not Collective 3273 3274 Input Parameters: 3275 + section - the section 3276 . field - the field number 3277 . numPoints - the number of points 3278 - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an 3279 arbitrary integer: its interpretation is up to sym. Orientations are used by DM: for their interpretation in that 3280 context, see DMPlexGetConeOrientation()). 3281 3282 Output Parameters: 3283 + perms - The permutations for the given orientations: set to NULL at conclusion 3284 - rots - The field rotations symmetries for the given orientations: set to NULL at conclusion 3285 3286 Level: developer 3287 3288 .seealso: PetscSectionRestorePointSyms(), petscSectionGetFieldPointSyms(), PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym() 3289 @*/ 3290 PetscErrorCode PetscSectionRestoreFieldPointSyms(PetscSection section, PetscInt field, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots) 3291 { 3292 PetscFunctionBegin; 3293 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3294 PetscCheckFalse(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); 3295 PetscCall(PetscSectionRestorePointSyms(section->field[field],numPoints,points,perms,rots)); 3296 PetscFunctionReturn(0); 3297 } 3298 3299 /*@ 3300 PetscSectionSymCopy - Copy the symmetries, assuming that the point structure is compatible 3301 3302 Not Collective 3303 3304 Input Parameter: 3305 . sym - the PetscSectionSym 3306 3307 Output Parameter: 3308 . nsym - the equivalent symmetries 3309 3310 Level: developer 3311 3312 .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms() 3313 @*/ 3314 PetscErrorCode PetscSectionSymCopy(PetscSectionSym sym, PetscSectionSym nsym) 3315 { 3316 PetscFunctionBegin; 3317 PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 1); 3318 PetscValidHeaderSpecific(nsym, PETSC_SECTION_SYM_CLASSID, 2); 3319 if (sym->ops->copy) PetscCall((*sym->ops->copy)(sym, nsym)); 3320 PetscFunctionReturn(0); 3321 } 3322 3323 /*@ 3324 PetscSectionSymDistribute - Distribute the symmetries in accordance with the input SF 3325 3326 Collective 3327 3328 Input Parameters: 3329 + sym - the PetscSectionSym 3330 - migrationSF - the distribution map from roots to leaves 3331 3332 Output Parameters: 3333 . dsym - the redistributed symmetries 3334 3335 Level: developer 3336 3337 .seealso: PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym(), PetscSectionSymLabelSetStratum(), PetscSectionGetPointSyms() 3338 @*/ 3339 PetscErrorCode PetscSectionSymDistribute(PetscSectionSym sym, PetscSF migrationSF, PetscSectionSym *dsym) 3340 { 3341 PetscFunctionBegin; 3342 PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID, 1); 3343 PetscValidHeaderSpecific(migrationSF, PETSCSF_CLASSID, 2); 3344 PetscValidPointer(dsym, 3); 3345 if (sym->ops->distribute) PetscCall((*sym->ops->distribute)(sym, migrationSF, dsym)); 3346 PetscFunctionReturn(0); 3347 } 3348 3349 /*@ 3350 PetscSectionGetUseFieldOffsets - Get the flag to use field offsets directly in a global section, rather than just the point offset 3351 3352 Not Collective 3353 3354 Input Parameter: 3355 . s - the global PetscSection 3356 3357 Output Parameters: 3358 . flg - the flag 3359 3360 Level: developer 3361 3362 .seealso: PetscSectionSetChart(), PetscSectionCreate() 3363 @*/ 3364 PetscErrorCode PetscSectionGetUseFieldOffsets(PetscSection s, PetscBool *flg) 3365 { 3366 PetscFunctionBegin; 3367 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 3368 *flg = s->useFieldOff; 3369 PetscFunctionReturn(0); 3370 } 3371 3372 /*@ 3373 PetscSectionSetUseFieldOffsets - Set the flag to use field offsets directly in a global section, rather than just the point offset 3374 3375 Not Collective 3376 3377 Input Parameters: 3378 + s - the global PetscSection 3379 - flg - the flag 3380 3381 Level: developer 3382 3383 .seealso: PetscSectionGetUseFieldOffsets(), PetscSectionSetChart(), PetscSectionCreate() 3384 @*/ 3385 PetscErrorCode PetscSectionSetUseFieldOffsets(PetscSection s, PetscBool flg) 3386 { 3387 PetscFunctionBegin; 3388 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 3389 s->useFieldOff = flg; 3390 PetscFunctionReturn(0); 3391 } 3392 3393 #define PetscSectionExpandPoints_Loop(TYPE) \ 3394 { \ 3395 PetscInt i, n, o0, o1, size; \ 3396 TYPE *a0 = (TYPE*)origArray, *a1; \ 3397 PetscCall(PetscSectionGetStorageSize(s, &size)); \ 3398 PetscCall(PetscMalloc1(size, &a1)); \ 3399 for (i=0; i<npoints; i++) { \ 3400 PetscCall(PetscSectionGetOffset(origSection, points_[i], &o0)); \ 3401 PetscCall(PetscSectionGetOffset(s, i, &o1)); \ 3402 PetscCall(PetscSectionGetDof(s, i, &n)); \ 3403 PetscCall(PetscMemcpy(&a1[o1], &a0[o0], n*unitsize)); \ 3404 } \ 3405 *newArray = (void*)a1; \ 3406 } 3407 3408 /*@ 3409 PetscSectionExtractDofsFromArray - Extracts elements of an array corresponding to DOFs of specified points. 3410 3411 Not Collective 3412 3413 Input Parameters: 3414 + origSection - the PetscSection describing the layout of the array 3415 . dataType - MPI_Datatype describing the data type of the array (currently only MPIU_INT, MPIU_SCALAR, MPIU_REAL) 3416 . origArray - the array; its size must be equal to the storage size of origSection 3417 - points - IS with points to extract; its indices must lie in the chart of origSection 3418 3419 Output Parameters: 3420 + newSection - the new PetscSection desribing the layout of the new array (with points renumbered 0,1,... but preserving numbers of DOFs) 3421 - newArray - the array of the extracted DOFs; its size is the storage size of newSection 3422 3423 Level: developer 3424 3425 .seealso: PetscSectionGetChart(), PetscSectionGetDof(), PetscSectionGetStorageSize(), PetscSectionCreate() 3426 @*/ 3427 PetscErrorCode PetscSectionExtractDofsFromArray(PetscSection origSection, MPI_Datatype dataType, const void *origArray, IS points, PetscSection *newSection, void *newArray[]) 3428 { 3429 PetscSection s; 3430 const PetscInt *points_; 3431 PetscInt i, n, npoints, pStart, pEnd; 3432 PetscMPIInt unitsize; 3433 3434 PetscFunctionBegin; 3435 PetscValidHeaderSpecific(origSection, PETSC_SECTION_CLASSID, 1); 3436 PetscValidPointer(origArray, 3); 3437 PetscValidHeaderSpecific(points, IS_CLASSID, 4); 3438 if (newSection) PetscValidPointer(newSection, 5); 3439 if (newArray) PetscValidPointer(newArray, 6); 3440 PetscCallMPI(MPI_Type_size(dataType, &unitsize)); 3441 PetscCall(ISGetLocalSize(points, &npoints)); 3442 PetscCall(ISGetIndices(points, &points_)); 3443 PetscCall(PetscSectionGetChart(origSection, &pStart, &pEnd)); 3444 PetscCall(PetscSectionCreate(PETSC_COMM_SELF, &s)); 3445 PetscCall(PetscSectionSetChart(s, 0, npoints)); 3446 for (i=0; i<npoints; i++) { 3447 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); 3448 PetscCall(PetscSectionGetDof(origSection, points_[i], &n)); 3449 PetscCall(PetscSectionSetDof(s, i, n)); 3450 } 3451 PetscCall(PetscSectionSetUp(s)); 3452 if (newArray) { 3453 if (dataType == MPIU_INT) {PetscSectionExpandPoints_Loop(PetscInt);} 3454 else if (dataType == MPIU_SCALAR) {PetscSectionExpandPoints_Loop(PetscScalar);} 3455 else if (dataType == MPIU_REAL) {PetscSectionExpandPoints_Loop(PetscReal);} 3456 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "not implemented for this MPI_Datatype"); 3457 } 3458 if (newSection) { 3459 *newSection = s; 3460 } else { 3461 PetscCall(PetscSectionDestroy(&s)); 3462 } 3463 PetscCall(ISRestoreIndices(points, &points_)); 3464 PetscFunctionReturn(0); 3465 } 3466