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