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