xref: /petsc/src/dm/impls/plex/plextree.c (revision 0dd1b1fe4d0133f02e75266bd5a0efecbf538cde)
1 #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2 #include <../src/sys/utils/hash.h>
3 #include <petsc/private/isimpl.h>
4 #include <petsc/private/petscfeimpl.h>
5 #include <petscsf.h>
6 #include <petscds.h>
7 
8 /** hierarchy routines */
9 
10 #undef __FUNCT__
11 #define __FUNCT__ "DMPlexSetReferenceTree"
12 /*@
13   DMPlexSetReferenceTree - set the reference tree for hierarchically non-conforming meshes.
14 
15   Not collective
16 
17   Input Parameters:
18 + dm - The DMPlex object
19 - ref - The reference tree DMPlex object
20 
21   Level: intermediate
22 
23 .seealso: DMPlexGetReferenceTree(), DMPlexCreateDefaultReferenceTree()
24 @*/
25 PetscErrorCode DMPlexSetReferenceTree(DM dm, DM ref)
26 {
27   DM_Plex        *mesh = (DM_Plex *)dm->data;
28   PetscErrorCode  ierr;
29 
30   PetscFunctionBegin;
31   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
32   if (ref) {PetscValidHeaderSpecific(ref, DM_CLASSID, 2);}
33   ierr = PetscObjectReference((PetscObject)ref);CHKERRQ(ierr);
34   ierr = DMDestroy(&mesh->referenceTree);CHKERRQ(ierr);
35   mesh->referenceTree = ref;
36   PetscFunctionReturn(0);
37 }
38 
39 #undef __FUNCT__
40 #define __FUNCT__ "DMPlexGetReferenceTree"
41 /*@
42   DMPlexGetReferenceTree - get the reference tree for hierarchically non-conforming meshes.
43 
44   Not collective
45 
46   Input Parameters:
47 . dm - The DMPlex object
48 
49   Output Parameters
50 . ref - The reference tree DMPlex object
51 
52   Level: intermediate
53 
54 .seealso: DMPlexSetReferenceTree(), DMPlexCreateDefaultReferenceTree()
55 @*/
56 PetscErrorCode DMPlexGetReferenceTree(DM dm, DM *ref)
57 {
58   DM_Plex        *mesh = (DM_Plex *)dm->data;
59 
60   PetscFunctionBegin;
61   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
62   PetscValidPointer(ref,2);
63   *ref = mesh->referenceTree;
64   PetscFunctionReturn(0);
65 }
66 
67 #undef __FUNCT__
68 #define __FUNCT__ "DMPlexReferenceTreeGetChildSymmetry_Default"
69 static PetscErrorCode DMPlexReferenceTreeGetChildSymmetry_Default(DM dm, PetscInt parent, PetscInt parentOrientA, PetscInt childOrientA, PetscInt childA, PetscInt parentOrientB, PetscInt *childOrientB, PetscInt *childB)
70 {
71   PetscInt       coneSize, dStart, dEnd, dim, ABswap, oAvert, oBvert, ABswapVert;
72   PetscErrorCode ierr;
73 
74   PetscFunctionBegin;
75   if (parentOrientA == parentOrientB) {
76     if (childOrientB) *childOrientB = childOrientA;
77     if (childB) *childB = childA;
78     PetscFunctionReturn(0);
79   }
80   for (dim = 0; dim < 3; dim++) {
81     ierr = DMPlexGetDepthStratum(dm,dim,&dStart,&dEnd);CHKERRQ(ierr);
82     if (parent >= dStart && parent <= dEnd) {
83       break;
84     }
85   }
86   if (dim > 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot perform child symmetry for %d-cells",dim);
87   if (!dim) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"A vertex has no children");
88   if (childA < dStart || childA >= dEnd) {
89     /* this is a lower-dimensional child: bootstrap */
90     PetscInt size, i, sA = -1, sB, sOrientB, sConeSize;
91     const PetscInt *supp, *coneA, *coneB, *oA, *oB;
92 
93     ierr = DMPlexGetSupportSize(dm,childA,&size);CHKERRQ(ierr);
94     ierr = DMPlexGetSupport(dm,childA,&supp);CHKERRQ(ierr);
95 
96     /* find a point sA in supp(childA) that has the same parent */
97     for (i = 0; i < size; i++) {
98       PetscInt sParent;
99 
100       sA   = supp[i];
101       if (sA == parent) continue;
102       ierr = DMPlexGetTreeParent(dm,sA,&sParent,NULL);CHKERRQ(ierr);
103       if (sParent == parent) {
104         break;
105       }
106     }
107     if (i == size) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"could not find support in children");
108     /* find out which point sB is in an equivalent position to sA under
109      * parentOrientB */
110     ierr = DMPlexReferenceTreeGetChildSymmetry_Default(dm,parent,parentOrientA,0,sA,parentOrientB,&sOrientB,&sB);CHKERRQ(ierr);
111     ierr = DMPlexGetConeSize(dm,sA,&sConeSize);CHKERRQ(ierr);
112     ierr = DMPlexGetCone(dm,sA,&coneA);CHKERRQ(ierr);
113     ierr = DMPlexGetCone(dm,sB,&coneB);CHKERRQ(ierr);
114     ierr = DMPlexGetConeOrientation(dm,sA,&oA);CHKERRQ(ierr);
115     ierr = DMPlexGetConeOrientation(dm,sB,&oB);CHKERRQ(ierr);
116     /* step through the cone of sA in natural order */
117     for (i = 0; i < sConeSize; i++) {
118       if (coneA[i] == childA) {
119         /* if childA is at position i in coneA,
120          * then we want the point that is at sOrientB*i in coneB */
121         PetscInt j = (sOrientB >= 0) ? ((sOrientB + i) % sConeSize) : ((sConeSize -(sOrientB+1) - i) % sConeSize);
122         if (childB) *childB = coneB[j];
123         if (childOrientB) {
124           PetscInt oBtrue;
125 
126           ierr          = DMPlexGetConeSize(dm,childA,&coneSize);CHKERRQ(ierr);
127           /* compose sOrientB and oB[j] */
128           if (coneSize != 0 && coneSize != 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expected a vertex or an edge");
129           /* we may have to flip an edge */
130           oBtrue        = coneSize ? ((sOrientB >= 0) ? oB[j] : -(oB[j] + 2)) : 0;
131           ABswap        = DihedralSwap(coneSize,oA[i],oBtrue);CHKERRQ(ierr);
132           *childOrientB = DihedralCompose(coneSize,childOrientA,ABswap);
133         }
134         break;
135       }
136     }
137     if (i == sConeSize) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"support cone mismatch");
138     PetscFunctionReturn(0);
139   }
140   /* get the cone size and symmetry swap */
141   ierr   = DMPlexGetConeSize(dm,parent,&coneSize);CHKERRQ(ierr);
142   ABswap = DihedralSwap(coneSize, parentOrientA, parentOrientB);
143   if (dim == 2) {
144     /* orientations refer to cones: we want them to refer to vertices:
145      * if it's a rotation, they are the same, but if the order is reversed, a
146      * permutation that puts side i first does *not* put vertex i first */
147     oAvert     = (parentOrientA >= 0) ? parentOrientA : -((-parentOrientA % coneSize) + 1);
148     oBvert     = (parentOrientB >= 0) ? parentOrientB : -((-parentOrientB % coneSize) + 1);
149     ABswapVert = DihedralSwap(coneSize, oAvert, oBvert);
150   }
151   else {
152     oAvert     = parentOrientA;
153     oBvert     = parentOrientB;
154     ABswapVert = ABswap;
155   }
156   if (childB) {
157     /* assume that each child corresponds to a vertex, in the same order */
158     PetscInt p, posA = -1, numChildren, i;
159     const PetscInt *children;
160 
161     /* count which position the child is in */
162     ierr = DMPlexGetTreeChildren(dm,parent,&numChildren,&children);CHKERRQ(ierr);
163     for (i = 0; i < numChildren; i++) {
164       p = children[i];
165       if (p == childA) {
166         posA = i;
167         break;
168       }
169     }
170     if (posA >= coneSize) {
171       /* this is the triangle in the middle of a uniformly refined triangle: it
172        * is invariant */
173       if (dim != 2 || posA != 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Expected a middle triangle, got something else");
174       *childB = childA;
175     }
176     else {
177       /* figure out position B by applying ABswapVert */
178       PetscInt posB;
179 
180       posB = (ABswapVert >= 0) ? ((ABswapVert + posA) % coneSize) : ((coneSize -(ABswapVert + 1) - posA) % coneSize);
181       if (childB) *childB = children[posB];
182     }
183   }
184   if (childOrientB) *childOrientB = DihedralCompose(coneSize,childOrientA,ABswap);
185   PetscFunctionReturn(0);
186 }
187 
188 #undef __FUNCT__
189 #define __FUNCT__ "DMPlexReferenceTreeGetChildSymmetry"
190 /*@
191   DMPlexReferenceTreeGetChildSymmetry - Given a reference tree, transform a childid and orientation from one parent frame to another
192 
193   Input Parameters:
194 + dm - the reference tree DMPlex object
195 . parent - the parent point
196 . parentOrientA - the reference orientation for describing the parent
197 . childOrientA - the reference orientation for describing the child
198 . childA - the reference childID for describing the child
199 - parentOrientB - the new orientation for describing the parent
200 
201   Output Parameters:
202 + childOrientB - if not NULL, set to the new oreintation for describing the child
203 . childB - if not NULL, the new childID for describing the child
204 
205   Level: developer
206 
207 .seealso: DMPlexGetReferenceTree(), DMPlexSetReferenceTree(), DMPlexSetTree()
208 @*/
209 PetscErrorCode DMPlexReferenceTreeGetChildSymmetry(DM dm, PetscInt parent, PetscInt parentOrientA, PetscInt childOrientA, PetscInt childA, PetscInt parentOrientB, PetscInt *childOrientB, PetscInt *childB)
210 {
211   DM_Plex        *mesh = (DM_Plex *)dm->data;
212   PetscErrorCode ierr;
213 
214   PetscFunctionBegin;
215   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
216   if (!mesh->getchildsymmetry) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"DMPlexReferenceTreeGetChildSymmetry not implemented");
217   ierr = mesh->getchildsymmetry(dm,parent,parentOrientA,childOrientA,childA,parentOrientB,childOrientB,childB);CHKERRQ(ierr);
218   PetscFunctionReturn(0);
219 }
220 
221 static PetscErrorCode DMPlexSetTree_Internal(DM,PetscSection,PetscInt*,PetscInt*,PetscBool,PetscBool);
222 
223 #undef __FUNCT__
224 #define __FUNCT__ "DMPlexCreateDefaultReferenceTree"
225 /*@
226   DMPlexCreateDefaultReferenceTree - create a reference tree for isotropic hierarchical mesh refinement.
227 
228   Collective on comm
229 
230   Input Parameters:
231 + comm    - the MPI communicator
232 . dim     - the spatial dimension
233 - simplex - Flag for simplex, otherwise use a tensor-product cell
234 
235   Output Parameters:
236 . ref     - the reference tree DMPlex object
237 
238   Level: intermediate
239 
240 .keywords: reference cell
241 .seealso: DMPlexSetReferenceTree(), DMPlexGetReferenceTree()
242 @*/
243 PetscErrorCode DMPlexCreateDefaultReferenceTree(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *ref)
244 {
245   DM_Plex       *mesh;
246   DM             K, Kref;
247   PetscInt       p, pStart, pEnd, pRefStart, pRefEnd, d, offset, parentSize, *parents, *childIDs;
248   PetscInt      *permvals, *unionCones, *coneSizes, *unionOrientations, numUnionPoints, *numDimPoints, numCones, numVerts;
249   DMLabel        identity, identityRef;
250   PetscSection   unionSection, unionConeSection, parentSection;
251   PetscScalar   *unionCoords;
252   IS             perm;
253   PetscErrorCode ierr;
254 
255   PetscFunctionBegin;
256 #if 1
257   comm = PETSC_COMM_SELF;
258 #endif
259   /* create a reference element */
260   ierr = DMPlexCreateReferenceCell(comm, dim, simplex, &K);CHKERRQ(ierr);
261   ierr = DMPlexCreateLabel(K, "identity");CHKERRQ(ierr);
262   ierr = DMPlexGetLabel(K, "identity", &identity);CHKERRQ(ierr);
263   ierr = DMPlexGetChart(K, &pStart, &pEnd);CHKERRQ(ierr);
264   for (p = pStart; p < pEnd; p++) {
265     ierr = DMLabelSetValue(identity, p, p);CHKERRQ(ierr);
266   }
267   /* refine it */
268   ierr = DMRefine(K,comm,&Kref);CHKERRQ(ierr);
269 
270   /* the reference tree is the union of these two, without duplicating
271    * points that appear in both */
272   ierr = DMPlexGetLabel(Kref, "identity", &identityRef);CHKERRQ(ierr);
273   ierr = DMPlexGetChart(Kref, &pRefStart, &pRefEnd);CHKERRQ(ierr);
274   ierr = PetscSectionCreate(comm, &unionSection);CHKERRQ(ierr);
275   ierr = PetscSectionSetChart(unionSection, 0, (pEnd - pStart) + (pRefEnd - pRefStart));CHKERRQ(ierr);
276   /* count points that will go in the union */
277   for (p = pStart; p < pEnd; p++) {
278     ierr = PetscSectionSetDof(unionSection, p - pStart, 1);CHKERRQ(ierr);
279   }
280   for (p = pRefStart; p < pRefEnd; p++) {
281     PetscInt q, qSize;
282     ierr = DMLabelGetValue(identityRef, p, &q);CHKERRQ(ierr);
283     ierr = DMLabelGetStratumSize(identityRef, q, &qSize);CHKERRQ(ierr);
284     if (qSize > 1) {
285       ierr = PetscSectionSetDof(unionSection, p - pRefStart + (pEnd - pStart), 1);CHKERRQ(ierr);
286     }
287   }
288   ierr = PetscMalloc1(pEnd - pStart + pRefEnd - pRefStart,&permvals);CHKERRQ(ierr);
289   offset = 0;
290   /* stratify points in the union by topological dimension */
291   for (d = 0; d <= dim; d++) {
292     PetscInt cStart, cEnd, c;
293 
294     ierr = DMPlexGetHeightStratum(K, d, &cStart, &cEnd);CHKERRQ(ierr);
295     for (c = cStart; c < cEnd; c++) {
296       permvals[offset++] = c;
297     }
298 
299     ierr = DMPlexGetHeightStratum(Kref, d, &cStart, &cEnd);CHKERRQ(ierr);
300     for (c = cStart; c < cEnd; c++) {
301       permvals[offset++] = c + (pEnd - pStart);
302     }
303   }
304   ierr = ISCreateGeneral(comm, (pEnd - pStart) + (pRefEnd - pRefStart), permvals, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr);
305   ierr = PetscSectionSetPermutation(unionSection,perm);CHKERRQ(ierr);
306   ierr = PetscSectionSetUp(unionSection);CHKERRQ(ierr);
307   ierr = PetscSectionGetStorageSize(unionSection,&numUnionPoints);CHKERRQ(ierr);
308   ierr = PetscMalloc2(numUnionPoints,&coneSizes,dim+1,&numDimPoints);CHKERRQ(ierr);
309   /* count dimension points */
310   for (d = 0; d <= dim; d++) {
311     PetscInt cStart, cOff, cOff2;
312     ierr = DMPlexGetHeightStratum(K,d,&cStart,NULL);CHKERRQ(ierr);
313     ierr = PetscSectionGetOffset(unionSection,cStart-pStart,&cOff);CHKERRQ(ierr);
314     if (d < dim) {
315       ierr = DMPlexGetHeightStratum(K,d+1,&cStart,NULL);CHKERRQ(ierr);
316       ierr = PetscSectionGetOffset(unionSection,cStart-pStart,&cOff2);CHKERRQ(ierr);
317     }
318     else {
319       cOff2 = numUnionPoints;
320     }
321     numDimPoints[dim - d] = cOff2 - cOff;
322   }
323   ierr = PetscSectionCreate(comm, &unionConeSection);CHKERRQ(ierr);
324   ierr = PetscSectionSetChart(unionConeSection, 0, numUnionPoints);CHKERRQ(ierr);
325   /* count the cones in the union */
326   for (p = pStart; p < pEnd; p++) {
327     PetscInt dof, uOff;
328 
329     ierr = DMPlexGetConeSize(K, p, &dof);CHKERRQ(ierr);
330     ierr = PetscSectionGetOffset(unionSection, p - pStart,&uOff);CHKERRQ(ierr);
331     ierr = PetscSectionSetDof(unionConeSection, uOff, dof);CHKERRQ(ierr);
332     coneSizes[uOff] = dof;
333   }
334   for (p = pRefStart; p < pRefEnd; p++) {
335     PetscInt dof, uDof, uOff;
336 
337     ierr = DMPlexGetConeSize(Kref, p, &dof);CHKERRQ(ierr);
338     ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr);
339     ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr);
340     if (uDof) {
341       ierr = PetscSectionSetDof(unionConeSection, uOff, dof);CHKERRQ(ierr);
342       coneSizes[uOff] = dof;
343     }
344   }
345   ierr = PetscSectionSetUp(unionConeSection);CHKERRQ(ierr);
346   ierr = PetscSectionGetStorageSize(unionConeSection,&numCones);CHKERRQ(ierr);
347   ierr = PetscMalloc2(numCones,&unionCones,numCones,&unionOrientations);CHKERRQ(ierr);
348   /* write the cones in the union */
349   for (p = pStart; p < pEnd; p++) {
350     PetscInt dof, uOff, c, cOff;
351     const PetscInt *cone, *orientation;
352 
353     ierr = DMPlexGetConeSize(K, p, &dof);CHKERRQ(ierr);
354     ierr = DMPlexGetCone(K, p, &cone);CHKERRQ(ierr);
355     ierr = DMPlexGetConeOrientation(K, p, &orientation);CHKERRQ(ierr);
356     ierr = PetscSectionGetOffset(unionSection, p - pStart,&uOff);CHKERRQ(ierr);
357     ierr = PetscSectionGetOffset(unionConeSection,uOff,&cOff);CHKERRQ(ierr);
358     for (c = 0; c < dof; c++) {
359       PetscInt e, eOff;
360       e                           = cone[c];
361       ierr                        = PetscSectionGetOffset(unionSection, e - pStart, &eOff);CHKERRQ(ierr);
362       unionCones[cOff + c]        = eOff;
363       unionOrientations[cOff + c] = orientation[c];
364     }
365   }
366   for (p = pRefStart; p < pRefEnd; p++) {
367     PetscInt dof, uDof, uOff, c, cOff;
368     const PetscInt *cone, *orientation;
369 
370     ierr = DMPlexGetConeSize(Kref, p, &dof);CHKERRQ(ierr);
371     ierr = DMPlexGetCone(Kref, p, &cone);CHKERRQ(ierr);
372     ierr = DMPlexGetConeOrientation(Kref, p, &orientation);CHKERRQ(ierr);
373     ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr);
374     ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr);
375     if (uDof) {
376       ierr = PetscSectionGetOffset(unionConeSection,uOff,&cOff);CHKERRQ(ierr);
377       for (c = 0; c < dof; c++) {
378         PetscInt e, eOff, eDof;
379 
380         e    = cone[c];
381         ierr = PetscSectionGetDof(unionSection, e - pRefStart + (pEnd - pStart),&eDof);CHKERRQ(ierr);
382         if (eDof) {
383           ierr = PetscSectionGetOffset(unionSection, e - pRefStart + (pEnd - pStart), &eOff);CHKERRQ(ierr);
384         }
385         else {
386           ierr = DMLabelGetValue(identityRef, e, &e);CHKERRQ(ierr);
387           ierr = PetscSectionGetOffset(unionSection, e - pStart, &eOff);CHKERRQ(ierr);
388         }
389         unionCones[cOff + c]        = eOff;
390         unionOrientations[cOff + c] = orientation[c];
391       }
392     }
393   }
394   /* get the coordinates */
395   {
396     PetscInt vStart, vEnd, vRefStart, vRefEnd, v, vDof, vOff;
397     PetscSection KcoordsSec, KrefCoordsSec;
398     Vec      KcoordsVec, KrefCoordsVec;
399     PetscScalar *Kcoords;
400 
401     DMGetCoordinateSection(K, &KcoordsSec);CHKERRQ(ierr);
402     DMGetCoordinatesLocal(K, &KcoordsVec);CHKERRQ(ierr);
403     DMGetCoordinateSection(Kref, &KrefCoordsSec);CHKERRQ(ierr);
404     DMGetCoordinatesLocal(Kref, &KrefCoordsVec);CHKERRQ(ierr);
405 
406     numVerts = numDimPoints[0];
407     ierr     = PetscMalloc1(numVerts * dim,&unionCoords);CHKERRQ(ierr);
408     ierr     = DMPlexGetDepthStratum(K,0,&vStart,&vEnd);CHKERRQ(ierr);
409 
410     offset = 0;
411     for (v = vStart; v < vEnd; v++) {
412       ierr = PetscSectionGetOffset(unionSection,v - pStart,&vOff);CHKERRQ(ierr);
413       ierr = VecGetValuesSection(KcoordsVec, KcoordsSec, v, &Kcoords);CHKERRQ(ierr);
414       for (d = 0; d < dim; d++) {
415         unionCoords[offset * dim + d] = Kcoords[d];
416       }
417       offset++;
418     }
419     ierr = DMPlexGetDepthStratum(Kref,0,&vRefStart,&vRefEnd);CHKERRQ(ierr);
420     for (v = vRefStart; v < vRefEnd; v++) {
421       ierr = PetscSectionGetDof(unionSection,v - pRefStart + (pEnd - pStart),&vDof);CHKERRQ(ierr);
422       ierr = PetscSectionGetOffset(unionSection,v - pRefStart + (pEnd - pStart),&vOff);CHKERRQ(ierr);
423       ierr = VecGetValuesSection(KrefCoordsVec, KrefCoordsSec, v, &Kcoords);CHKERRQ(ierr);
424       if (vDof) {
425         for (d = 0; d < dim; d++) {
426           unionCoords[offset * dim + d] = Kcoords[d];
427         }
428         offset++;
429       }
430     }
431   }
432   ierr = DMCreate(comm,ref);CHKERRQ(ierr);
433   ierr = DMSetType(*ref,DMPLEX);CHKERRQ(ierr);
434   ierr = DMSetDimension(*ref,dim);CHKERRQ(ierr);
435   ierr = DMPlexCreateFromDAG(*ref,dim,numDimPoints,coneSizes,unionCones,unionOrientations,unionCoords);CHKERRQ(ierr);
436   /* set the tree */
437   ierr = PetscSectionCreate(comm,&parentSection);CHKERRQ(ierr);
438   ierr = PetscSectionSetChart(parentSection,0,numUnionPoints);CHKERRQ(ierr);
439   for (p = pRefStart; p < pRefEnd; p++) {
440     PetscInt uDof, uOff;
441 
442     ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr);
443     ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr);
444     if (uDof) {
445       PetscSectionSetDof(parentSection,uOff,1);CHKERRQ(ierr);
446     }
447   }
448   ierr = PetscSectionSetUp(parentSection);CHKERRQ(ierr);
449   ierr = PetscSectionGetStorageSize(parentSection,&parentSize);CHKERRQ(ierr);
450   ierr = PetscMalloc2(parentSize,&parents,parentSize,&childIDs);CHKERRQ(ierr);
451   for (p = pRefStart; p < pRefEnd; p++) {
452     PetscInt uDof, uOff;
453 
454     ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr);
455     ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr);
456     if (uDof) {
457       PetscInt pOff, parent, parentU;
458       PetscSectionGetOffset(parentSection,uOff,&pOff);CHKERRQ(ierr);
459       DMLabelGetValue(identityRef,p,&parent);CHKERRQ(ierr);
460       ierr = PetscSectionGetOffset(unionSection, parent - pStart,&parentU);CHKERRQ(ierr);
461       parents[pOff] = parentU;
462       childIDs[pOff] = uOff;
463     }
464   }
465   ierr = DMPlexSetTree_Internal(*ref,parentSection,parents,childIDs,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr);
466   mesh = (DM_Plex *) (*ref)->data;
467   mesh->getchildsymmetry = DMPlexReferenceTreeGetChildSymmetry_Default;
468   ierr = PetscSectionDestroy(&parentSection);CHKERRQ(ierr);
469   ierr = PetscFree2(parents,childIDs);CHKERRQ(ierr);
470 
471   /* clean up */
472   ierr = PetscSectionDestroy(&unionSection);CHKERRQ(ierr);
473   ierr = PetscSectionDestroy(&unionConeSection);CHKERRQ(ierr);
474   ierr = ISDestroy(&perm);CHKERRQ(ierr);
475   ierr = PetscFree(unionCoords);CHKERRQ(ierr);
476   ierr = PetscFree2(unionCones,unionOrientations);CHKERRQ(ierr);
477   ierr = PetscFree2(coneSizes,numDimPoints);CHKERRQ(ierr);
478   ierr = DMDestroy(&K);CHKERRQ(ierr);
479   ierr = DMDestroy(&Kref);CHKERRQ(ierr);
480   PetscFunctionReturn(0);
481 }
482 
483 #undef __FUNCT__
484 #define __FUNCT__ "DMPlexTreeSymmetrize"
485 static PetscErrorCode DMPlexTreeSymmetrize(DM dm)
486 {
487   DM_Plex        *mesh = (DM_Plex *)dm->data;
488   PetscSection   childSec, pSec;
489   PetscInt       p, pSize, cSize, parMax = PETSC_MIN_INT, parMin = PETSC_MAX_INT;
490   PetscInt       *offsets, *children, pStart, pEnd;
491   PetscErrorCode ierr;
492 
493   PetscFunctionBegin;
494   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
495   ierr = PetscSectionDestroy(&mesh->childSection);CHKERRQ(ierr);
496   ierr = PetscFree(mesh->children);CHKERRQ(ierr);
497   pSec = mesh->parentSection;
498   if (!pSec) PetscFunctionReturn(0);
499   ierr = PetscSectionGetStorageSize(pSec,&pSize);CHKERRQ(ierr);
500   for (p = 0; p < pSize; p++) {
501     PetscInt par = mesh->parents[p];
502 
503     parMax = PetscMax(parMax,par+1);
504     parMin = PetscMin(parMin,par);
505   }
506   if (parMin > parMax) {
507     parMin = -1;
508     parMax = -1;
509   }
510   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)pSec),&childSec);CHKERRQ(ierr);
511   ierr = PetscSectionSetChart(childSec,parMin,parMax);CHKERRQ(ierr);
512   for (p = 0; p < pSize; p++) {
513     PetscInt par = mesh->parents[p];
514 
515     ierr = PetscSectionAddDof(childSec,par,1);CHKERRQ(ierr);
516   }
517   ierr = PetscSectionSetUp(childSec);CHKERRQ(ierr);
518   ierr = PetscSectionGetStorageSize(childSec,&cSize);CHKERRQ(ierr);
519   ierr = PetscMalloc1(cSize,&children);CHKERRQ(ierr);
520   ierr = PetscCalloc1(parMax-parMin,&offsets);CHKERRQ(ierr);
521   ierr = PetscSectionGetChart(pSec,&pStart,&pEnd);CHKERRQ(ierr);
522   for (p = pStart; p < pEnd; p++) {
523     PetscInt dof, off, i;
524 
525     ierr = PetscSectionGetDof(pSec,p,&dof);CHKERRQ(ierr);
526     ierr = PetscSectionGetOffset(pSec,p,&off);CHKERRQ(ierr);
527     for (i = 0; i < dof; i++) {
528       PetscInt par = mesh->parents[off + i], cOff;
529 
530       ierr = PetscSectionGetOffset(childSec,par,&cOff);CHKERRQ(ierr);
531       children[cOff + offsets[par-parMin]++] = p;
532     }
533   }
534   mesh->childSection = childSec;
535   mesh->children = children;
536   ierr = PetscFree(offsets);CHKERRQ(ierr);
537   PetscFunctionReturn(0);
538 }
539 
540 #undef __FUNCT__
541 #define __FUNCT__ "AnchorsFlatten"
542 static PetscErrorCode AnchorsFlatten (PetscSection section, IS is, PetscSection *sectionNew, IS *isNew)
543 {
544   PetscInt       pStart, pEnd, size, sizeNew, i, p, *valsNew = NULL;
545   const PetscInt *vals;
546   PetscSection   secNew;
547   PetscBool      anyNew, globalAnyNew;
548   PetscBool      compress;
549   PetscErrorCode ierr;
550 
551   PetscFunctionBegin;
552   ierr = PetscSectionGetChart(section,&pStart,&pEnd);CHKERRQ(ierr);
553   ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr);
554   ierr = ISGetIndices(is,&vals);CHKERRQ(ierr);
555   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)section),&secNew);CHKERRQ(ierr);
556   ierr = PetscSectionSetChart(secNew,pStart,pEnd);CHKERRQ(ierr);
557   for (i = 0; i < size; i++) {
558     PetscInt dof;
559 
560     p = vals[i];
561     if (p < pStart || p >= pEnd) continue;
562     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
563     if (dof) break;
564   }
565   if (i == size) {
566     ierr     = PetscSectionSetUp(secNew);CHKERRQ(ierr);
567     anyNew   = PETSC_FALSE;
568     compress = PETSC_FALSE;
569     sizeNew  = 0;
570   }
571   else {
572     anyNew = PETSC_TRUE;
573     for (p = pStart; p < pEnd; p++) {
574       PetscInt dof, off;
575 
576       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
577       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
578       for (i = 0; i < dof; i++) {
579         PetscInt q = vals[off + i], qDof = 0;
580 
581         if (q >= pStart && q < pEnd) {
582           ierr = PetscSectionGetDof(section, q, &qDof);CHKERRQ(ierr);
583         }
584         if (qDof) {
585           ierr = PetscSectionAddDof(secNew, p, qDof);CHKERRQ(ierr);
586         }
587         else {
588           ierr = PetscSectionAddDof(secNew, p, 1);CHKERRQ(ierr);
589         }
590       }
591     }
592     ierr = PetscSectionSetUp(secNew);CHKERRQ(ierr);
593     ierr = PetscSectionGetStorageSize(secNew,&sizeNew);CHKERRQ(ierr);
594     ierr = PetscMalloc1(sizeNew,&valsNew);CHKERRQ(ierr);
595     compress = PETSC_FALSE;
596     for (p = pStart; p < pEnd; p++) {
597       PetscInt dof, off, count, offNew, dofNew;
598 
599       ierr  = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
600       ierr  = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
601       ierr  = PetscSectionGetDof(secNew, p, &dofNew);CHKERRQ(ierr);
602       ierr  = PetscSectionGetOffset(secNew, p, &offNew);CHKERRQ(ierr);
603       count = 0;
604       for (i = 0; i < dof; i++) {
605         PetscInt q = vals[off + i], qDof = 0, qOff = 0, j;
606 
607         if (q >= pStart && q < pEnd) {
608           ierr = PetscSectionGetDof(section, q, &qDof);CHKERRQ(ierr);
609           ierr = PetscSectionGetOffset(section, q, &qOff);CHKERRQ(ierr);
610         }
611         if (qDof) {
612           PetscInt oldCount = count;
613 
614           for (j = 0; j < qDof; j++) {
615             PetscInt k, r = vals[qOff + j];
616 
617             for (k = 0; k < oldCount; k++) {
618               if (valsNew[offNew + k] == r) {
619                 break;
620               }
621             }
622             if (k == oldCount) {
623               valsNew[offNew + count++] = r;
624             }
625           }
626         }
627         else {
628           PetscInt k, oldCount = count;
629 
630           for (k = 0; k < oldCount; k++) {
631             if (valsNew[offNew + k] == q) {
632               break;
633             }
634           }
635           if (k == oldCount) {
636             valsNew[offNew + count++] = q;
637           }
638         }
639       }
640       if (count < dofNew) {
641         ierr = PetscSectionSetDof(secNew, p, count);CHKERRQ(ierr);
642         compress = PETSC_TRUE;
643       }
644     }
645   }
646   ierr = ISRestoreIndices(is,&vals);CHKERRQ(ierr);
647   ierr = MPI_Allreduce(&anyNew,&globalAnyNew,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)secNew));CHKERRQ(ierr);
648   if (!globalAnyNew) {
649     ierr = PetscSectionDestroy(&secNew);CHKERRQ(ierr);
650     *sectionNew = NULL;
651     *isNew = NULL;
652   }
653   else {
654     PetscBool globalCompress;
655 
656     ierr = MPI_Allreduce(&compress,&globalCompress,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)secNew));CHKERRQ(ierr);
657     if (compress) {
658       PetscSection secComp;
659       PetscInt *valsComp = NULL;
660 
661       ierr = PetscSectionCreate(PetscObjectComm((PetscObject)section),&secComp);CHKERRQ(ierr);
662       ierr = PetscSectionSetChart(secComp,pStart,pEnd);CHKERRQ(ierr);
663       for (p = pStart; p < pEnd; p++) {
664         PetscInt dof;
665 
666         ierr = PetscSectionGetDof(secNew, p, &dof);CHKERRQ(ierr);
667         ierr = PetscSectionSetDof(secComp, p, dof);CHKERRQ(ierr);
668       }
669       ierr = PetscSectionSetUp(secComp);CHKERRQ(ierr);
670       ierr = PetscSectionGetStorageSize(secComp,&sizeNew);CHKERRQ(ierr);
671       ierr = PetscMalloc1(sizeNew,&valsComp);CHKERRQ(ierr);
672       for (p = pStart; p < pEnd; p++) {
673         PetscInt dof, off, offNew, j;
674 
675         ierr = PetscSectionGetDof(secNew, p, &dof);CHKERRQ(ierr);
676         ierr = PetscSectionGetOffset(secNew, p, &off);CHKERRQ(ierr);
677         ierr = PetscSectionGetOffset(secComp, p, &offNew);CHKERRQ(ierr);
678         for (j = 0; j < dof; j++) {
679           valsComp[offNew + j] = valsNew[off + j];
680         }
681       }
682       ierr    = PetscSectionDestroy(&secNew);CHKERRQ(ierr);
683       secNew  = secComp;
684       ierr    = PetscFree(valsNew);CHKERRQ(ierr);
685       valsNew = valsComp;
686     }
687     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),sizeNew,valsNew,PETSC_OWN_POINTER,isNew);CHKERRQ(ierr);
688   }
689   PetscFunctionReturn(0);
690 }
691 
692 #undef __FUNCT__
693 #define __FUNCT__ "DMPlexCreateAnchors_Tree"
694 static PetscErrorCode DMPlexCreateAnchors_Tree(DM dm)
695 {
696   PetscInt       p, pStart, pEnd, *anchors, size;
697   PetscInt       aMin = PETSC_MAX_INT, aMax = PETSC_MIN_INT;
698   PetscSection   aSec;
699   DMLabel        canonLabel;
700   IS             aIS;
701   PetscErrorCode ierr;
702 
703   PetscFunctionBegin;
704   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
705   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
706   ierr = DMPlexGetLabel(dm,"canonical",&canonLabel);CHKERRQ(ierr);
707   for (p = pStart; p < pEnd; p++) {
708     PetscInt parent;
709 
710     if (canonLabel) {
711       PetscInt canon;
712 
713       ierr = DMLabelGetValue(canonLabel,p,&canon);CHKERRQ(ierr);
714       if (p != canon) continue;
715     }
716     ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
717     if (parent != p) {
718       aMin = PetscMin(aMin,p);
719       aMax = PetscMax(aMax,p+1);
720     }
721   }
722   if (aMin > aMax) {
723     aMin = -1;
724     aMax = -1;
725   }
726   ierr = PetscSectionCreate(PETSC_COMM_SELF,&aSec);CHKERRQ(ierr);
727   ierr = PetscSectionSetChart(aSec,aMin,aMax);CHKERRQ(ierr);
728   for (p = aMin; p < aMax; p++) {
729     PetscInt parent, ancestor = p;
730 
731     if (canonLabel) {
732       PetscInt canon;
733 
734       ierr = DMLabelGetValue(canonLabel,p,&canon);CHKERRQ(ierr);
735       if (p != canon) continue;
736     }
737     ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
738     while (parent != ancestor) {
739       ancestor = parent;
740       ierr     = DMPlexGetTreeParent(dm,ancestor,&parent,NULL);CHKERRQ(ierr);
741     }
742     if (ancestor != p) {
743       PetscInt closureSize, *closure = NULL;
744 
745       ierr = DMPlexGetTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
746       ierr = PetscSectionSetDof(aSec,p,closureSize);CHKERRQ(ierr);
747       ierr = DMPlexRestoreTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
748     }
749   }
750   ierr = PetscSectionSetUp(aSec);CHKERRQ(ierr);
751   ierr = PetscSectionGetStorageSize(aSec,&size);CHKERRQ(ierr);
752   ierr = PetscMalloc1(size,&anchors);CHKERRQ(ierr);
753   for (p = aMin; p < aMax; p++) {
754     PetscInt parent, ancestor = p;
755 
756     if (canonLabel) {
757       PetscInt canon;
758 
759       ierr = DMLabelGetValue(canonLabel,p,&canon);CHKERRQ(ierr);
760       if (p != canon) continue;
761     }
762     ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
763     while (parent != ancestor) {
764       ancestor = parent;
765       ierr     = DMPlexGetTreeParent(dm,ancestor,&parent,NULL);CHKERRQ(ierr);
766     }
767     if (ancestor != p) {
768       PetscInt j, closureSize, *closure = NULL, aOff;
769 
770       ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr);
771 
772       ierr = DMPlexGetTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
773       for (j = 0; j < closureSize; j++) {
774         anchors[aOff + j] = closure[2*j];
775       }
776       ierr = DMPlexRestoreTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
777     }
778   }
779   ierr = ISCreateGeneral(PETSC_COMM_SELF,size,anchors,PETSC_OWN_POINTER,&aIS);CHKERRQ(ierr);
780   {
781     PetscSection aSecNew = aSec;
782     IS           aISNew  = aIS;
783 
784     ierr = PetscObjectReference((PetscObject)aSec);CHKERRQ(ierr);
785     ierr = PetscObjectReference((PetscObject)aIS);CHKERRQ(ierr);
786     while (aSecNew) {
787       ierr    = PetscSectionDestroy(&aSec);CHKERRQ(ierr);
788       ierr    = ISDestroy(&aIS);CHKERRQ(ierr);
789       aSec    = aSecNew;
790       aIS     = aISNew;
791       aSecNew = NULL;
792       aISNew  = NULL;
793       ierr    = AnchorsFlatten(aSec,aIS,&aSecNew,&aISNew);CHKERRQ(ierr);
794     }
795   }
796   ierr = DMPlexSetAnchors(dm,aSec,aIS);CHKERRQ(ierr);
797   ierr = PetscSectionDestroy(&aSec);CHKERRQ(ierr);
798   ierr = ISDestroy(&aIS);CHKERRQ(ierr);
799   PetscFunctionReturn(0);
800 }
801 
802 #undef __FUNCT__
803 #define __FUNCT__ "DMPlexTreeExchangeSupports"
804 static PetscErrorCode DMPlexTreeExchangeSupports(DM dm)
805 {
806   DM_Plex *mesh = (DM_Plex *)dm->data;
807   PetscSection newSupportSection;
808   PetscInt newSize, *newSupports, pStart, pEnd, p, d, depth;
809   PetscInt *offsets;
810   PetscErrorCode ierr;
811 
812   PetscFunctionBegin;
813   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
814   /* symmetrize the hierarchy */
815   ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr);
816   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)(mesh->supportSection)),&newSupportSection);CHKERRQ(ierr);
817   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
818   ierr = PetscSectionSetChart(newSupportSection,pStart,pEnd);CHKERRQ(ierr);
819   ierr = PetscCalloc1(pEnd,&offsets);CHKERRQ(ierr);
820   /* if a point is in the support of q, it should be in the support of
821    * parent(q) */
822   for (d = 0; d <= depth; d++) {
823     ierr = DMPlexGetHeightStratum(dm,d,&pStart,&pEnd);CHKERRQ(ierr);
824     for (p = pStart; p < pEnd; ++p) {
825       PetscInt dof, q, qdof, parent;
826 
827       ierr = PetscSectionGetDof(mesh->supportSection, p, &dof);CHKERRQ(ierr);
828       ierr = PetscSectionAddDof(newSupportSection, p, dof);CHKERRQ(ierr);
829       q    = p;
830       ierr = DMPlexGetTreeParent(dm,q,&parent,NULL);CHKERRQ(ierr);
831       while (parent != q && parent >= pStart && parent < pEnd) {
832         q = parent;
833 
834         ierr = PetscSectionGetDof(mesh->supportSection, q, &qdof);CHKERRQ(ierr);
835         ierr = PetscSectionAddDof(newSupportSection,p,qdof);CHKERRQ(ierr);
836         ierr = PetscSectionAddDof(newSupportSection,q,dof);CHKERRQ(ierr);
837         ierr = DMPlexGetTreeParent(dm,q,&parent,NULL);CHKERRQ(ierr);
838       }
839     }
840   }
841   ierr = PetscSectionSetUp(newSupportSection);CHKERRQ(ierr);
842   ierr = PetscSectionGetStorageSize(newSupportSection,&newSize);CHKERRQ(ierr);
843   ierr = PetscMalloc1(newSize,&newSupports);CHKERRQ(ierr);
844   for (d = 0; d <= depth; d++) {
845     ierr = DMPlexGetHeightStratum(dm,d,&pStart,&pEnd);CHKERRQ(ierr);
846     for (p = pStart; p < pEnd; p++) {
847       PetscInt dof, off, q, qdof, qoff, newDof, newOff, newqOff, i, parent;
848 
849       ierr = PetscSectionGetDof(mesh->supportSection, p, &dof);CHKERRQ(ierr);
850       ierr = PetscSectionGetOffset(mesh->supportSection, p, &off);CHKERRQ(ierr);
851       ierr = PetscSectionGetDof(newSupportSection, p, &newDof);CHKERRQ(ierr);
852       ierr = PetscSectionGetOffset(newSupportSection, p, &newOff);CHKERRQ(ierr);
853       for (i = 0; i < dof; i++) {
854         newSupports[newOff+offsets[p]++] = mesh->supports[off + i];
855       }
856       mesh->maxSupportSize = PetscMax(mesh->maxSupportSize,newDof);
857 
858       q    = p;
859       ierr = DMPlexGetTreeParent(dm,q,&parent,NULL);CHKERRQ(ierr);
860       while (parent != q && parent >= pStart && parent < pEnd) {
861         q = parent;
862         ierr = PetscSectionGetDof(mesh->supportSection, q, &qdof);CHKERRQ(ierr);
863         ierr = PetscSectionGetOffset(mesh->supportSection, q, &qoff);CHKERRQ(ierr);
864         ierr = PetscSectionGetOffset(newSupportSection, q, &newqOff);CHKERRQ(ierr);
865         for (i = 0; i < qdof; i++) {
866           newSupports[newOff+offsets[p]++] = mesh->supports[qoff + i];
867         }
868         for (i = 0; i < dof; i++) {
869           newSupports[newqOff+offsets[q]++] = mesh->supports[off + i];
870         }
871         ierr = DMPlexGetTreeParent(dm,q,&parent,NULL);CHKERRQ(ierr);
872       }
873     }
874   }
875   ierr = PetscSectionDestroy(&mesh->supportSection);CHKERRQ(ierr);
876   mesh->supportSection = newSupportSection;
877   ierr = PetscFree(mesh->supports);CHKERRQ(ierr);
878   mesh->supports = newSupports;
879   ierr = PetscFree(offsets);CHKERRQ(ierr);
880 
881   PetscFunctionReturn(0);
882 }
883 
884 static PetscErrorCode DMPlexComputeAnchorMatrix_Tree_Direct(DM,PetscSection,PetscSection,Mat);
885 static PetscErrorCode DMPlexComputeAnchorMatrix_Tree_FromReference(DM,PetscSection,PetscSection,Mat);
886 
887 #undef __FUNCT__
888 #define __FUNCT__ "DMPlexSetTree_Internal"
889 static PetscErrorCode DMPlexSetTree_Internal(DM dm, PetscSection parentSection, PetscInt *parents, PetscInt *childIDs, PetscBool computeCanonical, PetscBool exchangeSupports)
890 {
891   DM_Plex       *mesh = (DM_Plex *)dm->data;
892   DM             refTree;
893   PetscInt       size;
894   PetscErrorCode ierr;
895 
896   PetscFunctionBegin;
897   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
898   PetscValidHeaderSpecific(parentSection, PETSC_SECTION_CLASSID, 2);
899   ierr = PetscObjectReference((PetscObject)parentSection);CHKERRQ(ierr);
900   ierr = PetscSectionDestroy(&mesh->parentSection);CHKERRQ(ierr);
901   mesh->parentSection = parentSection;
902   ierr = PetscSectionGetStorageSize(parentSection,&size);CHKERRQ(ierr);
903   if (parents != mesh->parents) {
904     ierr = PetscFree(mesh->parents);CHKERRQ(ierr);
905     ierr = PetscMalloc1(size,&mesh->parents);CHKERRQ(ierr);
906     ierr = PetscMemcpy(mesh->parents, parents, size * sizeof(*parents));CHKERRQ(ierr);
907   }
908   if (childIDs != mesh->childIDs) {
909     ierr = PetscFree(mesh->childIDs);CHKERRQ(ierr);
910     ierr = PetscMalloc1(size,&mesh->childIDs);CHKERRQ(ierr);
911     ierr = PetscMemcpy(mesh->childIDs, childIDs, size * sizeof(*childIDs));CHKERRQ(ierr);
912   }
913   ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr);
914   if (refTree) {
915     DMLabel canonLabel;
916 
917     ierr = DMPlexGetLabel(refTree,"canonical",&canonLabel);CHKERRQ(ierr);
918     if (canonLabel) {
919       PetscInt i;
920 
921       for (i = 0; i < size; i++) {
922         PetscInt canon;
923         ierr = DMLabelGetValue(canonLabel, mesh->childIDs[i], &canon);CHKERRQ(ierr);
924         if (canon >= 0) {
925           mesh->childIDs[i] = canon;
926         }
927       }
928     }
929     mesh->computeanchormatrix = DMPlexComputeAnchorMatrix_Tree_FromReference;
930   }
931   else {
932     mesh->computeanchormatrix = DMPlexComputeAnchorMatrix_Tree_Direct;
933   }
934   ierr = DMPlexTreeSymmetrize(dm);CHKERRQ(ierr);
935   if (computeCanonical) {
936     PetscInt d, dim;
937 
938     /* add the canonical label */
939     ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
940     ierr = DMPlexCreateLabel(dm,"canonical");CHKERRQ(ierr);
941     for (d = 0; d <= dim; d++) {
942       PetscInt p, dStart, dEnd, canon = -1, cNumChildren;
943       const PetscInt *cChildren;
944 
945       ierr = DMPlexGetDepthStratum(dm,d,&dStart,&dEnd);CHKERRQ(ierr);
946       for (p = dStart; p < dEnd; p++) {
947         ierr = DMPlexGetTreeChildren(dm,p,&cNumChildren,&cChildren);CHKERRQ(ierr);
948         if (cNumChildren) {
949           canon = p;
950           break;
951         }
952       }
953       if (canon == -1) continue;
954       for (p = dStart; p < dEnd; p++) {
955         PetscInt numChildren, i;
956         const PetscInt *children;
957 
958         ierr = DMPlexGetTreeChildren(dm,p,&numChildren,&children);CHKERRQ(ierr);
959         if (numChildren) {
960           if (numChildren != cNumChildren) SETERRQ2(PetscObjectComm((PetscObject)dm),PETSC_ERR_PLIB,"All parent points in a stratum should have the same number of children: %d != %d", numChildren, cNumChildren);
961           ierr = DMPlexSetLabelValue(dm,"canonical",p,canon);CHKERRQ(ierr);
962           for (i = 0; i < numChildren; i++) {
963             ierr = DMPlexSetLabelValue(dm,"canonical",children[i],cChildren[i]);CHKERRQ(ierr);
964           }
965         }
966       }
967     }
968   }
969   if (exchangeSupports) {
970     ierr = DMPlexTreeExchangeSupports(dm);CHKERRQ(ierr);
971   }
972   mesh->createanchors = DMPlexCreateAnchors_Tree;
973   /* reset anchors */
974   ierr = DMPlexSetAnchors(dm,NULL,NULL);CHKERRQ(ierr);
975   PetscFunctionReturn(0);
976 }
977 
978 #undef __FUNCT__
979 #define __FUNCT__ "DMPlexSetTree"
980 /*@
981   DMPlexSetTree - set the tree that describes the hierarchy of non-conforming mesh points.  This routine also creates
982   the point-to-point constraints determined by the tree: a point is constained to the points in the closure of its
983   tree root.
984 
985   Collective on dm
986 
987   Input Parameters:
988 + dm - the DMPlex object
989 . parentSection - a section describing the tree: a point has a parent if it has 1 dof in the section; the section
990                   offset indexes the parent and childID list; the reference count of parentSection is incremented
991 . parents - a list of the point parents; copied, can be destroyed
992 - childIDs - identifies the relationship of the child point to the parent point; if there is a reference tree, then
993              the child corresponds to the point in the reference tree with index childIDs; copied, can be destroyed
994 
995   Level: intermediate
996 
997 .seealso: DMPlexGetTree(), DMPlexSetReferenceTree(), DMPlexSetAnchors(), DMPlexGetTreeParent(), DMPlexGetTreeChildren()
998 @*/
999 PetscErrorCode DMPlexSetTree(DM dm, PetscSection parentSection, PetscInt parents[], PetscInt childIDs[])
1000 {
1001   PetscErrorCode ierr;
1002 
1003   PetscFunctionBegin;
1004   ierr = DMPlexSetTree_Internal(dm,parentSection,parents,childIDs,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr);
1005   PetscFunctionReturn(0);
1006 }
1007 
1008 #undef __FUNCT__
1009 #define __FUNCT__ "DMPlexGetTree"
1010 /*@
1011   DMPlexGetTree - get the tree that describes the hierarchy of non-conforming mesh points.
1012   Collective on dm
1013 
1014   Input Parameters:
1015 . dm - the DMPlex object
1016 
1017   Output Parameters:
1018 + parentSection - a section describing the tree: a point has a parent if it has 1 dof in the section; the section
1019                   offset indexes the parent and childID list
1020 . parents - a list of the point parents
1021 . childIDs - identifies the relationship of the child point to the parent point; if there is a reference tree, then
1022              the child corresponds to the point in the reference tree with index childID
1023 . childSection - the inverse of the parent section
1024 - children - a list of the point children
1025 
1026   Level: intermediate
1027 
1028 .seealso: DMPlexSetTree(), DMPlexSetReferenceTree(), DMPlexSetAnchors(), DMPlexGetTreeParent(), DMPlexGetTreeChildren()
1029 @*/
1030 PetscErrorCode DMPlexGetTree(DM dm, PetscSection *parentSection, PetscInt *parents[], PetscInt *childIDs[], PetscSection *childSection, PetscInt *children[])
1031 {
1032   DM_Plex        *mesh = (DM_Plex *)dm->data;
1033 
1034   PetscFunctionBegin;
1035   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1036   if (parentSection) *parentSection = mesh->parentSection;
1037   if (parents)       *parents       = mesh->parents;
1038   if (childIDs)      *childIDs      = mesh->childIDs;
1039   if (childSection)  *childSection  = mesh->childSection;
1040   if (children)      *children      = mesh->children;
1041   PetscFunctionReturn(0);
1042 }
1043 
1044 #undef __FUNCT__
1045 #define __FUNCT__ "DMPlexGetTreeParent"
1046 /*@
1047   DMPlexGetTreeParent - get the parent of a point in the tree describing the point hierarchy (not the Sieve DAG)
1048 
1049   Input Parameters:
1050 + dm - the DMPlex object
1051 - point - the query point
1052 
1053   Output Parameters:
1054 + parent - if not NULL, set to the parent of the point, or the point itself if the point does not have a parent
1055 - childID - if not NULL, set to the child ID of the point with respect to its parent, or 0 if the point
1056             does not have a parent
1057 
1058   Level: intermediate
1059 
1060 .seealso: DMPlexSetTree(), DMPlexGetTree(), DMPlexGetTreeChildren()
1061 @*/
1062 PetscErrorCode DMPlexGetTreeParent(DM dm, PetscInt point, PetscInt *parent, PetscInt *childID)
1063 {
1064   DM_Plex       *mesh = (DM_Plex *)dm->data;
1065   PetscSection   pSec;
1066   PetscErrorCode ierr;
1067 
1068   PetscFunctionBegin;
1069   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1070   pSec = mesh->parentSection;
1071   if (pSec && point >= pSec->pStart && point < pSec->pEnd) {
1072     PetscInt dof;
1073 
1074     ierr = PetscSectionGetDof (pSec, point, &dof);CHKERRQ(ierr);
1075     if (dof) {
1076       PetscInt off;
1077 
1078       ierr = PetscSectionGetOffset (pSec, point, &off);CHKERRQ(ierr);
1079       if (parent)  *parent = mesh->parents[off];
1080       if (childID) *childID = mesh->childIDs[off];
1081       PetscFunctionReturn(0);
1082     }
1083   }
1084   if (parent) {
1085     *parent = point;
1086   }
1087   if (childID) {
1088     *childID = 0;
1089   }
1090   PetscFunctionReturn(0);
1091 }
1092 
1093 #undef __FUNCT__
1094 #define __FUNCT__ "DMPlexGetTreeChildren"
1095 /*@C
1096   DMPlexGetTreeChildren - get the children of a point in the tree describing the point hierarchy (not the Sieve DAG)
1097 
1098   Input Parameters:
1099 + dm - the DMPlex object
1100 - point - the query point
1101 
1102   Output Parameters:
1103 + numChildren - if not NULL, set to the number of children
1104 - children - if not NULL, set to a list children, or set to NULL if the point has no children
1105 
1106   Level: intermediate
1107 
1108   Fortran Notes:
1109   Since it returns an array, this routine is only available in Fortran 90, and you must
1110   include petsc.h90 in your code.
1111 
1112 .seealso: DMPlexSetTree(), DMPlexGetTree(), DMPlexGetTreeParent()
1113 @*/
1114 PetscErrorCode DMPlexGetTreeChildren(DM dm, PetscInt point, PetscInt *numChildren, const PetscInt *children[])
1115 {
1116   DM_Plex       *mesh = (DM_Plex *)dm->data;
1117   PetscSection   childSec;
1118   PetscInt       dof = 0;
1119   PetscErrorCode ierr;
1120 
1121   PetscFunctionBegin;
1122   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1123   childSec = mesh->childSection;
1124   if (childSec && point >= childSec->pStart && point < childSec->pEnd) {
1125     ierr = PetscSectionGetDof (childSec, point, &dof);CHKERRQ(ierr);
1126   }
1127   if (numChildren) *numChildren = dof;
1128   if (children) {
1129     if (dof) {
1130       PetscInt off;
1131 
1132       ierr = PetscSectionGetOffset (childSec, point, &off);CHKERRQ(ierr);
1133       *children = &mesh->children[off];
1134     }
1135     else {
1136       *children = NULL;
1137     }
1138   }
1139   PetscFunctionReturn(0);
1140 }
1141 
1142 #undef __FUNCT__
1143 #define __FUNCT__ "DMPlexComputeAnchorMatrix_Tree_Direct"
1144 static PetscErrorCode DMPlexComputeAnchorMatrix_Tree_Direct(DM dm, PetscSection section, PetscSection cSec, Mat cMat)
1145 {
1146   PetscDS        ds;
1147   PetscInt       spdim;
1148   PetscInt       numFields, f, c, cStart, cEnd, pStart, pEnd, conStart, conEnd;
1149   const PetscInt *anchors;
1150   PetscSection   aSec;
1151   PetscReal      *v0, *v0parent, *vtmp, *J, *Jparent, *invJparent, detJ, detJparent;
1152   IS             aIS;
1153   PetscErrorCode ierr;
1154 
1155   PetscFunctionBegin;
1156   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
1157   ierr = DMGetDS(dm,&ds);CHKERRQ(ierr);
1158   ierr = PetscDSGetNumFields(ds,&numFields);CHKERRQ(ierr);
1159   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
1160   ierr = DMPlexGetAnchors(dm,&aSec,&aIS);CHKERRQ(ierr);
1161   ierr = ISGetIndices(aIS,&anchors);CHKERRQ(ierr);
1162   ierr = PetscSectionGetChart(cSec,&conStart,&conEnd);CHKERRQ(ierr);
1163   ierr = DMGetDimension(dm,&spdim);CHKERRQ(ierr);
1164   ierr = PetscMalloc6(spdim,&v0,spdim,&v0parent,spdim,&vtmp,spdim*spdim,&J,spdim*spdim,&Jparent,spdim*spdim,&invJparent);CHKERRQ(ierr);
1165 
1166   for (f = 0; f < numFields; f++) {
1167     PetscObject disc;
1168     PetscFE fe = NULL;
1169     PetscFV fv = NULL;
1170     PetscClassId id;
1171     PetscDualSpace space;
1172     PetscInt i, j, k, nPoints, offset;
1173     PetscInt fSize, fComp;
1174     PetscReal *B = NULL;
1175     PetscReal *weights, *pointsRef, *pointsReal;
1176     Mat Amat, Bmat, Xmat;
1177 
1178     ierr = PetscDSGetDiscretization(ds,f,&disc);CHKERRQ(ierr);
1179     ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
1180     if (id == PETSCFE_CLASSID) {
1181       fe = (PetscFE) disc;
1182       ierr = PetscFEGetDualSpace(fe,&space);CHKERRQ(ierr);
1183       ierr = PetscDualSpaceGetDimension(space,&fSize);CHKERRQ(ierr);
1184       ierr = PetscFEGetNumComponents(fe,&fComp);CHKERRQ(ierr);
1185     }
1186     else if (id == PETSCFV_CLASSID) {
1187       fv = (PetscFV) disc;
1188       ierr = PetscFVGetDualSpace(fv,&space);CHKERRQ(ierr);
1189       ierr = PetscDualSpaceGetDimension(space,&fSize);CHKERRQ(ierr);
1190       ierr = PetscFVGetNumComponents(fv,&fComp);CHKERRQ(ierr);
1191     }
1192     else SETERRQ1(PetscObjectComm(disc),PETSC_ERR_ARG_UNKNOWN_TYPE, "PetscDS discretization id %d not recognized.", id);
1193 
1194     ierr = MatCreate(PETSC_COMM_SELF,&Amat);CHKERRQ(ierr);
1195     ierr = MatSetSizes(Amat,fSize,fSize,fSize,fSize);CHKERRQ(ierr);
1196     ierr = MatSetType(Amat,MATSEQDENSE);CHKERRQ(ierr);
1197     ierr = MatSetUp(Amat);CHKERRQ(ierr);
1198     ierr = MatDuplicate(Amat,MAT_DO_NOT_COPY_VALUES,&Bmat);CHKERRQ(ierr);
1199     ierr = MatDuplicate(Amat,MAT_DO_NOT_COPY_VALUES,&Xmat);CHKERRQ(ierr);
1200     nPoints = 0;
1201     for (i = 0; i < fSize; i++) {
1202       PetscInt        qPoints;
1203       PetscQuadrature quad;
1204 
1205       ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr);
1206       ierr = PetscQuadratureGetData(quad,NULL,&qPoints,NULL,NULL);CHKERRQ(ierr);
1207       nPoints += qPoints;
1208     }
1209     ierr = PetscMalloc3(nPoints,&weights,spdim*nPoints,&pointsRef,spdim*nPoints,&pointsReal);CHKERRQ(ierr);
1210     offset = 0;
1211     for (i = 0; i < fSize; i++) {
1212       PetscInt        qPoints;
1213       const PetscReal    *p, *w;
1214       PetscQuadrature quad;
1215 
1216       ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr);
1217       ierr = PetscQuadratureGetData(quad,NULL,&qPoints,&p,&w);CHKERRQ(ierr);
1218       ierr = PetscMemcpy(weights+offset,w,qPoints*sizeof(*w));CHKERRQ(ierr);
1219       ierr = PetscMemcpy(pointsRef+spdim*offset,p,spdim*qPoints*sizeof(*p));CHKERRQ(ierr);
1220       offset += qPoints;
1221     }
1222     if (id == PETSCFE_CLASSID) {
1223       ierr = PetscFEGetTabulation(fe,nPoints,pointsRef,&B,NULL,NULL);CHKERRQ(ierr);
1224     }
1225     else {
1226       ierr = PetscFVGetTabulation(fv,nPoints,pointsRef,&B,NULL,NULL);CHKERRQ(ierr);
1227     }
1228     offset = 0;
1229     for (i = 0; i < fSize; i++) {
1230       PetscInt        qPoints;
1231       PetscQuadrature quad;
1232 
1233       ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr);
1234       ierr = PetscQuadratureGetData(quad,NULL,&qPoints,NULL,NULL);CHKERRQ(ierr);
1235       for (j = 0; j < fSize; j++) {
1236         PetscScalar val = 0.;
1237 
1238         for (k = 0; k < qPoints; k++) {
1239           val += B[((offset + k) * fSize + j) * fComp] * weights[k];
1240         }
1241         ierr = MatSetValue(Amat,i,j,val,INSERT_VALUES);CHKERRQ(ierr);
1242       }
1243       offset += qPoints;
1244     }
1245     if (id == PETSCFE_CLASSID) {
1246       ierr = PetscFERestoreTabulation(fe,nPoints,pointsRef,&B,NULL,NULL);CHKERRQ(ierr);
1247     }
1248     else {
1249       ierr = PetscFVRestoreTabulation(fv,nPoints,pointsRef,&B,NULL,NULL);CHKERRQ(ierr);
1250     }
1251     ierr = MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1252     ierr = MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1253     ierr = MatLUFactor(Amat,NULL,NULL,NULL);CHKERRQ(ierr);
1254     for (c = cStart; c < cEnd; c++) {
1255       PetscInt parent;
1256       PetscInt closureSize, closureSizeP, *closure = NULL, *closureP = NULL;
1257       PetscInt *childOffsets, *parentOffsets;
1258 
1259       ierr = DMPlexGetTreeParent(dm,c,&parent,NULL);CHKERRQ(ierr);
1260       if (parent == c) continue;
1261       ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
1262       for (i = 0; i < closureSize; i++) {
1263         PetscInt p = closure[2*i];
1264         PetscInt conDof;
1265 
1266         if (p < conStart || p >= conEnd) continue;
1267         if (numFields > 1) {
1268           ierr = PetscSectionGetFieldDof(cSec,p,f,&conDof);CHKERRQ(ierr);
1269         }
1270         else {
1271           ierr = PetscSectionGetDof(cSec,p,&conDof);CHKERRQ(ierr);
1272         }
1273         if (conDof) break;
1274       }
1275       if (i == closureSize) {
1276         ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
1277         continue;
1278       }
1279 
1280       ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, NULL, &detJ);CHKERRQ(ierr);
1281       ierr = DMPlexComputeCellGeometryFEM(dm, parent, NULL, v0parent, Jparent, invJparent, &detJparent);CHKERRQ(ierr);
1282       for (i = 0; i < nPoints; i++) {
1283         CoordinatesRefToReal(spdim, spdim, v0, J, &pointsRef[i*spdim],vtmp);CHKERRQ(ierr);
1284         CoordinatesRealToRef(spdim, spdim, v0parent, invJparent, vtmp, &pointsReal[i*spdim]);CHKERRQ(ierr);
1285       }
1286       if (id == PETSCFE_CLASSID) {
1287         ierr = PetscFEGetTabulation(fe,nPoints,pointsReal,&B,NULL,NULL);CHKERRQ(ierr);
1288       }
1289       else {
1290         ierr = PetscFVGetTabulation(fv,nPoints,pointsReal,&B,NULL,NULL);CHKERRQ(ierr);
1291       }
1292       offset = 0;
1293       for (i = 0; i < fSize; i++) {
1294         PetscInt        qPoints;
1295         PetscQuadrature quad;
1296 
1297         ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr);
1298         ierr = PetscQuadratureGetData(quad,NULL,&qPoints,NULL,NULL);CHKERRQ(ierr);
1299         for (j = 0; j < fSize; j++) {
1300           PetscScalar val = 0.;
1301 
1302           for (k = 0; k < qPoints; k++) {
1303             val += B[((offset + k) * fSize + j) * fComp] * weights[k];
1304           }
1305           MatSetValue(Bmat,i,j,val,INSERT_VALUES);CHKERRQ(ierr);
1306         }
1307         offset += qPoints;
1308       }
1309       if (id == PETSCFE_CLASSID) {
1310         ierr = PetscFERestoreTabulation(fe,nPoints,pointsReal,&B,NULL,NULL);CHKERRQ(ierr);
1311       }
1312       else {
1313         ierr = PetscFVRestoreTabulation(fv,nPoints,pointsReal,&B,NULL,NULL);CHKERRQ(ierr);
1314       }
1315       ierr = MatAssemblyBegin(Bmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1316       ierr = MatAssemblyEnd(Bmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1317       ierr = MatMatSolve(Amat,Bmat,Xmat);CHKERRQ(ierr);
1318       ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSizeP,&closureP);CHKERRQ(ierr);
1319       ierr = PetscMalloc2(closureSize+1,&childOffsets,closureSizeP+1,&parentOffsets);CHKERRQ(ierr);
1320       childOffsets[0] = 0;
1321       for (i = 0; i < closureSize; i++) {
1322         PetscInt p = closure[2*i];
1323         PetscInt dof;
1324 
1325         if (numFields > 1) {
1326           ierr = PetscSectionGetFieldDof(section,p,f,&dof);CHKERRQ(ierr);
1327         }
1328         else {
1329           ierr = PetscSectionGetDof(section,p,&dof);CHKERRQ(ierr);
1330         }
1331         childOffsets[i+1]=childOffsets[i]+dof / fComp;
1332       }
1333       parentOffsets[0] = 0;
1334       for (i = 0; i < closureSizeP; i++) {
1335         PetscInt p = closureP[2*i];
1336         PetscInt dof;
1337 
1338         if (numFields > 1) {
1339           ierr = PetscSectionGetFieldDof(section,p,f,&dof);CHKERRQ(ierr);
1340         }
1341         else {
1342           ierr = PetscSectionGetDof(section,p,&dof);CHKERRQ(ierr);
1343         }
1344         parentOffsets[i+1]=parentOffsets[i]+dof / fComp;
1345       }
1346       for (i = 0; i < closureSize; i++) {
1347         PetscInt conDof, conOff, aDof, aOff;
1348         PetscInt p = closure[2*i];
1349         PetscInt o = closure[2*i+1];
1350 
1351         if (p < conStart || p >= conEnd) continue;
1352         if (numFields > 1) {
1353           ierr = PetscSectionGetFieldDof(cSec,p,f,&conDof);CHKERRQ(ierr);
1354           ierr = PetscSectionGetFieldOffset(cSec,p,f,&conOff);CHKERRQ(ierr);
1355         }
1356         else {
1357           ierr = PetscSectionGetDof(cSec,p,&conDof);CHKERRQ(ierr);
1358           ierr = PetscSectionGetOffset(cSec,p,&conOff);CHKERRQ(ierr);
1359         }
1360         if (!conDof) continue;
1361         ierr = PetscSectionGetDof(aSec,p,&aDof);CHKERRQ(ierr);
1362         ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr);
1363         for (k = 0; k < aDof; k++) {
1364           PetscInt a = anchors[aOff + k];
1365           PetscInt aSecDof, aSecOff;
1366 
1367           if (numFields > 1) {
1368             ierr = PetscSectionGetFieldDof(section,a,f,&aSecDof);CHKERRQ(ierr);
1369             ierr = PetscSectionGetFieldOffset(section,a,f,&aSecOff);CHKERRQ(ierr);
1370           }
1371           else {
1372             ierr = PetscSectionGetDof(section,a,&aSecDof);CHKERRQ(ierr);
1373             ierr = PetscSectionGetOffset(section,a,&aSecOff);CHKERRQ(ierr);
1374           }
1375           if (!aSecDof) continue;
1376 
1377           for (j = 0; j < closureSizeP; j++) {
1378             PetscInt q = closureP[2*j];
1379             PetscInt oq = closureP[2*j+1];
1380 
1381             if (q == a) {
1382               PetscInt r, s, t;
1383 
1384               for (r = childOffsets[i]; r < childOffsets[i+1]; r++) {
1385                 for (s = parentOffsets[j]; s < parentOffsets[j+1]; s++) {
1386                   PetscScalar val;
1387                   PetscInt insertCol, insertRow;
1388 
1389                   ierr = MatGetValue(Xmat,r,s,&val);CHKERRQ(ierr);
1390                   if (o >= 0) {
1391                     insertRow = conOff + fComp * (r - childOffsets[i]);
1392                   }
1393                   else {
1394                     insertRow = conOff + fComp * (childOffsets[i + 1] - 1 - r);
1395                   }
1396                   if (oq >= 0) {
1397                     insertCol = aSecOff + fComp * (s - parentOffsets[j]);
1398                   }
1399                   else {
1400                     insertCol = aSecOff + fComp * (parentOffsets[j + 1] - 1 - s);
1401                   }
1402                   for (t = 0; t < fComp; t++) {
1403                     ierr = MatSetValue(cMat,insertRow + t,insertCol + t,val,INSERT_VALUES);CHKERRQ(ierr);
1404                   }
1405                 }
1406               }
1407             }
1408           }
1409         }
1410       }
1411       ierr = PetscFree2(childOffsets,parentOffsets);CHKERRQ(ierr);
1412       ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
1413       ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSizeP,&closureP);CHKERRQ(ierr);
1414     }
1415     ierr = MatDestroy(&Amat);CHKERRQ(ierr);
1416     ierr = MatDestroy(&Bmat);CHKERRQ(ierr);
1417     ierr = MatDestroy(&Xmat);CHKERRQ(ierr);
1418     ierr = PetscFree3(weights,pointsRef,pointsReal);CHKERRQ(ierr);
1419   }
1420   ierr = MatAssemblyBegin(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1421   ierr = MatAssemblyEnd(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1422   ierr = PetscFree6(v0,v0parent,vtmp,J,Jparent,invJparent);CHKERRQ(ierr);
1423   ierr = ISRestoreIndices(aIS,&anchors);CHKERRQ(ierr);
1424 
1425   PetscFunctionReturn(0);
1426 }
1427 
1428 #undef __FUNCT__
1429 #define __FUNCT__ "DMPlexComputeAnchorMatrix_Tree_FromReference"
1430 static PetscErrorCode DMPlexComputeAnchorMatrix_Tree_FromReference(DM dm, PetscSection section, PetscSection conSec, Mat cMat)
1431 {
1432   DM             refTree;
1433   PetscDS        ds;
1434   Mat            refCmat;
1435   PetscInt       numFields, f, pRefStart, pRefEnd, p, *rows, *cols, maxDof, maxAnDof, *perm, *iperm, pStart, pEnd, conStart, conEnd, **refPointFieldN;
1436   PetscScalar ***refPointFieldMats, *pointWork;
1437   PetscSection   refConSec, refAnSec, anSec, refSection;
1438   IS             refAnIS, anIS;
1439   const PetscInt *refAnchors, *anchors;
1440   PetscErrorCode ierr;
1441 
1442   PetscFunctionBegin;
1443   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1444   ierr = DMGetDS(dm,&ds);CHKERRQ(ierr);
1445   ierr = PetscDSGetNumFields(ds,&numFields);CHKERRQ(ierr);
1446   ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr);
1447   ierr = DMSetDS(refTree,ds);CHKERRQ(ierr);
1448   ierr = DMGetDefaultConstraints(refTree,&refConSec,&refCmat);CHKERRQ(ierr);
1449   ierr = DMPlexGetChart(refTree,&pRefStart,&pRefEnd);CHKERRQ(ierr);
1450   ierr = DMPlexGetAnchors(refTree,&refAnSec,&refAnIS);CHKERRQ(ierr);
1451   ierr = DMPlexGetAnchors(dm,&anSec,&anIS);CHKERRQ(ierr);
1452   ierr = ISGetIndices(refAnIS,&refAnchors);CHKERRQ(ierr);
1453   ierr = ISGetIndices(anIS,&anchors);CHKERRQ(ierr);
1454   ierr = DMGetDefaultSection(refTree,&refSection);CHKERRQ(ierr);
1455   ierr = PetscSectionGetChart(refConSec,&pRefStart,&pRefEnd);CHKERRQ(ierr);
1456   ierr = PetscSectionGetChart(conSec,&conStart,&conEnd);CHKERRQ(ierr);
1457   ierr = PetscMalloc1(pRefEnd-pRefStart,&refPointFieldMats);CHKERRQ(ierr);
1458   ierr = PetscMalloc1(pRefEnd-pRefStart,&refPointFieldN);CHKERRQ(ierr);
1459   ierr = PetscSectionGetMaxDof(refConSec,&maxDof);CHKERRQ(ierr);
1460   ierr = PetscSectionGetMaxDof(refAnSec,&maxAnDof);CHKERRQ(ierr);
1461   ierr = PetscMalloc1(maxDof,&rows);CHKERRQ(ierr);
1462   ierr = PetscMalloc1(maxDof*maxAnDof,&cols);CHKERRQ(ierr);
1463   ierr = PetscMalloc1(maxDof*maxDof*maxAnDof,&pointWork);CHKERRQ(ierr);
1464 
1465   /* step 1: get submats for every constrained point in the reference tree */
1466   for (p = pRefStart; p < pRefEnd; p++) {
1467     PetscInt parent, closureSize, *closure = NULL, pDof;
1468 
1469     ierr = DMPlexGetTreeParent(refTree,p,&parent,NULL);CHKERRQ(ierr);
1470     ierr = PetscSectionGetDof(refConSec,p,&pDof);CHKERRQ(ierr);
1471     if (!pDof || parent == p) continue;
1472 
1473     ierr = PetscMalloc1(numFields,&refPointFieldMats[p-pRefStart]);CHKERRQ(ierr);
1474     ierr = PetscCalloc1(numFields,&refPointFieldN[p-pRefStart]);CHKERRQ(ierr);
1475     ierr = DMPlexGetTransitiveClosure(refTree,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
1476     for (f = 0; f < numFields; f++) {
1477       PetscInt cDof, cOff, numCols, r, i, fComp;
1478       PetscObject disc;
1479       PetscClassId id;
1480       PetscFE fe = NULL;
1481       PetscFV fv = NULL;
1482 
1483       ierr = PetscDSGetDiscretization(ds,f,&disc);CHKERRQ(ierr);
1484       ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
1485       if (id == PETSCFE_CLASSID) {
1486         fe = (PetscFE) disc;
1487         ierr = PetscFEGetNumComponents(fe,&fComp);CHKERRQ(ierr);
1488       }
1489       else if (id == PETSCFV_CLASSID) {
1490         fv = (PetscFV) disc;
1491         ierr = PetscFVGetNumComponents(fv,&fComp);CHKERRQ(ierr);
1492       }
1493       else SETERRQ1(PetscObjectComm(disc),PETSC_ERR_ARG_UNKNOWN_TYPE, "PetscDS discretization id %d not recognized.", id);
1494 
1495       if (numFields > 1) {
1496         ierr = PetscSectionGetFieldDof(refConSec,p,f,&cDof);CHKERRQ(ierr);
1497         ierr = PetscSectionGetFieldOffset(refConSec,p,f,&cOff);CHKERRQ(ierr);
1498       }
1499       else {
1500         ierr = PetscSectionGetDof(refConSec,p,&cDof);CHKERRQ(ierr);
1501         ierr = PetscSectionGetOffset(refConSec,p,&cOff);CHKERRQ(ierr);
1502       }
1503 
1504       if (!cDof) continue;
1505       for (r = 0; r < cDof; r++) {
1506         rows[r] = cOff + r;
1507       }
1508       numCols = 0;
1509       for (i = 0; i < closureSize; i++) {
1510         PetscInt q = closure[2*i];
1511         PetscInt o = closure[2*i+1];
1512         PetscInt aDof, aOff, j;
1513 
1514         if (numFields > 1) {
1515           ierr = PetscSectionGetFieldDof(refSection,q,f,&aDof);CHKERRQ(ierr);
1516           ierr = PetscSectionGetFieldOffset(refSection,q,f,&aOff);CHKERRQ(ierr);
1517         }
1518         else {
1519           ierr = PetscSectionGetDof(refSection,q,&aDof);CHKERRQ(ierr);
1520           ierr = PetscSectionGetOffset(refSection,q,&aOff);CHKERRQ(ierr);
1521         }
1522 
1523         for (j = 0; j < aDof; j++) {
1524           PetscInt node = (o >= 0) ? (j / fComp) : ((aDof - 1 - j) / fComp);
1525           PetscInt comp = (j % fComp);
1526 
1527           cols[numCols++] = aOff + node * fComp + comp;
1528         }
1529       }
1530       refPointFieldN[p-pRefStart][f] = numCols;
1531       ierr = PetscMalloc1(cDof*numCols,&refPointFieldMats[p-pRefStart][f]);CHKERRQ(ierr);
1532       ierr = MatGetValues(refCmat,cDof,rows,numCols,cols,refPointFieldMats[p-pRefStart][f]);CHKERRQ(ierr);
1533     }
1534     ierr = DMPlexRestoreTransitiveClosure(refTree,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
1535   }
1536 
1537   /* step 2: compute the preorder */
1538   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
1539   ierr = PetscMalloc2(pEnd-pStart,&perm,pEnd-pStart,&iperm);CHKERRQ(ierr);
1540   for (p = pStart; p < pEnd; p++) {
1541     perm[p - pStart] = p;
1542     iperm[p - pStart] = p-pStart;
1543   }
1544   for (p = 0; p < pEnd - pStart;) {
1545     PetscInt point = perm[p];
1546     PetscInt parent;
1547 
1548     ierr = DMPlexGetTreeParent(dm,point,&parent,NULL);CHKERRQ(ierr);
1549     if (parent == point) {
1550       p++;
1551     }
1552     else {
1553       PetscInt size, closureSize, *closure = NULL, i;
1554 
1555       ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
1556       for (i = 0; i < closureSize; i++) {
1557         PetscInt q = closure[2*i];
1558         if (iperm[q-pStart] > iperm[point-pStart]) {
1559           /* swap */
1560           perm[p]               = q;
1561           perm[iperm[q-pStart]] = point;
1562           iperm[point-pStart]   = iperm[q-pStart];
1563           iperm[q-pStart]       = p;
1564           break;
1565         }
1566       }
1567       size = closureSize;
1568       ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
1569       if (i == size) {
1570         p++;
1571       }
1572     }
1573   }
1574 
1575   /* step 3: fill the constraint matrix */
1576   /* we are going to use a preorder progressive fill strategy.  Mat doesn't
1577    * allow progressive fill without assembly, so we are going to set up the
1578    * values outside of the Mat first.
1579    */
1580   {
1581     PetscInt nRows, row, nnz;
1582     PetscBool done;
1583     const PetscInt *ia, *ja;
1584     PetscScalar *vals;
1585 
1586     ierr = MatGetRowIJ(cMat,0,PETSC_FALSE,PETSC_FALSE,&nRows,&ia,&ja,&done);CHKERRQ(ierr);
1587     if (!done) SETERRQ(PetscObjectComm((PetscObject)cMat),PETSC_ERR_PLIB,"Could not get RowIJ of constraint matrix");
1588     nnz  = ia[nRows];
1589     /* malloc and then zero rows right before we fill them: this way valgrind
1590      * can tell if we are doing progressive fill in the wrong order */
1591     ierr = PetscMalloc1(nnz,&vals);CHKERRQ(ierr);
1592     for (p = 0; p < pEnd - pStart; p++) {
1593       PetscInt parent, childid, closureSize, *closure = NULL;
1594       PetscInt point = perm[p], pointDof;
1595 
1596       ierr = DMPlexGetTreeParent(dm,point,&parent,&childid);CHKERRQ(ierr);
1597       if ((point < conStart) || (point >= conEnd) || (parent == point)) continue;
1598       ierr = PetscSectionGetDof(conSec,point,&pointDof);CHKERRQ(ierr);
1599       if (!pointDof) continue;
1600       ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
1601       for (f = 0; f < numFields; f++) {
1602         PetscInt cDof, cOff, numCols, numFillCols, i, r, fComp, matOffset, offset;
1603         PetscScalar *pointMat;
1604         PetscObject disc;
1605         PetscClassId id;
1606         PetscFE fe = NULL;
1607         PetscFV fv = NULL;
1608 
1609         ierr = PetscDSGetDiscretization(ds,f,&disc);CHKERRQ(ierr);
1610         ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr);
1611         if (id == PETSCFE_CLASSID) {
1612           fe = (PetscFE) disc;
1613           ierr = PetscFEGetNumComponents(fe,&fComp);CHKERRQ(ierr);
1614         }
1615         else if (id == PETSCFV_CLASSID) {
1616           fv = (PetscFV) disc;
1617           ierr = PetscFVGetNumComponents(fv,&fComp);CHKERRQ(ierr);
1618         }
1619 
1620         if (numFields > 1) {
1621           ierr = PetscSectionGetFieldDof(conSec,point,f,&cDof);CHKERRQ(ierr);
1622           ierr = PetscSectionGetFieldOffset(conSec,point,f,&cOff);CHKERRQ(ierr);
1623         }
1624         else {
1625           ierr = PetscSectionGetDof(conSec,point,&cDof);CHKERRQ(ierr);
1626           ierr = PetscSectionGetOffset(conSec,point,&cOff);CHKERRQ(ierr);
1627         }
1628         if (!cDof) continue;
1629 
1630         /* make sure that every row for this point is the same size */
1631 #if defined(PETSC_USE_DEBUG)
1632         for (r = 0; r < cDof; r++) {
1633           if (cDof > 1 && r) {
1634             if ((ia[cOff+r+1]-ia[cOff+r]) != (ia[cOff+r]-ia[cOff+r-1])) {
1635               SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Two point rows have different nnz: %d vs. %d", (ia[rows[r]+1]-ia[rows[r]]), (ia[rows[r]]-ia[rows[r]-1]));
1636             }
1637           }
1638         }
1639 #endif
1640         /* zero rows */
1641         for (i = ia[cOff] ; i< ia[cOff+cDof];i++) {
1642           vals[i] = 0.;
1643         }
1644         matOffset = ia[cOff];
1645         numFillCols = ia[cOff+1] - matOffset;
1646         pointMat = refPointFieldMats[childid-pRefStart][f];
1647         numCols = refPointFieldN[childid-pRefStart][f];
1648         offset = 0;
1649         for (i = 0; i < closureSize; i++) {
1650           PetscInt q = closure[2*i];
1651           PetscInt o = closure[2*i+1];
1652           PetscInt aDof, aOff, j, k, qConDof, qConOff;
1653 
1654           qConDof = qConOff = 0;
1655           if (numFields > 1) {
1656             ierr = PetscSectionGetFieldDof(section,q,f,&aDof);CHKERRQ(ierr);
1657             ierr = PetscSectionGetFieldOffset(section,q,f,&aOff);CHKERRQ(ierr);
1658             if (q >= conStart && q < conEnd) {
1659               ierr = PetscSectionGetFieldDof(conSec,q,f,&qConDof);CHKERRQ(ierr);
1660               ierr = PetscSectionGetFieldOffset(conSec,q,f,&qConOff);CHKERRQ(ierr);
1661             }
1662           }
1663           else {
1664             ierr = PetscSectionGetDof(section,q,&aDof);CHKERRQ(ierr);
1665             ierr = PetscSectionGetOffset(section,q,&aOff);CHKERRQ(ierr);
1666             if (q >= conStart && q < conEnd) {
1667               ierr = PetscSectionGetDof(conSec,q,&qConDof);CHKERRQ(ierr);
1668               ierr = PetscSectionGetOffset(conSec,q,&qConOff);CHKERRQ(ierr);
1669             }
1670           }
1671           if (!aDof) continue;
1672           if (qConDof) {
1673             /* this point has anchors: its rows of the matrix should already
1674              * be filled, thanks to preordering */
1675             /* first multiply into pointWork, then set in matrix */
1676             PetscInt aMatOffset = ia[qConOff];
1677             PetscInt aNumFillCols = ia[qConOff + 1] - aMatOffset;
1678             for (r = 0; r < cDof; r++) {
1679               for (j = 0; j < aNumFillCols; j++) {
1680                 PetscScalar inVal = 0;
1681                 for (k = 0; k < aDof; k++) {
1682                   PetscInt node = (o >= 0) ? (k / fComp) : ((aDof - 1 - k) / fComp);
1683                   PetscInt comp = (k % fComp);
1684                   PetscInt col  = node * fComp + comp;
1685 
1686                   inVal += pointMat[r * numCols + offset + col] * vals[aMatOffset + aNumFillCols * k + j];
1687                 }
1688                 pointWork[r * aNumFillCols + j] = inVal;
1689               }
1690             }
1691             /* assume that the columns are sorted, spend less time searching */
1692             for (j = 0, k = 0; j < aNumFillCols; j++) {
1693               PetscInt col = ja[aMatOffset + j];
1694               for (;k < numFillCols; k++) {
1695                 if (ja[matOffset + k] == col) {
1696                   break;
1697                 }
1698               }
1699               if (k == numFillCols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"No nonzero space for (%d, %d)", cOff, col);
1700               for (r = 0; r < cDof; r++) {
1701                 vals[matOffset + numFillCols * r + k] = pointWork[r * aNumFillCols + j];
1702               }
1703             }
1704           }
1705           else {
1706             /* find where to put this portion of pointMat into the matrix */
1707             for (k = 0; k < numFillCols; k++) {
1708               if (ja[matOffset + k] == aOff) {
1709                 break;
1710               }
1711             }
1712             if (k == numFillCols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"No nonzero space for (%d, %d)", cOff, aOff);
1713             for (j = 0; j < aDof; j++) {
1714               PetscInt node = (o >= 0) ? (j / fComp) : ((aDof - 1 - j) / fComp);
1715               PetscInt comp = (j % fComp);
1716               PetscInt col  = node * fComp + comp;
1717               for (r = 0; r < cDof; r++) {
1718                 vals[matOffset + numFillCols * r + k + col] += pointMat[r * numCols + offset + col];
1719               }
1720             }
1721           }
1722           offset += aDof;
1723         }
1724       }
1725       ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
1726     }
1727     for (row = 0; row < nRows; row++) {
1728       ierr = MatSetValues(cMat,1,&row,ia[row+1]-ia[row],&ja[ia[row]],&vals[ia[row]],INSERT_VALUES);CHKERRQ(ierr);
1729     }
1730     ierr = MatRestoreRowIJ(cMat,0,PETSC_FALSE,PETSC_FALSE,&nRows,&ia,&ja,&done);CHKERRQ(ierr);
1731     if (!done) SETERRQ(PetscObjectComm((PetscObject)cMat),PETSC_ERR_PLIB,"Could not restore RowIJ of constraint matrix");
1732     ierr = MatAssemblyBegin(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1733     ierr = MatAssemblyEnd(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1734     ierr = PetscFree(vals);CHKERRQ(ierr);
1735   }
1736 
1737   /* clean up */
1738   ierr = ISRestoreIndices(refAnIS,&refAnchors);CHKERRQ(ierr);
1739   ierr = ISRestoreIndices(anIS,&anchors);CHKERRQ(ierr);
1740   ierr = PetscFree2(perm,iperm);CHKERRQ(ierr);
1741   ierr = PetscFree(rows);CHKERRQ(ierr);
1742   ierr = PetscFree(cols);CHKERRQ(ierr);
1743   ierr = PetscFree(pointWork);CHKERRQ(ierr);
1744   for (p = pRefStart; p < pRefEnd; p++) {
1745     PetscInt parent, pDof;
1746 
1747     ierr = DMPlexGetTreeParent(refTree,p,&parent,NULL);CHKERRQ(ierr);
1748     ierr = PetscSectionGetDof(refConSec,p,&pDof);CHKERRQ(ierr);
1749     if (!pDof || parent == p) continue;
1750 
1751     for (f = 0; f < numFields; f++) {
1752       PetscInt cDof;
1753 
1754       if (numFields > 1) {
1755         ierr = PetscSectionGetFieldDof(refConSec,p,f,&cDof);CHKERRQ(ierr);
1756       }
1757       else {
1758         ierr = PetscSectionGetDof(refConSec,p,&cDof);CHKERRQ(ierr);
1759       }
1760 
1761       if (!cDof) continue;
1762       ierr = PetscFree(refPointFieldMats[p - pRefStart][f]);CHKERRQ(ierr);
1763     }
1764     ierr = PetscFree(refPointFieldMats[p - pRefStart]);CHKERRQ(ierr);
1765     ierr = PetscFree(refPointFieldN[p - pRefStart]);CHKERRQ(ierr);
1766   }
1767   ierr = PetscFree(refPointFieldMats);CHKERRQ(ierr);
1768   ierr = PetscFree(refPointFieldN);CHKERRQ(ierr);
1769   PetscFunctionReturn(0);
1770 }
1771 
1772 #undef __FUNCT__
1773 #define __FUNCT__ "DMPlexTreeRefineCell"
1774 /* refine a single cell on rank 0: this is not intended to provide good local refinement, only to create an example of
1775  * a non-conforming mesh.  Local refinement comes later */
1776 PetscErrorCode DMPlexTreeRefineCell (DM dm, PetscInt cell, DM *ncdm)
1777 {
1778   DM K;
1779   PetscMPIInt rank;
1780   PetscInt dim, *pNewStart, *pNewEnd, *pNewCount, *pOldStart, *pOldEnd, offset, d, pStart, pEnd;
1781   PetscInt numNewCones, *newConeSizes, *newCones, *newOrientations;
1782   PetscInt *Kembedding;
1783   PetscInt *cellClosure=NULL, nc;
1784   PetscScalar *newVertexCoords;
1785   PetscInt numPointsWithParents, *parents, *childIDs, *perm, *iperm, *preOrient, pOffset;
1786   PetscSection parentSection;
1787   PetscErrorCode ierr;
1788 
1789   PetscFunctionBegin;
1790   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
1791   ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr);
1792   ierr = DMPlexCreate(PetscObjectComm((PetscObject)dm), ncdm);CHKERRQ(ierr);
1793   ierr = DMSetDimension(*ncdm,dim);CHKERRQ(ierr);
1794 
1795   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
1796   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&parentSection);CHKERRQ(ierr);
1797   ierr = DMPlexGetReferenceTree(dm,&K);CHKERRQ(ierr);
1798   if (!rank) {
1799     /* compute the new charts */
1800     ierr = PetscMalloc5(dim+1,&pNewCount,dim+1,&pNewStart,dim+1,&pNewEnd,dim+1,&pOldStart,dim+1,&pOldEnd);CHKERRQ(ierr);
1801     offset = 0;
1802     for (d = 0; d <= dim; d++) {
1803       PetscInt pOldCount, kStart, kEnd, k;
1804 
1805       pNewStart[d] = offset;
1806       ierr = DMPlexGetHeightStratum(dm,d,&pOldStart[d],&pOldEnd[d]);CHKERRQ(ierr);
1807       ierr = DMPlexGetHeightStratum(K,d,&kStart,&kEnd);CHKERRQ(ierr);
1808       pOldCount = pOldEnd[d] - pOldStart[d];
1809       /* adding the new points */
1810       pNewCount[d] = pOldCount + kEnd - kStart;
1811       if (!d) {
1812         /* removing the cell */
1813         pNewCount[d]--;
1814       }
1815       for (k = kStart; k < kEnd; k++) {
1816         PetscInt parent;
1817         ierr = DMPlexGetTreeParent(K,k,&parent,NULL);CHKERRQ(ierr);
1818         if (parent == k) {
1819           /* avoid double counting points that won't actually be new */
1820           pNewCount[d]--;
1821         }
1822       }
1823       pNewEnd[d] = pNewStart[d] + pNewCount[d];
1824       offset = pNewEnd[d];
1825 
1826     }
1827     if (cell < pOldStart[0] || cell >= pOldEnd[0]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"%d not in cell range [%d, %d)", cell, pOldStart[0], pOldEnd[0]);
1828     /* get the current closure of the cell that we are removing */
1829     ierr = DMPlexGetTransitiveClosure(dm,cell,PETSC_TRUE,&nc,&cellClosure);CHKERRQ(ierr);
1830 
1831     ierr = PetscMalloc1(pNewEnd[dim],&newConeSizes);CHKERRQ(ierr);
1832     {
1833       PetscInt kStart, kEnd, k, closureSizeK, *closureK = NULL, j;
1834 
1835       ierr = DMPlexGetChart(K,&kStart,&kEnd);CHKERRQ(ierr);
1836       ierr = PetscMalloc4(kEnd-kStart,&Kembedding,kEnd-kStart,&perm,kEnd-kStart,&iperm,kEnd-kStart,&preOrient);CHKERRQ(ierr);
1837 
1838       for (k = kStart; k < kEnd; k++) {
1839         perm[k - kStart] = k;
1840         iperm [k - kStart] = k - kStart;
1841         preOrient[k - kStart] = 0;
1842       }
1843 
1844       ierr = DMPlexGetTransitiveClosure(K,0,PETSC_TRUE,&closureSizeK,&closureK);CHKERRQ(ierr);
1845       for (j = 1; j < closureSizeK; j++) {
1846         PetscInt parentOrientA = closureK[2*j+1];
1847         PetscInt parentOrientB = cellClosure[2*j+1];
1848         PetscInt p, q;
1849 
1850         p = closureK[2*j];
1851         q = cellClosure[2*j];
1852         for (d = 0; d <= dim; d++) {
1853           if (q >= pOldStart[d] && q < pOldEnd[d]) {
1854             Kembedding[p] = (q - pOldStart[d]) + pNewStart[d];
1855           }
1856         }
1857         if (parentOrientA != parentOrientB) {
1858           PetscInt numChildren, i;
1859           const PetscInt *children;
1860 
1861           ierr = DMPlexGetTreeChildren(K,p,&numChildren,&children);CHKERRQ(ierr);
1862           for (i = 0; i < numChildren; i++) {
1863             PetscInt kPerm, oPerm;
1864 
1865             k    = children[i];
1866             ierr = DMPlexReferenceTreeGetChildSymmetry(K,p,parentOrientA,0,k,parentOrientB,&oPerm,&kPerm);CHKERRQ(ierr);
1867             /* perm = what refTree position I'm in */
1868             perm[kPerm-kStart]      = k;
1869             /* iperm = who is at this position */
1870             iperm[k-kStart]         = kPerm-kStart;
1871             preOrient[kPerm-kStart] = oPerm;
1872           }
1873         }
1874       }
1875       ierr = DMPlexRestoreTransitiveClosure(K,0,PETSC_TRUE,&closureSizeK,&closureK);CHKERRQ(ierr);
1876     }
1877     ierr = PetscSectionSetChart(parentSection,0,pNewEnd[dim]);CHKERRQ(ierr);
1878     offset = 0;
1879     numNewCones = 0;
1880     for (d = 0; d <= dim; d++) {
1881       PetscInt kStart, kEnd, k;
1882       PetscInt p;
1883       PetscInt size;
1884 
1885       for (p = pOldStart[d]; p < pOldEnd[d]; p++) {
1886         /* skip cell 0 */
1887         if (p == cell) continue;
1888         /* old cones to new cones */
1889         ierr = DMPlexGetConeSize(dm,p,&size);CHKERRQ(ierr);
1890         newConeSizes[offset++] = size;
1891         numNewCones += size;
1892       }
1893 
1894       ierr = DMPlexGetHeightStratum(K,d,&kStart,&kEnd);CHKERRQ(ierr);
1895       for (k = kStart; k < kEnd; k++) {
1896         PetscInt kParent;
1897 
1898         ierr = DMPlexGetTreeParent(K,k,&kParent,NULL);CHKERRQ(ierr);
1899         if (kParent != k) {
1900           Kembedding[k] = offset;
1901           ierr = DMPlexGetConeSize(K,k,&size);CHKERRQ(ierr);
1902           newConeSizes[offset++] = size;
1903           numNewCones += size;
1904           if (kParent != 0) {
1905             ierr = PetscSectionSetDof(parentSection,Kembedding[k],1);CHKERRQ(ierr);
1906           }
1907         }
1908       }
1909     }
1910 
1911     ierr = PetscSectionSetUp(parentSection);CHKERRQ(ierr);
1912     ierr = PetscSectionGetStorageSize(parentSection,&numPointsWithParents);CHKERRQ(ierr);
1913     ierr = PetscMalloc2(numNewCones,&newCones,numNewCones,&newOrientations);CHKERRQ(ierr);
1914     ierr = PetscMalloc2(numPointsWithParents,&parents,numPointsWithParents,&childIDs);CHKERRQ(ierr);
1915 
1916     /* fill new cones */
1917     offset = 0;
1918     for (d = 0; d <= dim; d++) {
1919       PetscInt kStart, kEnd, k, l;
1920       PetscInt p;
1921       PetscInt size;
1922       const PetscInt *cone, *orientation;
1923 
1924       for (p = pOldStart[d]; p < pOldEnd[d]; p++) {
1925         /* skip cell 0 */
1926         if (p == cell) continue;
1927         /* old cones to new cones */
1928         ierr = DMPlexGetConeSize(dm,p,&size);CHKERRQ(ierr);
1929         ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr);
1930         ierr = DMPlexGetConeOrientation(dm,p,&orientation);CHKERRQ(ierr);
1931         for (l = 0; l < size; l++) {
1932           newCones[offset]          = (cone[l] - pOldStart[d + 1]) + pNewStart[d + 1];
1933           newOrientations[offset++] = orientation[l];
1934         }
1935       }
1936 
1937       ierr = DMPlexGetHeightStratum(K,d,&kStart,&kEnd);CHKERRQ(ierr);
1938       for (k = kStart; k < kEnd; k++) {
1939         PetscInt kPerm = perm[k], kParent;
1940         PetscInt preO  = preOrient[k];
1941 
1942         ierr = DMPlexGetTreeParent(K,k,&kParent,NULL);CHKERRQ(ierr);
1943         if (kParent != k) {
1944           /* embed new cones */
1945           ierr = DMPlexGetConeSize(K,k,&size);CHKERRQ(ierr);
1946           ierr = DMPlexGetCone(K,kPerm,&cone);CHKERRQ(ierr);
1947           ierr = DMPlexGetConeOrientation(K,kPerm,&orientation);CHKERRQ(ierr);
1948           for (l = 0; l < size; l++) {
1949             PetscInt q, m = (preO >= 0) ? ((preO + l) % size) : ((size -(preO + 1) - l) % size);
1950             PetscInt newO, lSize, oTrue;
1951 
1952             q                         = iperm[cone[m]];
1953             newCones[offset]          = Kembedding[q];
1954             ierr                      = DMPlexGetConeSize(K,q,&lSize);CHKERRQ(ierr);
1955             oTrue                     = orientation[m];
1956             oTrue                     = ((!lSize) || (preOrient[k] >= 0)) ? oTrue : -(oTrue + 2);
1957             newO                      = DihedralCompose(lSize,oTrue,preOrient[q]);
1958             newOrientations[offset++] = newO;
1959           }
1960           if (kParent != 0) {
1961             PetscInt newPoint = Kembedding[kParent];
1962             ierr              = PetscSectionGetOffset(parentSection,Kembedding[k],&pOffset);CHKERRQ(ierr);
1963             parents[pOffset]  = newPoint;
1964             childIDs[pOffset] = k;
1965           }
1966         }
1967       }
1968     }
1969 
1970     ierr = PetscMalloc1(dim*(pNewEnd[dim]-pNewStart[dim]),&newVertexCoords);CHKERRQ(ierr);
1971 
1972     /* fill coordinates */
1973     offset = 0;
1974     {
1975       PetscInt kStart, kEnd, l;
1976       PetscSection vSection;
1977       PetscInt v;
1978       Vec coords;
1979       PetscScalar *coordvals;
1980       PetscInt dof, off;
1981       PetscReal v0[3], J[9], detJ;
1982 
1983 #if defined(PETSC_USE_DEBUG)
1984       {
1985         PetscInt k;
1986         ierr = DMPlexGetHeightStratum(K,0,&kStart,&kEnd);CHKERRQ(ierr);
1987         for (k = kStart; k < kEnd; k++) {
1988           ierr = DMPlexComputeCellGeometryFEM(K, k, NULL, v0, J, NULL, &detJ);CHKERRQ(ierr);
1989           if (detJ <= 0.) SETERRQ1 (PETSC_COMM_SELF,PETSC_ERR_PLIB,"reference tree cell %d has bad determinant",k);
1990         }
1991       }
1992 #endif
1993       ierr = DMPlexComputeCellGeometryFEM(dm, cell, NULL, v0, J, NULL, &detJ);CHKERRQ(ierr);
1994       ierr = DMGetCoordinateSection(dm,&vSection);CHKERRQ(ierr);
1995       ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
1996       ierr = VecGetArray(coords,&coordvals);CHKERRQ(ierr);
1997       for (v = pOldStart[dim]; v < pOldEnd[dim]; v++) {
1998 
1999         ierr = PetscSectionGetDof(vSection,v,&dof);CHKERRQ(ierr);
2000         ierr = PetscSectionGetOffset(vSection,v,&off);CHKERRQ(ierr);
2001         for (l = 0; l < dof; l++) {
2002           newVertexCoords[offset++] = coordvals[off + l];
2003         }
2004       }
2005       ierr = VecRestoreArray(coords,&coordvals);CHKERRQ(ierr);
2006 
2007       ierr = DMGetCoordinateSection(K,&vSection);CHKERRQ(ierr);
2008       ierr = DMGetCoordinatesLocal(K,&coords);CHKERRQ(ierr);
2009       ierr = VecGetArray(coords,&coordvals);CHKERRQ(ierr);
2010       ierr = DMPlexGetDepthStratum(K,0,&kStart,&kEnd);CHKERRQ(ierr);
2011       for (v = kStart; v < kEnd; v++) {
2012         PetscReal coord[3], newCoord[3];
2013         PetscInt  vPerm = perm[v];
2014         PetscInt  kParent;
2015 
2016         ierr = DMPlexGetTreeParent(K,v,&kParent,NULL);CHKERRQ(ierr);
2017         if (kParent != v) {
2018           /* this is a new vertex */
2019           ierr = PetscSectionGetOffset(vSection,vPerm,&off);CHKERRQ(ierr);
2020           for (l = 0; l < dim; ++l) coord[l] = PetscRealPart(coordvals[off+l]);
2021           CoordinatesRefToReal(dim, dim, v0, J, coord, newCoord);CHKERRQ(ierr);
2022           for (l = 0; l < dim; ++l) newVertexCoords[offset+l] = newCoord[l];
2023           offset += dim;
2024         }
2025       }
2026       ierr = VecRestoreArray(coords,&coordvals);CHKERRQ(ierr);
2027     }
2028 
2029     /* need to reverse the order of pNewCount: vertices first, cells last */
2030     for (d = 0; d < (dim + 1) / 2; d++) {
2031       PetscInt tmp;
2032 
2033       tmp = pNewCount[d];
2034       pNewCount[d] = pNewCount[dim - d];
2035       pNewCount[dim - d] = tmp;
2036     }
2037 
2038     ierr = DMPlexCreateFromDAG(*ncdm,dim,pNewCount,newConeSizes,newCones,newOrientations,newVertexCoords);CHKERRQ(ierr);
2039     ierr = DMPlexSetReferenceTree(*ncdm,K);CHKERRQ(ierr);
2040     ierr = DMPlexSetTree(*ncdm,parentSection,parents,childIDs);CHKERRQ(ierr);
2041 
2042     /* clean up */
2043     ierr = DMPlexRestoreTransitiveClosure(dm,cell,PETSC_TRUE,&nc,&cellClosure);CHKERRQ(ierr);
2044     ierr = PetscFree5(pNewCount,pNewStart,pNewEnd,pOldStart,pOldEnd);CHKERRQ(ierr);
2045     ierr = PetscFree(newConeSizes);CHKERRQ(ierr);
2046     ierr = PetscFree2(newCones,newOrientations);CHKERRQ(ierr);
2047     ierr = PetscFree(newVertexCoords);CHKERRQ(ierr);
2048     ierr = PetscFree2(parents,childIDs);CHKERRQ(ierr);
2049     ierr = PetscFree4(Kembedding,perm,iperm,preOrient);CHKERRQ(ierr);
2050   }
2051   else {
2052     PetscInt    p, counts[4];
2053     PetscInt    *coneSizes, *cones, *orientations;
2054     Vec         coordVec;
2055     PetscScalar *coords;
2056 
2057     for (d = 0; d <= dim; d++) {
2058       PetscInt dStart, dEnd;
2059 
2060       ierr = DMPlexGetDepthStratum(dm,d,&dStart,&dEnd);CHKERRQ(ierr);
2061       counts[d] = dEnd - dStart;
2062     }
2063     ierr = PetscMalloc1(pEnd-pStart,&coneSizes);CHKERRQ(ierr);
2064     for (p = pStart; p < pEnd; p++) {
2065       ierr = DMPlexGetConeSize(dm,p,&coneSizes[p-pStart]);CHKERRQ(ierr);
2066     }
2067     ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
2068     ierr = DMPlexGetConeOrientations(dm, &orientations);CHKERRQ(ierr);
2069     ierr = DMGetCoordinatesLocal(dm,&coordVec);CHKERRQ(ierr);
2070     ierr = VecGetArray(coordVec,&coords);CHKERRQ(ierr);
2071 
2072     ierr = PetscSectionSetChart(parentSection,pStart,pEnd);CHKERRQ(ierr);
2073     ierr = PetscSectionSetUp(parentSection);CHKERRQ(ierr);
2074     ierr = DMPlexCreateFromDAG(*ncdm,dim,counts,coneSizes,cones,orientations,NULL);CHKERRQ(ierr);
2075     ierr = DMPlexSetReferenceTree(*ncdm,K);CHKERRQ(ierr);
2076     ierr = DMPlexSetTree(*ncdm,parentSection,NULL,NULL);CHKERRQ(ierr);
2077     ierr = VecRestoreArray(coordVec,&coords);CHKERRQ(ierr);
2078   }
2079   ierr = PetscSectionDestroy(&parentSection);CHKERRQ(ierr);
2080 
2081   PetscFunctionReturn(0);
2082 }
2083