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 <petscsection.h> 7 #include <petscsf.h> 8 #include <petscviewer.h> 9 10 PetscClassId PETSC_SECTION_CLASSID; 11 12 /*@ 13 PetscSectionCreate - Allocates PetscSection space and sets the map contents to the default. 14 15 Collective 16 17 Input Parameters: 18 + comm - the MPI communicator 19 - s - pointer to the section 20 21 Level: beginner 22 23 Notes: 24 Typical calling sequence 25 $ PetscSectionCreate(MPI_Comm,PetscSection *); 26 $ PetscSectionSetNumFields(PetscSection, numFields); 27 $ PetscSectionSetChart(PetscSection,low,high); 28 $ PetscSectionSetDof(PetscSection,point,numdof); 29 $ PetscSectionSetUp(PetscSection); 30 $ PetscSectionGetOffset(PetscSection,point,PetscInt *); 31 $ PetscSectionDestroy(PetscSection); 32 33 The PetscSection object and methods are intended to be used in the PETSc Vec and Mat implementions; it is 34 recommended they not be used in user codes unless you really gain something in their use. 35 36 .seealso: PetscSection, PetscSectionDestroy() 37 @*/ 38 PetscErrorCode PetscSectionCreate(MPI_Comm comm, PetscSection *s) 39 { 40 PetscErrorCode ierr; 41 42 PetscFunctionBegin; 43 PetscValidPointer(s,2); 44 ierr = ISInitializePackage();CHKERRQ(ierr); 45 46 ierr = PetscHeaderCreate(*s,PETSC_SECTION_CLASSID,"PetscSection","Section","IS",comm,PetscSectionDestroy,PetscSectionView);CHKERRQ(ierr); 47 48 (*s)->pStart = -1; 49 (*s)->pEnd = -1; 50 (*s)->perm = NULL; 51 (*s)->pointMajor = PETSC_TRUE; 52 (*s)->includesConstraints = PETSC_TRUE; 53 (*s)->maxDof = 0; 54 (*s)->atlasDof = NULL; 55 (*s)->atlasOff = NULL; 56 (*s)->bc = NULL; 57 (*s)->bcIndices = NULL; 58 (*s)->setup = PETSC_FALSE; 59 (*s)->numFields = 0; 60 (*s)->fieldNames = NULL; 61 (*s)->field = NULL; 62 (*s)->useFieldOff = PETSC_FALSE; 63 (*s)->compNames = NULL; 64 (*s)->clObj = NULL; 65 (*s)->clHash = NULL; 66 (*s)->clSection = NULL; 67 (*s)->clPoints = NULL; 68 PetscFunctionReturn(0); 69 } 70 71 /*@ 72 PetscSectionCopy - Creates a shallow (if possible) copy of the PetscSection 73 74 Collective 75 76 Input Parameter: 77 . section - the PetscSection 78 79 Output Parameter: 80 . newSection - the copy 81 82 Level: intermediate 83 84 .seealso: PetscSection, PetscSectionCreate(), PetscSectionDestroy() 85 @*/ 86 PetscErrorCode PetscSectionCopy(PetscSection section, PetscSection newSection) 87 { 88 PetscSectionSym sym; 89 IS perm; 90 PetscInt numFields, f, c, pStart, pEnd, p; 91 PetscErrorCode ierr; 92 93 PetscFunctionBegin; 94 PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1); 95 PetscValidHeaderSpecific(newSection, PETSC_SECTION_CLASSID, 2); 96 ierr = PetscSectionReset(newSection);CHKERRQ(ierr); 97 ierr = PetscSectionGetNumFields(section, &numFields);CHKERRQ(ierr); 98 if (numFields) {ierr = PetscSectionSetNumFields(newSection, numFields);CHKERRQ(ierr);} 99 for (f = 0; f < numFields; ++f) { 100 const char *fieldName = NULL, *compName = NULL; 101 PetscInt numComp = 0; 102 103 ierr = PetscSectionGetFieldName(section, f, &fieldName);CHKERRQ(ierr); 104 ierr = PetscSectionSetFieldName(newSection, f, fieldName);CHKERRQ(ierr); 105 ierr = PetscSectionGetFieldComponents(section, f, &numComp);CHKERRQ(ierr); 106 ierr = PetscSectionSetFieldComponents(newSection, f, numComp);CHKERRQ(ierr); 107 for (c = 0; c < numComp; ++c) { 108 ierr = PetscSectionGetComponentName(section, f, c, &compName);CHKERRQ(ierr); 109 ierr = PetscSectionSetComponentName(newSection, f, c, compName);CHKERRQ(ierr); 110 } 111 ierr = PetscSectionGetFieldSym(section, f, &sym);CHKERRQ(ierr); 112 ierr = PetscSectionSetFieldSym(newSection, f, sym);CHKERRQ(ierr); 113 } 114 ierr = PetscSectionGetChart(section, &pStart, &pEnd);CHKERRQ(ierr); 115 ierr = PetscSectionSetChart(newSection, pStart, pEnd);CHKERRQ(ierr); 116 ierr = PetscSectionGetPermutation(section, &perm);CHKERRQ(ierr); 117 ierr = PetscSectionSetPermutation(newSection, perm);CHKERRQ(ierr); 118 ierr = PetscSectionGetSym(section, &sym);CHKERRQ(ierr); 119 ierr = PetscSectionSetSym(newSection, sym);CHKERRQ(ierr); 120 for (p = pStart; p < pEnd; ++p) { 121 PetscInt dof, cdof, fcdof = 0; 122 123 ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 124 ierr = PetscSectionSetDof(newSection, p, dof);CHKERRQ(ierr); 125 ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); 126 if (cdof) {ierr = PetscSectionSetConstraintDof(newSection, p, cdof);CHKERRQ(ierr);} 127 for (f = 0; f < numFields; ++f) { 128 ierr = PetscSectionGetFieldDof(section, p, f, &dof);CHKERRQ(ierr); 129 ierr = PetscSectionSetFieldDof(newSection, p, f, dof);CHKERRQ(ierr); 130 if (cdof) { 131 ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr); 132 if (fcdof) {ierr = PetscSectionSetFieldConstraintDof(newSection, p, f, fcdof);CHKERRQ(ierr);} 133 } 134 } 135 } 136 ierr = PetscSectionSetUp(newSection);CHKERRQ(ierr); 137 for (p = pStart; p < pEnd; ++p) { 138 PetscInt off, cdof, fcdof = 0; 139 const PetscInt *cInd; 140 141 /* Must set offsets in case they do not agree with the prefix sums */ 142 ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 143 ierr = PetscSectionSetOffset(newSection, p, off);CHKERRQ(ierr); 144 ierr = PetscSectionGetConstraintDof(section, p, &cdof);CHKERRQ(ierr); 145 if (cdof) { 146 ierr = PetscSectionGetConstraintIndices(section, p, &cInd);CHKERRQ(ierr); 147 ierr = PetscSectionSetConstraintIndices(newSection, p, cInd);CHKERRQ(ierr); 148 for (f = 0; f < numFields; ++f) { 149 ierr = PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);CHKERRQ(ierr); 150 if (fcdof) { 151 ierr = PetscSectionGetFieldConstraintIndices(section, p, f, &cInd);CHKERRQ(ierr); 152 ierr = PetscSectionSetFieldConstraintIndices(newSection, p, f, cInd);CHKERRQ(ierr); 153 } 154 } 155 } 156 } 157 PetscFunctionReturn(0); 158 } 159 160 /*@ 161 PetscSectionClone - Creates a shallow (if possible) copy of the PetscSection 162 163 Collective 164 165 Input Parameter: 166 . section - the PetscSection 167 168 Output Parameter: 169 . newSection - the copy 170 171 Level: beginner 172 173 .seealso: PetscSection, PetscSectionCreate(), PetscSectionDestroy() 174 @*/ 175 PetscErrorCode PetscSectionClone(PetscSection section, PetscSection *newSection) 176 { 177 PetscErrorCode ierr; 178 179 PetscFunctionBegin; 180 PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1); 181 PetscValidPointer(newSection, 2); 182 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) section), newSection);CHKERRQ(ierr); 183 ierr = PetscSectionCopy(section, *newSection);CHKERRQ(ierr); 184 PetscFunctionReturn(0); 185 } 186 187 /*@ 188 PetscSectionSetFromOptions - sets parameters in a PetscSection from the options database 189 190 Collective on PetscSection 191 192 Input Parameter: 193 . section - the PetscSection 194 195 Options Database: 196 . -petscsection_point_major the dof order 197 198 Level: intermediate 199 200 .seealso: PetscSection, PetscSectionCreate(), PetscSectionDestroy() 201 @*/ 202 PetscErrorCode PetscSectionSetFromOptions(PetscSection s) 203 { 204 PetscErrorCode ierr; 205 206 PetscFunctionBegin; 207 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 208 ierr = PetscObjectOptionsBegin((PetscObject) s);CHKERRQ(ierr); 209 ierr = PetscOptionsBool("-petscsection_point_major", "The for ordering, either point major or field major", "PetscSectionSetPointMajor", s->pointMajor, &s->pointMajor, NULL);CHKERRQ(ierr); 210 /* process any options handlers added with PetscObjectAddOptionsHandler() */ 211 ierr = PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject) s);CHKERRQ(ierr); 212 ierr = PetscOptionsEnd();CHKERRQ(ierr); 213 ierr = PetscObjectViewFromOptions((PetscObject) s, NULL, "-petscsection_view");CHKERRQ(ierr); 214 PetscFunctionReturn(0); 215 } 216 217 /*@ 218 PetscSectionCompare - Compares two sections 219 220 Collective on PetscSection 221 222 Input Parameters: 223 + s1 - the first PetscSection 224 - s2 - the second PetscSection 225 226 Output Parameter: 227 . congruent - PETSC_TRUE if the two sections are congruent, PETSC_FALSE otherwise 228 229 Level: intermediate 230 231 Notes: 232 Field names are disregarded. 233 234 .seealso: PetscSection, PetscSectionCreate(), PetscSectionCopy(), PetscSectionClone() 235 @*/ 236 PetscErrorCode PetscSectionCompare(PetscSection s1, PetscSection s2, PetscBool *congruent) 237 { 238 PetscInt pStart, pEnd, nfields, ncdof, nfcdof, p, f, n1, n2; 239 const PetscInt *idx1, *idx2; 240 IS perm1, perm2; 241 PetscBool flg; 242 PetscMPIInt mflg; 243 PetscErrorCode ierr; 244 245 PetscFunctionBegin; 246 PetscValidHeaderSpecific(s1, PETSC_SECTION_CLASSID, 1); 247 PetscValidHeaderSpecific(s2, PETSC_SECTION_CLASSID, 2); 248 PetscValidBoolPointer(congruent,3); 249 flg = PETSC_FALSE; 250 251 ierr = MPI_Comm_compare(PetscObjectComm((PetscObject)s1),PetscObjectComm((PetscObject)s2),&mflg);CHKERRMPI(ierr); 252 if (mflg != MPI_CONGRUENT && mflg != MPI_IDENT) { 253 *congruent = PETSC_FALSE; 254 PetscFunctionReturn(0); 255 } 256 257 ierr = PetscSectionGetChart(s1, &pStart, &pEnd);CHKERRQ(ierr); 258 ierr = PetscSectionGetChart(s2, &n1, &n2);CHKERRQ(ierr); 259 if (pStart != n1 || pEnd != n2) goto not_congruent; 260 261 ierr = PetscSectionGetPermutation(s1, &perm1);CHKERRQ(ierr); 262 ierr = PetscSectionGetPermutation(s2, &perm2);CHKERRQ(ierr); 263 if (perm1 && perm2) { 264 ierr = ISEqual(perm1, perm2, congruent);CHKERRQ(ierr); 265 if (!(*congruent)) goto not_congruent; 266 } else if (perm1 != perm2) goto not_congruent; 267 268 for (p = pStart; p < pEnd; ++p) { 269 ierr = PetscSectionGetOffset(s1, p, &n1);CHKERRQ(ierr); 270 ierr = PetscSectionGetOffset(s2, p, &n2);CHKERRQ(ierr); 271 if (n1 != n2) goto not_congruent; 272 273 ierr = PetscSectionGetDof(s1, p, &n1);CHKERRQ(ierr); 274 ierr = PetscSectionGetDof(s2, p, &n2);CHKERRQ(ierr); 275 if (n1 != n2) goto not_congruent; 276 277 ierr = PetscSectionGetConstraintDof(s1, p, &ncdof);CHKERRQ(ierr); 278 ierr = PetscSectionGetConstraintDof(s2, p, &n2);CHKERRQ(ierr); 279 if (ncdof != n2) goto not_congruent; 280 281 ierr = PetscSectionGetConstraintIndices(s1, p, &idx1);CHKERRQ(ierr); 282 ierr = PetscSectionGetConstraintIndices(s2, p, &idx2);CHKERRQ(ierr); 283 ierr = PetscArraycmp(idx1, idx2, ncdof, congruent);CHKERRQ(ierr); 284 if (!(*congruent)) goto not_congruent; 285 } 286 287 ierr = PetscSectionGetNumFields(s1, &nfields);CHKERRQ(ierr); 288 ierr = PetscSectionGetNumFields(s2, &n2);CHKERRQ(ierr); 289 if (nfields != n2) goto not_congruent; 290 291 for (f = 0; f < nfields; ++f) { 292 ierr = PetscSectionGetFieldComponents(s1, f, &n1);CHKERRQ(ierr); 293 ierr = PetscSectionGetFieldComponents(s2, f, &n2);CHKERRQ(ierr); 294 if (n1 != n2) goto not_congruent; 295 296 for (p = pStart; p < pEnd; ++p) { 297 ierr = PetscSectionGetFieldOffset(s1, p, f, &n1);CHKERRQ(ierr); 298 ierr = PetscSectionGetFieldOffset(s2, p, f, &n2);CHKERRQ(ierr); 299 if (n1 != n2) goto not_congruent; 300 301 ierr = PetscSectionGetFieldDof(s1, p, f, &n1);CHKERRQ(ierr); 302 ierr = PetscSectionGetFieldDof(s2, p, f, &n2);CHKERRQ(ierr); 303 if (n1 != n2) goto not_congruent; 304 305 ierr = PetscSectionGetFieldConstraintDof(s1, p, f, &nfcdof);CHKERRQ(ierr); 306 ierr = PetscSectionGetFieldConstraintDof(s2, p, f, &n2);CHKERRQ(ierr); 307 if (nfcdof != n2) goto not_congruent; 308 309 ierr = PetscSectionGetFieldConstraintIndices(s1, p, f, &idx1);CHKERRQ(ierr); 310 ierr = PetscSectionGetFieldConstraintIndices(s2, p, f, &idx2);CHKERRQ(ierr); 311 ierr = PetscArraycmp(idx1, idx2, nfcdof, congruent);CHKERRQ(ierr); 312 if (!(*congruent)) goto not_congruent; 313 } 314 } 315 316 flg = PETSC_TRUE; 317 not_congruent: 318 ierr = MPIU_Allreduce(&flg,congruent,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)s1));CHKERRMPI(ierr); 319 PetscFunctionReturn(0); 320 } 321 322 /*@ 323 PetscSectionGetNumFields - Returns the number of fields, or 0 if no fields were defined. 324 325 Not collective 326 327 Input Parameter: 328 . s - the PetscSection 329 330 Output Parameter: 331 . numFields - the number of fields defined, or 0 if none were defined 332 333 Level: intermediate 334 335 .seealso: PetscSectionSetNumFields() 336 @*/ 337 PetscErrorCode PetscSectionGetNumFields(PetscSection s, PetscInt *numFields) 338 { 339 PetscFunctionBegin; 340 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 341 PetscValidPointer(numFields,2); 342 *numFields = s->numFields; 343 PetscFunctionReturn(0); 344 } 345 346 /*@ 347 PetscSectionSetNumFields - Sets the number of fields. 348 349 Not collective 350 351 Input Parameters: 352 + s - the PetscSection 353 - numFields - the number of fields 354 355 Level: intermediate 356 357 .seealso: PetscSectionGetNumFields() 358 @*/ 359 PetscErrorCode PetscSectionSetNumFields(PetscSection s, PetscInt numFields) 360 { 361 PetscInt f; 362 PetscErrorCode ierr; 363 364 PetscFunctionBegin; 365 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 366 if (numFields <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "The number of fields %D must be positive", numFields); 367 ierr = PetscSectionReset(s);CHKERRQ(ierr); 368 369 s->numFields = numFields; 370 ierr = PetscMalloc1(s->numFields, &s->numFieldComponents);CHKERRQ(ierr); 371 ierr = PetscMalloc1(s->numFields, &s->fieldNames);CHKERRQ(ierr); 372 ierr = PetscMalloc1(s->numFields, &s->compNames);CHKERRQ(ierr); 373 ierr = PetscMalloc1(s->numFields, &s->field);CHKERRQ(ierr); 374 for (f = 0; f < s->numFields; ++f) { 375 char name[64]; 376 377 s->numFieldComponents[f] = 1; 378 379 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), &s->field[f]);CHKERRQ(ierr); 380 ierr = PetscSNPrintf(name, 64, "Field_%D", f);CHKERRQ(ierr); 381 ierr = PetscStrallocpy(name, (char **) &s->fieldNames[f]);CHKERRQ(ierr); 382 ierr = PetscSNPrintf(name, 64, "Component_0");CHKERRQ(ierr); 383 ierr = PetscMalloc1(s->numFieldComponents[f], &s->compNames[f]);CHKERRQ(ierr); 384 ierr = PetscStrallocpy(name, (char **) &s->compNames[f][0]);CHKERRQ(ierr); 385 } 386 PetscFunctionReturn(0); 387 } 388 389 /*@C 390 PetscSectionGetFieldName - Returns the name of a field in the PetscSection 391 392 Not Collective 393 394 Input Parameters: 395 + s - the PetscSection 396 - field - the field number 397 398 Output Parameter: 399 . fieldName - the field name 400 401 Level: intermediate 402 403 .seealso: PetscSectionSetFieldName() 404 @*/ 405 PetscErrorCode PetscSectionGetFieldName(PetscSection s, PetscInt field, const char *fieldName[]) 406 { 407 PetscFunctionBegin; 408 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 409 PetscValidPointer(fieldName, 3); 410 if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields); 411 *fieldName = s->fieldNames[field]; 412 PetscFunctionReturn(0); 413 } 414 415 /*@C 416 PetscSectionSetFieldName - Sets the name of a field in the PetscSection 417 418 Not Collective 419 420 Input Parameters: 421 + s - the PetscSection 422 . field - the field number 423 - fieldName - the field name 424 425 Level: intermediate 426 427 .seealso: PetscSectionGetFieldName() 428 @*/ 429 PetscErrorCode PetscSectionSetFieldName(PetscSection s, PetscInt field, const char fieldName[]) 430 { 431 PetscErrorCode ierr; 432 433 PetscFunctionBegin; 434 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 435 if (fieldName) PetscValidCharPointer(fieldName, 3); 436 if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields); 437 ierr = PetscFree(s->fieldNames[field]);CHKERRQ(ierr); 438 ierr = PetscStrallocpy(fieldName, (char**) &s->fieldNames[field]);CHKERRQ(ierr); 439 PetscFunctionReturn(0); 440 } 441 442 /*@C 443 PetscSectionGetComponentName - Gets the name of a field component in the PetscSection 444 445 Not Collective 446 447 Input Parameters: 448 + s - the PetscSection 449 . field - the field number 450 . comp - the component number 451 - compName - the component name 452 453 Level: intermediate 454 455 .seealso: PetscSectionSetComponentName() 456 @*/ 457 PetscErrorCode PetscSectionGetComponentName(PetscSection s, PetscInt field, PetscInt comp, const char *compName[]) 458 { 459 PetscFunctionBegin; 460 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 461 PetscValidPointer(compName, 4); 462 if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields); 463 if ((comp < 0) || (comp >= s->numFieldComponents[field])) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field component %D should be in [%D, %D)", comp, 0, s->numFieldComponents[field]); 464 *compName = s->compNames[field][comp]; 465 PetscFunctionReturn(0); 466 } 467 468 /*@C 469 PetscSectionSetComponentName - Sets the name of a field component in the PetscSection 470 471 Not Collective 472 473 Input Parameters: 474 + s - the PetscSection 475 . field - the field number 476 . comp - the component number 477 - compName - the component name 478 479 Level: intermediate 480 481 .seealso: PetscSectionGetComponentName() 482 @*/ 483 PetscErrorCode PetscSectionSetComponentName(PetscSection s, PetscInt field, PetscInt comp, const char compName[]) 484 { 485 PetscErrorCode ierr; 486 487 PetscFunctionBegin; 488 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 489 if (compName) PetscValidCharPointer(compName, 4); 490 if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields); 491 if ((comp < 0) || (comp >= s->numFieldComponents[field])) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field component %D should be in [%D, %D)", comp, 0, s->numFieldComponents[field]); 492 ierr = PetscFree(s->compNames[field][comp]);CHKERRQ(ierr); 493 ierr = PetscStrallocpy(compName, (char**) &s->compNames[field][comp]);CHKERRQ(ierr); 494 PetscFunctionReturn(0); 495 } 496 497 /*@ 498 PetscSectionGetFieldComponents - Returns the number of field components for the given field. 499 500 Not collective 501 502 Input Parameters: 503 + s - the PetscSection 504 - field - the field number 505 506 Output Parameter: 507 . numComp - the number of field components 508 509 Level: intermediate 510 511 .seealso: PetscSectionSetFieldComponents(), PetscSectionGetNumFields() 512 @*/ 513 PetscErrorCode PetscSectionGetFieldComponents(PetscSection s, PetscInt field, PetscInt *numComp) 514 { 515 PetscFunctionBegin; 516 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 517 PetscValidPointer(numComp, 3); 518 if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields); 519 *numComp = s->numFieldComponents[field]; 520 PetscFunctionReturn(0); 521 } 522 523 /*@ 524 PetscSectionSetFieldComponents - Sets the number of field components for the given field. 525 526 Not collective 527 528 Input Parameters: 529 + s - the PetscSection 530 . field - the field number 531 - numComp - the number of field components 532 533 Level: intermediate 534 535 .seealso: PetscSectionGetFieldComponents(), PetscSectionGetNumFields() 536 @*/ 537 PetscErrorCode PetscSectionSetFieldComponents(PetscSection s, PetscInt field, PetscInt numComp) 538 { 539 PetscErrorCode ierr; 540 PetscInt c; 541 char name[64]; 542 543 PetscFunctionBegin; 544 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 545 if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields); 546 if (s->compNames) { 547 for (c = 0; c < s->numFieldComponents[field]; ++c) { 548 ierr = PetscFree(s->compNames[field][c]);CHKERRQ(ierr); 549 } 550 ierr = PetscFree(s->compNames[field]);CHKERRQ(ierr); 551 } 552 553 s->numFieldComponents[field] = numComp; 554 if (numComp) { 555 ierr = PetscMalloc1(numComp, (char ***) &s->compNames[field]);CHKERRQ(ierr); 556 for (c = 0; c < numComp; ++c) { 557 ierr = PetscSNPrintf(name, 64, "%D", c);CHKERRQ(ierr); 558 ierr = PetscStrallocpy(name, (char **) &s->compNames[field][c]);CHKERRQ(ierr); 559 } 560 } 561 PetscFunctionReturn(0); 562 } 563 564 static PetscErrorCode PetscSectionCheckConstraints_Static(PetscSection s) 565 { 566 PetscErrorCode ierr; 567 568 PetscFunctionBegin; 569 if (!s->bc) { 570 ierr = PetscSectionCreate(PETSC_COMM_SELF, &s->bc);CHKERRQ(ierr); 571 ierr = PetscSectionSetChart(s->bc, s->pStart, s->pEnd);CHKERRQ(ierr); 572 } 573 PetscFunctionReturn(0); 574 } 575 576 /*@ 577 PetscSectionGetChart - Returns the range [pStart, pEnd) in which points lie. 578 579 Not collective 580 581 Input Parameter: 582 . s - the PetscSection 583 584 Output Parameters: 585 + pStart - the first point 586 - pEnd - one past the last point 587 588 Level: intermediate 589 590 .seealso: PetscSectionSetChart(), PetscSectionCreate() 591 @*/ 592 PetscErrorCode PetscSectionGetChart(PetscSection s, PetscInt *pStart, PetscInt *pEnd) 593 { 594 PetscFunctionBegin; 595 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 596 if (pStart) *pStart = s->pStart; 597 if (pEnd) *pEnd = s->pEnd; 598 PetscFunctionReturn(0); 599 } 600 601 /*@ 602 PetscSectionSetChart - Sets the range [pStart, pEnd) in which points lie. 603 604 Not collective 605 606 Input Parameters: 607 + s - the PetscSection 608 . pStart - the first point 609 - pEnd - one past the last point 610 611 Level: intermediate 612 613 .seealso: PetscSectionGetChart(), PetscSectionCreate() 614 @*/ 615 PetscErrorCode PetscSectionSetChart(PetscSection s, PetscInt pStart, PetscInt pEnd) 616 { 617 PetscInt f; 618 PetscErrorCode ierr; 619 620 PetscFunctionBegin; 621 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 622 if (pStart == s->pStart && pEnd == s->pEnd) PetscFunctionReturn(0); 623 /* Cannot Reset() because it destroys field information */ 624 s->setup = PETSC_FALSE; 625 ierr = PetscSectionDestroy(&s->bc);CHKERRQ(ierr); 626 ierr = PetscFree(s->bcIndices);CHKERRQ(ierr); 627 ierr = PetscFree2(s->atlasDof, s->atlasOff);CHKERRQ(ierr); 628 629 s->pStart = pStart; 630 s->pEnd = pEnd; 631 ierr = PetscMalloc2((pEnd - pStart), &s->atlasDof, (pEnd - pStart), &s->atlasOff);CHKERRQ(ierr); 632 ierr = PetscArrayzero(s->atlasDof, pEnd - pStart);CHKERRQ(ierr); 633 for (f = 0; f < s->numFields; ++f) { 634 ierr = PetscSectionSetChart(s->field[f], pStart, pEnd);CHKERRQ(ierr); 635 } 636 PetscFunctionReturn(0); 637 } 638 639 /*@ 640 PetscSectionGetPermutation - Returns the permutation of [0, pEnd-pStart) or NULL 641 642 Not collective 643 644 Input Parameter: 645 . s - the PetscSection 646 647 Output Parameters: 648 . perm - The permutation as an IS 649 650 Level: intermediate 651 652 .seealso: PetscSectionSetPermutation(), PetscSectionCreate() 653 @*/ 654 PetscErrorCode PetscSectionGetPermutation(PetscSection s, IS *perm) 655 { 656 PetscFunctionBegin; 657 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 658 if (perm) {PetscValidPointer(perm, 2); *perm = s->perm;} 659 PetscFunctionReturn(0); 660 } 661 662 /*@ 663 PetscSectionSetPermutation - Sets the permutation for [0, pEnd-pStart) 664 665 Not collective 666 667 Input Parameters: 668 + s - the PetscSection 669 - perm - the permutation of points 670 671 Level: intermediate 672 673 .seealso: PetscSectionGetPermutation(), PetscSectionCreate() 674 @*/ 675 PetscErrorCode PetscSectionSetPermutation(PetscSection s, IS perm) 676 { 677 PetscErrorCode ierr; 678 679 PetscFunctionBegin; 680 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 681 if (perm) PetscValidHeaderSpecific(perm, IS_CLASSID, 2); 682 if (s->setup) SETERRQ(PetscObjectComm((PetscObject) s), PETSC_ERR_ARG_WRONGSTATE, "Cannot set a permutation after the section is setup"); 683 if (s->perm != perm) { 684 ierr = ISDestroy(&s->perm);CHKERRQ(ierr); 685 if (perm) { 686 s->perm = perm; 687 ierr = PetscObjectReference((PetscObject) s->perm);CHKERRQ(ierr); 688 } 689 } 690 PetscFunctionReturn(0); 691 } 692 693 /*@ 694 PetscSectionGetPointMajor - Returns the flag for dof ordering, true if it is point major, otherwise field major 695 696 Not collective 697 698 Input Parameter: 699 . s - the PetscSection 700 701 Output Parameter: 702 . pm - the flag for point major ordering 703 704 Level: intermediate 705 706 .seealso: PetscSectionSetPointMajor() 707 @*/ 708 PetscErrorCode PetscSectionGetPointMajor(PetscSection s, PetscBool *pm) 709 { 710 PetscFunctionBegin; 711 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 712 PetscValidPointer(pm,2); 713 *pm = s->pointMajor; 714 PetscFunctionReturn(0); 715 } 716 717 /*@ 718 PetscSectionSetPointMajor - Sets the flag for dof ordering, true if it is point major, otherwise field major 719 720 Not collective 721 722 Input Parameters: 723 + s - the PetscSection 724 - pm - the flag for point major ordering 725 726 Not collective 727 728 Level: intermediate 729 730 .seealso: PetscSectionGetPointMajor() 731 @*/ 732 PetscErrorCode PetscSectionSetPointMajor(PetscSection s, PetscBool pm) 733 { 734 PetscFunctionBegin; 735 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 736 if (s->setup) SETERRQ(PetscObjectComm((PetscObject) s), PETSC_ERR_ARG_WRONGSTATE, "Cannot set the dof ordering after the section is setup"); 737 s->pointMajor = pm; 738 PetscFunctionReturn(0); 739 } 740 741 /*@ 742 PetscSectionGetIncludesConstraints - Returns the flag indicating if constrained dofs were included when computing offsets 743 744 Not collective 745 746 Input Parameter: 747 . s - the PetscSection 748 749 Output Parameter: 750 . includesConstraints - the flag indicating if constrained dofs were included when computing offsets 751 752 Level: intermediate 753 754 .seealso: PetscSectionSetIncludesConstraints() 755 @*/ 756 PetscErrorCode PetscSectionGetIncludesConstraints(PetscSection s, PetscBool *includesConstraints) 757 { 758 PetscFunctionBegin; 759 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 760 PetscValidBoolPointer(includesConstraints,2); 761 *includesConstraints = s->includesConstraints; 762 PetscFunctionReturn(0); 763 } 764 765 /*@ 766 PetscSectionSetIncludesConstraints - Sets the flag indicating if constrained dofs are to be included when computing offsets 767 768 Not collective 769 770 Input Parameters: 771 + s - the PetscSection 772 - includesConstraints - the flag indicating if constrained dofs are to be included when computing offsets 773 774 Not collective 775 776 Level: intermediate 777 778 .seealso: PetscSectionGetIncludesConstraints() 779 @*/ 780 PetscErrorCode PetscSectionSetIncludesConstraints(PetscSection s, PetscBool includesConstraints) 781 { 782 PetscFunctionBegin; 783 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 784 if (s->setup) SETERRQ(PetscObjectComm((PetscObject) s), PETSC_ERR_ARG_WRONGSTATE, "Cannot set includesConstraints after the section is set up"); 785 s->includesConstraints = includesConstraints; 786 PetscFunctionReturn(0); 787 } 788 789 /*@ 790 PetscSectionGetDof - Return the number of degrees of freedom associated with a given point. 791 792 Not collective 793 794 Input Parameters: 795 + s - the PetscSection 796 - point - the point 797 798 Output Parameter: 799 . numDof - the number of dof 800 801 Level: intermediate 802 803 .seealso: PetscSectionSetDof(), PetscSectionCreate() 804 @*/ 805 PetscErrorCode PetscSectionGetDof(PetscSection s, PetscInt point, PetscInt *numDof) 806 { 807 PetscFunctionBeginHot; 808 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 809 PetscValidPointer(numDof, 3); 810 if (PetscDefined(USE_DEBUG)) { 811 if ((point < s->pStart) || (point >= s->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %D should be in [%D, %D)", point, s->pStart, s->pEnd); 812 } 813 *numDof = s->atlasDof[point - s->pStart]; 814 PetscFunctionReturn(0); 815 } 816 817 /*@ 818 PetscSectionSetDof - Sets the number of degrees of freedom associated with a given point. 819 820 Not collective 821 822 Input Parameters: 823 + s - the PetscSection 824 . point - the point 825 - numDof - the number of dof 826 827 Level: intermediate 828 829 .seealso: PetscSectionGetDof(), PetscSectionAddDof(), PetscSectionCreate() 830 @*/ 831 PetscErrorCode PetscSectionSetDof(PetscSection s, PetscInt point, PetscInt numDof) 832 { 833 PetscFunctionBegin; 834 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 835 if ((point < s->pStart) || (point >= s->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %D should be in [%D, %D)", point, s->pStart, s->pEnd); 836 s->atlasDof[point - s->pStart] = numDof; 837 PetscFunctionReturn(0); 838 } 839 840 /*@ 841 PetscSectionAddDof - Adds to the number of degrees of freedom associated with a given point. 842 843 Not collective 844 845 Input Parameters: 846 + s - the PetscSection 847 . point - the point 848 - numDof - the number of additional dof 849 850 Level: intermediate 851 852 .seealso: PetscSectionGetDof(), PetscSectionSetDof(), PetscSectionCreate() 853 @*/ 854 PetscErrorCode PetscSectionAddDof(PetscSection s, PetscInt point, PetscInt numDof) 855 { 856 PetscFunctionBeginHot; 857 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 858 if (PetscDefined(USE_DEBUG)) { 859 if ((point < s->pStart) || (point >= s->pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %D should be in [%D, %D)", point, s->pStart, s->pEnd); 860 } 861 s->atlasDof[point - s->pStart] += numDof; 862 PetscFunctionReturn(0); 863 } 864 865 /*@ 866 PetscSectionGetFieldDof - Return the number of degrees of freedom associated with a field on a given point. 867 868 Not collective 869 870 Input Parameters: 871 + s - the PetscSection 872 . point - the point 873 - field - the field 874 875 Output Parameter: 876 . numDof - the number of dof 877 878 Level: intermediate 879 880 .seealso: PetscSectionSetFieldDof(), PetscSectionCreate() 881 @*/ 882 PetscErrorCode PetscSectionGetFieldDof(PetscSection s, PetscInt point, PetscInt field, PetscInt *numDof) 883 { 884 PetscErrorCode ierr; 885 886 PetscFunctionBegin; 887 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 888 if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, 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 if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, 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 if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, 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 PetscValidPointer(numDof, 4); 1052 if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, 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 if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, 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 if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, 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) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "SF roots %D < pEnd %D", nroots, pEnd); 1331 if (maxleaf >= nroots) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Max local leaf %D >= nroots %D", 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) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Global dof %D for point %D is not the unconstrained %D", -(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) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %D < %D 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)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %D should be in [%D, %D)", 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)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %D should be in [%D, %D)", 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 PetscValidPointer(offset, 4); 1668 if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, 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 if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, 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 PetscValidPointer(offset, 4); 1729 if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, 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 Parameter: 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) SETERRQ2(PetscObjectComm((PetscObject) s), PETSC_ERR_ARG_WRONG, "Number of requested fields %D greater than number of fields %D", 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, " (%4D) dim %2D offset %3D 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, " %D", s->bcIndices[s->bc->atlasOff[p]+b]);CHKERRQ(ierr); 2109 } 2110 } 2111 ierr = PetscViewerASCIISynchronizedPrintf(viewer, "\n");CHKERRQ(ierr); 2112 } else { 2113 ierr = PetscViewerASCIISynchronizedPrintf(viewer, " (%4D) dim %2D offset %3D\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() 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, "%D fields\n", s->numFields);CHKERRQ(ierr); 2186 for (f = 0; f < s->numFields; ++f) { 2187 ierr = PetscViewerASCIIPrintf(viewer, " field %D with %D 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 static PetscErrorCode PetscSectionResetClosurePermutation(PetscSection section) 2204 { 2205 PetscErrorCode ierr; 2206 PetscSectionClosurePermVal clVal; 2207 2208 PetscFunctionBegin; 2209 if (!section->clHash) PetscFunctionReturn(0); 2210 kh_foreach_value(section->clHash, clVal, { 2211 ierr = PetscFree(clVal.perm);CHKERRQ(ierr); 2212 ierr = PetscFree(clVal.invPerm);CHKERRQ(ierr); 2213 }); 2214 kh_destroy(ClPerm, section->clHash); 2215 section->clHash = NULL; 2216 PetscFunctionReturn(0); 2217 } 2218 2219 /*@ 2220 PetscSectionReset - Frees all section data. 2221 2222 Not collective 2223 2224 Input Parameters: 2225 . s - the PetscSection 2226 2227 Level: beginner 2228 2229 .seealso: PetscSection, PetscSectionCreate() 2230 @*/ 2231 PetscErrorCode PetscSectionReset(PetscSection s) 2232 { 2233 PetscInt f, c; 2234 PetscErrorCode ierr; 2235 2236 PetscFunctionBegin; 2237 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2238 for (f = 0; f < s->numFields; ++f) { 2239 ierr = PetscSectionDestroy(&s->field[f]);CHKERRQ(ierr); 2240 ierr = PetscFree(s->fieldNames[f]);CHKERRQ(ierr); 2241 for (c = 0; c < s->numFieldComponents[f]; ++c) { 2242 ierr = PetscFree(s->compNames[f][c]);CHKERRQ(ierr); 2243 } 2244 ierr = PetscFree(s->compNames[f]);CHKERRQ(ierr); 2245 } 2246 ierr = PetscFree(s->numFieldComponents);CHKERRQ(ierr); 2247 ierr = PetscFree(s->fieldNames);CHKERRQ(ierr); 2248 ierr = PetscFree(s->compNames);CHKERRQ(ierr); 2249 ierr = PetscFree(s->field);CHKERRQ(ierr); 2250 ierr = PetscSectionDestroy(&s->bc);CHKERRQ(ierr); 2251 ierr = PetscFree(s->bcIndices);CHKERRQ(ierr); 2252 ierr = PetscFree2(s->atlasDof, s->atlasOff);CHKERRQ(ierr); 2253 ierr = PetscSectionDestroy(&s->clSection);CHKERRQ(ierr); 2254 ierr = ISDestroy(&s->clPoints);CHKERRQ(ierr); 2255 ierr = ISDestroy(&s->perm);CHKERRQ(ierr); 2256 ierr = PetscSectionResetClosurePermutation(s);CHKERRQ(ierr); 2257 ierr = PetscSectionSymDestroy(&s->sym);CHKERRQ(ierr); 2258 ierr = PetscSectionDestroy(&s->clSection);CHKERRQ(ierr); 2259 ierr = ISDestroy(&s->clPoints);CHKERRQ(ierr); 2260 2261 s->pStart = -1; 2262 s->pEnd = -1; 2263 s->maxDof = 0; 2264 s->setup = PETSC_FALSE; 2265 s->numFields = 0; 2266 s->clObj = NULL; 2267 PetscFunctionReturn(0); 2268 } 2269 2270 /*@ 2271 PetscSectionDestroy - Frees a section object and frees its range if that exists. 2272 2273 Not collective 2274 2275 Input Parameters: 2276 . s - the PetscSection 2277 2278 Level: beginner 2279 2280 .seealso: PetscSection, PetscSectionCreate() 2281 @*/ 2282 PetscErrorCode PetscSectionDestroy(PetscSection *s) 2283 { 2284 PetscErrorCode ierr; 2285 2286 PetscFunctionBegin; 2287 if (!*s) PetscFunctionReturn(0); 2288 PetscValidHeaderSpecific(*s,PETSC_SECTION_CLASSID,1); 2289 if (--((PetscObject)(*s))->refct > 0) { 2290 *s = NULL; 2291 PetscFunctionReturn(0); 2292 } 2293 ierr = PetscSectionReset(*s);CHKERRQ(ierr); 2294 ierr = PetscHeaderDestroy(s);CHKERRQ(ierr); 2295 PetscFunctionReturn(0); 2296 } 2297 2298 PetscErrorCode VecIntGetValuesSection(PetscInt *baseArray, PetscSection s, PetscInt point, const PetscInt **values) 2299 { 2300 const PetscInt p = point - s->pStart; 2301 2302 PetscFunctionBegin; 2303 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 2); 2304 *values = &baseArray[s->atlasOff[p]]; 2305 PetscFunctionReturn(0); 2306 } 2307 2308 PetscErrorCode VecIntSetValuesSection(PetscInt *baseArray, PetscSection s, PetscInt point, const PetscInt values[], InsertMode mode) 2309 { 2310 PetscInt *array; 2311 const PetscInt p = point - s->pStart; 2312 const PetscInt orientation = 0; /* Needs to be included for use in closure operations */ 2313 PetscInt cDim = 0; 2314 PetscErrorCode ierr; 2315 2316 PetscFunctionBegin; 2317 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 2); 2318 ierr = PetscSectionGetConstraintDof(s, p, &cDim);CHKERRQ(ierr); 2319 array = &baseArray[s->atlasOff[p]]; 2320 if (!cDim) { 2321 if (orientation >= 0) { 2322 const PetscInt dim = s->atlasDof[p]; 2323 PetscInt i; 2324 2325 if (mode == INSERT_VALUES) { 2326 for (i = 0; i < dim; ++i) array[i] = values[i]; 2327 } else { 2328 for (i = 0; i < dim; ++i) array[i] += values[i]; 2329 } 2330 } else { 2331 PetscInt offset = 0; 2332 PetscInt j = -1, field, i; 2333 2334 for (field = 0; field < s->numFields; ++field) { 2335 const PetscInt dim = s->field[field]->atlasDof[p]; 2336 2337 for (i = dim-1; i >= 0; --i) array[++j] = values[i+offset]; 2338 offset += dim; 2339 } 2340 } 2341 } else { 2342 if (orientation >= 0) { 2343 const PetscInt dim = s->atlasDof[p]; 2344 PetscInt cInd = 0, i; 2345 const PetscInt *cDof; 2346 2347 ierr = PetscSectionGetConstraintIndices(s, point, &cDof);CHKERRQ(ierr); 2348 if (mode == INSERT_VALUES) { 2349 for (i = 0; i < dim; ++i) { 2350 if ((cInd < cDim) && (i == cDof[cInd])) {++cInd; continue;} 2351 array[i] = values[i]; 2352 } 2353 } else { 2354 for (i = 0; i < dim; ++i) { 2355 if ((cInd < cDim) && (i == cDof[cInd])) {++cInd; continue;} 2356 array[i] += values[i]; 2357 } 2358 } 2359 } else { 2360 const PetscInt *cDof; 2361 PetscInt offset = 0; 2362 PetscInt cOffset = 0; 2363 PetscInt j = 0, field; 2364 2365 ierr = PetscSectionGetConstraintIndices(s, point, &cDof);CHKERRQ(ierr); 2366 for (field = 0; field < s->numFields; ++field) { 2367 const PetscInt dim = s->field[field]->atlasDof[p]; /* PetscSectionGetFieldDof() */ 2368 const PetscInt tDim = s->field[field]->bc->atlasDof[p]; /* PetscSectionGetFieldConstraintDof() */ 2369 const PetscInt sDim = dim - tDim; 2370 PetscInt cInd = 0, i,k; 2371 2372 for (i = 0, k = dim+offset-1; i < dim; ++i, ++j, --k) { 2373 if ((cInd < sDim) && (j == cDof[cInd+cOffset])) {++cInd; continue;} 2374 array[j] = values[k]; 2375 } 2376 offset += dim; 2377 cOffset += dim - tDim; 2378 } 2379 } 2380 } 2381 PetscFunctionReturn(0); 2382 } 2383 2384 /*@C 2385 PetscSectionHasConstraints - Determine whether a section has constrained dofs 2386 2387 Not collective 2388 2389 Input Parameter: 2390 . s - The PetscSection 2391 2392 Output Parameter: 2393 . hasConstraints - flag indicating that the section has constrained dofs 2394 2395 Level: intermediate 2396 2397 .seealso: PetscSectionSetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection 2398 @*/ 2399 PetscErrorCode PetscSectionHasConstraints(PetscSection s, PetscBool *hasConstraints) 2400 { 2401 PetscFunctionBegin; 2402 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2403 PetscValidPointer(hasConstraints, 2); 2404 *hasConstraints = s->bc ? PETSC_TRUE : PETSC_FALSE; 2405 PetscFunctionReturn(0); 2406 } 2407 2408 /*@C 2409 PetscSectionGetConstraintIndices - Get the point dof numbers, in [0, dof), which are constrained 2410 2411 Not collective 2412 2413 Input Parameters: 2414 + s - The PetscSection 2415 - point - The point 2416 2417 Output Parameter: 2418 . indices - The constrained dofs 2419 2420 Note: In Fortran, you call PetscSectionGetConstraintIndicesF90() and PetscSectionRestoreConstraintIndicesF90() 2421 2422 Level: intermediate 2423 2424 .seealso: PetscSectionSetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection 2425 @*/ 2426 PetscErrorCode PetscSectionGetConstraintIndices(PetscSection s, PetscInt point, const PetscInt **indices) 2427 { 2428 PetscErrorCode ierr; 2429 2430 PetscFunctionBegin; 2431 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2432 if (s->bc) { 2433 ierr = VecIntGetValuesSection(s->bcIndices, s->bc, point, indices);CHKERRQ(ierr); 2434 } else *indices = NULL; 2435 PetscFunctionReturn(0); 2436 } 2437 2438 /*@C 2439 PetscSectionSetConstraintIndices - Set the point dof numbers, in [0, dof), which are constrained 2440 2441 Not collective 2442 2443 Input Parameters: 2444 + s - The PetscSection 2445 . point - The point 2446 - indices - The constrained dofs 2447 2448 Note: The Fortran is PetscSectionSetConstraintIndicesF90() 2449 2450 Level: intermediate 2451 2452 .seealso: PetscSectionGetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection 2453 @*/ 2454 PetscErrorCode PetscSectionSetConstraintIndices(PetscSection s, PetscInt point, const PetscInt indices[]) 2455 { 2456 PetscErrorCode ierr; 2457 2458 PetscFunctionBegin; 2459 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2460 if (s->bc) { 2461 const PetscInt dof = s->atlasDof[point]; 2462 const PetscInt cdof = s->bc->atlasDof[point]; 2463 PetscInt d; 2464 2465 for (d = 0; d < cdof; ++d) { 2466 if (indices[d] >= dof) SETERRQ4(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Point %D dof %D, invalid constraint index[%D]: %D", point, dof, d, indices[d]); 2467 } 2468 ierr = VecIntSetValuesSection(s->bcIndices, s->bc, point, indices, INSERT_VALUES);CHKERRQ(ierr); 2469 } 2470 PetscFunctionReturn(0); 2471 } 2472 2473 /*@C 2474 PetscSectionGetFieldConstraintIndices - Get the field dof numbers, in [0, fdof), which are constrained 2475 2476 Not collective 2477 2478 Input Parameters: 2479 + s - The PetscSection 2480 . field - The field number 2481 - point - The point 2482 2483 Output Parameter: 2484 . indices - The constrained dofs sorted in ascending order 2485 2486 Notes: 2487 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(). 2488 2489 Fortran Note: 2490 In Fortran, you call PetscSectionGetFieldConstraintIndicesF90() and PetscSectionRestoreFieldConstraintIndicesF90() 2491 2492 Level: intermediate 2493 2494 .seealso: PetscSectionSetFieldConstraintIndices(), PetscSectionGetConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection 2495 @*/ 2496 PetscErrorCode PetscSectionGetFieldConstraintIndices(PetscSection s, PetscInt point, PetscInt field, const PetscInt **indices) 2497 { 2498 PetscErrorCode ierr; 2499 2500 PetscFunctionBegin; 2501 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2502 if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields); 2503 ierr = PetscSectionGetConstraintIndices(s->field[field], point, indices);CHKERRQ(ierr); 2504 PetscFunctionReturn(0); 2505 } 2506 2507 /*@C 2508 PetscSectionSetFieldConstraintIndices - Set the field dof numbers, in [0, fdof), which are constrained 2509 2510 Not collective 2511 2512 Input Parameters: 2513 + s - The PetscSection 2514 . point - The point 2515 . field - The field number 2516 - indices - The constrained dofs 2517 2518 Note: The Fortran is PetscSectionSetFieldConstraintIndicesF90() 2519 2520 Level: intermediate 2521 2522 .seealso: PetscSectionSetConstraintIndices(), PetscSectionGetFieldConstraintIndices(), PetscSectionGetConstraintDof(), PetscSection 2523 @*/ 2524 PetscErrorCode PetscSectionSetFieldConstraintIndices(PetscSection s, PetscInt point, PetscInt field, const PetscInt indices[]) 2525 { 2526 PetscErrorCode ierr; 2527 2528 PetscFunctionBegin; 2529 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 2530 if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields); 2531 ierr = PetscSectionSetConstraintIndices(s->field[field], point, indices);CHKERRQ(ierr); 2532 PetscFunctionReturn(0); 2533 } 2534 2535 /*@ 2536 PetscSectionPermute - Reorder the section according to the input point permutation 2537 2538 Collective on PetscSection 2539 2540 Input Parameter: 2541 + section - The PetscSection object 2542 - perm - The point permutation, old point p becomes new point perm[p] 2543 2544 Output Parameter: 2545 . sectionNew - The permuted PetscSection 2546 2547 Level: intermediate 2548 2549 .seealso: MatPermute() 2550 @*/ 2551 PetscErrorCode PetscSectionPermute(PetscSection section, IS permutation, PetscSection *sectionNew) 2552 { 2553 PetscSection s = section, sNew; 2554 const PetscInt *perm; 2555 PetscInt numFields, f, c, numPoints, pStart, pEnd, p; 2556 PetscErrorCode ierr; 2557 2558 PetscFunctionBegin; 2559 PetscValidHeaderSpecific(section, PETSC_SECTION_CLASSID, 1); 2560 PetscValidHeaderSpecific(permutation, IS_CLASSID, 2); 2561 PetscValidPointer(sectionNew, 3); 2562 ierr = PetscSectionCreate(PetscObjectComm((PetscObject) s), &sNew);CHKERRQ(ierr); 2563 ierr = PetscSectionGetNumFields(s, &numFields);CHKERRQ(ierr); 2564 if (numFields) {ierr = PetscSectionSetNumFields(sNew, numFields);CHKERRQ(ierr);} 2565 for (f = 0; f < numFields; ++f) { 2566 const char *name; 2567 PetscInt numComp; 2568 2569 ierr = PetscSectionGetFieldName(s, f, &name);CHKERRQ(ierr); 2570 ierr = PetscSectionSetFieldName(sNew, f, name);CHKERRQ(ierr); 2571 ierr = PetscSectionGetFieldComponents(s, f, &numComp);CHKERRQ(ierr); 2572 ierr = PetscSectionSetFieldComponents(sNew, f, numComp);CHKERRQ(ierr); 2573 for (c = 0; c < s->numFieldComponents[f]; ++c) { 2574 ierr = PetscSectionGetComponentName(s, f, c, &name);CHKERRQ(ierr); 2575 ierr = PetscSectionSetComponentName(sNew, f, c, name);CHKERRQ(ierr); 2576 } 2577 } 2578 ierr = ISGetLocalSize(permutation, &numPoints);CHKERRQ(ierr); 2579 ierr = ISGetIndices(permutation, &perm);CHKERRQ(ierr); 2580 ierr = PetscSectionGetChart(s, &pStart, &pEnd);CHKERRQ(ierr); 2581 ierr = PetscSectionSetChart(sNew, pStart, pEnd);CHKERRQ(ierr); 2582 if (numPoints < pEnd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Permutation size %D is less than largest Section point %D", numPoints, pEnd); 2583 for (p = pStart; p < pEnd; ++p) { 2584 PetscInt dof, cdof; 2585 2586 ierr = PetscSectionGetDof(s, p, &dof);CHKERRQ(ierr); 2587 ierr = PetscSectionSetDof(sNew, perm[p], dof);CHKERRQ(ierr); 2588 ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 2589 if (cdof) {ierr = PetscSectionSetConstraintDof(sNew, perm[p], cdof);CHKERRQ(ierr);} 2590 for (f = 0; f < numFields; ++f) { 2591 ierr = PetscSectionGetFieldDof(s, p, f, &dof);CHKERRQ(ierr); 2592 ierr = PetscSectionSetFieldDof(sNew, perm[p], f, dof);CHKERRQ(ierr); 2593 ierr = PetscSectionGetFieldConstraintDof(s, p, f, &cdof);CHKERRQ(ierr); 2594 if (cdof) {ierr = PetscSectionSetFieldConstraintDof(sNew, perm[p], f, cdof);CHKERRQ(ierr);} 2595 } 2596 } 2597 ierr = PetscSectionSetUp(sNew);CHKERRQ(ierr); 2598 for (p = pStart; p < pEnd; ++p) { 2599 const PetscInt *cind; 2600 PetscInt cdof; 2601 2602 ierr = PetscSectionGetConstraintDof(s, p, &cdof);CHKERRQ(ierr); 2603 if (cdof) { 2604 ierr = PetscSectionGetConstraintIndices(s, p, &cind);CHKERRQ(ierr); 2605 ierr = PetscSectionSetConstraintIndices(sNew, perm[p], cind);CHKERRQ(ierr); 2606 } 2607 for (f = 0; f < numFields; ++f) { 2608 ierr = PetscSectionGetFieldConstraintDof(s, p, f, &cdof);CHKERRQ(ierr); 2609 if (cdof) { 2610 ierr = PetscSectionGetFieldConstraintIndices(s, p, f, &cind);CHKERRQ(ierr); 2611 ierr = PetscSectionSetFieldConstraintIndices(sNew, perm[p], f, cind);CHKERRQ(ierr); 2612 } 2613 } 2614 } 2615 ierr = ISRestoreIndices(permutation, &perm);CHKERRQ(ierr); 2616 *sectionNew = sNew; 2617 PetscFunctionReturn(0); 2618 } 2619 2620 /*@ 2621 PetscSectionSetClosureIndex - Set a cache of points in the closure of each point in the section 2622 2623 Collective on section 2624 2625 Input Parameters: 2626 + section - The PetscSection 2627 . obj - A PetscObject which serves as the key for this index 2628 . clSection - Section giving the size of the closure of each point 2629 - clPoints - IS giving the points in each closure 2630 2631 Note: We compress out closure points with no dofs in this section 2632 2633 Level: advanced 2634 2635 .seealso: PetscSectionGetClosureIndex(), DMPlexCreateClosureIndex() 2636 @*/ 2637 PetscErrorCode PetscSectionSetClosureIndex(PetscSection section, PetscObject obj, PetscSection clSection, IS clPoints) 2638 { 2639 PetscErrorCode ierr; 2640 2641 PetscFunctionBegin; 2642 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 2643 PetscValidHeaderSpecific(clSection,PETSC_SECTION_CLASSID,3); 2644 PetscValidHeaderSpecific(clPoints,IS_CLASSID,4); 2645 if (section->clObj != obj) {ierr = PetscSectionResetClosurePermutation(section);CHKERRQ(ierr);} 2646 section->clObj = obj; 2647 ierr = PetscObjectReference((PetscObject)clSection);CHKERRQ(ierr); 2648 ierr = PetscObjectReference((PetscObject)clPoints);CHKERRQ(ierr); 2649 ierr = PetscSectionDestroy(§ion->clSection);CHKERRQ(ierr); 2650 ierr = ISDestroy(§ion->clPoints);CHKERRQ(ierr); 2651 section->clSection = clSection; 2652 section->clPoints = clPoints; 2653 PetscFunctionReturn(0); 2654 } 2655 2656 /*@ 2657 PetscSectionGetClosureIndex - Get the cache of points in the closure of each point in the section 2658 2659 Collective on section 2660 2661 Input Parameters: 2662 + section - The PetscSection 2663 - obj - A PetscObject which serves as the key for this index 2664 2665 Output Parameters: 2666 + clSection - Section giving the size of the closure of each point 2667 - clPoints - IS giving the points in each closure 2668 2669 Note: We compress out closure points with no dofs in this section 2670 2671 Level: advanced 2672 2673 .seealso: PetscSectionSetClosureIndex(), DMPlexCreateClosureIndex() 2674 @*/ 2675 PetscErrorCode PetscSectionGetClosureIndex(PetscSection section, PetscObject obj, PetscSection *clSection, IS *clPoints) 2676 { 2677 PetscFunctionBegin; 2678 if (section->clObj == obj) { 2679 if (clSection) *clSection = section->clSection; 2680 if (clPoints) *clPoints = section->clPoints; 2681 } else { 2682 if (clSection) *clSection = NULL; 2683 if (clPoints) *clPoints = NULL; 2684 } 2685 PetscFunctionReturn(0); 2686 } 2687 2688 PetscErrorCode PetscSectionSetClosurePermutation_Internal(PetscSection section, PetscObject obj, PetscInt depth, PetscInt clSize, PetscCopyMode mode, PetscInt *clPerm) 2689 { 2690 PetscInt i; 2691 khiter_t iter; 2692 int new_entry; 2693 PetscSectionClosurePermKey key = {depth, clSize}; 2694 PetscSectionClosurePermVal *val; 2695 PetscErrorCode ierr; 2696 2697 PetscFunctionBegin; 2698 if (section->clObj != obj) { 2699 ierr = PetscSectionDestroy(§ion->clSection);CHKERRQ(ierr); 2700 ierr = ISDestroy(§ion->clPoints);CHKERRQ(ierr); 2701 } 2702 section->clObj = obj; 2703 if (!section->clHash) {ierr = PetscClPermCreate(§ion->clHash);CHKERRQ(ierr);} 2704 iter = kh_put(ClPerm, section->clHash, key, &new_entry); 2705 val = &kh_val(section->clHash, iter); 2706 if (!new_entry) { 2707 ierr = PetscFree(val->perm);CHKERRQ(ierr); 2708 ierr = PetscFree(val->invPerm);CHKERRQ(ierr); 2709 } 2710 if (mode == PETSC_COPY_VALUES) { 2711 ierr = PetscMalloc1(clSize, &val->perm);CHKERRQ(ierr); 2712 ierr = PetscLogObjectMemory((PetscObject) obj, clSize*sizeof(PetscInt));CHKERRQ(ierr); 2713 ierr = PetscArraycpy(val->perm, clPerm, clSize);CHKERRQ(ierr); 2714 } else if (mode == PETSC_OWN_POINTER) { 2715 val->perm = clPerm; 2716 } else SETERRQ(PetscObjectComm(obj), PETSC_ERR_SUP, "Do not support borrowed arrays"); 2717 ierr = PetscMalloc1(clSize, &val->invPerm);CHKERRQ(ierr); 2718 for (i = 0; i < clSize; ++i) val->invPerm[clPerm[i]] = i; 2719 PetscFunctionReturn(0); 2720 } 2721 2722 /*@ 2723 PetscSectionSetClosurePermutation - Set the dof permutation for the closure of each cell in the section, meaning clPerm[newIndex] = oldIndex. 2724 2725 Not Collective 2726 2727 Input Parameters: 2728 + section - The PetscSection 2729 . obj - A PetscObject which serves as the key for this index (usually a DM) 2730 . depth - Depth of points on which to apply the given permutation 2731 - perm - Permutation of the cell dof closure 2732 2733 Note: 2734 The specified permutation will only be applied to points at depth whose closure size matches the length of perm. In a 2735 mixed-topology or variable-degree finite element space, this function can be called multiple times at each depth for 2736 each topology and degree. 2737 2738 This approach assumes that (depth, len(perm)) uniquely identifies the desired permutation; this might not be true for 2739 exotic/enriched spaces on mixed topology meshes. 2740 2741 Level: intermediate 2742 2743 .seealso: PetscSectionGetClosurePermutation(), PetscSectionGetClosureIndex(), DMPlexCreateClosureIndex(), PetscCopyMode 2744 @*/ 2745 PetscErrorCode PetscSectionSetClosurePermutation(PetscSection section, PetscObject obj, PetscInt depth, IS perm) 2746 { 2747 const PetscInt *clPerm = NULL; 2748 PetscInt clSize = 0; 2749 PetscErrorCode ierr; 2750 2751 PetscFunctionBegin; 2752 if (perm) { 2753 ierr = ISGetLocalSize(perm, &clSize);CHKERRQ(ierr); 2754 ierr = ISGetIndices(perm, &clPerm);CHKERRQ(ierr); 2755 } 2756 ierr = PetscSectionSetClosurePermutation_Internal(section, obj, depth, clSize, PETSC_COPY_VALUES, (PetscInt *) clPerm);CHKERRQ(ierr); 2757 if (perm) {ierr = ISRestoreIndices(perm, &clPerm);CHKERRQ(ierr);} 2758 PetscFunctionReturn(0); 2759 } 2760 2761 PetscErrorCode PetscSectionGetClosurePermutation_Internal(PetscSection section, PetscObject obj, PetscInt depth, PetscInt size, const PetscInt *perm[]) 2762 { 2763 PetscErrorCode ierr; 2764 2765 PetscFunctionBegin; 2766 if (section->clObj == obj) { 2767 PetscSectionClosurePermKey k = {depth, size}; 2768 PetscSectionClosurePermVal v; 2769 ierr = PetscClPermGet(section->clHash, k, &v);CHKERRQ(ierr); 2770 if (perm) *perm = v.perm; 2771 } else { 2772 if (perm) *perm = NULL; 2773 } 2774 PetscFunctionReturn(0); 2775 } 2776 2777 /*@ 2778 PetscSectionGetClosurePermutation - Get the dof permutation for the closure of each cell in the section, meaning clPerm[newIndex] = oldIndex. 2779 2780 Not collective 2781 2782 Input Parameters: 2783 + section - The PetscSection 2784 . obj - A PetscObject which serves as the key for this index (usually a DM) 2785 . depth - Depth stratum on which to obtain closure permutation 2786 - clSize - Closure size to be permuted (e.g., may vary with element topology and degree) 2787 2788 Output Parameter: 2789 . perm - The dof closure permutation 2790 2791 Note: 2792 The user must destroy the IS that is returned. 2793 2794 Level: intermediate 2795 2796 .seealso: PetscSectionSetClosurePermutation(), PetscSectionGetClosureInversePermutation(), PetscSectionGetClosureIndex(), PetscSectionSetClosureIndex(), DMPlexCreateClosureIndex() 2797 @*/ 2798 PetscErrorCode PetscSectionGetClosurePermutation(PetscSection section, PetscObject obj, PetscInt depth, PetscInt clSize, IS *perm) 2799 { 2800 const PetscInt *clPerm; 2801 PetscErrorCode ierr; 2802 2803 PetscFunctionBegin; 2804 ierr = PetscSectionGetClosurePermutation_Internal(section, obj, depth, clSize, &clPerm);CHKERRQ(ierr); 2805 ierr = ISCreateGeneral(PETSC_COMM_SELF, clSize, clPerm, PETSC_USE_POINTER, perm);CHKERRQ(ierr); 2806 PetscFunctionReturn(0); 2807 } 2808 2809 PetscErrorCode PetscSectionGetClosureInversePermutation_Internal(PetscSection section, PetscObject obj, PetscInt depth, PetscInt size, const PetscInt *perm[]) 2810 { 2811 PetscErrorCode ierr; 2812 2813 PetscFunctionBegin; 2814 if (section->clObj == obj && section->clHash) { 2815 PetscSectionClosurePermKey k = {depth, size}; 2816 PetscSectionClosurePermVal v; 2817 ierr = PetscClPermGet(section->clHash, k, &v);CHKERRQ(ierr); 2818 if (perm) *perm = v.invPerm; 2819 } else { 2820 if (perm) *perm = NULL; 2821 } 2822 PetscFunctionReturn(0); 2823 } 2824 2825 /*@ 2826 PetscSectionGetClosureInversePermutation - Get the inverse dof permutation for the closure of each cell in the section, meaning clPerm[oldIndex] = newIndex. 2827 2828 Not collective 2829 2830 Input Parameters: 2831 + section - The PetscSection 2832 . obj - A PetscObject which serves as the key for this index (usually a DM) 2833 . depth - Depth stratum on which to obtain closure permutation 2834 - clSize - Closure size to be permuted (e.g., may vary with element topology and degree) 2835 2836 Output Parameters: 2837 . perm - The dof closure permutation 2838 2839 Note: 2840 The user must destroy the IS that is returned. 2841 2842 Level: intermediate 2843 2844 .seealso: PetscSectionSetClosurePermutation(), PetscSectionGetClosureIndex(), PetscSectionSetClosureIndex(), DMPlexCreateClosureIndex() 2845 @*/ 2846 PetscErrorCode PetscSectionGetClosureInversePermutation(PetscSection section, PetscObject obj, PetscInt depth, PetscInt clSize, IS *perm) 2847 { 2848 const PetscInt *clPerm; 2849 PetscErrorCode ierr; 2850 2851 PetscFunctionBegin; 2852 ierr = PetscSectionGetClosureInversePermutation_Internal(section, obj, depth, clSize, &clPerm);CHKERRQ(ierr); 2853 ierr = ISCreateGeneral(PETSC_COMM_SELF, clSize, clPerm, PETSC_USE_POINTER, perm);CHKERRQ(ierr); 2854 PetscFunctionReturn(0); 2855 } 2856 2857 /*@ 2858 PetscSectionGetField - Get the subsection associated with a single field 2859 2860 Input Parameters: 2861 + s - The PetscSection 2862 - field - The field number 2863 2864 Output Parameter: 2865 . subs - The subsection for the given field 2866 2867 Level: intermediate 2868 2869 .seealso: PetscSectionSetNumFields() 2870 @*/ 2871 PetscErrorCode PetscSectionGetField(PetscSection s, PetscInt field, PetscSection *subs) 2872 { 2873 PetscFunctionBegin; 2874 PetscValidHeaderSpecific(s,PETSC_SECTION_CLASSID,1); 2875 PetscValidPointer(subs,3); 2876 if ((field < 0) || (field >= s->numFields)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section field %D should be in [%D, %D)", field, 0, s->numFields); 2877 *subs = s->field[field]; 2878 PetscFunctionReturn(0); 2879 } 2880 2881 PetscClassId PETSC_SECTION_SYM_CLASSID; 2882 PetscFunctionList PetscSectionSymList = NULL; 2883 2884 /*@ 2885 PetscSectionSymCreate - Creates an empty PetscSectionSym object. 2886 2887 Collective 2888 2889 Input Parameter: 2890 . comm - the MPI communicator 2891 2892 Output Parameter: 2893 . sym - pointer to the new set of symmetries 2894 2895 Level: developer 2896 2897 .seealso: PetscSectionSym, PetscSectionSymDestroy() 2898 @*/ 2899 PetscErrorCode PetscSectionSymCreate(MPI_Comm comm, PetscSectionSym *sym) 2900 { 2901 PetscErrorCode ierr; 2902 2903 PetscFunctionBegin; 2904 PetscValidPointer(sym,2); 2905 ierr = ISInitializePackage();CHKERRQ(ierr); 2906 ierr = PetscHeaderCreate(*sym,PETSC_SECTION_SYM_CLASSID,"PetscSectionSym","Section Symmetry","IS",comm,PetscSectionSymDestroy,PetscSectionSymView);CHKERRQ(ierr); 2907 PetscFunctionReturn(0); 2908 } 2909 2910 /*@C 2911 PetscSectionSymSetType - Builds a PetscSection symmetry, for a particular implementation. 2912 2913 Collective on PetscSectionSym 2914 2915 Input Parameters: 2916 + sym - The section symmetry object 2917 - method - The name of the section symmetry type 2918 2919 Level: developer 2920 2921 .seealso: PetscSectionSymGetType(), PetscSectionSymCreate() 2922 @*/ 2923 PetscErrorCode PetscSectionSymSetType(PetscSectionSym sym, PetscSectionSymType method) 2924 { 2925 PetscErrorCode (*r)(PetscSectionSym); 2926 PetscBool match; 2927 PetscErrorCode ierr; 2928 2929 PetscFunctionBegin; 2930 PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID,1); 2931 ierr = PetscObjectTypeCompare((PetscObject) sym, method, &match);CHKERRQ(ierr); 2932 if (match) PetscFunctionReturn(0); 2933 2934 ierr = PetscFunctionListFind(PetscSectionSymList,method,&r);CHKERRQ(ierr); 2935 if (!r) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE, "Unknown PetscSectionSym type: %s", method); 2936 if (sym->ops->destroy) { 2937 ierr = (*sym->ops->destroy)(sym);CHKERRQ(ierr); 2938 sym->ops->destroy = NULL; 2939 } 2940 ierr = (*r)(sym);CHKERRQ(ierr); 2941 ierr = PetscObjectChangeTypeName((PetscObject)sym,method);CHKERRQ(ierr); 2942 PetscFunctionReturn(0); 2943 } 2944 2945 /*@C 2946 PetscSectionSymGetType - Gets the section symmetry type name (as a string) from the PetscSectionSym. 2947 2948 Not Collective 2949 2950 Input Parameter: 2951 . sym - The section symmetry 2952 2953 Output Parameter: 2954 . type - The index set type name 2955 2956 Level: developer 2957 2958 .seealso: PetscSectionSymSetType(), PetscSectionSymCreate() 2959 @*/ 2960 PetscErrorCode PetscSectionSymGetType(PetscSectionSym sym, PetscSectionSymType *type) 2961 { 2962 PetscFunctionBegin; 2963 PetscValidHeaderSpecific(sym, PETSC_SECTION_SYM_CLASSID,1); 2964 PetscValidCharPointer(type,2); 2965 *type = ((PetscObject)sym)->type_name; 2966 PetscFunctionReturn(0); 2967 } 2968 2969 /*@C 2970 PetscSectionSymRegister - Adds a new section symmetry implementation 2971 2972 Not Collective 2973 2974 Input Parameters: 2975 + name - The name of a new user-defined creation routine 2976 - create_func - The creation routine itself 2977 2978 Notes: 2979 PetscSectionSymRegister() may be called multiple times to add several user-defined vectors 2980 2981 Level: developer 2982 2983 .seealso: PetscSectionSymCreate(), PetscSectionSymSetType() 2984 @*/ 2985 PetscErrorCode PetscSectionSymRegister(const char sname[], PetscErrorCode (*function)(PetscSectionSym)) 2986 { 2987 PetscErrorCode ierr; 2988 2989 PetscFunctionBegin; 2990 ierr = ISInitializePackage();CHKERRQ(ierr); 2991 ierr = PetscFunctionListAdd(&PetscSectionSymList,sname,function);CHKERRQ(ierr); 2992 PetscFunctionReturn(0); 2993 } 2994 2995 /*@ 2996 PetscSectionSymDestroy - Destroys a section symmetry. 2997 2998 Collective on PetscSectionSym 2999 3000 Input Parameters: 3001 . sym - the section symmetry 3002 3003 Level: developer 3004 3005 .seealso: PetscSectionSymCreate(), PetscSectionSymDestroy() 3006 @*/ 3007 PetscErrorCode PetscSectionSymDestroy(PetscSectionSym *sym) 3008 { 3009 SymWorkLink link,next; 3010 PetscErrorCode ierr; 3011 3012 PetscFunctionBegin; 3013 if (!*sym) PetscFunctionReturn(0); 3014 PetscValidHeaderSpecific((*sym),PETSC_SECTION_SYM_CLASSID,1); 3015 if (--((PetscObject)(*sym))->refct > 0) {*sym = NULL; PetscFunctionReturn(0);} 3016 if ((*sym)->ops->destroy) { 3017 ierr = (*(*sym)->ops->destroy)(*sym);CHKERRQ(ierr); 3018 } 3019 if ((*sym)->workout) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Work array still checked out"); 3020 for (link=(*sym)->workin; link; link=next) { 3021 next = link->next; 3022 ierr = PetscFree2(*(PetscInt***)&link->perms,*(PetscScalar***)&link->rots);CHKERRQ(ierr); 3023 ierr = PetscFree(link);CHKERRQ(ierr); 3024 } 3025 (*sym)->workin = NULL; 3026 ierr = PetscHeaderDestroy(sym);CHKERRQ(ierr); 3027 PetscFunctionReturn(0); 3028 } 3029 3030 /*@C 3031 PetscSectionSymView - Displays a section symmetry 3032 3033 Collective on PetscSectionSym 3034 3035 Input Parameters: 3036 + sym - the index set 3037 - viewer - viewer used to display the set, for example PETSC_VIEWER_STDOUT_SELF. 3038 3039 Level: developer 3040 3041 .seealso: PetscViewerASCIIOpen() 3042 @*/ 3043 PetscErrorCode PetscSectionSymView(PetscSectionSym sym,PetscViewer viewer) 3044 { 3045 PetscErrorCode ierr; 3046 3047 PetscFunctionBegin; 3048 PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,1); 3049 if (!viewer) { 3050 ierr = PetscViewerASCIIGetStdout(PetscObjectComm((PetscObject)sym),&viewer);CHKERRQ(ierr); 3051 } 3052 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,2); 3053 PetscCheckSameComm(sym,1,viewer,2); 3054 ierr = PetscObjectPrintClassNamePrefixType((PetscObject)sym,viewer);CHKERRQ(ierr); 3055 if (sym->ops->view) { 3056 ierr = (*sym->ops->view)(sym,viewer);CHKERRQ(ierr); 3057 } 3058 PetscFunctionReturn(0); 3059 } 3060 3061 /*@ 3062 PetscSectionSetSym - Set the symmetries for the data referred to by the section 3063 3064 Collective on PetscSection 3065 3066 Input Parameters: 3067 + section - the section describing data layout 3068 - sym - the symmetry describing the affect of orientation on the access of the data 3069 3070 Level: developer 3071 3072 .seealso: PetscSectionGetSym(), PetscSectionSymCreate() 3073 @*/ 3074 PetscErrorCode PetscSectionSetSym(PetscSection section, PetscSectionSym sym) 3075 { 3076 PetscErrorCode ierr; 3077 3078 PetscFunctionBegin; 3079 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3080 ierr = PetscSectionSymDestroy(&(section->sym));CHKERRQ(ierr); 3081 if (sym) { 3082 PetscValidHeaderSpecific(sym,PETSC_SECTION_SYM_CLASSID,2); 3083 PetscCheckSameComm(section,1,sym,2); 3084 ierr = PetscObjectReference((PetscObject) sym);CHKERRQ(ierr); 3085 } 3086 section->sym = sym; 3087 PetscFunctionReturn(0); 3088 } 3089 3090 /*@ 3091 PetscSectionGetSym - Get the symmetries for the data referred to by the section 3092 3093 Not collective 3094 3095 Input Parameters: 3096 . section - the section describing data layout 3097 3098 Output Parameters: 3099 . sym - the symmetry describing the affect of orientation on the access of the data 3100 3101 Level: developer 3102 3103 .seealso: PetscSectionSetSym(), PetscSectionSymCreate() 3104 @*/ 3105 PetscErrorCode PetscSectionGetSym(PetscSection section, PetscSectionSym *sym) 3106 { 3107 PetscFunctionBegin; 3108 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3109 *sym = section->sym; 3110 PetscFunctionReturn(0); 3111 } 3112 3113 /*@ 3114 PetscSectionSetFieldSym - Set the symmetries for the data referred to by a field of the section 3115 3116 Collective on PetscSection 3117 3118 Input Parameters: 3119 + section - the section describing data layout 3120 . field - the field number 3121 - sym - the symmetry describing the affect of orientation on the access of the data 3122 3123 Level: developer 3124 3125 .seealso: PetscSectionGetFieldSym(), PetscSectionSymCreate() 3126 @*/ 3127 PetscErrorCode PetscSectionSetFieldSym(PetscSection section, PetscInt field, PetscSectionSym sym) 3128 { 3129 PetscErrorCode ierr; 3130 3131 PetscFunctionBegin; 3132 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3133 if (field < 0 || field >= section->numFields) SETERRQ2(PetscObjectComm((PetscObject)section),PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number %D (not in [0,%D)", field, section->numFields); 3134 ierr = PetscSectionSetSym(section->field[field],sym);CHKERRQ(ierr); 3135 PetscFunctionReturn(0); 3136 } 3137 3138 /*@ 3139 PetscSectionGetFieldSym - Get the symmetries for the data referred to by a field of the section 3140 3141 Collective on PetscSection 3142 3143 Input Parameters: 3144 + section - the section describing data layout 3145 - field - the field number 3146 3147 Output Parameters: 3148 . sym - the symmetry describing the affect of orientation on the access of the data 3149 3150 Level: developer 3151 3152 .seealso: PetscSectionSetFieldSym(), PetscSectionSymCreate() 3153 @*/ 3154 PetscErrorCode PetscSectionGetFieldSym(PetscSection section, PetscInt field, PetscSectionSym *sym) 3155 { 3156 PetscFunctionBegin; 3157 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3158 if (field < 0 || field >= section->numFields) SETERRQ2(PetscObjectComm((PetscObject)section),PETSC_ERR_ARG_OUTOFRANGE,"Invalid field number %D (not in [0,%D)", field, section->numFields); 3159 *sym = section->field[field]->sym; 3160 PetscFunctionReturn(0); 3161 } 3162 3163 /*@C 3164 PetscSectionGetPointSyms - Get the symmetries for a set of points in a PetscSection under specific orientations. 3165 3166 Not collective 3167 3168 Input Parameters: 3169 + section - the section 3170 . numPoints - the number of points 3171 - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an 3172 arbitrary integer: its interpretation is up to sym. Orientations are used by DM: for their interpretation in that 3173 context, see DMPlexGetConeOrientation()). 3174 3175 Output Parameter: 3176 + perms - The permutations for the given orientations (or NULL if there is no symmetry or the permutation is the identity). 3177 - rots - The field rotations symmetries for the given orientations (or NULL if there is no symmetry or the rotations are all 3178 identity). 3179 3180 Example of usage, gathering dofs into a local array (lArray) from a section array (sArray): 3181 .vb 3182 const PetscInt **perms; 3183 const PetscScalar **rots; 3184 PetscInt lOffset; 3185 3186 PetscSectionGetPointSyms(section,numPoints,points,&perms,&rots); 3187 for (i = 0, lOffset = 0; i < numPoints; i++) { 3188 PetscInt point = points[2*i], dof, sOffset; 3189 const PetscInt *perm = perms ? perms[i] : NULL; 3190 const PetscScalar *rot = rots ? rots[i] : NULL; 3191 3192 PetscSectionGetDof(section,point,&dof); 3193 PetscSectionGetOffset(section,point,&sOffset); 3194 3195 if (perm) {for (j = 0; j < dof; j++) {lArray[lOffset + perm[j]] = sArray[sOffset + j];}} 3196 else {for (j = 0; j < dof; j++) {lArray[lOffset + j ] = sArray[sOffset + j];}} 3197 if (rot) {for (j = 0; j < dof; j++) {lArray[lOffset + j ] *= rot[j]; }} 3198 lOffset += dof; 3199 } 3200 PetscSectionRestorePointSyms(section,numPoints,points,&perms,&rots); 3201 .ve 3202 3203 Example of usage, adding dofs into a section array (sArray) from a local array (lArray): 3204 .vb 3205 const PetscInt **perms; 3206 const PetscScalar **rots; 3207 PetscInt lOffset; 3208 3209 PetscSectionGetPointSyms(section,numPoints,points,&perms,&rots); 3210 for (i = 0, lOffset = 0; i < numPoints; i++) { 3211 PetscInt point = points[2*i], dof, sOffset; 3212 const PetscInt *perm = perms ? perms[i] : NULL; 3213 const PetscScalar *rot = rots ? rots[i] : NULL; 3214 3215 PetscSectionGetDof(section,point,&dof); 3216 PetscSectionGetOffset(section,point,&sOff); 3217 3218 if (perm) {for (j = 0; j < dof; j++) {sArray[sOffset + j] += lArray[lOffset + perm[j]] * (rot ? PetscConj(rot[perm[j]]) : 1.);}} 3219 else {for (j = 0; j < dof; j++) {sArray[sOffset + j] += lArray[lOffset + j ] * (rot ? PetscConj(rot[ j ]) : 1.);}} 3220 offset += dof; 3221 } 3222 PetscSectionRestorePointSyms(section,numPoints,points,&perms,&rots); 3223 .ve 3224 3225 Level: developer 3226 3227 .seealso: PetscSectionRestorePointSyms(), PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym() 3228 @*/ 3229 PetscErrorCode PetscSectionGetPointSyms(PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots) 3230 { 3231 PetscSectionSym sym; 3232 PetscErrorCode ierr; 3233 3234 PetscFunctionBegin; 3235 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3236 if (numPoints) PetscValidIntPointer(points,3); 3237 if (perms) *perms = NULL; 3238 if (rots) *rots = NULL; 3239 sym = section->sym; 3240 if (sym && (perms || rots)) { 3241 SymWorkLink link; 3242 3243 if (sym->workin) { 3244 link = sym->workin; 3245 sym->workin = sym->workin->next; 3246 } else { 3247 ierr = PetscNewLog(sym,&link);CHKERRQ(ierr); 3248 } 3249 if (numPoints > link->numPoints) { 3250 ierr = PetscFree2(*(PetscInt***)&link->perms,*(PetscInt***)&link->rots);CHKERRQ(ierr); 3251 ierr = PetscMalloc2(numPoints,(PetscInt***)&link->perms,numPoints,(PetscScalar***)&link->rots);CHKERRQ(ierr); 3252 link->numPoints = numPoints; 3253 } 3254 link->next = sym->workout; 3255 sym->workout = link; 3256 ierr = PetscArrayzero((PetscInt**)link->perms,numPoints);CHKERRQ(ierr); 3257 ierr = PetscArrayzero((PetscInt**)link->rots,numPoints);CHKERRQ(ierr); 3258 ierr = (*sym->ops->getpoints) (sym, section, numPoints, points, link->perms, link->rots);CHKERRQ(ierr); 3259 if (perms) *perms = link->perms; 3260 if (rots) *rots = link->rots; 3261 } 3262 PetscFunctionReturn(0); 3263 } 3264 3265 /*@C 3266 PetscSectionRestorePointSyms - Restore the symmetries returned by PetscSectionGetPointSyms() 3267 3268 Not collective 3269 3270 Input Parameters: 3271 + section - the section 3272 . numPoints - the number of points 3273 - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an 3274 arbitrary integer: its interpretation is up to sym. Orientations are used by DM: for their interpretation in that 3275 context, see DMPlexGetConeOrientation()). 3276 3277 Output Parameter: 3278 + perms - The permutations for the given orientations: set to NULL at conclusion 3279 - rots - The field rotations symmetries for the given orientations: set to NULL at conclusion 3280 3281 Level: developer 3282 3283 .seealso: PetscSectionGetPointSyms(), PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym() 3284 @*/ 3285 PetscErrorCode PetscSectionRestorePointSyms(PetscSection section, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots) 3286 { 3287 PetscSectionSym sym; 3288 3289 PetscFunctionBegin; 3290 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3291 sym = section->sym; 3292 if (sym && (perms || rots)) { 3293 SymWorkLink *p,link; 3294 3295 for (p=&sym->workout; (link=*p); p=&link->next) { 3296 if ((perms && link->perms == *perms) || (rots && link->rots == *rots)) { 3297 *p = link->next; 3298 link->next = sym->workin; 3299 sym->workin = link; 3300 if (perms) *perms = NULL; 3301 if (rots) *rots = NULL; 3302 PetscFunctionReturn(0); 3303 } 3304 } 3305 SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Array was not checked out"); 3306 } 3307 PetscFunctionReturn(0); 3308 } 3309 3310 /*@C 3311 PetscSectionGetFieldPointSyms - Get the symmetries for a set of points in a field of a PetscSection under specific orientations. 3312 3313 Not collective 3314 3315 Input Parameters: 3316 + section - the section 3317 . field - the field of the section 3318 . numPoints - the number of points 3319 - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an 3320 arbitrary integer: its interpretation is up to sym. Orientations are used by DM: for their interpretation in that 3321 context, see DMPlexGetConeOrientation()). 3322 3323 Output Parameter: 3324 + perms - The permutations for the given orientations (or NULL if there is no symmetry or the permutation is the identity). 3325 - rots - The field rotations symmetries for the given orientations (or NULL if there is no symmetry or the rotations are all 3326 identity). 3327 3328 Level: developer 3329 3330 .seealso: PetscSectionGetPointSyms(), PetscSectionRestoreFieldPointSyms() 3331 @*/ 3332 PetscErrorCode PetscSectionGetFieldPointSyms(PetscSection section, PetscInt field, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots) 3333 { 3334 PetscErrorCode ierr; 3335 3336 PetscFunctionBegin; 3337 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3338 if (field > section->numFields) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"field %D greater than number of fields (%D) in section",field,section->numFields); 3339 ierr = PetscSectionGetPointSyms(section->field[field],numPoints,points,perms,rots);CHKERRQ(ierr); 3340 PetscFunctionReturn(0); 3341 } 3342 3343 /*@C 3344 PetscSectionRestoreFieldPointSyms - Restore the symmetries returned by PetscSectionGetFieldPointSyms() 3345 3346 Not collective 3347 3348 Input Parameters: 3349 + section - the section 3350 . field - the field number 3351 . numPoints - the number of points 3352 - points - an array of size 2 * numPoints, containing a list of (point, orientation) pairs. (An orientation is an 3353 arbitrary integer: its interpretation is up to sym. Orientations are used by DM: for their interpretation in that 3354 context, see DMPlexGetConeOrientation()). 3355 3356 Output Parameter: 3357 + perms - The permutations for the given orientations: set to NULL at conclusion 3358 - rots - The field rotations symmetries for the given orientations: set to NULL at conclusion 3359 3360 Level: developer 3361 3362 .seealso: PetscSectionRestorePointSyms(), petscSectionGetFieldPointSyms(), PetscSectionSymCreate(), PetscSectionSetSym(), PetscSectionGetSym() 3363 @*/ 3364 PetscErrorCode PetscSectionRestoreFieldPointSyms(PetscSection section, PetscInt field, PetscInt numPoints, const PetscInt *points, const PetscInt ***perms, const PetscScalar ***rots) 3365 { 3366 PetscErrorCode ierr; 3367 3368 PetscFunctionBegin; 3369 PetscValidHeaderSpecific(section,PETSC_SECTION_CLASSID,1); 3370 if (field > section->numFields) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"field %D greater than number of fields (%D) in section",field,section->numFields); 3371 ierr = PetscSectionRestorePointSyms(section->field[field],numPoints,points,perms,rots);CHKERRQ(ierr); 3372 PetscFunctionReturn(0); 3373 } 3374 3375 /*@ 3376 PetscSectionGetUseFieldOffsets - Get the flag to use field offsets directly in a global section, rather than just the point offset 3377 3378 Not collective 3379 3380 Input Parameter: 3381 . s - the global PetscSection 3382 3383 Output Parameters: 3384 . flg - the flag 3385 3386 Level: developer 3387 3388 .seealso: PetscSectionSetChart(), PetscSectionCreate() 3389 @*/ 3390 PetscErrorCode PetscSectionGetUseFieldOffsets(PetscSection s, PetscBool *flg) 3391 { 3392 PetscFunctionBegin; 3393 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 3394 *flg = s->useFieldOff; 3395 PetscFunctionReturn(0); 3396 } 3397 3398 /*@ 3399 PetscSectionSetUseFieldOffsets - Set the flag to use field offsets directly in a global section, rather than just the point offset 3400 3401 Not collective 3402 3403 Input Parameters: 3404 + s - the global PetscSection 3405 - flg - the flag 3406 3407 Level: developer 3408 3409 .seealso: PetscSectionGetUseFieldOffsets(), PetscSectionSetChart(), PetscSectionCreate() 3410 @*/ 3411 PetscErrorCode PetscSectionSetUseFieldOffsets(PetscSection s, PetscBool flg) 3412 { 3413 PetscFunctionBegin; 3414 PetscValidHeaderSpecific(s, PETSC_SECTION_CLASSID, 1); 3415 s->useFieldOff = flg; 3416 PetscFunctionReturn(0); 3417 } 3418 3419 #define PetscSectionExpandPoints_Loop(TYPE) \ 3420 { \ 3421 PetscInt i, n, o0, o1, size; \ 3422 TYPE *a0 = (TYPE*)origArray, *a1; \ 3423 ierr = PetscSectionGetStorageSize(s, &size);CHKERRQ(ierr); \ 3424 ierr = PetscMalloc1(size, &a1);CHKERRQ(ierr); \ 3425 for (i=0; i<npoints; i++) { \ 3426 ierr = PetscSectionGetOffset(origSection, points_[i], &o0);CHKERRQ(ierr); \ 3427 ierr = PetscSectionGetOffset(s, i, &o1);CHKERRQ(ierr); \ 3428 ierr = PetscSectionGetDof(s, i, &n);CHKERRQ(ierr); \ 3429 ierr = PetscMemcpy(&a1[o1], &a0[o0], n*unitsize);CHKERRQ(ierr); \ 3430 } \ 3431 *newArray = (void*)a1; \ 3432 } 3433 3434 /*@ 3435 PetscSectionExtractDofsFromArray - Extracts elements of an array corresponding to DOFs of specified points. 3436 3437 Not collective 3438 3439 Input Parameters: 3440 + origSection - the PetscSection describing the layout of the array 3441 . dataType - MPI_Datatype describing the data type of the array (currently only MPIU_INT, MPIU_SCALAR, MPIU_REAL) 3442 . origArray - the array; its size must be equal to the storage size of origSection 3443 - points - IS with points to extract; its indices must lie in the chart of origSection 3444 3445 Output Parameters: 3446 + newSection - the new PetscSection desribing the layout of the new array (with points renumbered 0,1,... but preserving numbers of DOFs) 3447 - newArray - the array of the extracted DOFs; its size is the storage size of newSection 3448 3449 Level: developer 3450 3451 .seealso: PetscSectionGetChart(), PetscSectionGetDof(), PetscSectionGetStorageSize(), PetscSectionCreate() 3452 @*/ 3453 PetscErrorCode PetscSectionExtractDofsFromArray(PetscSection origSection, MPI_Datatype dataType, const void *origArray, IS points, PetscSection *newSection, void *newArray[]) 3454 { 3455 PetscSection s; 3456 const PetscInt *points_; 3457 PetscInt i, n, npoints, pStart, pEnd; 3458 PetscMPIInt unitsize; 3459 PetscErrorCode ierr; 3460 3461 PetscFunctionBegin; 3462 PetscValidHeaderSpecific(origSection, PETSC_SECTION_CLASSID, 1); 3463 PetscValidPointer(origArray, 3); 3464 PetscValidHeaderSpecific(points, IS_CLASSID, 4); 3465 if (newSection) PetscValidPointer(newSection, 5); 3466 if (newArray) PetscValidPointer(newArray, 6); 3467 ierr = MPI_Type_size(dataType, &unitsize);CHKERRMPI(ierr); 3468 ierr = ISGetLocalSize(points, &npoints);CHKERRQ(ierr); 3469 ierr = ISGetIndices(points, &points_);CHKERRQ(ierr); 3470 ierr = PetscSectionGetChart(origSection, &pStart, &pEnd);CHKERRQ(ierr); 3471 ierr = PetscSectionCreate(PETSC_COMM_SELF, &s);CHKERRQ(ierr); 3472 ierr = PetscSectionSetChart(s, 0, npoints);CHKERRQ(ierr); 3473 for (i=0; i<npoints; i++) { 3474 if (PetscUnlikely(points_[i] < pStart || points_[i] >= pEnd)) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "point %d (index %d) in input IS out of input section's chart", points_[i], i); 3475 ierr = PetscSectionGetDof(origSection, points_[i], &n);CHKERRQ(ierr); 3476 ierr = PetscSectionSetDof(s, i, n);CHKERRQ(ierr); 3477 } 3478 ierr = PetscSectionSetUp(s);CHKERRQ(ierr); 3479 if (newArray) { 3480 if (dataType == MPIU_INT) {PetscSectionExpandPoints_Loop(PetscInt);} 3481 else if (dataType == MPIU_SCALAR) {PetscSectionExpandPoints_Loop(PetscScalar);} 3482 else if (dataType == MPIU_REAL) {PetscSectionExpandPoints_Loop(PetscReal);} 3483 else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "not implemented for this MPI_Datatype"); 3484 } 3485 if (newSection) { 3486 *newSection = s; 3487 } else { 3488 ierr = PetscSectionDestroy(&s);CHKERRQ(ierr); 3489 } 3490 PetscFunctionReturn(0); 3491 } 3492