xref: /petsc/src/dm/impls/plex/plextree.c (revision 66af876c15e99dd06ca9b64a8060c74d7bf4c8c6)
1d6a7ad0dSToby Isaac #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2d6a7ad0dSToby Isaac #include <../src/sys/utils/hash.h>
3d6a7ad0dSToby Isaac #include <petsc-private/isimpl.h>
4d6a7ad0dSToby Isaac #include <petscsf.h>
5d6a7ad0dSToby Isaac 
6d6a7ad0dSToby Isaac /** hierarchy routines */
7d6a7ad0dSToby Isaac 
8d6a7ad0dSToby Isaac #undef __FUNCT__
9d6a7ad0dSToby Isaac #define __FUNCT__ "DMPlexSetReferenceTree"
10d6a7ad0dSToby Isaac /*@
11d6a7ad0dSToby Isaac   DMPlexSetReferenceTree - set the reference tree for hierarchically non-conforming meshes.
12d6a7ad0dSToby Isaac 
13d6a7ad0dSToby Isaac   Not collective
14d6a7ad0dSToby Isaac 
15d6a7ad0dSToby Isaac   Input Parameters:
16d6a7ad0dSToby Isaac + dm - The DMPlex object
17d6a7ad0dSToby Isaac - ref - The reference tree DMPlex object
18d6a7ad0dSToby Isaac 
190b7167a0SToby Isaac   Level: intermediate
20d6a7ad0dSToby Isaac 
21da43764aSToby Isaac .seealso: DMPlexGetReferenceTree(), DMPlexCreateDefaultReferenceTree()
22d6a7ad0dSToby Isaac @*/
23d6a7ad0dSToby Isaac PetscErrorCode DMPlexSetReferenceTree(DM dm, DM ref)
24d6a7ad0dSToby Isaac {
25d6a7ad0dSToby Isaac   DM_Plex        *mesh = (DM_Plex *)dm->data;
26d6a7ad0dSToby Isaac   PetscErrorCode  ierr;
27d6a7ad0dSToby Isaac 
28d6a7ad0dSToby Isaac   PetscFunctionBegin;
29d6a7ad0dSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
30d6a7ad0dSToby Isaac   PetscValidHeaderSpecific(ref, DM_CLASSID, 2);
31d6a7ad0dSToby Isaac   ierr = PetscObjectReference((PetscObject)ref);CHKERRQ(ierr);
32d6a7ad0dSToby Isaac   ierr = DMDestroy(&mesh->referenceTree);CHKERRQ(ierr);
33d6a7ad0dSToby Isaac   mesh->referenceTree = ref;
34d6a7ad0dSToby Isaac   PetscFunctionReturn(0);
35d6a7ad0dSToby Isaac }
36d6a7ad0dSToby Isaac 
37d6a7ad0dSToby Isaac #undef __FUNCT__
38d6a7ad0dSToby Isaac #define __FUNCT__ "DMPlexGetReferenceTree"
39d6a7ad0dSToby Isaac /*@
40d6a7ad0dSToby Isaac   DMPlexGetReferenceTree - get the reference tree for hierarchically non-conforming meshes.
41d6a7ad0dSToby Isaac 
42d6a7ad0dSToby Isaac   Not collective
43d6a7ad0dSToby Isaac 
44d6a7ad0dSToby Isaac   Input Parameters:
45d6a7ad0dSToby Isaac . dm - The DMPlex object
46d6a7ad0dSToby Isaac 
47d6a7ad0dSToby Isaac   Output Parameters
48d6a7ad0dSToby Isaac . ref - The reference tree DMPlex object
49d6a7ad0dSToby Isaac 
500b7167a0SToby Isaac   Level: intermediate
51d6a7ad0dSToby Isaac 
52da43764aSToby Isaac .seealso: DMPlexSetReferenceTree(), DMPlexCreateDefaultReferenceTree()
53d6a7ad0dSToby Isaac @*/
54d6a7ad0dSToby Isaac PetscErrorCode DMPlexGetReferenceTree(DM dm, DM *ref)
55d6a7ad0dSToby Isaac {
56d6a7ad0dSToby Isaac   DM_Plex        *mesh = (DM_Plex *)dm->data;
57d6a7ad0dSToby Isaac 
58d6a7ad0dSToby Isaac   PetscFunctionBegin;
59d6a7ad0dSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
60d6a7ad0dSToby Isaac   PetscValidPointer(ref,2);
61d6a7ad0dSToby Isaac   *ref = mesh->referenceTree;
62d6a7ad0dSToby Isaac   PetscFunctionReturn(0);
63d6a7ad0dSToby Isaac }
64d6a7ad0dSToby Isaac 
65da43764aSToby Isaac #undef __FUNCT__
66da43764aSToby Isaac #define __FUNCT__ "DMPlexCreateDefaultReferenceTree"
67da43764aSToby Isaac /*@
68da43764aSToby Isaac   DMPlexCreateDefaultReferenceTree - create a reference tree for isotropic hierarchical mesh refinement.
69da43764aSToby Isaac 
70da43764aSToby Isaac   Collective on comm
71da43764aSToby Isaac 
72da43764aSToby Isaac   Input Parameters:
73da43764aSToby Isaac + comm    - the MPI communicator
74da43764aSToby Isaac . dim     - the spatial dimension
75da43764aSToby Isaac - simplex - Flag for simplex, otherwise use a tensor-product cell
76da43764aSToby Isaac 
77da43764aSToby Isaac   Output Parameters:
78da43764aSToby Isaac . ref     - the reference tree DMPlex object
79da43764aSToby Isaac 
800b7167a0SToby Isaac   Level: intermediate
81da43764aSToby Isaac 
82da43764aSToby Isaac .keywords: reference cell
83da43764aSToby Isaac .seealso: DMPlexSetReferenceTree(), DMPlexGetReferenceTree()
84da43764aSToby Isaac @*/
85da43764aSToby Isaac PetscErrorCode DMPlexCreateDefaultReferenceTree(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *ref)
86da43764aSToby Isaac {
87da43764aSToby Isaac   DM             K, Kref;
8810f7e118SToby Isaac   PetscInt       p, pStart, pEnd, pRefStart, pRefEnd, d, offset, parentSize, *parents, *childIDs;
89da43764aSToby Isaac   PetscInt      *permvals, *unionCones, *coneSizes, *unionOrientations, numUnionPoints, *numDimPoints, numCones, numVerts;
90da43764aSToby Isaac   DMLabel        identity, identityRef;
9110f7e118SToby Isaac   PetscSection   unionSection, unionConeSection, parentSection;
92da43764aSToby Isaac   PetscScalar   *unionCoords;
93da43764aSToby Isaac   IS             perm;
94da43764aSToby Isaac   PetscErrorCode ierr;
95da43764aSToby Isaac 
96da43764aSToby Isaac   PetscFunctionBegin;
97da43764aSToby Isaac   /* create a reference element */
98da43764aSToby Isaac   ierr = DMPlexCreateReferenceCell(comm, dim, simplex, &K);CHKERRQ(ierr);
99da43764aSToby Isaac   ierr = DMPlexCreateLabel(K, "identity");CHKERRQ(ierr);
100da43764aSToby Isaac   ierr = DMPlexGetLabel(K, "identity", &identity);CHKERRQ(ierr);
101da43764aSToby Isaac   ierr = DMPlexGetChart(K, &pStart, &pEnd);CHKERRQ(ierr);
102da43764aSToby Isaac   for (p = pStart; p < pEnd; p++) {
103da43764aSToby Isaac     ierr = DMLabelSetValue(identity, p, p);CHKERRQ(ierr);
104da43764aSToby Isaac   }
105da43764aSToby Isaac   /* refine it */
106da43764aSToby Isaac   ierr = DMRefine(K,comm,&Kref);CHKERRQ(ierr);
107da43764aSToby Isaac 
108da43764aSToby Isaac   /* the reference tree is the union of these two, without duplicating
109da43764aSToby Isaac    * points that appear in both */
110da43764aSToby Isaac   ierr = DMPlexGetLabel(Kref, "identity", &identityRef);CHKERRQ(ierr);
111da43764aSToby Isaac   ierr = DMPlexGetChart(Kref, &pRefStart, &pRefEnd);CHKERRQ(ierr);
112da43764aSToby Isaac   ierr = PetscSectionCreate(comm, &unionSection);CHKERRQ(ierr);
113da43764aSToby Isaac   ierr = PetscSectionSetChart(unionSection, 0, (pEnd - pStart) + (pRefEnd - pRefStart));CHKERRQ(ierr);
114da43764aSToby Isaac   /* count points that will go in the union */
115da43764aSToby Isaac   for (p = pStart; p < pEnd; p++) {
116da43764aSToby Isaac     ierr = PetscSectionSetDof(unionSection, p - pStart, 1);CHKERRQ(ierr);
117da43764aSToby Isaac   }
118da43764aSToby Isaac   for (p = pRefStart; p < pRefEnd; p++) {
119da43764aSToby Isaac     PetscInt q, qSize;
120da43764aSToby Isaac     ierr = DMLabelGetValue(identityRef, p, &q);CHKERRQ(ierr);
121da43764aSToby Isaac     ierr = DMLabelGetStratumSize(identityRef, q, &qSize);CHKERRQ(ierr);
122da43764aSToby Isaac     if (qSize > 1) {
123da43764aSToby Isaac       ierr = PetscSectionSetDof(unionSection, p - pRefStart + (pEnd - pStart), 1);CHKERRQ(ierr);
124da43764aSToby Isaac     }
125da43764aSToby Isaac   }
126da43764aSToby Isaac   ierr = PetscMalloc1((pEnd - pStart) + (pRefEnd - pRefStart),&permvals);CHKERRQ(ierr);
127da43764aSToby Isaac   offset = 0;
128da43764aSToby Isaac   /* stratify points in the union by topological dimension */
129da43764aSToby Isaac   for (d = 0; d <= dim; d++) {
130da43764aSToby Isaac     PetscInt cStart, cEnd, c;
131da43764aSToby Isaac 
132da43764aSToby Isaac     ierr = DMPlexGetHeightStratum(K, d, &cStart, &cEnd);CHKERRQ(ierr);
133da43764aSToby Isaac     for (c = cStart; c < cEnd; c++) {
134da43764aSToby Isaac       permvals[offset++] = c;
135da43764aSToby Isaac     }
136da43764aSToby Isaac 
137da43764aSToby Isaac     ierr = DMPlexGetHeightStratum(Kref, d, &cStart, &cEnd);CHKERRQ(ierr);
138da43764aSToby Isaac     for (c = cStart; c < cEnd; c++) {
139da43764aSToby Isaac       permvals[offset++] = c + (pEnd - pStart);
140da43764aSToby Isaac     }
141da43764aSToby Isaac   }
142da43764aSToby Isaac   ierr = ISCreateGeneral(comm, (pEnd - pStart) + (pRefEnd - pRefStart), permvals, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr);
143da43764aSToby Isaac   ierr = PetscSectionSetPermutation(unionSection,perm);CHKERRQ(ierr);
144da43764aSToby Isaac   ierr = PetscSectionSetUp(unionSection);CHKERRQ(ierr);
145da43764aSToby Isaac   ierr = PetscSectionGetStorageSize(unionSection,&numUnionPoints);CHKERRQ(ierr);
146da43764aSToby Isaac   ierr = PetscMalloc2(numUnionPoints,&coneSizes,dim+1,&numDimPoints);CHKERRQ(ierr);
147da43764aSToby Isaac   /* count dimension points */
148da43764aSToby Isaac   for (d = 0; d <= dim; d++) {
149da43764aSToby Isaac     PetscInt cStart, cOff, cOff2;
150da43764aSToby Isaac     ierr = DMPlexGetHeightStratum(K,d,&cStart,NULL);CHKERRQ(ierr);
151da43764aSToby Isaac     ierr = PetscSectionGetOffset(unionSection,cStart-pStart,&cOff);CHKERRQ(ierr);
152da43764aSToby Isaac     if (d < dim) {
153da43764aSToby Isaac       ierr = DMPlexGetHeightStratum(K,d+1,&cStart,NULL);CHKERRQ(ierr);
154da43764aSToby Isaac       ierr = PetscSectionGetOffset(unionSection,cStart-pStart,&cOff2);CHKERRQ(ierr);
155da43764aSToby Isaac     }
156da43764aSToby Isaac     else {
157da43764aSToby Isaac       cOff2 = numUnionPoints;
158da43764aSToby Isaac     }
159da43764aSToby Isaac     numDimPoints[dim - d] = cOff2 - cOff;
160da43764aSToby Isaac   }
161da43764aSToby Isaac   ierr = PetscSectionCreate(comm, &unionConeSection);CHKERRQ(ierr);
162da43764aSToby Isaac   ierr = PetscSectionSetChart(unionConeSection, 0, numUnionPoints);CHKERRQ(ierr);
163da43764aSToby Isaac   /* count the cones in the union */
164da43764aSToby Isaac   for (p = pStart; p < pEnd; p++) {
165da43764aSToby Isaac     PetscInt dof, uOff;
166da43764aSToby Isaac 
167da43764aSToby Isaac     ierr = DMPlexGetConeSize(K, p, &dof);CHKERRQ(ierr);
168da43764aSToby Isaac     ierr = PetscSectionGetOffset(unionSection, p - pStart,&uOff);CHKERRQ(ierr);
169da43764aSToby Isaac     ierr = PetscSectionSetDof(unionConeSection, uOff, dof);CHKERRQ(ierr);
170da43764aSToby Isaac     coneSizes[uOff] = dof;
171da43764aSToby Isaac   }
172da43764aSToby Isaac   for (p = pRefStart; p < pRefEnd; p++) {
173da43764aSToby Isaac     PetscInt dof, uDof, uOff;
174da43764aSToby Isaac 
175da43764aSToby Isaac     ierr = DMPlexGetConeSize(Kref, p, &dof);CHKERRQ(ierr);
176da43764aSToby Isaac     ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr);
177da43764aSToby Isaac     ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr);
178da43764aSToby Isaac     if (uDof) {
179da43764aSToby Isaac       ierr = PetscSectionSetDof(unionConeSection, uOff, dof);CHKERRQ(ierr);
180da43764aSToby Isaac       coneSizes[uOff] = dof;
181da43764aSToby Isaac     }
182da43764aSToby Isaac   }
183da43764aSToby Isaac   ierr = PetscSectionSetUp(unionConeSection);CHKERRQ(ierr);
184da43764aSToby Isaac   ierr = PetscSectionGetStorageSize(unionConeSection,&numCones);CHKERRQ(ierr);
185da43764aSToby Isaac   ierr = PetscMalloc2(numCones,&unionCones,numCones,&unionOrientations);CHKERRQ(ierr);
186da43764aSToby Isaac   /* write the cones in the union */
187da43764aSToby Isaac   for (p = pStart; p < pEnd; p++) {
188da43764aSToby Isaac     PetscInt dof, uOff, c, cOff;
189da43764aSToby Isaac     const PetscInt *cone, *orientation;
190da43764aSToby Isaac 
191da43764aSToby Isaac     ierr = DMPlexGetConeSize(K, p, &dof);CHKERRQ(ierr);
192da43764aSToby Isaac     ierr = DMPlexGetCone(K, p, &cone);CHKERRQ(ierr);
193da43764aSToby Isaac     ierr = DMPlexGetConeOrientation(K, p, &orientation);CHKERRQ(ierr);
194da43764aSToby Isaac     ierr = PetscSectionGetOffset(unionSection, p - pStart,&uOff);CHKERRQ(ierr);
195da43764aSToby Isaac     ierr = PetscSectionGetOffset(unionConeSection,uOff,&cOff);CHKERRQ(ierr);
196da43764aSToby Isaac     for (c = 0; c < dof; c++) {
197da43764aSToby Isaac       PetscInt e, eOff;
198da43764aSToby Isaac       e                           = cone[c];
199da43764aSToby Isaac       ierr                        = PetscSectionGetOffset(unionSection, e - pStart, &eOff);CHKERRQ(ierr);
200da43764aSToby Isaac       unionCones[cOff + c]        = eOff;
201da43764aSToby Isaac       unionOrientations[cOff + c] = orientation[c];
202da43764aSToby Isaac     }
203da43764aSToby Isaac   }
204da43764aSToby Isaac   for (p = pRefStart; p < pRefEnd; p++) {
205da43764aSToby Isaac     PetscInt dof, uDof, uOff, c, cOff;
206da43764aSToby Isaac     const PetscInt *cone, *orientation;
207da43764aSToby Isaac 
208da43764aSToby Isaac     ierr = DMPlexGetConeSize(Kref, p, &dof);CHKERRQ(ierr);
209da43764aSToby Isaac     ierr = DMPlexGetCone(Kref, p, &cone);CHKERRQ(ierr);
210da43764aSToby Isaac     ierr = DMPlexGetConeOrientation(Kref, p, &orientation);CHKERRQ(ierr);
211da43764aSToby Isaac     ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr);
212da43764aSToby Isaac     ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr);
213da43764aSToby Isaac     if (uDof) {
214da43764aSToby Isaac       ierr = PetscSectionGetOffset(unionConeSection,uOff,&cOff);CHKERRQ(ierr);
215da43764aSToby Isaac       for (c = 0; c < dof; c++) {
216da43764aSToby Isaac         PetscInt e, eOff, eDof;
217da43764aSToby Isaac 
218da43764aSToby Isaac         e    = cone[c];
219da43764aSToby Isaac         ierr = PetscSectionGetDof(unionSection, e - pRefStart + (pEnd - pStart),&eDof);CHKERRQ(ierr);
220da43764aSToby Isaac         if (eDof) {
221da43764aSToby Isaac           ierr = PetscSectionGetOffset(unionSection, e - pRefStart + (pEnd - pStart), &eOff);CHKERRQ(ierr);
222da43764aSToby Isaac         }
223da43764aSToby Isaac         else {
224da43764aSToby Isaac           ierr = DMLabelGetValue(identityRef, e, &e);CHKERRQ(ierr);
225da43764aSToby Isaac           ierr = PetscSectionGetOffset(unionSection, e - pStart, &eOff);CHKERRQ(ierr);
226da43764aSToby Isaac         }
227da43764aSToby Isaac         unionCones[cOff + c]        = eOff;
228da43764aSToby Isaac         unionOrientations[cOff + c] = orientation[c];
229da43764aSToby Isaac       }
230da43764aSToby Isaac     }
231da43764aSToby Isaac   }
232da43764aSToby Isaac   /* get the coordinates */
233da43764aSToby Isaac   {
234da43764aSToby Isaac     PetscInt vStart, vEnd, vRefStart, vRefEnd, v, vDof, vOff;
235da43764aSToby Isaac     PetscSection KcoordsSec, KrefCoordsSec;
236da43764aSToby Isaac     Vec      KcoordsVec, KrefCoordsVec;
237da43764aSToby Isaac     PetscScalar *Kcoords;
238da43764aSToby Isaac 
239da43764aSToby Isaac     DMGetCoordinateSection(K, &KcoordsSec);CHKERRQ(ierr);
240da43764aSToby Isaac     DMGetCoordinatesLocal(K, &KcoordsVec);CHKERRQ(ierr);
241da43764aSToby Isaac     DMGetCoordinateSection(Kref, &KrefCoordsSec);CHKERRQ(ierr);
242da43764aSToby Isaac     DMGetCoordinatesLocal(Kref, &KrefCoordsVec);CHKERRQ(ierr);
243da43764aSToby Isaac 
244da43764aSToby Isaac     numVerts = numDimPoints[0];
245da43764aSToby Isaac     ierr     = PetscMalloc1(numVerts * dim,&unionCoords);CHKERRQ(ierr);
246da43764aSToby Isaac     ierr     = DMPlexGetDepthStratum(K,0,&vStart,&vEnd);CHKERRQ(ierr);
247da43764aSToby Isaac 
248da43764aSToby Isaac     offset = 0;
249da43764aSToby Isaac     for (v = vStart; v < vEnd; v++) {
250da43764aSToby Isaac       ierr = PetscSectionGetOffset(unionSection,v - pStart,&vOff);CHKERRQ(ierr);
251da43764aSToby Isaac       ierr = VecGetValuesSection(KcoordsVec, KcoordsSec, v, &Kcoords);CHKERRQ(ierr);
252da43764aSToby Isaac       for (d = 0; d < dim; d++) {
253da43764aSToby Isaac         unionCoords[offset * dim + d] = Kcoords[d];
254da43764aSToby Isaac       }
255da43764aSToby Isaac       offset++;
256da43764aSToby Isaac     }
257da43764aSToby Isaac     ierr = DMPlexGetDepthStratum(Kref,0,&vRefStart,&vRefEnd);CHKERRQ(ierr);
258da43764aSToby Isaac     for (v = vRefStart; v < vRefEnd; v++) {
259da43764aSToby Isaac       ierr = PetscSectionGetDof(unionSection,v - pRefStart + (pEnd - pStart),&vDof);CHKERRQ(ierr);
260da43764aSToby Isaac       ierr = PetscSectionGetOffset(unionSection,v - pRefStart + (pEnd - pStart),&vOff);CHKERRQ(ierr);
261da43764aSToby Isaac       ierr = VecGetValuesSection(KrefCoordsVec, KrefCoordsSec, v, &Kcoords);CHKERRQ(ierr);
262da43764aSToby Isaac       if (vDof) {
263da43764aSToby Isaac         for (d = 0; d < dim; d++) {
264da43764aSToby Isaac           unionCoords[offset * dim + d] = Kcoords[d];
265da43764aSToby Isaac         }
266da43764aSToby Isaac         offset++;
267da43764aSToby Isaac       }
268da43764aSToby Isaac     }
269da43764aSToby Isaac   }
270da43764aSToby Isaac   ierr = DMCreate(comm,ref);CHKERRQ(ierr);
271da43764aSToby Isaac   ierr = DMSetType(*ref,DMPLEX);CHKERRQ(ierr);
272da43764aSToby Isaac   ierr = DMPlexSetDimension(*ref,dim);CHKERRQ(ierr);
273da43764aSToby Isaac   ierr = DMPlexCreateFromDAG(*ref,dim,numDimPoints,coneSizes,unionCones,unionOrientations,unionCoords);CHKERRQ(ierr);
27410f7e118SToby Isaac   /* set the tree */
27510f7e118SToby Isaac   ierr = PetscSectionCreate(comm,&parentSection);CHKERRQ(ierr);
27610f7e118SToby Isaac   ierr = PetscSectionSetChart(parentSection,0,numUnionPoints);CHKERRQ(ierr);
27710f7e118SToby Isaac   for (p = pRefStart; p < pRefEnd; p++) {
27810f7e118SToby Isaac     PetscInt uDof, uOff;
27910f7e118SToby Isaac 
28010f7e118SToby Isaac     ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr);
28110f7e118SToby Isaac     ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr);
28210f7e118SToby Isaac     if (uDof) {
28310f7e118SToby Isaac       PetscSectionSetDof(parentSection,uOff,1);CHKERRQ(ierr);
28410f7e118SToby Isaac     }
28510f7e118SToby Isaac   }
28610f7e118SToby Isaac   ierr = PetscSectionSetUp(parentSection);CHKERRQ(ierr);
28710f7e118SToby Isaac   ierr = PetscSectionGetStorageSize(parentSection,&parentSize);CHKERRQ(ierr);
28810f7e118SToby Isaac   ierr = PetscMalloc2(parentSize,&parents,parentSize,&childIDs);CHKERRQ(ierr);
28910f7e118SToby Isaac   for (p = pRefStart; p < pRefEnd; p++) {
29010f7e118SToby Isaac     PetscInt uDof, uOff;
29110f7e118SToby Isaac 
29210f7e118SToby Isaac     ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr);
29310f7e118SToby Isaac     ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr);
29410f7e118SToby Isaac     if (uDof) {
29510f7e118SToby Isaac       PetscInt pOff, parent, parentU;
29610f7e118SToby Isaac       PetscSectionGetOffset(parentSection,uOff,&pOff);CHKERRQ(ierr);
29710f7e118SToby Isaac       DMLabelGetValue(identityRef,p,&parent);CHKERRQ(ierr);
29810f7e118SToby Isaac       ierr = PetscSectionGetOffset(unionSection, parent - pStart,&parentU);CHKERRQ(ierr);
29910f7e118SToby Isaac       parents[pOff] = parentU;
30010f7e118SToby Isaac       childIDs[pOff] = uOff;
30110f7e118SToby Isaac     }
30210f7e118SToby Isaac   }
30310f7e118SToby Isaac   ierr = DMPlexSetTree(*ref,parentSection,parents,childIDs);CHKERRQ(ierr);
30410f7e118SToby Isaac   ierr = PetscSectionDestroy(&parentSection);CHKERRQ(ierr);
30510f7e118SToby Isaac   ierr = PetscFree2(parents,childIDs);CHKERRQ(ierr);
30610f7e118SToby Isaac 
307da43764aSToby Isaac   /* clean up */
308da43764aSToby Isaac   ierr = PetscSectionDestroy(&unionSection);CHKERRQ(ierr);
309da43764aSToby Isaac   ierr = PetscSectionDestroy(&unionConeSection);CHKERRQ(ierr);
310da43764aSToby Isaac   ierr = ISDestroy(&perm);CHKERRQ(ierr);
311da43764aSToby Isaac   ierr = PetscFree(unionCoords);CHKERRQ(ierr);
312da43764aSToby Isaac   ierr = PetscFree2(unionCones,unionOrientations);CHKERRQ(ierr);
313da43764aSToby Isaac   ierr = PetscFree2(coneSizes,numDimPoints);CHKERRQ(ierr);
314da43764aSToby Isaac   ierr = DMDestroy(&K);CHKERRQ(ierr);
315da43764aSToby Isaac   ierr = DMDestroy(&Kref);CHKERRQ(ierr);
316da43764aSToby Isaac   PetscFunctionReturn(0);
317da43764aSToby Isaac }
318da43764aSToby Isaac 
319d961a43aSToby Isaac #undef __FUNCT__
320878b19aaSToby Isaac #define __FUNCT__ "DMPlexTreeSymmetrize"
321878b19aaSToby Isaac static PetscErrorCode DMPlexTreeSymmetrize(DM dm)
322878b19aaSToby Isaac {
323878b19aaSToby Isaac   DM_Plex        *mesh = (DM_Plex *)dm->data;
324878b19aaSToby Isaac   PetscSection   childSec, pSec;
325878b19aaSToby Isaac   PetscInt       p, pSize, cSize, parMax = PETSC_MIN_INT, parMin = PETSC_MAX_INT;
326878b19aaSToby Isaac   PetscInt       *offsets, *children, pStart, pEnd;
327878b19aaSToby Isaac   PetscErrorCode ierr;
328878b19aaSToby Isaac 
329878b19aaSToby Isaac   PetscFunctionBegin;
330878b19aaSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
331878b19aaSToby Isaac   ierr = PetscSectionDestroy(&mesh->childSection);CHKERRQ(ierr);
332878b19aaSToby Isaac   ierr = PetscFree(mesh->children);CHKERRQ(ierr);
333878b19aaSToby Isaac   pSec = mesh->parentSection;
334878b19aaSToby Isaac   if (!pSec) PetscFunctionReturn(0);
335878b19aaSToby Isaac   ierr = PetscSectionGetStorageSize(pSec,&pSize);CHKERRQ(ierr);
336878b19aaSToby Isaac   for (p = 0; p < pSize; p++) {
337878b19aaSToby Isaac     PetscInt par = mesh->parents[p];
338878b19aaSToby Isaac 
339878b19aaSToby Isaac     parMax = PetscMax(parMax,par+1);
340878b19aaSToby Isaac     parMin = PetscMin(parMin,par);
341878b19aaSToby Isaac   }
342878b19aaSToby Isaac   if (parMin > parMax) {
343878b19aaSToby Isaac     parMin = -1;
344878b19aaSToby Isaac     parMax = -1;
345878b19aaSToby Isaac   }
346878b19aaSToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)pSec),&childSec);CHKERRQ(ierr);
347878b19aaSToby Isaac   ierr = PetscSectionSetChart(childSec,parMin,parMax);CHKERRQ(ierr);
348878b19aaSToby Isaac   for (p = 0; p < pSize; p++) {
349878b19aaSToby Isaac     PetscInt par = mesh->parents[p];
350878b19aaSToby Isaac 
351878b19aaSToby Isaac     ierr = PetscSectionAddDof(childSec,par,1);CHKERRQ(ierr);
352878b19aaSToby Isaac   }
353878b19aaSToby Isaac   ierr = PetscSectionSetUp(childSec);CHKERRQ(ierr);
354878b19aaSToby Isaac   ierr = PetscSectionGetStorageSize(childSec,&cSize);CHKERRQ(ierr);
355878b19aaSToby Isaac   ierr = PetscMalloc1(cSize,&children);CHKERRQ(ierr);
356878b19aaSToby Isaac   ierr = PetscCalloc1(parMax-parMin,&offsets);CHKERRQ(ierr);
357878b19aaSToby Isaac   ierr = PetscSectionGetChart(pSec,&pStart,&pEnd);CHKERRQ(ierr);
358878b19aaSToby Isaac   for (p = pStart; p < pEnd; p++) {
359878b19aaSToby Isaac     PetscInt dof, off, i;
360878b19aaSToby Isaac 
361878b19aaSToby Isaac     ierr = PetscSectionGetDof(pSec,p,&dof);CHKERRQ(ierr);
362878b19aaSToby Isaac     ierr = PetscSectionGetOffset(pSec,p,&off);CHKERRQ(ierr);
363878b19aaSToby Isaac     for (i = 0; i < dof; i++) {
364878b19aaSToby Isaac       PetscInt par = mesh->parents[off + i], cOff;
365878b19aaSToby Isaac 
366878b19aaSToby Isaac       ierr = PetscSectionGetOffset(childSec,par,&cOff);CHKERRQ(ierr);
367878b19aaSToby Isaac       children[cOff + offsets[par-parMin]++] = p;
368878b19aaSToby Isaac     }
369878b19aaSToby Isaac   }
370878b19aaSToby Isaac   mesh->childSection = childSec;
371878b19aaSToby Isaac   mesh->children = children;
372878b19aaSToby Isaac   ierr = PetscFree(offsets);CHKERRQ(ierr);
373878b19aaSToby Isaac   PetscFunctionReturn(0);
374878b19aaSToby Isaac }
375878b19aaSToby Isaac 
376878b19aaSToby Isaac #undef __FUNCT__
377*66af876cSToby Isaac #define __FUNCT__ "DMPlexComputeConstraints_Tree"
378*66af876cSToby Isaac static PetscErrorCode DMPlexComputeConstraints_Tree(DM dm)
379*66af876cSToby Isaac {
380*66af876cSToby Isaac   PetscInt       p, pStart, pEnd, *anchors, size;
381*66af876cSToby Isaac   PetscInt       aMin = PETSC_MAX_INT, aMax = PETSC_MIN_INT;
382*66af876cSToby Isaac   PetscSection   aSec;
383*66af876cSToby Isaac   IS             aIS;
384*66af876cSToby Isaac   PetscErrorCode ierr;
385*66af876cSToby Isaac 
386*66af876cSToby Isaac   PetscFunctionBegin;
387*66af876cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
388*66af876cSToby Isaac   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
389*66af876cSToby Isaac   for (p = pStart; p < pEnd; p++) {
390*66af876cSToby Isaac     PetscInt parent;
391*66af876cSToby Isaac 
392*66af876cSToby Isaac     ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
393*66af876cSToby Isaac     if (parent != p) {
394*66af876cSToby Isaac       aMin = PetscMin(aMin,p);
395*66af876cSToby Isaac       aMax = PetscMax(aMax,p+1);
396*66af876cSToby Isaac     }
397*66af876cSToby Isaac   }
398*66af876cSToby Isaac   if (aMin > aMax) {
399*66af876cSToby Isaac     aMin = -1;
400*66af876cSToby Isaac     aMax = -1;
401*66af876cSToby Isaac   }
402*66af876cSToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&aSec);CHKERRQ(ierr);
403*66af876cSToby Isaac   ierr = PetscSectionSetChart(aSec,aMin,aMax);CHKERRQ(ierr);
404*66af876cSToby Isaac   for (p = aMin; p < aMax; p++) {
405*66af876cSToby Isaac     PetscInt parent, ancestor = p;
406*66af876cSToby Isaac 
407*66af876cSToby Isaac     ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
408*66af876cSToby Isaac     while (parent != ancestor) {
409*66af876cSToby Isaac       ancestor = parent;
410*66af876cSToby Isaac       ierr     = DMPlexGetTreeParent(dm,ancestor,&parent,NULL);CHKERRQ(ierr);
411*66af876cSToby Isaac     }
412*66af876cSToby Isaac     if (ancestor != p) {
413*66af876cSToby Isaac       PetscInt closureSize, *closure = NULL;
414*66af876cSToby Isaac 
415*66af876cSToby Isaac       ierr = DMPlexGetTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
416*66af876cSToby Isaac       ierr = PetscSectionSetDof(aSec,p,closureSize);CHKERRQ(ierr);
417*66af876cSToby Isaac       ierr = DMPlexRestoreTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
418*66af876cSToby Isaac     }
419*66af876cSToby Isaac   }
420*66af876cSToby Isaac   ierr = PetscSectionSetUp(aSec);CHKERRQ(ierr);
421*66af876cSToby Isaac   ierr = PetscSectionGetStorageSize(aSec,&size);CHKERRQ(ierr);
422*66af876cSToby Isaac   ierr = PetscMalloc1(size,&anchors);CHKERRQ(ierr);
423*66af876cSToby Isaac   for (p = aMin; p < aMax; p++) {
424*66af876cSToby Isaac     PetscInt parent, ancestor = p;
425*66af876cSToby Isaac 
426*66af876cSToby Isaac     ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
427*66af876cSToby Isaac     while (parent != ancestor) {
428*66af876cSToby Isaac       ancestor = parent;
429*66af876cSToby Isaac       ierr     = DMPlexGetTreeParent(dm,ancestor,&parent,NULL);CHKERRQ(ierr);
430*66af876cSToby Isaac     }
431*66af876cSToby Isaac     if (ancestor != p) {
432*66af876cSToby Isaac       PetscInt j, closureSize, *closure = NULL, aOff;
433*66af876cSToby Isaac 
434*66af876cSToby Isaac       ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr);
435*66af876cSToby Isaac 
436*66af876cSToby Isaac       ierr = DMPlexGetTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
437*66af876cSToby Isaac       for (j = 0; j < closureSize; j++) {
438*66af876cSToby Isaac         anchors[aOff + j] = closure[2*j];
439*66af876cSToby Isaac       }
440*66af876cSToby Isaac       ierr = DMPlexRestoreTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
441*66af876cSToby Isaac     }
442*66af876cSToby Isaac   }
443*66af876cSToby Isaac   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm),size,anchors,PETSC_OWN_POINTER,&aIS);CHKERRQ(ierr);
444*66af876cSToby Isaac   ierr = DMPlexSetConstraints(dm,aSec,aIS);CHKERRQ(ierr);
445*66af876cSToby Isaac   ierr = PetscSectionDestroy(&aSec);CHKERRQ(ierr);
446*66af876cSToby Isaac   ierr = ISDestroy(&aIS);CHKERRQ(ierr);
447*66af876cSToby Isaac   PetscFunctionReturn(0);
448*66af876cSToby Isaac }
449*66af876cSToby Isaac 
450*66af876cSToby Isaac #undef __FUNCT__
4510b7167a0SToby Isaac #define __FUNCT__ "DMPlexSetTree"
4520b7167a0SToby Isaac /*@
4530b7167a0SToby Isaac   DMPlexSetTree - set the tree that describes the hierarchy of non-conforming mesh points.  This routine also creates
4540b7167a0SToby Isaac   the point-to-point constraints determined by the tree: a point is constained to the points in the closure of its
4550b7167a0SToby Isaac   tree root.
4560b7167a0SToby Isaac 
4570b7167a0SToby Isaac   Collective on dm
4580b7167a0SToby Isaac 
4590b7167a0SToby Isaac   Input Parameters:
4600b7167a0SToby Isaac + dm - the DMPlex object
4610b7167a0SToby Isaac . parentSection - a section describing the tree: a point has a parent if it has 1 dof in the section; the section
4620b7167a0SToby Isaac                   offset indexes the parent and childID list; the reference count of parentSection is incremented
4630b7167a0SToby Isaac . parents - a list of the point parents; copied, can be destroyed
4640b7167a0SToby Isaac - childIDs - identifies the relationship of the child point to the parent point; if there is a reference tree, then
4650b7167a0SToby Isaac              the child corresponds to the point in the reference tree with index childIDs; copied, can be destroyed
4660b7167a0SToby Isaac 
4670b7167a0SToby Isaac   Level: intermediate
4680b7167a0SToby Isaac 
469b2f41788SToby Isaac .seealso: DMPlexGetTree(), DMPlexSetReferenceTree(), DMPlexSetConstraints(), DMPlexGetTreeParent(), DMPlexGetTreeChildren()
4700b7167a0SToby Isaac @*/
471b2f41788SToby Isaac PetscErrorCode DMPlexSetTree(DM dm, PetscSection parentSection, PetscInt parents[], PetscInt childIDs[])
4720b7167a0SToby Isaac {
4730b7167a0SToby Isaac   DM_Plex       *mesh = (DM_Plex *)dm->data;
4740b7167a0SToby Isaac   PetscInt       size;
4750b7167a0SToby Isaac   PetscErrorCode ierr;
4760b7167a0SToby Isaac 
4770b7167a0SToby Isaac   PetscFunctionBegin;
4780b7167a0SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
4790b7167a0SToby Isaac   PetscValidHeaderSpecific(parentSection, PETSC_SECTION_CLASSID, 2);
4800b7167a0SToby Isaac   ierr = PetscObjectReference((PetscObject)parentSection);CHKERRQ(ierr);
4810b7167a0SToby Isaac   ierr = PetscSectionDestroy(&mesh->parentSection);CHKERRQ(ierr);
4820b7167a0SToby Isaac   mesh->parentSection = parentSection;
4830b7167a0SToby Isaac   ierr = PetscSectionGetStorageSize(parentSection,&size);CHKERRQ(ierr);
4840b7167a0SToby Isaac   if (parents != mesh->parents) {
4850b7167a0SToby Isaac     ierr = PetscFree(mesh->parents);CHKERRQ(ierr);
4860b7167a0SToby Isaac     ierr = PetscMalloc1(size,&mesh->parents);CHKERRQ(ierr);
4870b7167a0SToby Isaac     ierr = PetscMemcpy(mesh->parents, parents, size * sizeof(*parents));CHKERRQ(ierr);
4880b7167a0SToby Isaac   }
4890b7167a0SToby Isaac   if (childIDs != mesh->childIDs) {
4900b7167a0SToby Isaac     ierr = PetscFree(mesh->childIDs);CHKERRQ(ierr);
4910b7167a0SToby Isaac     ierr = PetscMalloc1(size,&mesh->childIDs);CHKERRQ(ierr);
4920b7167a0SToby Isaac     ierr = PetscMemcpy(mesh->childIDs, childIDs, size * sizeof(*childIDs));CHKERRQ(ierr);
4930b7167a0SToby Isaac   }
494878b19aaSToby Isaac   ierr = DMPlexTreeSymmetrize(dm);CHKERRQ(ierr);
495*66af876cSToby Isaac   ierr = DMPlexComputeConstraints_Tree(dm);CHKERRQ(ierr);
4960b7167a0SToby Isaac   PetscFunctionReturn(0);
4970b7167a0SToby Isaac }
4980b7167a0SToby Isaac 
4990b7167a0SToby Isaac #undef __FUNCT__
500b2f41788SToby Isaac #define __FUNCT__ "DMPlexGetTree"
501b2f41788SToby Isaac /*@
502b2f41788SToby Isaac   DMPlexGetTree - get the tree that describes the hierarchy of non-conforming mesh points.
503b2f41788SToby Isaac   Collective on dm
504b2f41788SToby Isaac 
505b2f41788SToby Isaac   Input Parameters:
506b2f41788SToby Isaac . dm - the DMPlex object
507b2f41788SToby Isaac 
508b2f41788SToby Isaac   Output Parameters:
509b2f41788SToby Isaac + parentSection - a section describing the tree: a point has a parent if it has 1 dof in the section; the section
510b2f41788SToby Isaac                   offset indexes the parent and childID list
511b2f41788SToby Isaac . parents - a list of the point parents
512b2f41788SToby Isaac . childIDs - identifies the relationship of the child point to the parent point; if there is a reference tree, then
513b2f41788SToby Isaac              the child corresponds to the point in the reference tree with index childID
514b2f41788SToby Isaac . childSection - the inverse of the parent section
515b2f41788SToby Isaac - children - a list of the point children
516b2f41788SToby Isaac 
517b2f41788SToby Isaac   Level: intermediate
518b2f41788SToby Isaac 
519b2f41788SToby Isaac .seealso: DMPlexSetTree(), DMPlexSetReferenceTree(), DMPlexSetConstraints(), DMPlexGetTreeParent(), DMPlexGetTreeChildren()
520b2f41788SToby Isaac @*/
521b2f41788SToby Isaac PetscErrorCode DMPlexGetTree(DM dm, PetscSection *parentSection, PetscInt *parents[], PetscInt *childIDs[], PetscSection *childSection, PetscInt *children[])
522b2f41788SToby Isaac {
523b2f41788SToby Isaac   DM_Plex        *mesh = (DM_Plex *)dm->data;
524b2f41788SToby Isaac 
525b2f41788SToby Isaac   PetscFunctionBegin;
526b2f41788SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
527b2f41788SToby Isaac   if (parentSection) *parentSection = mesh->parentSection;
528b2f41788SToby Isaac   if (parents)       *parents       = mesh->parents;
529b2f41788SToby Isaac   if (childIDs)      *childIDs      = mesh->childIDs;
530b2f41788SToby Isaac   if (childSection)  *childSection  = mesh->childSection;
531b2f41788SToby Isaac   if (children)      *children      = mesh->children;
532b2f41788SToby Isaac   PetscFunctionReturn(0);
533b2f41788SToby Isaac }
534b2f41788SToby Isaac 
535b2f41788SToby Isaac #undef __FUNCT__
536d961a43aSToby Isaac #define __FUNCT__ "DMPlexGetTreeParent"
537d961a43aSToby Isaac /*@
538d961a43aSToby Isaac   DMPlexGetTreeParent - get the parent of a point in the tree describing the point hierarchy (not the Sieve DAG)
539d961a43aSToby Isaac 
540d961a43aSToby Isaac   Input Parameters:
541d961a43aSToby Isaac + dm - the DMPlex object
542d961a43aSToby Isaac - point - the query point
543d961a43aSToby Isaac 
544d961a43aSToby Isaac   Output Parameters:
545d961a43aSToby Isaac + parent - if not NULL, set to the parent of the point, or the point itself if the point does not have a parent
546d961a43aSToby Isaac - childID - if not NULL, set to the child ID of the point with respect to its parent, or 0 if the point
547d961a43aSToby Isaac             does not have a parent
548d961a43aSToby Isaac 
549d961a43aSToby Isaac   Level: intermediate
550d961a43aSToby Isaac 
551d961a43aSToby Isaac .seealso: DMPlexSetTree(), DMPlexGetTree(), DMPlexGetTreeChildren()
552d961a43aSToby Isaac @*/
553d961a43aSToby Isaac PetscErrorCode DMPlexGetTreeParent(DM dm, PetscInt point, PetscInt *parent, PetscInt *childID)
554d961a43aSToby Isaac {
555d961a43aSToby Isaac   DM_Plex       *mesh = (DM_Plex *)dm->data;
556d961a43aSToby Isaac   PetscSection   pSec;
557d961a43aSToby Isaac   PetscErrorCode ierr;
558d961a43aSToby Isaac 
559d961a43aSToby Isaac   PetscFunctionBegin;
560d961a43aSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
561d961a43aSToby Isaac   pSec = mesh->parentSection;
562d961a43aSToby Isaac   if (pSec && point >= pSec->pStart && point < pSec->pEnd) {
563d961a43aSToby Isaac     PetscInt dof;
564d961a43aSToby Isaac 
565d961a43aSToby Isaac     ierr = PetscSectionGetDof (pSec, point, &dof);CHKERRQ(ierr);
566d961a43aSToby Isaac     if (dof) {
567d961a43aSToby Isaac       PetscInt off;
568d961a43aSToby Isaac 
569d961a43aSToby Isaac       ierr = PetscSectionGetOffset (pSec, point, &off);CHKERRQ(ierr);
570d961a43aSToby Isaac       if (parent)  *parent = mesh->parents[off];
571d961a43aSToby Isaac       if (childID) *childID = mesh->childIDs[off];
572d961a43aSToby Isaac       PetscFunctionReturn(0);
573d961a43aSToby Isaac     }
574d961a43aSToby Isaac   }
575d961a43aSToby Isaac   if (parent) {
576d961a43aSToby Isaac     *parent = point;
577d961a43aSToby Isaac   }
578d961a43aSToby Isaac   if (childID) {
579d961a43aSToby Isaac     *childID = 0;
580d961a43aSToby Isaac   }
581d961a43aSToby Isaac   PetscFunctionReturn(0);
582d961a43aSToby Isaac }
583d961a43aSToby Isaac 
584d961a43aSToby Isaac #undef __FUNCT__
585d961a43aSToby Isaac #define __FUNCT__ "DMPlexGetTreeChildren"
586d961a43aSToby Isaac /*@C
587d961a43aSToby Isaac   DMPlexGetTreeChildren - get the children of a point in the tree describing the point hierarchy (not the Sieve DAG)
588d961a43aSToby Isaac 
589d961a43aSToby Isaac   Input Parameters:
590d961a43aSToby Isaac + dm - the DMPlex object
591d961a43aSToby Isaac - point - the query point
592d961a43aSToby Isaac 
593d961a43aSToby Isaac   Output Parameters:
594d961a43aSToby Isaac + numChildren - if not NULL, set to the number of children
595d961a43aSToby Isaac - children - if not NULL, set to a list children, or set to NULL if the point has no children
596d961a43aSToby Isaac 
597d961a43aSToby Isaac   Level: intermediate
598d961a43aSToby Isaac 
599d961a43aSToby Isaac   Fortran Notes:
600d961a43aSToby Isaac   Since it returns an array, this routine is only available in Fortran 90, and you must
601d961a43aSToby Isaac   include petsc.h90 in your code.
602d961a43aSToby Isaac 
603d961a43aSToby Isaac .seealso: DMPlexSetTree(), DMPlexGetTree(), DMPlexGetTreeParent()
604d961a43aSToby Isaac @*/
605d961a43aSToby Isaac PetscErrorCode DMPlexGetTreeChildren(DM dm, PetscInt point, PetscInt *numChildren, const PetscInt *children[])
606d961a43aSToby Isaac {
607d961a43aSToby Isaac   DM_Plex       *mesh = (DM_Plex *)dm->data;
608d961a43aSToby Isaac   PetscSection   childSec;
609d961a43aSToby Isaac   PetscInt       dof = 0;
610d961a43aSToby Isaac   PetscErrorCode ierr;
611d961a43aSToby Isaac 
612d961a43aSToby Isaac   PetscFunctionBegin;
613d961a43aSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
614d961a43aSToby Isaac   childSec = mesh->childSection;
615d961a43aSToby Isaac   if (childSec && point >= childSec->pStart && point < childSec->pEnd) {
616d961a43aSToby Isaac     ierr = PetscSectionGetDof (childSec, point, &dof);CHKERRQ(ierr);
617d961a43aSToby Isaac   }
618d961a43aSToby Isaac   if (numChildren) *numChildren = dof;
619d961a43aSToby Isaac   if (children) {
620d961a43aSToby Isaac     if (dof) {
621d961a43aSToby Isaac       PetscInt off;
622d961a43aSToby Isaac 
623d961a43aSToby Isaac       ierr = PetscSectionGetOffset (childSec, point, &off);CHKERRQ(ierr);
624d961a43aSToby Isaac       *children = &mesh->children[off];
625d961a43aSToby Isaac     }
626d961a43aSToby Isaac     else {
627d961a43aSToby Isaac       *children = NULL;
628d961a43aSToby Isaac     }
629d961a43aSToby Isaac   }
630d961a43aSToby Isaac   PetscFunctionReturn(0);
631d961a43aSToby Isaac }
632