Actual source code: vsectionis.c
petsc-dev 2014-02-02
1: /*
2: This file contains routines for basic section object implementation.
3: */
5: #include <petsc-private/isimpl.h> /*I "petscvec.h" I*/
6: #include <petscsf.h>
7: #include <petscviewer.h>
9: PetscClassId PETSC_SECTION_CLASSID;
13: /*@
14: PetscSectionCreate - Allocates PetscSection space and sets the map contents to the default.
16: Collective on MPI_Comm
18: Input Parameters:
19: + comm - the MPI communicator
20: - s - pointer to the section
22: Level: developer
24: Notes: Typical calling sequence
25: PetscSectionCreate(MPI_Comm,PetscSection *);
26: PetscSectionSetChart(PetscSection,low,high);
27: PetscSectionSetDof(PetscSection,point,numdof);
28: PetscSectionSetUp(PetscSection);
29: PetscSectionGetOffset(PetscSection,point,PetscInt *);
30: PetscSectionDestroy(PetscSection);
32: The PetscSection object and methods are intended to be used in the PETSc Vec and Mat implementions; it is
33: recommended they not be used in user codes unless you really gain something in their use.
35: .seealso: PetscSection, PetscSectionDestroy()
36: @*/
37: PetscErrorCode PetscSectionCreate(MPI_Comm comm, PetscSection *s)
38: {
43: ISInitializePackage();
45: PetscHeaderCreate(*s,_p_PetscSection,int,PETSC_SECTION_CLASSID,"PetscSection","Section","IS",comm,PetscSectionDestroy,PetscSectionView);
47: (*s)->atlasLayout.comm = comm;
48: (*s)->atlasLayout.pStart = -1;
49: (*s)->atlasLayout.pEnd = -1;
50: (*s)->atlasLayout.numDof = 1;
51: (*s)->maxDof = 0;
52: (*s)->atlasDof = NULL;
53: (*s)->atlasOff = NULL;
54: (*s)->bc = NULL;
55: (*s)->bcIndices = NULL;
56: (*s)->setup = PETSC_FALSE;
57: (*s)->numFields = 0;
58: (*s)->fieldNames = NULL;
59: (*s)->field = NULL;
60: (*s)->clObj = NULL;
61: (*s)->clSection = NULL;
62: (*s)->clPoints = NULL;
63: return(0);
64: }
68: /*@
69: PetscSectionClone - Creates a shallow (if possible) copy of the PetscSection
71: Collective on MPI_Comm
73: Input Parameter:
74: . section - the PetscSection
76: Output Parameter:
77: . newSection - the copy
79: Level: developer
81: .seealso: PetscSection, PetscSectionCreate(), PetscSectionDestroy()
82: @*/
83: PetscErrorCode PetscSectionClone(PetscSection section, PetscSection *newSection)
84: {
85: PetscInt numFields, f, pStart, pEnd, p;
89: PetscSectionCreate(section->atlasLayout.comm, newSection);
90: PetscSectionGetNumFields(section, &numFields);
91: if (numFields) {PetscSectionSetNumFields(*newSection, numFields);}
92: for (f = 0; f < numFields; ++f) {
93: const char *name = NULL;
94: PetscInt numComp = 0;
96: PetscSectionGetFieldName(section, f, &name);
97: PetscSectionSetFieldName(*newSection, f, name);
98: PetscSectionGetFieldComponents(section, f, &numComp);
99: PetscSectionSetFieldComponents(*newSection, f, numComp);
100: }
101: PetscSectionGetChart(section, &pStart, &pEnd);
102: PetscSectionSetChart(*newSection, pStart, pEnd);
103: for (p = pStart; p < pEnd; ++p) {
104: PetscInt dof, cdof, fcdof = 0;
106: PetscSectionGetDof(section, p, &dof);
107: PetscSectionSetDof(*newSection, p, dof);
108: PetscSectionGetConstraintDof(section, p, &cdof);
109: if (cdof) {PetscSectionSetConstraintDof(*newSection, p, cdof);}
110: for (f = 0; f < numFields; ++f) {
111: PetscSectionGetFieldDof(section, p, f, &dof);
112: PetscSectionSetFieldDof(*newSection, p, f, dof);
113: if (cdof) {
114: PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);
115: if (fcdof) {PetscSectionSetFieldConstraintDof(*newSection, p, f, fcdof);}
116: }
117: }
118: }
119: PetscSectionSetUp(*newSection);
120: for (p = pStart; p < pEnd; ++p) {
121: PetscInt off, cdof, fcdof = 0;
122: const PetscInt *cInd;
124: /* Must set offsets in case they do not agree with the prefix sums */
125: PetscSectionGetOffset(section, p, &off);
126: PetscSectionSetOffset(*newSection, p, off);
127: PetscSectionGetConstraintDof(section, p, &cdof);
128: if (cdof) {
129: PetscSectionGetConstraintIndices(section, p, &cInd);
130: PetscSectionSetConstraintIndices(*newSection, p, cInd);
131: for (f = 0; f < numFields; ++f) {
132: PetscSectionGetFieldConstraintDof(section, p, f, &fcdof);
133: if (fcdof) {
134: PetscSectionGetFieldConstraintIndices(section, p, f, &cInd);
135: PetscSectionSetFieldConstraintIndices(*newSection, p, f, cInd);
136: }
137: }
138: }
139: }
140: return(0);
141: }
145: /*@
146: PetscSectionGetNumFields - Returns the number of fields, or 0 if no fields were defined.
148: Not collective
150: Input Parameter:
151: . s - the PetscSection
153: Output Parameter:
154: . numFields - the number of fields defined, or 0 if none were defined
156: Level: intermediate
158: .seealso: PetscSectionSetNumFields()
159: @*/
160: PetscErrorCode PetscSectionGetNumFields(PetscSection s, PetscInt *numFields)
161: {
164: *numFields = s->numFields;
165: return(0);
166: }
170: /*@
171: PetscSectionSetNumFields - Sets the number of fields.
173: Not collective
175: Input Parameters:
176: + s - the PetscSection
177: - numFields - the number of fields
179: Level: intermediate
181: .seealso: PetscSectionGetNumFields()
182: @*/
183: PetscErrorCode PetscSectionSetNumFields(PetscSection s, PetscInt numFields)
184: {
185: PetscInt f;
189: if (numFields <= 0) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "The number of fields %d must be positive", numFields);
190: PetscSectionReset(s);
192: s->numFields = numFields;
193: PetscMalloc1(s->numFields, &s->numFieldComponents);
194: PetscMalloc1(s->numFields, &s->fieldNames);
195: PetscMalloc1(s->numFields, &s->field);
196: for (f = 0; f < s->numFields; ++f) {
197: char name[64];
199: s->numFieldComponents[f] = 1;
201: PetscSectionCreate(s->atlasLayout.comm, &s->field[f]);
202: PetscSNPrintf(name, 64, "Field_%D", f);
203: PetscStrallocpy(name, (char **) &s->fieldNames[f]);
204: }
205: return(0);
206: }
210: /*@C
211: PetscSectionGetFieldName - Returns the name of a field in the PetscSection
213: Not Collective
215: Input Parameters:
216: + s - the PetscSection
217: - field - the field number
219: Output Parameter:
220: . fieldName - the field name
222: Level: developer
224: .seealso: PetscSectionSetFieldName()
225: @*/
226: PetscErrorCode PetscSectionGetFieldName(PetscSection s, PetscInt field, const char *fieldName[])
227: {
230: 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);
231: *fieldName = s->fieldNames[field];
232: return(0);
233: }
237: /*@C
238: PetscSectionSetFieldName - Sets the name of a field in the PetscSection
240: Not Collective
242: Input Parameters:
243: + s - the PetscSection
244: . field - the field number
245: - fieldName - the field name
247: Level: developer
249: .seealso: PetscSectionGetFieldName()
250: @*/
251: PetscErrorCode PetscSectionSetFieldName(PetscSection s, PetscInt field, const char fieldName[])
252: {
257: 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);
258: PetscFree(s->fieldNames[field]);
259: PetscStrallocpy(fieldName, (char**) &s->fieldNames[field]);
260: return(0);
261: }
265: /*@
266: PetscSectionGetFieldComponents - Returns the number of field components for the given field.
268: Not collective
270: Input Parameters:
271: + s - the PetscSection
272: - field - the field number
274: Output Parameter:
275: . numComp - the number of field components
277: Level: intermediate
279: .seealso: PetscSectionSetNumFieldComponents(), PetscSectionGetNumFields()
280: @*/
281: PetscErrorCode PetscSectionGetFieldComponents(PetscSection s, PetscInt field, PetscInt *numComp)
282: {
285: 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);
286: *numComp = s->numFieldComponents[field];
287: return(0);
288: }
292: /*@
293: PetscSectionSetFieldComponents - Sets the number of field components for the given field.
295: Not collective
297: Input Parameters:
298: + s - the PetscSection
299: . field - the field number
300: - numComp - the number of field components
302: Level: intermediate
304: .seealso: PetscSectionGetNumFieldComponents(), PetscSectionGetNumFields()
305: @*/
306: PetscErrorCode PetscSectionSetFieldComponents(PetscSection s, PetscInt field, PetscInt numComp)
307: {
309: 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);
310: s->numFieldComponents[field] = numComp;
311: return(0);
312: }
316: PetscErrorCode PetscSectionCheckConstraints(PetscSection s)
317: {
321: if (!s->bc) {
322: PetscSectionCreate(PETSC_COMM_SELF, &s->bc);
323: PetscSectionSetChart(s->bc, s->atlasLayout.pStart, s->atlasLayout.pEnd);
324: }
325: return(0);
326: }
330: /*@
331: PetscSectionGetChart - Returns the range [pStart, pEnd) in which points in the lie.
333: Not collective
335: Input Parameter:
336: . s - the PetscSection
338: Output Parameters:
339: + pStart - the first point
340: - pEnd - one past the last point
342: Level: intermediate
344: .seealso: PetscSectionSetChart(), PetscSectionCreate()
345: @*/
346: PetscErrorCode PetscSectionGetChart(PetscSection s, PetscInt *pStart, PetscInt *pEnd)
347: {
349: if (pStart) *pStart = s->atlasLayout.pStart;
350: if (pEnd) *pEnd = s->atlasLayout.pEnd;
351: return(0);
352: }
356: /*@
357: PetscSectionSetChart - Sets the range [pStart, pEnd) in which points in the lie.
359: Not collective
361: Input Parameters:
362: + s - the PetscSection
363: . pStart - the first point
364: - pEnd - one past the last point
366: Level: intermediate
368: .seealso: PetscSectionSetChart(), PetscSectionCreate()
369: @*/
370: PetscErrorCode PetscSectionSetChart(PetscSection s, PetscInt pStart, PetscInt pEnd)
371: {
372: PetscInt f;
376: /* Cannot Reset() because it destroys field information */
377: s->setup = PETSC_FALSE;
378: PetscSectionDestroy(&s->bc);
379: PetscFree(s->bcIndices);
380: PetscFree2(s->atlasDof, s->atlasOff);
382: s->atlasLayout.pStart = pStart;
383: s->atlasLayout.pEnd = pEnd;
384: PetscMalloc2((pEnd - pStart), &s->atlasDof, (pEnd - pStart), &s->atlasOff);
385: PetscMemzero(s->atlasDof, (pEnd - pStart)*sizeof(PetscInt));
386: for (f = 0; f < s->numFields; ++f) {
387: PetscSectionSetChart(s->field[f], pStart, pEnd);
388: }
389: return(0);
390: }
394: /*@
395: PetscSectionGetDof - Return the number of degrees of freedom associated with a given point.
397: Not collective
399: Input Parameters:
400: + s - the PetscSection
401: - point - the point
403: Output Parameter:
404: . numDof - the number of dof
406: Level: intermediate
408: .seealso: PetscSectionSetDof(), PetscSectionCreate()
409: @*/
410: PetscErrorCode PetscSectionGetDof(PetscSection s, PetscInt point, PetscInt *numDof)
411: {
413: #if defined(PETSC_USE_DEBUG)
414: if ((point < s->atlasLayout.pStart) || (point >= s->atlasLayout.pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %d should be in [%d, %d)", point, s->atlasLayout.pStart, s->atlasLayout.pEnd);
415: #endif
416: *numDof = s->atlasDof[point - s->atlasLayout.pStart];
417: return(0);
418: }
422: /*@
423: PetscSectionSetDof - Sets the number of degrees of freedom associated with a given point.
425: Not collective
427: Input Parameters:
428: + s - the PetscSection
429: . point - the point
430: - numDof - the number of dof
432: Level: intermediate
434: .seealso: PetscSectionGetDof(), PetscSectionAddDof(), PetscSectionCreate()
435: @*/
436: PetscErrorCode PetscSectionSetDof(PetscSection s, PetscInt point, PetscInt numDof)
437: {
439: if ((point < s->atlasLayout.pStart) || (point >= s->atlasLayout.pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %d should be in [%d, %d)", point, s->atlasLayout.pStart, s->atlasLayout.pEnd);
440: s->atlasDof[point - s->atlasLayout.pStart] = numDof;
441: return(0);
442: }
446: /*@
447: PetscSectionAddDof - Adds to the number of degrees of freedom associated with a given point.
449: Not collective
451: Input Parameters:
452: + s - the PetscSection
453: . point - the point
454: - numDof - the number of additional dof
456: Level: intermediate
458: .seealso: PetscSectionGetDof(), PetscSectionSetDof(), PetscSectionCreate()
459: @*/
460: PetscErrorCode PetscSectionAddDof(PetscSection s, PetscInt point, PetscInt numDof)
461: {
463: if ((point < s->atlasLayout.pStart) || (point >= s->atlasLayout.pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %d should be in [%d, %d)", point, s->atlasLayout.pStart, s->atlasLayout.pEnd);
464: s->atlasDof[point - s->atlasLayout.pStart] += numDof;
465: return(0);
466: }
470: /*@
471: PetscSectionGetFieldDof - Return the number of degrees of freedom associated with a field on a given point.
473: Not collective
475: Input Parameters:
476: + s - the PetscSection
477: . point - the point
478: - field - the field
480: Output Parameter:
481: . numDof - the number of dof
483: Level: intermediate
485: .seealso: PetscSectionSetFieldDof(), PetscSectionCreate()
486: @*/
487: PetscErrorCode PetscSectionGetFieldDof(PetscSection s, PetscInt point, PetscInt field, PetscInt *numDof)
488: {
492: 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);
493: PetscSectionGetDof(s->field[field], point, numDof);
494: return(0);
495: }
499: /*@
500: PetscSectionSetFieldDof - Sets the number of degrees of freedom associated with a field on a given point.
502: Not collective
504: Input Parameters:
505: + s - the PetscSection
506: . point - the point
507: . field - the field
508: - numDof - the number of dof
510: Level: intermediate
512: .seealso: PetscSectionGetFieldDof(), PetscSectionCreate()
513: @*/
514: PetscErrorCode PetscSectionSetFieldDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof)
515: {
519: 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);
520: PetscSectionSetDof(s->field[field], point, numDof);
521: return(0);
522: }
526: /*@
527: PetscSectionAddFieldDof - Adds a number of degrees of freedom associated with a field on a given point.
529: Not collective
531: Input Parameters:
532: + s - the PetscSection
533: . point - the point
534: . field - the field
535: - numDof - the number of dof
537: Level: intermediate
539: .seealso: PetscSectionSetFieldDof(), PetscSectionGetFieldDof(), PetscSectionCreate()
540: @*/
541: PetscErrorCode PetscSectionAddFieldDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof)
542: {
546: 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);
547: PetscSectionAddDof(s->field[field], point, numDof);
548: return(0);
549: }
553: /*@
554: PetscSectionGetConstraintDof - Return the number of constrained degrees of freedom associated with a given point.
556: Not collective
558: Input Parameters:
559: + s - the PetscSection
560: - point - the point
562: Output Parameter:
563: . numDof - the number of dof which are fixed by constraints
565: Level: intermediate
567: .seealso: PetscSectionGetDof(), PetscSectionSetConstraintDof(), PetscSectionCreate()
568: @*/
569: PetscErrorCode PetscSectionGetConstraintDof(PetscSection s, PetscInt point, PetscInt *numDof)
570: {
574: if (s->bc) {
575: PetscSectionGetDof(s->bc, point, numDof);
576: } else *numDof = 0;
577: return(0);
578: }
582: /*@
583: PetscSectionSetConstraintDof - Set the number of constrained degrees of freedom associated with a given point.
585: Not collective
587: Input Parameters:
588: + s - the PetscSection
589: . point - the point
590: - numDof - the number of dof which are fixed by constraints
592: Level: intermediate
594: .seealso: PetscSectionSetDof(), PetscSectionGetConstraintDof(), PetscSectionCreate()
595: @*/
596: PetscErrorCode PetscSectionSetConstraintDof(PetscSection s, PetscInt point, PetscInt numDof)
597: {
601: if (numDof) {
602: PetscSectionCheckConstraints(s);
603: PetscSectionSetDof(s->bc, point, numDof);
604: }
605: return(0);
606: }
610: /*@
611: PetscSectionAddConstraintDof - Increment the number of constrained degrees of freedom associated with a given point.
613: Not collective
615: Input Parameters:
616: + s - the PetscSection
617: . point - the point
618: - numDof - the number of additional dof which are fixed by constraints
620: Level: intermediate
622: .seealso: PetscSectionAddDof(), PetscSectionGetConstraintDof(), PetscSectionCreate()
623: @*/
624: PetscErrorCode PetscSectionAddConstraintDof(PetscSection s, PetscInt point, PetscInt numDof)
625: {
629: if (numDof) {
630: PetscSectionCheckConstraints(s);
631: PetscSectionAddDof(s->bc, point, numDof);
632: }
633: return(0);
634: }
638: /*@
639: PetscSectionGetFieldConstraintDof - Return the number of constrained degrees of freedom associated with a given field on a point.
641: Not collective
643: Input Parameters:
644: + s - the PetscSection
645: . point - the point
646: - field - the field
648: Output Parameter:
649: . numDof - the number of dof which are fixed by constraints
651: Level: intermediate
653: .seealso: PetscSectionGetDof(), PetscSectionSetFieldConstraintDof(), PetscSectionCreate()
654: @*/
655: PetscErrorCode PetscSectionGetFieldConstraintDof(PetscSection s, PetscInt point, PetscInt field, PetscInt *numDof)
656: {
660: 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);
661: PetscSectionGetConstraintDof(s->field[field], point, numDof);
662: return(0);
663: }
667: /*@
668: PetscSectionSetFieldConstraintDof - Set the number of constrained degrees of freedom associated with a given field on a point.
670: Not collective
672: Input Parameters:
673: + s - the PetscSection
674: . point - the point
675: . field - the field
676: - numDof - the number of dof which are fixed by constraints
678: Level: intermediate
680: .seealso: PetscSectionSetDof(), PetscSectionGetFieldConstraintDof(), PetscSectionCreate()
681: @*/
682: PetscErrorCode PetscSectionSetFieldConstraintDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof)
683: {
687: 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);
688: PetscSectionSetConstraintDof(s->field[field], point, numDof);
689: return(0);
690: }
694: /*@
695: PetscSectionAddFieldConstraintDof - Increment the number of constrained degrees of freedom associated with a given field on a point.
697: Not collective
699: Input Parameters:
700: + s - the PetscSection
701: . point - the point
702: . field - the field
703: - numDof - the number of additional dof which are fixed by constraints
705: Level: intermediate
707: .seealso: PetscSectionAddDof(), PetscSectionGetFieldConstraintDof(), PetscSectionCreate()
708: @*/
709: PetscErrorCode PetscSectionAddFieldConstraintDof(PetscSection s, PetscInt point, PetscInt field, PetscInt numDof)
710: {
714: 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);
715: PetscSectionAddConstraintDof(s->field[field], point, numDof);
716: return(0);
717: }
721: PetscErrorCode PetscSectionSetUpBC(PetscSection s)
722: {
726: if (s->bc) {
727: const PetscInt last = (s->bc->atlasLayout.pEnd-s->bc->atlasLayout.pStart) - 1;
729: PetscSectionSetUp(s->bc);
730: PetscMalloc1((s->bc->atlasOff[last] + s->bc->atlasDof[last]), &s->bcIndices);
731: }
732: return(0);
733: }
737: /*@
738: PetscSectionSetUp - Calculate offsets based upon the number of degrees of freedom for each point.
740: Not collective
742: Input Parameter:
743: . s - the PetscSection
745: Level: intermediate
747: .seealso: PetscSectionCreate()
748: @*/
749: PetscErrorCode PetscSectionSetUp(PetscSection s)
750: {
751: PetscInt offset = 0, p, f;
755: if (s->setup) return(0);
756: s->setup = PETSC_TRUE;
757: for (p = 0; p < s->atlasLayout.pEnd - s->atlasLayout.pStart; ++p) {
758: s->atlasOff[p] = offset;
759: offset += s->atlasDof[p];
760: s->maxDof = PetscMax(s->maxDof, s->atlasDof[p]);
761: }
762: PetscSectionSetUpBC(s);
763: /* Assume that all fields have the same chart */
764: for (p = 0; p < s->atlasLayout.pEnd - s->atlasLayout.pStart; ++p) {
765: PetscInt off = s->atlasOff[p];
767: for (f = 0; f < s->numFields; ++f) {
768: PetscSection sf = s->field[f];
770: sf->atlasOff[p] = off;
771: off += sf->atlasDof[p];
772: }
773: }
774: for (f = 0; f < s->numFields; ++f) {
775: PetscSectionSetUpBC(s->field[f]);
776: }
777: return(0);
778: }
782: /*@
783: PetscSectionGetMaxDof - Return the maximum number of degrees of freedom on any point in the chart
785: Not collective
787: Input Parameters:
788: . s - the PetscSection
790: Output Parameter:
791: . maxDof - the maximum dof
793: Level: intermediate
795: .seealso: PetscSectionGetDof(), PetscSectionSetDof(), PetscSectionCreate()
796: @*/
797: PetscErrorCode PetscSectionGetMaxDof(PetscSection s, PetscInt *maxDof)
798: {
800: *maxDof = s->maxDof;
801: return(0);
802: }
806: /*@
807: PetscSectionGetStorageSize - Return the size of an arary or local Vec capable of holding all the degrees of freedom.
809: Not collective
811: Input Parameters:
812: + s - the PetscSection
813: - point - the point
815: Output Parameter:
816: . size - the size of an array which can hold all the dofs
818: Level: intermediate
820: .seealso: PetscSectionGetOffset(), PetscSectionCreate()
821: @*/
822: PetscErrorCode PetscSectionGetStorageSize(PetscSection s, PetscInt *size)
823: {
824: PetscInt p, n = 0;
827: for (p = 0; p < s->atlasLayout.pEnd - s->atlasLayout.pStart; ++p) n += s->atlasDof[p] > 0 ? s->atlasDof[p] : 0;
828: *size = n;
829: return(0);
830: }
834: PetscErrorCode PetscSectionGetConstrainedStorageSize(PetscSection s, PetscInt *size)
835: {
836: PetscInt p, n = 0;
839: for (p = 0; p < s->atlasLayout.pEnd - s->atlasLayout.pStart; ++p) {
840: const PetscInt cdof = s->bc ? s->bc->atlasDof[p] : 0;
841: n += s->atlasDof[p] > 0 ? s->atlasDof[p] - cdof : 0;
842: }
843: *size = n;
844: return(0);
845: }
849: /*@
850: PetscSectionCreateGlobalSection - Create a section describing the global field layout using
851: the local section and an SF describing the section point overlap.
853: Input Parameters:
854: + s - The PetscSection for the local field layout
855: . sf - The SF describing parallel layout of the section points (leaves are unowned local points)
856: - includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
858: Output Parameter:
859: . gsection - The PetscSection for the global field layout
861: Note: This gives negative sizes and offsets to points not owned by this process
863: Level: developer
865: .seealso: PetscSectionCreate()
866: @*/
867: PetscErrorCode PetscSectionCreateGlobalSection(PetscSection s, PetscSF sf, PetscBool includeConstraints, PetscSection *gsection)
868: {
869: PetscInt *neg = NULL, *recv = NULL;
870: PetscInt pStart, pEnd, p, dof, cdof, off, globalOff = 0, nroots, nlocal;
874: PetscSectionCreate(s->atlasLayout.comm, gsection);
875: PetscSectionGetChart(s, &pStart, &pEnd);
876: PetscSectionSetChart(*gsection, pStart, pEnd);
877: PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);
878: nlocal = nroots; /* The local/leaf space matches global/root space */
879: /* Must allocate for all points visible to SF, which may be more than this section */
880: if (nroots >= 0) { /* nroots < 0 means that the graph has not been set, only happens in serial */
881: PetscMalloc2(nroots,&neg,nlocal,&recv);
882: PetscMemzero(neg,nroots*sizeof(PetscInt));
883: }
884: /* Mark all local points with negative dof */
885: for (p = pStart; p < pEnd; ++p) {
886: PetscSectionGetDof(s, p, &dof);
887: PetscSectionSetDof(*gsection, p, dof);
888: PetscSectionGetConstraintDof(s, p, &cdof);
889: if (!includeConstraints && cdof > 0) {PetscSectionSetConstraintDof(*gsection, p, cdof);}
890: if (neg) neg[p] = -(dof+1);
891: }
892: PetscSectionSetUpBC(*gsection);
893: if (nroots >= 0) {
894: PetscMemzero(recv,nlocal*sizeof(PetscInt));
895: PetscSFBcastBegin(sf, MPIU_INT, neg, recv);
896: PetscSFBcastEnd(sf, MPIU_INT, neg, recv);
897: for (p = pStart; p < pEnd; ++p) {
898: if (recv[p] < 0) (*gsection)->atlasDof[p-pStart] = recv[p];
899: }
900: }
901: /* Calculate new sizes, get proccess offset, and calculate point offsets */
902: for (p = 0, off = 0; p < pEnd-pStart; ++p) {
903: cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
904: (*gsection)->atlasOff[p] = off;
905: off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0;
906: }
907: MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, s->atlasLayout.comm);
908: globalOff -= off;
909: for (p = pStart, off = 0; p < pEnd; ++p) {
910: (*gsection)->atlasOff[p-pStart] += globalOff;
911: if (neg) neg[p] = -((*gsection)->atlasOff[p-pStart]+1);
912: }
913: /* Put in negative offsets for ghost points */
914: if (nroots >= 0) {
915: PetscMemzero(recv,nlocal*sizeof(PetscInt));
916: PetscSFBcastBegin(sf, MPIU_INT, neg, recv);
917: PetscSFBcastEnd(sf, MPIU_INT, neg, recv);
918: for (p = pStart; p < pEnd; ++p) {
919: if (recv[p] < 0) (*gsection)->atlasOff[p-pStart] = recv[p];
920: }
921: }
922: PetscFree2(neg,recv);
923: return(0);
924: }
928: /*@
929: PetscSectionCreateGlobalSectionCensored - Create a section describing the global field layout using
930: the local section and an SF describing the section point overlap.
932: Input Parameters:
933: + s - The PetscSection for the local field layout
934: . sf - The SF describing parallel layout of the section points
935: . includeConstraints - By default this is PETSC_FALSE, meaning that the global field vector will not possess constrained dofs
936: . numExcludes - The number of exclusion ranges
937: - excludes - An array [start_0, end_0, start_1, end_1, ...] where there are numExcludes pairs
939: Output Parameter:
940: . gsection - The PetscSection for the global field layout
942: Note: This gives negative sizes and offsets to points not owned by this process
944: Level: developer
946: .seealso: PetscSectionCreate()
947: @*/
948: PetscErrorCode PetscSectionCreateGlobalSectionCensored(PetscSection s, PetscSF sf, PetscBool includeConstraints, PetscInt numExcludes, const PetscInt excludes[], PetscSection *gsection)
949: {
950: PetscInt *neg = NULL, *tmpOff = NULL;
951: PetscInt pStart, pEnd, p, e, dof, cdof, off, globalOff = 0, nroots;
955: PetscSectionCreate(s->atlasLayout.comm, gsection);
956: PetscSectionGetChart(s, &pStart, &pEnd);
957: PetscSectionSetChart(*gsection, pStart, pEnd);
958: PetscSFGetGraph(sf, &nroots, NULL, NULL, NULL);
959: if (nroots >= 0) {
960: if (nroots < pEnd-pStart) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "PetscSF nroots %d < %d section size", nroots, pEnd-pStart);
961: PetscCalloc1(nroots, &neg);
962: if (nroots > pEnd-pStart) {
963: PetscCalloc1(nroots, &tmpOff);
964: } else {
965: tmpOff = &(*gsection)->atlasDof[-pStart];
966: }
967: }
968: /* Mark ghost points with negative dof */
969: for (p = pStart; p < pEnd; ++p) {
970: for (e = 0; e < numExcludes; ++e) {
971: if ((p >= excludes[e*2+0]) && (p < excludes[e*2+1])) {
972: PetscSectionSetDof(*gsection, p, 0);
973: break;
974: }
975: }
976: if (e < numExcludes) continue;
977: PetscSectionGetDof(s, p, &dof);
978: PetscSectionSetDof(*gsection, p, dof);
979: PetscSectionGetConstraintDof(s, p, &cdof);
980: if (!includeConstraints && cdof > 0) {PetscSectionSetConstraintDof(*gsection, p, cdof);}
981: if (neg) neg[p] = -(dof+1);
982: }
983: PetscSectionSetUpBC(*gsection);
984: if (nroots >= 0) {
985: PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);
986: PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);
987: if (nroots > pEnd - pStart) {
988: for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasDof[p-pStart] = tmpOff[p];}
989: }
990: }
991: /* Calculate new sizes, get proccess offset, and calculate point offsets */
992: for (p = 0, off = 0; p < pEnd-pStart; ++p) {
993: cdof = (!includeConstraints && s->bc) ? s->bc->atlasDof[p] : 0;
994: (*gsection)->atlasOff[p] = off;
995: off += (*gsection)->atlasDof[p] > 0 ? (*gsection)->atlasDof[p]-cdof : 0;
996: }
997: MPI_Scan(&off, &globalOff, 1, MPIU_INT, MPI_SUM, s->atlasLayout.comm);
998: globalOff -= off;
999: for (p = 0, off = 0; p < pEnd-pStart; ++p) {
1000: (*gsection)->atlasOff[p] += globalOff;
1001: if (neg) neg[p+pStart] = -((*gsection)->atlasOff[p]+1);
1002: }
1003: /* Put in negative offsets for ghost points */
1004: if (nroots >= 0) {
1005: PetscSFBcastBegin(sf, MPIU_INT, neg, tmpOff);
1006: PetscSFBcastEnd(sf, MPIU_INT, neg, tmpOff);
1007: if (nroots > pEnd - pStart) {
1008: for (p = pStart; p < pEnd; ++p) {if (tmpOff[p] < 0) (*gsection)->atlasOff[p-pStart] = tmpOff[p];}
1009: }
1010: }
1011: if (nroots >= 0 && nroots > pEnd-pStart) {PetscFree(tmpOff);}
1012: PetscFree(neg);
1013: return(0);
1014: }
1018: PetscErrorCode PetscSectionGetPointLayout(MPI_Comm comm, PetscSection s, PetscLayout *layout)
1019: {
1020: PetscInt pStart, pEnd, p, localSize = 0;
1024: PetscSectionGetChart(s, &pStart, &pEnd);
1025: for (p = pStart; p < pEnd; ++p) {
1026: PetscInt dof;
1028: PetscSectionGetDof(s, p, &dof);
1029: if (dof > 0) ++localSize;
1030: }
1031: PetscLayoutCreate(comm, layout);
1032: PetscLayoutSetLocalSize(*layout, localSize);
1033: PetscLayoutSetBlockSize(*layout, 1);
1034: PetscLayoutSetUp(*layout);
1035: return(0);
1036: }
1040: PetscErrorCode PetscSectionGetValueLayout(MPI_Comm comm, PetscSection s, PetscLayout *layout)
1041: {
1042: PetscInt pStart, pEnd, p, localSize = 0;
1046: PetscSectionGetChart(s, &pStart, &pEnd);
1047: for (p = pStart; p < pEnd; ++p) {
1048: PetscInt dof,cdof;
1050: PetscSectionGetDof(s, p, &dof);
1051: PetscSectionGetConstraintDof(s, p, &cdof);
1052: if (dof-cdof > 0) localSize += dof-cdof;
1053: }
1054: PetscLayoutCreate(comm, layout);
1055: PetscLayoutSetLocalSize(*layout, localSize);
1056: PetscLayoutSetBlockSize(*layout, 1);
1057: PetscLayoutSetUp(*layout);
1058: return(0);
1059: }
1063: /*@
1064: PetscSectionGetOffset - Return the offset into an array or local Vec for the dof associated with the given point.
1066: Not collective
1068: Input Parameters:
1069: + s - the PetscSection
1070: - point - the point
1072: Output Parameter:
1073: . offset - the offset
1075: Level: intermediate
1077: .seealso: PetscSectionGetFieldOffset(), PetscSectionCreate()
1078: @*/
1079: PetscErrorCode PetscSectionGetOffset(PetscSection s, PetscInt point, PetscInt *offset)
1080: {
1082: #if defined(PETSC_USE_DEBUG)
1083: if ((point < s->atlasLayout.pStart) || (point >= s->atlasLayout.pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %d should be in [%d, %d)", point, s->atlasLayout.pStart, s->atlasLayout.pEnd);
1084: #endif
1085: *offset = s->atlasOff[point - s->atlasLayout.pStart];
1086: return(0);
1087: }
1091: /*@
1092: PetscSectionSetOffset - Set the offset into an array or local Vec for the dof associated with the given point.
1094: Not collective
1096: Input Parameters:
1097: + s - the PetscSection
1098: . point - the point
1099: - offset - the offset
1101: Note: The user usually does not call this function, but uses PetscSectionSetUp()
1103: Level: intermediate
1105: .seealso: PetscSectionGetFieldOffset(), PetscSectionCreate(), PetscSectionSetUp()
1106: @*/
1107: PetscErrorCode PetscSectionSetOffset(PetscSection s, PetscInt point, PetscInt offset)
1108: {
1110: if ((point < s->atlasLayout.pStart) || (point >= s->atlasLayout.pEnd)) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Section point %d should be in [%d, %d)", point, s->atlasLayout.pStart, s->atlasLayout.pEnd);
1111: s->atlasOff[point - s->atlasLayout.pStart] = offset;
1112: return(0);
1113: }
1117: /*@
1118: PetscSectionGetFieldOffset - Return the offset into an array or local Vec for the dof associated with the given point.
1120: Not collective
1122: Input Parameters:
1123: + s - the PetscSection
1124: . point - the point
1125: - field - the field
1127: Output Parameter:
1128: . offset - the offset
1130: Level: intermediate
1132: .seealso: PetscSectionGetOffset(), PetscSectionCreate()
1133: @*/
1134: PetscErrorCode PetscSectionGetFieldOffset(PetscSection s, PetscInt point, PetscInt field, PetscInt *offset)
1135: {
1139: 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);
1140: PetscSectionGetOffset(s->field[field], point, offset);
1141: return(0);
1142: }
1146: /*@
1147: PetscSectionSetFieldOffset - Set the offset into an array or local Vec for the dof associated with the given point.
1149: Not collective
1151: Input Parameters:
1152: + s - the PetscSection
1153: . point - the point
1154: . field - the field
1155: - offset - the offset
1157: Note: The user usually does not call this function, but uses PetscSectionSetUp()
1159: Level: intermediate
1161: .seealso: PetscSectionGetOffset(), PetscSectionCreate(), PetscSectionSetUp()
1162: @*/
1163: PetscErrorCode PetscSectionSetFieldOffset(PetscSection s, PetscInt point, PetscInt field, PetscInt offset)
1164: {
1168: 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);
1169: PetscSectionSetOffset(s->field[field], point, offset);
1170: return(0);
1171: }
1175: /* This gives the offset on a point of the field, ignoring constraints */
1176: PetscErrorCode PetscSectionGetFieldPointOffset(PetscSection s, PetscInt point, PetscInt field, PetscInt *offset)
1177: {
1178: PetscInt off, foff;
1182: 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);
1183: PetscSectionGetOffset(s, point, &off);
1184: PetscSectionGetOffset(s->field[field], point, &foff);
1185: *offset = foff - off;
1186: return(0);
1187: }
1191: /*@
1192: PetscSectionGetOffsetRange - Return the full range of offsets [start, end)
1194: Not collective
1196: Input Parameter:
1197: . s - the PetscSection
1199: Output Parameters:
1200: + start - the minimum offset
1201: - end - one more than the maximum offset
1203: Level: intermediate
1205: .seealso: PetscSectionGetOffset(), PetscSectionCreate()
1206: @*/
1207: PetscErrorCode PetscSectionGetOffsetRange(PetscSection s, PetscInt *start, PetscInt *end)
1208: {
1209: PetscInt os = 0, oe = 0, pStart, pEnd, p;
1213: if (s->atlasOff) {os = s->atlasOff[0]; oe = s->atlasOff[0];}
1214: PetscSectionGetChart(s, &pStart, &pEnd);
1215: for (p = 0; p < pEnd-pStart; ++p) {
1216: PetscInt dof = s->atlasDof[p], off = s->atlasOff[p];
1218: if (off >= 0) {
1219: os = PetscMin(os, off);
1220: oe = PetscMax(oe, off+dof);
1221: }
1222: }
1223: if (start) *start = os;
1224: if (end) *end = oe;
1225: return(0);
1226: }
1230: PetscErrorCode PetscSectionCreateSubsection(PetscSection s, PetscInt numFields, PetscInt fields[], PetscSection *subs)
1231: {
1232: PetscInt nF, f, pStart, pEnd, p, maxCdof = 0;
1236: if (!numFields) return(0);
1237: PetscSectionGetNumFields(s, &nF);
1238: if (numFields > nF) SETERRQ2(s->atlasLayout.comm, PETSC_ERR_ARG_WRONG, "Number of requested fields %d greater than number of fields %d", numFields, nF);
1239: PetscSectionCreate(s->atlasLayout.comm, subs);
1240: PetscSectionSetNumFields(*subs, numFields);
1241: for (f = 0; f < numFields; ++f) {
1242: const char *name = NULL;
1243: PetscInt numComp = 0;
1245: PetscSectionGetFieldName(s, fields[f], &name);
1246: PetscSectionSetFieldName(*subs, f, name);
1247: PetscSectionGetFieldComponents(s, fields[f], &numComp);
1248: PetscSectionSetFieldComponents(*subs, f, numComp);
1249: }
1250: PetscSectionGetChart(s, &pStart, &pEnd);
1251: PetscSectionSetChart(*subs, pStart, pEnd);
1252: for (p = pStart; p < pEnd; ++p) {
1253: PetscInt dof = 0, cdof = 0, fdof = 0, cfdof = 0;
1255: for (f = 0; f < numFields; ++f) {
1256: PetscSectionGetFieldDof(s, p, fields[f], &fdof);
1257: PetscSectionSetFieldDof(*subs, p, f, fdof);
1258: PetscSectionGetFieldConstraintDof(s, p, fields[f], &cfdof);
1259: if (cfdof) {PetscSectionSetFieldConstraintDof(*subs, p, f, cfdof);}
1260: dof += fdof;
1261: cdof += cfdof;
1262: }
1263: PetscSectionSetDof(*subs, p, dof);
1264: if (cdof) {PetscSectionSetConstraintDof(*subs, p, cdof);}
1265: maxCdof = PetscMax(cdof, maxCdof);
1266: }
1267: PetscSectionSetUp(*subs);
1268: if (maxCdof) {
1269: PetscInt *indices;
1271: PetscMalloc1(maxCdof, &indices);
1272: for (p = pStart; p < pEnd; ++p) {
1273: PetscInt cdof;
1275: PetscSectionGetConstraintDof(*subs, p, &cdof);
1276: if (cdof) {
1277: const PetscInt *oldIndices = NULL;
1278: PetscInt fdof = 0, cfdof = 0, fc, numConst = 0, fOff = 0;
1280: for (f = 0; f < numFields; ++f) {
1281: PetscInt oldFoff = 0, oldf;
1283: PetscSectionGetFieldDof(s, p, fields[f], &fdof);
1284: PetscSectionGetFieldConstraintDof(s, p, fields[f], &cfdof);
1285: PetscSectionGetFieldConstraintIndices(s, p, fields[f], &oldIndices);
1286: /* This can be sped up if we assume sorted fields */
1287: for (oldf = 0; oldf < fields[f]; ++oldf) {
1288: PetscInt oldfdof = 0;
1289: PetscSectionGetFieldDof(s, p, oldf, &oldfdof);
1290: oldFoff += oldfdof;
1291: }
1292: for (fc = 0; fc < cfdof; ++fc) indices[numConst+fc] = oldIndices[fc] + (fOff - oldFoff);
1293: PetscSectionSetFieldConstraintIndices(*subs, p, f, &indices[numConst]);
1294: numConst += cfdof;
1295: fOff += fdof;
1296: }
1297: PetscSectionSetConstraintIndices(*subs, p, indices);
1298: }
1299: }
1300: PetscFree(indices);
1301: }
1302: return(0);
1303: }
1307: PetscErrorCode PetscSectionCreateSubmeshSection(PetscSection s, IS subpointMap, PetscSection *subs)
1308: {
1309: const PetscInt *points = NULL, *indices = NULL;
1310: PetscInt numFields, f, numSubpoints = 0, pStart, pEnd, p, subp;
1314: PetscSectionGetNumFields(s, &numFields);
1315: PetscSectionCreate(s->atlasLayout.comm, subs);
1316: if (numFields) {PetscSectionSetNumFields(*subs, numFields);}
1317: for (f = 0; f < numFields; ++f) {
1318: const char *name = NULL;
1319: PetscInt numComp = 0;
1321: PetscSectionGetFieldName(s, f, &name);
1322: PetscSectionSetFieldName(*subs, f, name);
1323: PetscSectionGetFieldComponents(s, f, &numComp);
1324: PetscSectionSetFieldComponents(*subs, f, numComp);
1325: }
1326: /* For right now, we do not try to squeeze the subchart */
1327: if (subpointMap) {
1328: ISGetSize(subpointMap, &numSubpoints);
1329: ISGetIndices(subpointMap, &points);
1330: }
1331: PetscSectionGetChart(s, &pStart, &pEnd);
1332: PetscSectionSetChart(*subs, 0, numSubpoints);
1333: for (p = pStart; p < pEnd; ++p) {
1334: PetscInt dof, cdof, fdof = 0, cfdof = 0;
1336: PetscFindInt(p, numSubpoints, points, &subp);
1337: if (subp < 0) continue;
1338: for (f = 0; f < numFields; ++f) {
1339: PetscSectionGetFieldDof(s, p, f, &fdof);
1340: PetscSectionSetFieldDof(*subs, subp, f, fdof);
1341: PetscSectionGetFieldConstraintDof(s, p, f, &cfdof);
1342: if (cfdof) {PetscSectionSetFieldConstraintDof(*subs, subp, f, cfdof);}
1343: }
1344: PetscSectionGetDof(s, p, &dof);
1345: PetscSectionSetDof(*subs, subp, dof);
1346: PetscSectionGetConstraintDof(s, p, &cdof);
1347: if (cdof) {PetscSectionSetConstraintDof(*subs, subp, cdof);}
1348: }
1349: PetscSectionSetUp(*subs);
1350: /* Change offsets to original offsets */
1351: for (p = pStart; p < pEnd; ++p) {
1352: PetscInt off, foff = 0;
1354: PetscFindInt(p, numSubpoints, points, &subp);
1355: if (subp < 0) continue;
1356: for (f = 0; f < numFields; ++f) {
1357: PetscSectionGetFieldOffset(s, p, f, &foff);
1358: PetscSectionSetFieldOffset(*subs, subp, f, foff);
1359: }
1360: PetscSectionGetOffset(s, p, &off);
1361: PetscSectionSetOffset(*subs, subp, off);
1362: }
1363: /* Copy constraint indices */
1364: for (subp = 0; subp < numSubpoints; ++subp) {
1365: PetscInt cdof;
1367: PetscSectionGetConstraintDof(*subs, subp, &cdof);
1368: if (cdof) {
1369: for (f = 0; f < numFields; ++f) {
1370: PetscSectionGetFieldConstraintIndices(s, points[subp], f, &indices);
1371: PetscSectionSetFieldConstraintIndices(*subs, subp, f, indices);
1372: }
1373: PetscSectionGetConstraintIndices(s, points[subp], &indices);
1374: PetscSectionSetConstraintIndices(*subs, subp, indices);
1375: }
1376: }
1377: if (subpointMap) {ISRestoreIndices(subpointMap, &points);}
1378: return(0);
1379: }
1383: PetscErrorCode PetscSectionView_ASCII(PetscSection s, PetscViewer viewer)
1384: {
1385: PetscInt p;
1386: PetscMPIInt rank;
1390: if (s->atlasLayout.numDof != 1) SETERRQ1(PETSC_COMM_SELF, PETSC_ERR_SUP, "Cannot handle %d dof in a uniform section", s->atlasLayout.numDof);
1391: MPI_Comm_rank(PetscObjectComm((PetscObject)viewer), &rank);
1392: PetscViewerASCIISynchronizedAllow(viewer, PETSC_TRUE);
1393: PetscViewerASCIISynchronizedPrintf(viewer, "Process %d:\n", rank);
1394: for (p = 0; p < s->atlasLayout.pEnd - s->atlasLayout.pStart; ++p) {
1395: if ((s->bc) && (s->bc->atlasDof[p] > 0)) {
1396: PetscInt b;
1398: PetscViewerASCIISynchronizedPrintf(viewer, " (%4d) dim %2d offset %3d constrained", p+s->atlasLayout.pStart, s->atlasDof[p], s->atlasOff[p]);
1399: for (b = 0; b < s->bc->atlasDof[p]; ++b) {
1400: PetscViewerASCIISynchronizedPrintf(viewer, " %d", s->bcIndices[s->bc->atlasOff[p]+b]);
1401: }
1402: PetscViewerASCIISynchronizedPrintf(viewer, "\n");
1403: } else {
1404: PetscViewerASCIISynchronizedPrintf(viewer, " (%4d) dim %2d offset %3d\n", p+s->atlasLayout.pStart, s->atlasDof[p], s->atlasOff[p]);
1405: }
1406: }
1407: PetscViewerFlush(viewer);
1408: return(0);
1409: }
1413: /*@C
1414: PetscSectionView - Views a PetscSection
1416: Collective on PetscSection
1418: Input Parameters:
1419: + s - the PetscSection object to view
1420: - v - the viewer
1422: Level: developer
1424: .seealso PetscSectionCreate(), PetscSectionDestroy()
1425: @*/
1426: PetscErrorCode PetscSectionView(PetscSection s, PetscViewer viewer)
1427: {
1428: PetscBool isascii;
1429: PetscInt f;
1433: if (!viewer) {PetscViewerASCIIGetStdout(PETSC_COMM_SELF, &viewer);}
1435: PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &isascii);
1436: if (isascii) {
1437: if (s->numFields) {
1438: PetscViewerASCIIPrintf(viewer, "PetscSection with %d fields\n", s->numFields);
1439: for (f = 0; f < s->numFields; ++f) {
1440: PetscViewerASCIIPrintf(viewer, " field %d with %d components\n", f, s->numFieldComponents[f]);
1441: PetscSectionView_ASCII(s->field[f], viewer);
1442: }
1443: } else {
1444: PetscViewerASCIIPrintf(viewer, "PetscSection\n");
1445: PetscSectionView_ASCII(s, viewer);
1446: }
1447: }
1448: return(0);
1449: }
1453: /*@
1454: PetscSectionReset - Frees all section data.
1456: Not collective
1458: Input Parameters:
1459: . s - the PetscSection
1461: Level: developer
1463: .seealso: PetscSection, PetscSectionCreate()
1464: @*/
1465: PetscErrorCode PetscSectionReset(PetscSection s)
1466: {
1467: PetscInt f;
1471: PetscFree(s->numFieldComponents);
1472: for (f = 0; f < s->numFields; ++f) {
1473: PetscSectionDestroy(&s->field[f]);
1474: PetscFree(s->fieldNames[f]);
1475: }
1476: PetscFree(s->fieldNames);
1477: PetscFree(s->field);
1478: PetscSectionDestroy(&s->bc);
1479: PetscFree(s->bcIndices);
1480: PetscFree2(s->atlasDof, s->atlasOff);
1481: PetscSectionDestroy(&s->clSection);
1482: ISDestroy(&s->clPoints);
1484: s->atlasLayout.pStart = -1;
1485: s->atlasLayout.pEnd = -1;
1486: s->atlasLayout.numDof = 1;
1487: s->maxDof = 0;
1488: s->setup = PETSC_FALSE;
1489: s->numFields = 0;
1490: s->clObj = NULL;
1491: return(0);
1492: }
1496: /*@
1497: PetscSectionDestroy - Frees a section object and frees its range if that exists.
1499: Not collective
1501: Input Parameters:
1502: . s - the PetscSection
1504: Level: developer
1506: The PetscSection object and methods are intended to be used in the PETSc Vec and Mat implementions; it is
1507: recommended they not be used in user codes unless you really gain something in their use.
1509: .seealso: PetscSection, PetscSectionCreate()
1510: @*/
1511: PetscErrorCode PetscSectionDestroy(PetscSection *s)
1512: {
1516: if (!*s) return(0);
1518: if (--((PetscObject)(*s))->refct > 0) {
1519: *s = NULL;
1520: return(0);
1521: }
1522: PetscSectionReset(*s);
1523: PetscHeaderDestroy(s);
1524: return(0);
1525: }
1529: PetscErrorCode VecIntGetValuesSection(PetscInt *baseArray, PetscSection s, PetscInt point, const PetscInt **values)
1530: {
1531: const PetscInt p = point - s->atlasLayout.pStart;
1534: *values = &baseArray[s->atlasOff[p]];
1535: return(0);
1536: }
1540: PetscErrorCode VecIntSetValuesSection(PetscInt *baseArray, PetscSection s, PetscInt point, const PetscInt values[], InsertMode mode)
1541: {
1542: PetscInt *array;
1543: const PetscInt p = point - s->atlasLayout.pStart;
1544: const PetscInt orientation = 0; /* Needs to be included for use in closure operations */
1545: PetscInt cDim = 0;
1549: PetscSectionGetConstraintDof(s, p, &cDim);
1550: array = &baseArray[s->atlasOff[p]];
1551: if (!cDim) {
1552: if (orientation >= 0) {
1553: const PetscInt dim = s->atlasDof[p];
1554: PetscInt i;
1556: if (mode == INSERT_VALUES) {
1557: for (i = 0; i < dim; ++i) array[i] = values[i];
1558: } else {
1559: for (i = 0; i < dim; ++i) array[i] += values[i];
1560: }
1561: } else {
1562: PetscInt offset = 0;
1563: PetscInt j = -1, field, i;
1565: for (field = 0; field < s->numFields; ++field) {
1566: const PetscInt dim = s->field[field]->atlasDof[p];
1568: for (i = dim-1; i >= 0; --i) array[++j] = values[i+offset];
1569: offset += dim;
1570: }
1571: }
1572: } else {
1573: if (orientation >= 0) {
1574: const PetscInt dim = s->atlasDof[p];
1575: PetscInt cInd = 0, i;
1576: const PetscInt *cDof;
1578: PetscSectionGetConstraintIndices(s, point, &cDof);
1579: if (mode == INSERT_VALUES) {
1580: for (i = 0; i < dim; ++i) {
1581: if ((cInd < cDim) && (i == cDof[cInd])) {++cInd; continue;}
1582: array[i] = values[i];
1583: }
1584: } else {
1585: for (i = 0; i < dim; ++i) {
1586: if ((cInd < cDim) && (i == cDof[cInd])) {++cInd; continue;}
1587: array[i] += values[i];
1588: }
1589: }
1590: } else {
1591: const PetscInt *cDof;
1592: PetscInt offset = 0;
1593: PetscInt cOffset = 0;
1594: PetscInt j = 0, field;
1596: PetscSectionGetConstraintIndices(s, point, &cDof);
1597: for (field = 0; field < s->numFields; ++field) {
1598: const PetscInt dim = s->field[field]->atlasDof[p]; /* PetscSectionGetFieldDof() */
1599: const PetscInt tDim = s->field[field]->bc->atlasDof[p]; /* PetscSectionGetFieldConstraintDof() */
1600: const PetscInt sDim = dim - tDim;
1601: PetscInt cInd = 0, i,k;
1603: for (i = 0, k = dim+offset-1; i < dim; ++i, ++j, --k) {
1604: if ((cInd < sDim) && (j == cDof[cInd+cOffset])) {++cInd; continue;}
1605: array[j] = values[k];
1606: }
1607: offset += dim;
1608: cOffset += dim - tDim;
1609: }
1610: }
1611: }
1612: return(0);
1613: }
1617: PetscErrorCode PetscSectionGetConstraintIndices(PetscSection s, PetscInt point, const PetscInt **indices)
1618: {
1622: if (s->bc) {
1623: VecIntGetValuesSection(s->bcIndices, s->bc, point, indices);
1624: } else *indices = NULL;
1625: return(0);
1626: }
1630: PetscErrorCode PetscSectionSetConstraintIndices(PetscSection s, PetscInt point, const PetscInt indices[])
1631: {
1635: if (s->bc) {
1636: VecIntSetValuesSection(s->bcIndices, s->bc, point, indices, INSERT_VALUES);
1637: }
1638: return(0);
1639: }
1643: PetscErrorCode PetscSectionGetFieldConstraintIndices(PetscSection s, PetscInt point, PetscInt field, const PetscInt **indices)
1644: {
1648: 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);
1649: PetscSectionGetConstraintIndices(s->field[field], point, indices);
1650: return(0);
1651: }
1655: PetscErrorCode PetscSectionSetFieldConstraintIndices(PetscSection s, PetscInt point, PetscInt field, const PetscInt indices[])
1656: {
1660: 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);
1661: PetscSectionSetConstraintIndices(s->field[field], point, indices);
1662: return(0);
1663: }
1667: /*@
1668: PetscSectionPermute - Reorder the section according to the input point permutation
1670: Collective on PetscSection
1672: Input Parameter:
1673: + section - The PetscSection object
1674: - perm - The point permutation, old point p becomes new point perm[p]
1676: Output Parameter:
1677: . sectionNew - The permuted PetscSection
1679: Level: intermediate
1681: .keywords: mesh
1682: .seealso: MatPermute()
1683: @*/
1684: PetscErrorCode PetscSectionPermute(PetscSection section, IS permutation, PetscSection *sectionNew)
1685: {
1686: PetscSection s = section, sNew;
1687: const PetscInt *perm;
1688: PetscInt numFields, f, numPoints, pStart, pEnd, p;
1689: PetscErrorCode ierr;
1695: PetscSectionCreate(PetscObjectComm((PetscObject) s), &sNew);
1696: PetscSectionGetNumFields(s, &numFields);
1697: if (numFields) {PetscSectionSetNumFields(sNew, numFields);}
1698: for (f = 0; f < numFields; ++f) {
1699: const char *name;
1700: PetscInt numComp;
1702: PetscSectionGetFieldName(s, f, &name);
1703: PetscSectionSetFieldName(sNew, f, name);
1704: PetscSectionGetFieldComponents(s, f, &numComp);
1705: PetscSectionSetFieldComponents(sNew, f, numComp);
1706: }
1707: ISGetLocalSize(permutation, &numPoints);
1708: ISGetIndices(permutation, &perm);
1709: PetscSectionGetChart(s, &pStart, &pEnd);
1710: PetscSectionSetChart(sNew, pStart, pEnd);
1711: if (numPoints < pEnd) SETERRQ2(PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "Permutation size %d is less than largest Section point %d", numPoints, pEnd);
1712: for (p = pStart; p < pEnd; ++p) {
1713: PetscInt dof, cdof;
1715: PetscSectionGetDof(s, p, &dof);
1716: PetscSectionSetDof(sNew, perm[p], dof);
1717: PetscSectionGetConstraintDof(s, p, &cdof);
1718: if (cdof) {PetscSectionSetConstraintDof(sNew, perm[p], cdof);}
1719: for (f = 0; f < numFields; ++f) {
1720: PetscSectionGetFieldDof(s, p, f, &dof);
1721: PetscSectionSetFieldDof(sNew, perm[p], f, dof);
1722: PetscSectionGetFieldConstraintDof(s, p, f, &cdof);
1723: if (cdof) {PetscSectionSetFieldConstraintDof(sNew, perm[p], f, cdof);}
1724: }
1725: }
1726: PetscSectionSetUp(sNew);
1727: for (p = pStart; p < pEnd; ++p) {
1728: const PetscInt *cind;
1729: PetscInt cdof;
1731: PetscSectionGetConstraintDof(s, p, &cdof);
1732: if (cdof) {
1733: PetscSectionGetConstraintIndices(s, p, &cind);
1734: PetscSectionSetConstraintIndices(sNew, perm[p], cind);
1735: }
1736: for (f = 0; f < numFields; ++f) {
1737: PetscSectionGetFieldConstraintDof(s, p, f, &cdof);
1738: if (cdof) {
1739: PetscSectionGetFieldConstraintIndices(s, p, f, &cind);
1740: PetscSectionSetFieldConstraintIndices(sNew, perm[p], f, cind);
1741: }
1742: }
1743: }
1744: ISRestoreIndices(permutation, &perm);
1745: *sectionNew = sNew;
1746: return(0);
1747: }
1749: /*
1750: I need extract and merge routines for section based on fields. This should be trivial except for updating the
1751: constraint indices.
1753: Then I need a new interface for DMCreateDecomposition that takes groups of fields and returns a real DMPlex
1754: that shares the mesh parts, and has the extracted section
1755: */
1759: PetscErrorCode PetscSFConvertPartition(PetscSF sfPart, PetscSection partSection, IS partition, ISLocalToGlobalMapping *renumbering, PetscSF *sf)
1760: {
1761: MPI_Comm comm;
1762: PetscSF sfPoints;
1763: PetscInt *partSizes,*partOffsets,p,i,numParts,numMyPoints,numPoints,count;
1764: const PetscInt *partArray;
1765: PetscSFNode *sendPoints;
1766: PetscMPIInt rank;
1770: PetscObjectGetComm((PetscObject)sfPart,&comm);
1771: MPI_Comm_rank(comm, &rank);
1773: /* Get the number of parts and sizes that I have to distribute */
1774: PetscSFGetGraph(sfPart,NULL,&numParts,NULL,NULL);
1775: PetscMalloc2(numParts,&partSizes,numParts,&partOffsets);
1776: for (p=0,numPoints=0; p<numParts; p++) {
1777: PetscSectionGetDof(partSection, p, &partSizes[p]);
1778: numPoints += partSizes[p];
1779: }
1780: numMyPoints = 0;
1781: PetscSFFetchAndOpBegin(sfPart,MPIU_INT,&numMyPoints,partSizes,partOffsets,MPIU_SUM);
1782: PetscSFFetchAndOpEnd(sfPart,MPIU_INT,&numMyPoints,partSizes,partOffsets,MPIU_SUM);
1783: /* I will receive numMyPoints. I will send a total of numPoints, to be placed on remote procs at partOffsets */
1785: /* Create SF mapping locations (addressed through partition, as indexed leaves) to new owners (roots) */
1786: PetscMalloc1(numPoints,&sendPoints);
1787: for (p=0,count=0; p<numParts; p++) {
1788: for (i=0; i<partSizes[p]; i++) {
1789: sendPoints[count].rank = p;
1790: sendPoints[count].index = partOffsets[p]+i;
1791: count++;
1792: }
1793: }
1794: if (count != numPoints) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Count %D should equal numPoints=%D",count,numPoints);
1795: PetscFree2(partSizes,partOffsets);
1796: ISGetIndices(partition,&partArray);
1797: PetscSFCreate(comm,&sfPoints);
1798: PetscSFSetFromOptions(sfPoints);
1799: PetscSFSetGraph(sfPoints,numMyPoints,numPoints,partArray,PETSC_USE_POINTER,sendPoints,PETSC_OWN_POINTER);
1801: /* Invert SF so that the new owners are leaves and the locations indexed through partition are the roots */
1802: PetscSFCreateInverseSF(sfPoints,sf);
1803: PetscSFDestroy(&sfPoints);
1804: ISRestoreIndices(partition,&partArray);
1806: /* Create the new local-to-global mapping */
1807: ISLocalToGlobalMappingCreateSF(*sf,0,renumbering);
1808: return(0);
1809: }
1813: /*@C
1814: PetscSFDistributeSection - Create a new PetscSection reorganized, moving from the root to the leaves of the SF
1816: Collective
1818: Input Parameters:
1819: + sf - The SF
1820: - rootSection - Section defined on root space
1822: Output Parameters:
1823: + remoteOffsets - root offsets in leaf storage, or NULL
1824: - leafSection - Section defined on the leaf space
1826: Level: intermediate
1828: .seealso: PetscSFCreate()
1829: @*/
1830: PetscErrorCode PetscSFDistributeSection(PetscSF sf, PetscSection rootSection, PetscInt **remoteOffsets, PetscSection leafSection)
1831: {
1832: PetscSF embedSF;
1833: const PetscInt *ilocal, *indices;
1834: IS selected;
1835: PetscInt numFields, nleaves, rpStart, rpEnd, lpStart = PETSC_MAX_INT, lpEnd = -1, i, f;
1839: PetscSectionGetNumFields(rootSection, &numFields);
1840: if (numFields) {PetscSectionSetNumFields(leafSection, numFields);}
1841: for (f = 0; f < numFields; ++f) {
1842: PetscInt numComp = 0;
1843: PetscSectionGetFieldComponents(rootSection, f, &numComp);
1844: PetscSectionSetFieldComponents(leafSection, f, numComp);
1845: }
1846: PetscSectionGetChart(rootSection, &rpStart, &rpEnd);
1847: ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected);
1848: ISGetIndices(selected, &indices);
1849: PetscSFCreateEmbeddedSF(sf, rpEnd - rpStart, indices, &embedSF);
1850: ISRestoreIndices(selected, &indices);
1851: ISDestroy(&selected);
1852: PetscSFGetGraph(embedSF, NULL, &nleaves, &ilocal, NULL);
1853: if (ilocal) {
1854: for (i = 0; i < nleaves; ++i) {
1855: lpStart = PetscMin(lpStart, ilocal[i]);
1856: lpEnd = PetscMax(lpEnd, ilocal[i]);
1857: }
1858: ++lpEnd;
1859: } else {
1860: lpStart = 0;
1861: lpEnd = nleaves;
1862: }
1863: PetscSectionSetChart(leafSection, lpStart, lpEnd);
1864: /* Could fuse these at the cost of a copy and extra allocation */
1865: PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart]);
1866: PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasDof[-rpStart], &leafSection->atlasDof[-lpStart]);
1867: for (f = 0; f < numFields; ++f) {
1868: PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart]);
1869: PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->field[f]->atlasDof[-rpStart], &leafSection->field[f]->atlasDof[-lpStart]);
1870: }
1871: if (remoteOffsets) {
1872: PetscMalloc1((lpEnd - lpStart), remoteOffsets);
1873: PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart]);
1874: PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart]);
1875: }
1876: PetscSFDestroy(&embedSF);
1877: PetscSectionSetUp(leafSection);
1878: return(0);
1879: }
1883: PetscErrorCode PetscSFCreateRemoteOffsets(PetscSF sf, PetscSection rootSection, PetscSection leafSection, PetscInt **remoteOffsets)
1884: {
1885: PetscSF embedSF;
1886: const PetscInt *indices;
1887: IS selected;
1888: PetscInt numRoots, rpStart = 0, rpEnd = 0, lpStart = 0, lpEnd = 0;
1889: PetscErrorCode ierr;
1892: *remoteOffsets = NULL;
1893: PetscSFGetGraph(sf, &numRoots, NULL, NULL, NULL);
1894: if (numRoots < 0) return(0);
1895: PetscSectionGetChart(rootSection, &rpStart, &rpEnd);
1896: PetscSectionGetChart(leafSection, &lpStart, &lpEnd);
1897: ISCreateStride(PETSC_COMM_SELF, rpEnd - rpStart, rpStart, 1, &selected);
1898: ISGetIndices(selected, &indices);
1899: PetscSFCreateEmbeddedSF(sf, rpEnd - rpStart, indices, &embedSF);
1900: ISRestoreIndices(selected, &indices);
1901: ISDestroy(&selected);
1902: PetscMalloc1((lpEnd - lpStart), remoteOffsets);
1903: PetscSFBcastBegin(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart]);
1904: PetscSFBcastEnd(embedSF, MPIU_INT, &rootSection->atlasOff[-rpStart], &(*remoteOffsets)[-lpStart]);
1905: PetscSFDestroy(&embedSF);
1906: return(0);
1907: }
1911: /*@C
1912: PetscSFCreateSectionSF - Create an expanded SF of dofs, assuming the input SF relates points
1914: Input Parameters:
1915: + sf - The SF
1916: . rootSection - Data layout of remote points for outgoing data (this is usually the serial section), or NULL
1917: - remoteOffsets - Offsets for point data on remote processes (these are offsets from the root section), or NULL
1919: Output Parameters:
1920: + leafSection - Data layout of local points for incoming data (this is the distributed section)
1921: - sectionSF - The new SF
1923: Note: Either rootSection or remoteOffsets can be specified
1925: Level: intermediate
1927: .seealso: PetscSFCreate()
1928: @*/
1929: PetscErrorCode PetscSFCreateSectionSF(PetscSF sf, PetscSection rootSection, PetscInt remoteOffsets[], PetscSection leafSection, PetscSF *sectionSF)
1930: {
1931: MPI_Comm comm;
1932: const PetscInt *localPoints;
1933: const PetscSFNode *remotePoints;
1934: PetscInt lpStart, lpEnd;
1935: PetscInt numRoots, numSectionRoots, numPoints, numIndices = 0;
1936: PetscInt *localIndices;
1937: PetscSFNode *remoteIndices;
1938: PetscInt i, ind;
1939: PetscErrorCode ierr;
1947: PetscObjectGetComm((PetscObject)sf,&comm);
1948: PetscSFCreate(comm, sectionSF);
1949: PetscSectionGetChart(leafSection, &lpStart, &lpEnd);
1950: PetscSectionGetStorageSize(rootSection, &numSectionRoots);
1951: PetscSFGetGraph(sf, &numRoots, &numPoints, &localPoints, &remotePoints);
1952: if (numRoots < 0) return(0);
1953: for (i = 0; i < numPoints; ++i) {
1954: PetscInt localPoint = localPoints ? localPoints[i] : i;
1955: PetscInt dof;
1957: if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
1958: PetscSectionGetDof(leafSection, localPoint, &dof);
1959: numIndices += dof;
1960: }
1961: }
1962: PetscMalloc1(numIndices, &localIndices);
1963: PetscMalloc1(numIndices, &remoteIndices);
1964: /* Create new index graph */
1965: for (i = 0, ind = 0; i < numPoints; ++i) {
1966: PetscInt localPoint = localPoints ? localPoints[i] : i;
1967: PetscInt rank = remotePoints[i].rank;
1969: if ((localPoint >= lpStart) && (localPoint < lpEnd)) {
1970: PetscInt remoteOffset = remoteOffsets[localPoint-lpStart];
1971: PetscInt loff, dof, d;
1973: PetscSectionGetOffset(leafSection, localPoint, &loff);
1974: PetscSectionGetDof(leafSection, localPoint, &dof);
1975: for (d = 0; d < dof; ++d, ++ind) {
1976: localIndices[ind] = loff+d;
1977: remoteIndices[ind].rank = rank;
1978: remoteIndices[ind].index = remoteOffset+d;
1979: }
1980: }
1981: }
1982: PetscFree(remoteOffsets);
1983: if (numIndices != ind) SETERRQ2(comm, PETSC_ERR_PLIB, "Inconsistency in indices, %d should be %d", ind, numIndices);
1984: PetscSFSetGraph(*sectionSF, numSectionRoots, numIndices, localIndices, PETSC_OWN_POINTER, remoteIndices, PETSC_OWN_POINTER);
1985: return(0);
1986: }
1990: /*@
1991: PetscSectionSetClosureIndex - Set a cache of points in the closure of each point in the section
1993: Input Parameters:
1994: + section - The PetscSection
1995: . obj - A PetscObject which serves as the key for this index
1996: . clSection - Section giving the size of the closure of each point
1997: - clPoints - IS giving the points in each closure
1999: Note: We compress out closure points with no dofs in this section
2001: Level: intermediate
2003: .seealso: PetscSectionGetClosureIndex(), DMPlexCreateClosureIndex()
2004: @*/
2005: PetscErrorCode PetscSectionSetClosureIndex(PetscSection section, PetscObject obj, PetscSection clSection, IS clPoints)
2006: {
2010: section->clObj = obj;
2011: PetscSectionDestroy(§ion->clSection);
2012: ISDestroy(§ion->clPoints);
2013: section->clSection = clSection;
2014: section->clPoints = clPoints;
2015: return(0);
2016: }
2020: /*@
2021: PetscSectionGetClosureIndex - Get the cache of points in the closure of each point in the section
2023: Input Parameters:
2024: + section - The PetscSection
2025: - obj - A PetscObject which serves as the key for this index
2027: Output Parameters:
2028: + clSection - Section giving the size of the closure of each point
2029: - clPoints - IS giving the points in each closure
2031: Note: We compress out closure points with no dofs in this section
2033: Level: intermediate
2035: .seealso: PetscSectionSetClosureIndex(), DMPlexCreateClosureIndex()
2036: @*/
2037: PetscErrorCode PetscSectionGetClosureIndex(PetscSection section, PetscObject obj, PetscSection *clSection, IS *clPoints)
2038: {
2040: if (section->clObj == obj) {
2041: if (clSection) *clSection = section->clSection;
2042: if (clPoints) *clPoints = section->clPoints;
2043: } else {
2044: if (clSection) *clSection = NULL;
2045: if (clPoints) *clPoints = NULL;
2046: }
2047: return(0);
2048: }
2052: /*@
2053: PetscSectionGetField - Get the subsection associated with a single field
2055: Input Parameters:
2056: + s - The PetscSection
2057: - field - The field number
2059: Output Parameter:
2060: . subs - The subsection for the given field
2062: Level: intermediate
2064: .seealso: PetscSectionSetNumFields()
2065: @*/
2066: PetscErrorCode PetscSectionGetField(PetscSection s, PetscInt field, PetscSection *subs)
2067: {
2073: 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);
2074: PetscObjectReference((PetscObject) s->field[field]);
2075: *subs = s->field[field];
2076: return(0);
2077: }