xref: /petsc/src/dm/impls/plex/plextree.c (revision 10f7e118563f850904d540fe4bb7f55e1ac84ee8)
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;
88*10f7e118SToby 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;
91*10f7e118SToby 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);
274*10f7e118SToby Isaac   /* set the tree */
275*10f7e118SToby Isaac   ierr = PetscSectionCreate(comm,&parentSection);CHKERRQ(ierr);
276*10f7e118SToby Isaac   ierr = PetscSectionSetChart(parentSection,0,numUnionPoints);CHKERRQ(ierr);
277*10f7e118SToby Isaac   for (p = pRefStart; p < pRefEnd; p++) {
278*10f7e118SToby Isaac     PetscInt uDof, uOff;
279*10f7e118SToby Isaac 
280*10f7e118SToby Isaac     ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr);
281*10f7e118SToby Isaac     ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr);
282*10f7e118SToby Isaac     if (uDof) {
283*10f7e118SToby Isaac       PetscSectionSetDof(parentSection,uOff,1);CHKERRQ(ierr);
284*10f7e118SToby Isaac     }
285*10f7e118SToby Isaac   }
286*10f7e118SToby Isaac   ierr = PetscSectionSetUp(parentSection);CHKERRQ(ierr);
287*10f7e118SToby Isaac   ierr = PetscSectionGetStorageSize(parentSection,&parentSize);CHKERRQ(ierr);
288*10f7e118SToby Isaac   ierr = PetscMalloc2(parentSize,&parents,parentSize,&childIDs);CHKERRQ(ierr);
289*10f7e118SToby Isaac   for (p = pRefStart; p < pRefEnd; p++) {
290*10f7e118SToby Isaac     PetscInt uDof, uOff;
291*10f7e118SToby Isaac 
292*10f7e118SToby Isaac     ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr);
293*10f7e118SToby Isaac     ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr);
294*10f7e118SToby Isaac     if (uDof) {
295*10f7e118SToby Isaac       PetscInt pOff, parent, parentU;
296*10f7e118SToby Isaac       PetscSectionGetOffset(parentSection,uOff,&pOff);CHKERRQ(ierr);
297*10f7e118SToby Isaac       DMLabelGetValue(identityRef,p,&parent);CHKERRQ(ierr);
298*10f7e118SToby Isaac       ierr = PetscSectionGetOffset(unionSection, parent - pStart,&parentU);CHKERRQ(ierr);
299*10f7e118SToby Isaac       parents[pOff] = parentU;
300*10f7e118SToby Isaac       childIDs[pOff] = uOff;
301*10f7e118SToby Isaac     }
302*10f7e118SToby Isaac   }
303*10f7e118SToby Isaac   ierr = DMPlexSetTree(*ref,parentSection,parents,childIDs);CHKERRQ(ierr);
304*10f7e118SToby Isaac   ierr = PetscSectionDestroy(&parentSection);CHKERRQ(ierr);
305*10f7e118SToby Isaac   ierr = PetscFree2(parents,childIDs);CHKERRQ(ierr);
306*10f7e118SToby 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__
3200b7167a0SToby Isaac #define __FUNCT__ "DMPlexSetTree"
3210b7167a0SToby Isaac /*@
3220b7167a0SToby Isaac   DMPlexSetTree - set the tree that describes the hierarchy of non-conforming mesh points.  This routine also creates
3230b7167a0SToby Isaac   the point-to-point constraints determined by the tree: a point is constained to the points in the closure of its
3240b7167a0SToby Isaac   tree root.
3250b7167a0SToby Isaac 
3260b7167a0SToby Isaac   Collective on dm
3270b7167a0SToby Isaac 
3280b7167a0SToby Isaac   Input Parameters:
3290b7167a0SToby Isaac + dm - the DMPlex object
3300b7167a0SToby Isaac . parentSection - a section describing the tree: a point has a parent if it has 1 dof in the section; the section
3310b7167a0SToby Isaac                   offset indexes the parent and childID list; the reference count of parentSection is incremented
3320b7167a0SToby Isaac . parents - a list of the point parents; copied, can be destroyed
3330b7167a0SToby Isaac - childIDs - identifies the relationship of the child point to the parent point; if there is a reference tree, then
3340b7167a0SToby Isaac              the child corresponds to the point in the reference tree with index childIDs; copied, can be destroyed
3350b7167a0SToby Isaac 
3360b7167a0SToby Isaac   Level: intermediate
3370b7167a0SToby Isaac 
3380b7167a0SToby Isaac .seealso: DMPlexSetReferenceTree(), DMPlexSetConstraints(), DMPlexGetTreeParent()
3390b7167a0SToby Isaac @*/
3400b7167a0SToby Isaac PetscErrorCode DMPlexSetTree(DM dm, PetscSection parentSection, PetscInt *parents, PetscInt *childIDs)
3410b7167a0SToby Isaac {
3420b7167a0SToby Isaac   DM_Plex       *mesh = (DM_Plex *)dm->data;
3430b7167a0SToby Isaac   PetscInt       size;
3440b7167a0SToby Isaac   PetscErrorCode ierr;
3450b7167a0SToby Isaac 
3460b7167a0SToby Isaac   PetscFunctionBegin;
3470b7167a0SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
3480b7167a0SToby Isaac   PetscValidHeaderSpecific(parentSection, PETSC_SECTION_CLASSID, 2);
3490b7167a0SToby Isaac   ierr = PetscObjectReference((PetscObject)parentSection);CHKERRQ(ierr);
3500b7167a0SToby Isaac   ierr = PetscSectionDestroy(&mesh->parentSection);CHKERRQ(ierr);
3510b7167a0SToby Isaac   mesh->parentSection = parentSection;
3520b7167a0SToby Isaac   ierr = PetscSectionGetStorageSize(parentSection,&size);CHKERRQ(ierr);
3530b7167a0SToby Isaac   if (parents != mesh->parents) {
3540b7167a0SToby Isaac     ierr = PetscFree(mesh->parents);CHKERRQ(ierr);
3550b7167a0SToby Isaac     ierr = PetscMalloc1(size,&mesh->parents);CHKERRQ(ierr);
3560b7167a0SToby Isaac     ierr = PetscMemcpy(mesh->parents, parents, size * sizeof(*parents));CHKERRQ(ierr);
3570b7167a0SToby Isaac   }
3580b7167a0SToby Isaac   if (childIDs != mesh->childIDs) {
3590b7167a0SToby Isaac     ierr = PetscFree(mesh->childIDs);CHKERRQ(ierr);
3600b7167a0SToby Isaac     ierr = PetscMalloc1(size,&mesh->childIDs);CHKERRQ(ierr);
3610b7167a0SToby Isaac     ierr = PetscMemcpy(mesh->childIDs, childIDs, size * sizeof(*childIDs));CHKERRQ(ierr);
3620b7167a0SToby Isaac   }
3630b7167a0SToby Isaac   PetscFunctionReturn(0);
3640b7167a0SToby Isaac }
3650b7167a0SToby Isaac 
3660b7167a0SToby Isaac #undef __FUNCT__
367d961a43aSToby Isaac #define __FUNCT__ "DMPlexGetTreeParent"
368d961a43aSToby Isaac /*@
369d961a43aSToby Isaac   DMPlexGetTreeParent - get the parent of a point in the tree describing the point hierarchy (not the Sieve DAG)
370d961a43aSToby Isaac 
371d961a43aSToby Isaac   Input Parameters:
372d961a43aSToby Isaac + dm - the DMPlex object
373d961a43aSToby Isaac - point - the query point
374d961a43aSToby Isaac 
375d961a43aSToby Isaac   Output Parameters:
376d961a43aSToby Isaac + parent - if not NULL, set to the parent of the point, or the point itself if the point does not have a parent
377d961a43aSToby Isaac - childID - if not NULL, set to the child ID of the point with respect to its parent, or 0 if the point
378d961a43aSToby Isaac             does not have a parent
379d961a43aSToby Isaac 
380d961a43aSToby Isaac   Level: intermediate
381d961a43aSToby Isaac 
382d961a43aSToby Isaac .seealso: DMPlexSetTree(), DMPlexGetTree(), DMPlexGetTreeChildren()
383d961a43aSToby Isaac @*/
384d961a43aSToby Isaac PetscErrorCode DMPlexGetTreeParent(DM dm, PetscInt point, PetscInt *parent, PetscInt *childID)
385d961a43aSToby Isaac {
386d961a43aSToby Isaac   DM_Plex       *mesh = (DM_Plex *)dm->data;
387d961a43aSToby Isaac   PetscSection   pSec;
388d961a43aSToby Isaac   PetscErrorCode ierr;
389d961a43aSToby Isaac 
390d961a43aSToby Isaac   PetscFunctionBegin;
391d961a43aSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
392d961a43aSToby Isaac   pSec = mesh->parentSection;
393d961a43aSToby Isaac   if (pSec && point >= pSec->pStart && point < pSec->pEnd) {
394d961a43aSToby Isaac     PetscInt dof;
395d961a43aSToby Isaac 
396d961a43aSToby Isaac     ierr = PetscSectionGetDof (pSec, point, &dof);CHKERRQ(ierr);
397d961a43aSToby Isaac     if (dof) {
398d961a43aSToby Isaac       PetscInt off;
399d961a43aSToby Isaac 
400d961a43aSToby Isaac       ierr = PetscSectionGetOffset (pSec, point, &off);CHKERRQ(ierr);
401d961a43aSToby Isaac       if (parent)  *parent = mesh->parents[off];
402d961a43aSToby Isaac       if (childID) *childID = mesh->childIDs[off];
403d961a43aSToby Isaac       PetscFunctionReturn(0);
404d961a43aSToby Isaac     }
405d961a43aSToby Isaac   }
406d961a43aSToby Isaac   if (parent) {
407d961a43aSToby Isaac     *parent = point;
408d961a43aSToby Isaac   }
409d961a43aSToby Isaac   if (childID) {
410d961a43aSToby Isaac     *childID = 0;
411d961a43aSToby Isaac   }
412d961a43aSToby Isaac   PetscFunctionReturn(0);
413d961a43aSToby Isaac }
414d961a43aSToby Isaac 
415d961a43aSToby Isaac #undef __FUNCT__
416d961a43aSToby Isaac #define __FUNCT__ "DMPlexGetTreeChildren"
417d961a43aSToby Isaac /*@C
418d961a43aSToby Isaac   DMPlexGetTreeChildren - get the children of a point in the tree describing the point hierarchy (not the Sieve DAG)
419d961a43aSToby Isaac 
420d961a43aSToby Isaac   Input Parameters:
421d961a43aSToby Isaac + dm - the DMPlex object
422d961a43aSToby Isaac - point - the query point
423d961a43aSToby Isaac 
424d961a43aSToby Isaac   Output Parameters:
425d961a43aSToby Isaac + numChildren - if not NULL, set to the number of children
426d961a43aSToby Isaac - children - if not NULL, set to a list children, or set to NULL if the point has no children
427d961a43aSToby Isaac 
428d961a43aSToby Isaac   Level: intermediate
429d961a43aSToby Isaac 
430d961a43aSToby Isaac   Fortran Notes:
431d961a43aSToby Isaac   Since it returns an array, this routine is only available in Fortran 90, and you must
432d961a43aSToby Isaac   include petsc.h90 in your code.
433d961a43aSToby Isaac 
434d961a43aSToby Isaac .seealso: DMPlexSetTree(), DMPlexGetTree(), DMPlexGetTreeParent()
435d961a43aSToby Isaac @*/
436d961a43aSToby Isaac PetscErrorCode DMPlexGetTreeChildren(DM dm, PetscInt point, PetscInt *numChildren, const PetscInt *children[])
437d961a43aSToby Isaac {
438d961a43aSToby Isaac   DM_Plex       *mesh = (DM_Plex *)dm->data;
439d961a43aSToby Isaac   PetscSection   childSec;
440d961a43aSToby Isaac   PetscInt       dof = 0;
441d961a43aSToby Isaac   PetscErrorCode ierr;
442d961a43aSToby Isaac 
443d961a43aSToby Isaac   PetscFunctionBegin;
444d961a43aSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
445d961a43aSToby Isaac   childSec = mesh->childSection;
446d961a43aSToby Isaac   if (childSec && point >= childSec->pStart && point < childSec->pEnd) {
447d961a43aSToby Isaac     ierr = PetscSectionGetDof (childSec, point, &dof);CHKERRQ(ierr);
448d961a43aSToby Isaac   }
449d961a43aSToby Isaac   if (numChildren) *numChildren = dof;
450d961a43aSToby Isaac   if (children) {
451d961a43aSToby Isaac     if (dof) {
452d961a43aSToby Isaac       PetscInt off;
453d961a43aSToby Isaac 
454d961a43aSToby Isaac       ierr = PetscSectionGetOffset (childSec, point, &off);CHKERRQ(ierr);
455d961a43aSToby Isaac       *children = &mesh->children[off];
456d961a43aSToby Isaac     }
457d961a43aSToby Isaac     else {
458d961a43aSToby Isaac       *children = NULL;
459d961a43aSToby Isaac     }
460d961a43aSToby Isaac   }
461d961a43aSToby Isaac   PetscFunctionReturn(0);
462d961a43aSToby Isaac }
463