xref: /petsc/src/dm/impls/plex/plextree.c (revision ccffa8f8cbde339a35c113c5d2513090583a52b7)
1d6a7ad0dSToby Isaac #include <petsc-private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
2d6a7ad0dSToby Isaac #include <../src/sys/utils/hash.h>
3d6a7ad0dSToby Isaac #include <petsc-private/isimpl.h>
40c37af3bSToby Isaac #include <petsc-private/petscfeimpl.h>
5d6a7ad0dSToby Isaac #include <petscsf.h>
60c37af3bSToby Isaac #include <petscds.h>
7d6a7ad0dSToby Isaac 
8d6a7ad0dSToby Isaac /** hierarchy routines */
9d6a7ad0dSToby Isaac 
10d6a7ad0dSToby Isaac #undef __FUNCT__
11d6a7ad0dSToby Isaac #define __FUNCT__ "DMPlexSetReferenceTree"
12d6a7ad0dSToby Isaac /*@
13d6a7ad0dSToby Isaac   DMPlexSetReferenceTree - set the reference tree for hierarchically non-conforming meshes.
14d6a7ad0dSToby Isaac 
15d6a7ad0dSToby Isaac   Not collective
16d6a7ad0dSToby Isaac 
17d6a7ad0dSToby Isaac   Input Parameters:
18d6a7ad0dSToby Isaac + dm - The DMPlex object
19d6a7ad0dSToby Isaac - ref - The reference tree DMPlex object
20d6a7ad0dSToby Isaac 
210b7167a0SToby Isaac   Level: intermediate
22d6a7ad0dSToby Isaac 
23da43764aSToby Isaac .seealso: DMPlexGetReferenceTree(), DMPlexCreateDefaultReferenceTree()
24d6a7ad0dSToby Isaac @*/
25d6a7ad0dSToby Isaac PetscErrorCode DMPlexSetReferenceTree(DM dm, DM ref)
26d6a7ad0dSToby Isaac {
27d6a7ad0dSToby Isaac   DM_Plex        *mesh = (DM_Plex *)dm->data;
28d6a7ad0dSToby Isaac   PetscErrorCode  ierr;
29d6a7ad0dSToby Isaac 
30d6a7ad0dSToby Isaac   PetscFunctionBegin;
31d6a7ad0dSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
32d6a7ad0dSToby Isaac   PetscValidHeaderSpecific(ref, DM_CLASSID, 2);
33d6a7ad0dSToby Isaac   ierr = PetscObjectReference((PetscObject)ref);CHKERRQ(ierr);
34d6a7ad0dSToby Isaac   ierr = DMDestroy(&mesh->referenceTree);CHKERRQ(ierr);
35d6a7ad0dSToby Isaac   mesh->referenceTree = ref;
36d6a7ad0dSToby Isaac   PetscFunctionReturn(0);
37d6a7ad0dSToby Isaac }
38d6a7ad0dSToby Isaac 
39d6a7ad0dSToby Isaac #undef __FUNCT__
40d6a7ad0dSToby Isaac #define __FUNCT__ "DMPlexGetReferenceTree"
41d6a7ad0dSToby Isaac /*@
42d6a7ad0dSToby Isaac   DMPlexGetReferenceTree - get the reference tree for hierarchically non-conforming meshes.
43d6a7ad0dSToby Isaac 
44d6a7ad0dSToby Isaac   Not collective
45d6a7ad0dSToby Isaac 
46d6a7ad0dSToby Isaac   Input Parameters:
47d6a7ad0dSToby Isaac . dm - The DMPlex object
48d6a7ad0dSToby Isaac 
49d6a7ad0dSToby Isaac   Output Parameters
50d6a7ad0dSToby Isaac . ref - The reference tree DMPlex object
51d6a7ad0dSToby Isaac 
520b7167a0SToby Isaac   Level: intermediate
53d6a7ad0dSToby Isaac 
54da43764aSToby Isaac .seealso: DMPlexSetReferenceTree(), DMPlexCreateDefaultReferenceTree()
55d6a7ad0dSToby Isaac @*/
56d6a7ad0dSToby Isaac PetscErrorCode DMPlexGetReferenceTree(DM dm, DM *ref)
57d6a7ad0dSToby Isaac {
58d6a7ad0dSToby Isaac   DM_Plex        *mesh = (DM_Plex *)dm->data;
59d6a7ad0dSToby Isaac 
60d6a7ad0dSToby Isaac   PetscFunctionBegin;
61d6a7ad0dSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
62d6a7ad0dSToby Isaac   PetscValidPointer(ref,2);
63d6a7ad0dSToby Isaac   *ref = mesh->referenceTree;
64d6a7ad0dSToby Isaac   PetscFunctionReturn(0);
65d6a7ad0dSToby Isaac }
66d6a7ad0dSToby Isaac 
67dcbd3bf7SToby Isaac #undef __FUNCT__
68dcbd3bf7SToby Isaac #define __FUNCT__ "DMPlexReferenceTreeGetChildSymmetry_Default"
69dcbd3bf7SToby Isaac static PetscErrorCode DMPlexReferenceTreeGetChildSymmetry_Default(DM dm, PetscInt parent, PetscInt parentOrientA, PetscInt childOrientA, PetscInt childA, PetscInt parentOrientB, PetscInt *childOrientB, PetscInt *childB)
70dcbd3bf7SToby Isaac {
71dcbd3bf7SToby Isaac   PetscInt       coneSize, dStart, dEnd, dim, ABswap, oAvert, oBvert, ABswapVert;
72dcbd3bf7SToby Isaac   PetscErrorCode ierr;
73dcbd3bf7SToby Isaac 
74dcbd3bf7SToby Isaac   PetscFunctionBegin;
75dcbd3bf7SToby Isaac   if (parentOrientA == parentOrientB) {
76dcbd3bf7SToby Isaac     if (childOrientB) *childOrientB = childOrientA;
77dcbd3bf7SToby Isaac     if (childB) *childB = childA;
78dcbd3bf7SToby Isaac     PetscFunctionReturn(0);
79dcbd3bf7SToby Isaac   }
80dcbd3bf7SToby Isaac   for (dim = 0; dim < 3; dim++) {
81dcbd3bf7SToby Isaac     ierr = DMPlexGetDepthStratum(dm,dim,&dStart,&dEnd);CHKERRQ(ierr);
82dcbd3bf7SToby Isaac     if (parent >= dStart && parent <= dEnd) {
83dcbd3bf7SToby Isaac       break;
84dcbd3bf7SToby Isaac     }
85dcbd3bf7SToby Isaac   }
86dcbd3bf7SToby Isaac   if (dim > 2) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot perform child symmetry for %d-cells",dim);
87dcbd3bf7SToby Isaac   if (!dim) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"A vertex has no children");
88dcbd3bf7SToby Isaac   if (childA < dStart || childA >= dEnd) {
89dcbd3bf7SToby Isaac     /* this is a lower-dimensional child: bootstrap */
90dcbd3bf7SToby Isaac     PetscInt size, i, sA = -1, sB, sOrientB, sConeSize;
91dcbd3bf7SToby Isaac     const PetscInt *supp, *coneA, *coneB, *oA, *oB;
92dcbd3bf7SToby Isaac 
93dcbd3bf7SToby Isaac     ierr = DMPlexGetSupportSize(dm,childA,&size);CHKERRQ(ierr);
94dcbd3bf7SToby Isaac     ierr = DMPlexGetSupport(dm,childA,&supp);CHKERRQ(ierr);
95dcbd3bf7SToby Isaac 
96dcbd3bf7SToby Isaac     /* find a point sA in supp(childA) that has the same parent */
97dcbd3bf7SToby Isaac     for (i = 0; i < size; i++) {
98dcbd3bf7SToby Isaac       PetscInt sParent;
99dcbd3bf7SToby Isaac 
100dcbd3bf7SToby Isaac       sA   = supp[i];
101dcbd3bf7SToby Isaac       if (sA == parent) continue;
102dcbd3bf7SToby Isaac       ierr = DMPlexGetTreeParent(dm,sA,&sParent,NULL);CHKERRQ(ierr);
103dcbd3bf7SToby Isaac       if (sParent == parent) {
104dcbd3bf7SToby Isaac         break;
105dcbd3bf7SToby Isaac       }
106dcbd3bf7SToby Isaac     }
107dcbd3bf7SToby Isaac     if (i == size) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"could not find support in children");
108dcbd3bf7SToby Isaac     /* find out which point sB is in an equivalent position to sA under
109dcbd3bf7SToby Isaac      * parentOrientB */
110dcbd3bf7SToby Isaac     ierr = DMPlexReferenceTreeGetChildSymmetry_Default(dm,parent,parentOrientA,0,sA,parentOrientB,&sOrientB,&sB);CHKERRQ(ierr);
111dcbd3bf7SToby Isaac     ierr = DMPlexGetConeSize(dm,sA,&sConeSize);CHKERRQ(ierr);
112dcbd3bf7SToby Isaac     ierr = DMPlexGetCone(dm,sA,&coneA);CHKERRQ(ierr);
113dcbd3bf7SToby Isaac     ierr = DMPlexGetCone(dm,sB,&coneB);CHKERRQ(ierr);
114dcbd3bf7SToby Isaac     ierr = DMPlexGetConeOrientation(dm,sA,&oA);CHKERRQ(ierr);
115dcbd3bf7SToby Isaac     ierr = DMPlexGetConeOrientation(dm,sB,&oB);CHKERRQ(ierr);
116dcbd3bf7SToby Isaac     /* step through the cone of sA in natural order */
117dcbd3bf7SToby Isaac     for (i = 0; i < sConeSize; i++) {
118dcbd3bf7SToby Isaac       if (coneA[i] == childA) {
119dcbd3bf7SToby Isaac         /* if childA is at position i in coneA,
120dcbd3bf7SToby Isaac          * then we want the point that is at sOrientB*i in coneB */
121dcbd3bf7SToby Isaac         PetscInt j = (sOrientB >= 0) ? ((sOrientB + i) % sConeSize) : ((sConeSize -(sOrientB+1) - i) % sConeSize);
122dcbd3bf7SToby Isaac         if (childB) *childB = coneB[j];
123dcbd3bf7SToby Isaac         if (childOrientB) {
124dcbd3bf7SToby Isaac           PetscInt oBtrue;
125dcbd3bf7SToby Isaac 
126dcbd3bf7SToby Isaac           ierr          = DMPlexGetConeSize(dm,childA,&coneSize);CHKERRQ(ierr);
127dcbd3bf7SToby Isaac           /* compose sOrientB and oB[j] */
128dcbd3bf7SToby Isaac           if (coneSize != 0 && coneSize != 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Expected a vertex or an edge");
129dcbd3bf7SToby Isaac           /* we may have to flip an edge */
130dcbd3bf7SToby Isaac           oBtrue        = coneSize ? ((sOrientB >= 0) ? oB[j] : -(oB[j] + 2)) : 0;
131dcbd3bf7SToby Isaac           ABswap        = DihedralSwap(coneSize,oA[i],oBtrue);CHKERRQ(ierr);
132dcbd3bf7SToby Isaac           *childOrientB = DihedralCompose(coneSize,childOrientA,ABswap);
133dcbd3bf7SToby Isaac         }
134dcbd3bf7SToby Isaac         break;
135dcbd3bf7SToby Isaac       }
136dcbd3bf7SToby Isaac     }
137dcbd3bf7SToby Isaac     if (i == sConeSize) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"support cone mismatch");
138dcbd3bf7SToby Isaac     PetscFunctionReturn(0);
139dcbd3bf7SToby Isaac   }
140dcbd3bf7SToby Isaac   /* get the cone size and symmetry swap */
141dcbd3bf7SToby Isaac   ierr   = DMPlexGetConeSize(dm,parent,&coneSize);CHKERRQ(ierr);
142dcbd3bf7SToby Isaac   ABswap = DihedralSwap(coneSize, parentOrientA, parentOrientB);
143dcbd3bf7SToby Isaac   if (dim == 2) {
144dcbd3bf7SToby Isaac     /* orientations refer to cones: we want them to refer to vertices:
145dcbd3bf7SToby Isaac      * if it's a rotation, they are the same, but if the order is reversed, a
146dcbd3bf7SToby Isaac      * permutation that puts side i first does *not* put vertex i first */
147dcbd3bf7SToby Isaac     oAvert     = (parentOrientA >= 0) ? parentOrientA : -((-parentOrientA % coneSize) + 1);
148dcbd3bf7SToby Isaac     oBvert     = (parentOrientB >= 0) ? parentOrientB : -((-parentOrientB % coneSize) + 1);
149dcbd3bf7SToby Isaac     ABswapVert = DihedralSwap(coneSize, oAvert, oBvert);
150dcbd3bf7SToby Isaac   }
151dcbd3bf7SToby Isaac   else {
152dcbd3bf7SToby Isaac     oAvert     = parentOrientA;
153dcbd3bf7SToby Isaac     oBvert     = parentOrientB;
154dcbd3bf7SToby Isaac     ABswapVert = ABswap;
155dcbd3bf7SToby Isaac   }
156dcbd3bf7SToby Isaac   if (childB) {
157dcbd3bf7SToby Isaac     /* assume that each child corresponds to a vertex, in the same order */
158dcbd3bf7SToby Isaac     PetscInt p, posA = -1, numChildren, i;
159dcbd3bf7SToby Isaac     const PetscInt *children;
160dcbd3bf7SToby Isaac 
161dcbd3bf7SToby Isaac     /* count which position the child is in */
162dcbd3bf7SToby Isaac     ierr = DMPlexGetTreeChildren(dm,parent,&numChildren,&children);CHKERRQ(ierr);
163dcbd3bf7SToby Isaac     for (i = 0; i < numChildren; i++) {
164dcbd3bf7SToby Isaac       p = children[i];
165dcbd3bf7SToby Isaac       if (p == childA) {
166dcbd3bf7SToby Isaac         posA = i;
167dcbd3bf7SToby Isaac         break;
168dcbd3bf7SToby Isaac       }
169dcbd3bf7SToby Isaac     }
170dcbd3bf7SToby Isaac     if (posA >= coneSize) {
171dcbd3bf7SToby Isaac       /* this is the triangle in the middle of a uniformly refined triangle: it
172dcbd3bf7SToby Isaac        * is invariant */
173dcbd3bf7SToby Isaac       if (dim != 2 || posA != 3) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Expected a middle triangle, got something else");
174dcbd3bf7SToby Isaac       *childB = childA;
175dcbd3bf7SToby Isaac     }
176dcbd3bf7SToby Isaac     else {
177dcbd3bf7SToby Isaac       /* figure out position B by applying ABswapVert */
178dcbd3bf7SToby Isaac       PetscInt posB;
179dcbd3bf7SToby Isaac 
180dcbd3bf7SToby Isaac       posB = (ABswapVert >= 0) ? ((ABswapVert + posA) % coneSize) : ((coneSize -(ABswapVert + 1) - posA) % coneSize);
181dcbd3bf7SToby Isaac       if (childB) *childB = children[posB];
182dcbd3bf7SToby Isaac     }
183dcbd3bf7SToby Isaac   }
184dcbd3bf7SToby Isaac   if (childOrientB) *childOrientB = DihedralCompose(coneSize,childOrientA,ABswap);
185dcbd3bf7SToby Isaac   PetscFunctionReturn(0);
186dcbd3bf7SToby Isaac }
187dcbd3bf7SToby Isaac 
188dcbd3bf7SToby Isaac #undef __FUNCT__
189dcbd3bf7SToby Isaac #define __FUNCT__ "DMPlexReferenceTreeGetChildSymmetry"
190dcbd3bf7SToby Isaac /*@
191dcbd3bf7SToby Isaac   DMPlexReferenceTreeGetChildSymmetry - Given a reference tree, transform a childid and orientation from one parent frame to another
192dcbd3bf7SToby Isaac 
193dcbd3bf7SToby Isaac   Input Parameters:
194dcbd3bf7SToby Isaac + dm - the reference tree DMPlex object
195dcbd3bf7SToby Isaac . parent - the parent point
196dcbd3bf7SToby Isaac . parentOrientA - the reference orientation for describing the parent
197dcbd3bf7SToby Isaac . childOrientA - the reference orientation for describing the child
198dcbd3bf7SToby Isaac . childA - the reference childID for describing the child
199dcbd3bf7SToby Isaac - parentOrientB - the new orientation for describing the parent
200dcbd3bf7SToby Isaac 
201dcbd3bf7SToby Isaac   Output Parameters:
202dcbd3bf7SToby Isaac + childOrientB - if not NULL, set to the new oreintation for describing the child
203dcbd3bf7SToby Isaac . childB - if not NULL, the new childID for describing the child
204dcbd3bf7SToby Isaac 
205dcbd3bf7SToby Isaac   Level: developer
206dcbd3bf7SToby Isaac 
207dcbd3bf7SToby Isaac .seealso: DMPlexGetReferenceTree(), DMPlexSetReferenceTree(), DMPlexSetTree()
208dcbd3bf7SToby Isaac @*/
209dcbd3bf7SToby Isaac PetscErrorCode DMPlexReferenceTreeGetChildSymmetry(DM dm, PetscInt parent, PetscInt parentOrientA, PetscInt childOrientA, PetscInt childA, PetscInt parentOrientB, PetscInt *childOrientB, PetscInt *childB)
210dcbd3bf7SToby Isaac {
211dcbd3bf7SToby Isaac   DM_Plex        *mesh = (DM_Plex *)dm->data;
212dcbd3bf7SToby Isaac   PetscErrorCode ierr;
213dcbd3bf7SToby Isaac 
214dcbd3bf7SToby Isaac   PetscFunctionBegin;
215dcbd3bf7SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
216dcbd3bf7SToby Isaac   if (!mesh->getchildsymmetry) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"DMPlexReferenceTreeGetChildSymmetry not implemented");
217dcbd3bf7SToby Isaac   ierr = mesh->getchildsymmetry(dm,parent,parentOrientA,childOrientA,childA,parentOrientB,childOrientB,childB);CHKERRQ(ierr);
218dcbd3bf7SToby Isaac   PetscFunctionReturn(0);
219dcbd3bf7SToby Isaac }
220dcbd3bf7SToby Isaac 
221f9f063d4SToby Isaac static PetscErrorCode DMPlexSetTree_Internal(DM,PetscSection,PetscInt*,PetscInt*,PetscBool);
222f9f063d4SToby Isaac 
223da43764aSToby Isaac #undef __FUNCT__
224da43764aSToby Isaac #define __FUNCT__ "DMPlexCreateDefaultReferenceTree"
225da43764aSToby Isaac /*@
226da43764aSToby Isaac   DMPlexCreateDefaultReferenceTree - create a reference tree for isotropic hierarchical mesh refinement.
227da43764aSToby Isaac 
228da43764aSToby Isaac   Collective on comm
229da43764aSToby Isaac 
230da43764aSToby Isaac   Input Parameters:
231da43764aSToby Isaac + comm    - the MPI communicator
232da43764aSToby Isaac . dim     - the spatial dimension
233da43764aSToby Isaac - simplex - Flag for simplex, otherwise use a tensor-product cell
234da43764aSToby Isaac 
235da43764aSToby Isaac   Output Parameters:
236da43764aSToby Isaac . ref     - the reference tree DMPlex object
237da43764aSToby Isaac 
2380b7167a0SToby Isaac   Level: intermediate
239da43764aSToby Isaac 
240da43764aSToby Isaac .keywords: reference cell
241da43764aSToby Isaac .seealso: DMPlexSetReferenceTree(), DMPlexGetReferenceTree()
242da43764aSToby Isaac @*/
243da43764aSToby Isaac PetscErrorCode DMPlexCreateDefaultReferenceTree(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *ref)
244da43764aSToby Isaac {
245dcbd3bf7SToby Isaac   DM_Plex       *mesh;
246da43764aSToby Isaac   DM             K, Kref;
24710f7e118SToby Isaac   PetscInt       p, pStart, pEnd, pRefStart, pRefEnd, d, offset, parentSize, *parents, *childIDs;
248da43764aSToby Isaac   PetscInt      *permvals, *unionCones, *coneSizes, *unionOrientations, numUnionPoints, *numDimPoints, numCones, numVerts;
249da43764aSToby Isaac   DMLabel        identity, identityRef;
25010f7e118SToby Isaac   PetscSection   unionSection, unionConeSection, parentSection;
251da43764aSToby Isaac   PetscScalar   *unionCoords;
252da43764aSToby Isaac   IS             perm;
253da43764aSToby Isaac   PetscErrorCode ierr;
254da43764aSToby Isaac 
255da43764aSToby Isaac   PetscFunctionBegin;
256da43764aSToby Isaac   /* create a reference element */
257da43764aSToby Isaac   ierr = DMPlexCreateReferenceCell(comm, dim, simplex, &K);CHKERRQ(ierr);
258da43764aSToby Isaac   ierr = DMPlexCreateLabel(K, "identity");CHKERRQ(ierr);
259da43764aSToby Isaac   ierr = DMPlexGetLabel(K, "identity", &identity);CHKERRQ(ierr);
260da43764aSToby Isaac   ierr = DMPlexGetChart(K, &pStart, &pEnd);CHKERRQ(ierr);
261da43764aSToby Isaac   for (p = pStart; p < pEnd; p++) {
262da43764aSToby Isaac     ierr = DMLabelSetValue(identity, p, p);CHKERRQ(ierr);
263da43764aSToby Isaac   }
264da43764aSToby Isaac   /* refine it */
265da43764aSToby Isaac   ierr = DMRefine(K,comm,&Kref);CHKERRQ(ierr);
266da43764aSToby Isaac 
267da43764aSToby Isaac   /* the reference tree is the union of these two, without duplicating
268da43764aSToby Isaac    * points that appear in both */
269da43764aSToby Isaac   ierr = DMPlexGetLabel(Kref, "identity", &identityRef);CHKERRQ(ierr);
270da43764aSToby Isaac   ierr = DMPlexGetChart(Kref, &pRefStart, &pRefEnd);CHKERRQ(ierr);
271da43764aSToby Isaac   ierr = PetscSectionCreate(comm, &unionSection);CHKERRQ(ierr);
272da43764aSToby Isaac   ierr = PetscSectionSetChart(unionSection, 0, (pEnd - pStart) + (pRefEnd - pRefStart));CHKERRQ(ierr);
273da43764aSToby Isaac   /* count points that will go in the union */
274da43764aSToby Isaac   for (p = pStart; p < pEnd; p++) {
275da43764aSToby Isaac     ierr = PetscSectionSetDof(unionSection, p - pStart, 1);CHKERRQ(ierr);
276da43764aSToby Isaac   }
277da43764aSToby Isaac   for (p = pRefStart; p < pRefEnd; p++) {
278da43764aSToby Isaac     PetscInt q, qSize;
279da43764aSToby Isaac     ierr = DMLabelGetValue(identityRef, p, &q);CHKERRQ(ierr);
280da43764aSToby Isaac     ierr = DMLabelGetStratumSize(identityRef, q, &qSize);CHKERRQ(ierr);
281da43764aSToby Isaac     if (qSize > 1) {
282da43764aSToby Isaac       ierr = PetscSectionSetDof(unionSection, p - pRefStart + (pEnd - pStart), 1);CHKERRQ(ierr);
283da43764aSToby Isaac     }
284da43764aSToby Isaac   }
285da43764aSToby Isaac   ierr = PetscMalloc1((pEnd - pStart) + (pRefEnd - pRefStart),&permvals);CHKERRQ(ierr);
286da43764aSToby Isaac   offset = 0;
287da43764aSToby Isaac   /* stratify points in the union by topological dimension */
288da43764aSToby Isaac   for (d = 0; d <= dim; d++) {
289da43764aSToby Isaac     PetscInt cStart, cEnd, c;
290da43764aSToby Isaac 
291da43764aSToby Isaac     ierr = DMPlexGetHeightStratum(K, d, &cStart, &cEnd);CHKERRQ(ierr);
292da43764aSToby Isaac     for (c = cStart; c < cEnd; c++) {
293da43764aSToby Isaac       permvals[offset++] = c;
294da43764aSToby Isaac     }
295da43764aSToby Isaac 
296da43764aSToby Isaac     ierr = DMPlexGetHeightStratum(Kref, d, &cStart, &cEnd);CHKERRQ(ierr);
297da43764aSToby Isaac     for (c = cStart; c < cEnd; c++) {
298da43764aSToby Isaac       permvals[offset++] = c + (pEnd - pStart);
299da43764aSToby Isaac     }
300da43764aSToby Isaac   }
301da43764aSToby Isaac   ierr = ISCreateGeneral(comm, (pEnd - pStart) + (pRefEnd - pRefStart), permvals, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr);
302da43764aSToby Isaac   ierr = PetscSectionSetPermutation(unionSection,perm);CHKERRQ(ierr);
303da43764aSToby Isaac   ierr = PetscSectionSetUp(unionSection);CHKERRQ(ierr);
304da43764aSToby Isaac   ierr = PetscSectionGetStorageSize(unionSection,&numUnionPoints);CHKERRQ(ierr);
305da43764aSToby Isaac   ierr = PetscMalloc2(numUnionPoints,&coneSizes,dim+1,&numDimPoints);CHKERRQ(ierr);
306da43764aSToby Isaac   /* count dimension points */
307da43764aSToby Isaac   for (d = 0; d <= dim; d++) {
308da43764aSToby Isaac     PetscInt cStart, cOff, cOff2;
309da43764aSToby Isaac     ierr = DMPlexGetHeightStratum(K,d,&cStart,NULL);CHKERRQ(ierr);
310da43764aSToby Isaac     ierr = PetscSectionGetOffset(unionSection,cStart-pStart,&cOff);CHKERRQ(ierr);
311da43764aSToby Isaac     if (d < dim) {
312da43764aSToby Isaac       ierr = DMPlexGetHeightStratum(K,d+1,&cStart,NULL);CHKERRQ(ierr);
313da43764aSToby Isaac       ierr = PetscSectionGetOffset(unionSection,cStart-pStart,&cOff2);CHKERRQ(ierr);
314da43764aSToby Isaac     }
315da43764aSToby Isaac     else {
316da43764aSToby Isaac       cOff2 = numUnionPoints;
317da43764aSToby Isaac     }
318da43764aSToby Isaac     numDimPoints[dim - d] = cOff2 - cOff;
319da43764aSToby Isaac   }
320da43764aSToby Isaac   ierr = PetscSectionCreate(comm, &unionConeSection);CHKERRQ(ierr);
321da43764aSToby Isaac   ierr = PetscSectionSetChart(unionConeSection, 0, numUnionPoints);CHKERRQ(ierr);
322da43764aSToby Isaac   /* count the cones in the union */
323da43764aSToby Isaac   for (p = pStart; p < pEnd; p++) {
324da43764aSToby Isaac     PetscInt dof, uOff;
325da43764aSToby Isaac 
326da43764aSToby Isaac     ierr = DMPlexGetConeSize(K, p, &dof);CHKERRQ(ierr);
327da43764aSToby Isaac     ierr = PetscSectionGetOffset(unionSection, p - pStart,&uOff);CHKERRQ(ierr);
328da43764aSToby Isaac     ierr = PetscSectionSetDof(unionConeSection, uOff, dof);CHKERRQ(ierr);
329da43764aSToby Isaac     coneSizes[uOff] = dof;
330da43764aSToby Isaac   }
331da43764aSToby Isaac   for (p = pRefStart; p < pRefEnd; p++) {
332da43764aSToby Isaac     PetscInt dof, uDof, uOff;
333da43764aSToby Isaac 
334da43764aSToby Isaac     ierr = DMPlexGetConeSize(Kref, p, &dof);CHKERRQ(ierr);
335da43764aSToby Isaac     ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr);
336da43764aSToby Isaac     ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr);
337da43764aSToby Isaac     if (uDof) {
338da43764aSToby Isaac       ierr = PetscSectionSetDof(unionConeSection, uOff, dof);CHKERRQ(ierr);
339da43764aSToby Isaac       coneSizes[uOff] = dof;
340da43764aSToby Isaac     }
341da43764aSToby Isaac   }
342da43764aSToby Isaac   ierr = PetscSectionSetUp(unionConeSection);CHKERRQ(ierr);
343da43764aSToby Isaac   ierr = PetscSectionGetStorageSize(unionConeSection,&numCones);CHKERRQ(ierr);
344da43764aSToby Isaac   ierr = PetscMalloc2(numCones,&unionCones,numCones,&unionOrientations);CHKERRQ(ierr);
345da43764aSToby Isaac   /* write the cones in the union */
346da43764aSToby Isaac   for (p = pStart; p < pEnd; p++) {
347da43764aSToby Isaac     PetscInt dof, uOff, c, cOff;
348da43764aSToby Isaac     const PetscInt *cone, *orientation;
349da43764aSToby Isaac 
350da43764aSToby Isaac     ierr = DMPlexGetConeSize(K, p, &dof);CHKERRQ(ierr);
351da43764aSToby Isaac     ierr = DMPlexGetCone(K, p, &cone);CHKERRQ(ierr);
352da43764aSToby Isaac     ierr = DMPlexGetConeOrientation(K, p, &orientation);CHKERRQ(ierr);
353da43764aSToby Isaac     ierr = PetscSectionGetOffset(unionSection, p - pStart,&uOff);CHKERRQ(ierr);
354da43764aSToby Isaac     ierr = PetscSectionGetOffset(unionConeSection,uOff,&cOff);CHKERRQ(ierr);
355da43764aSToby Isaac     for (c = 0; c < dof; c++) {
356da43764aSToby Isaac       PetscInt e, eOff;
357da43764aSToby Isaac       e                           = cone[c];
358da43764aSToby Isaac       ierr                        = PetscSectionGetOffset(unionSection, e - pStart, &eOff);CHKERRQ(ierr);
359da43764aSToby Isaac       unionCones[cOff + c]        = eOff;
360da43764aSToby Isaac       unionOrientations[cOff + c] = orientation[c];
361da43764aSToby Isaac     }
362da43764aSToby Isaac   }
363da43764aSToby Isaac   for (p = pRefStart; p < pRefEnd; p++) {
364da43764aSToby Isaac     PetscInt dof, uDof, uOff, c, cOff;
365da43764aSToby Isaac     const PetscInt *cone, *orientation;
366da43764aSToby Isaac 
367da43764aSToby Isaac     ierr = DMPlexGetConeSize(Kref, p, &dof);CHKERRQ(ierr);
368da43764aSToby Isaac     ierr = DMPlexGetCone(Kref, p, &cone);CHKERRQ(ierr);
369da43764aSToby Isaac     ierr = DMPlexGetConeOrientation(Kref, p, &orientation);CHKERRQ(ierr);
370da43764aSToby Isaac     ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr);
371da43764aSToby Isaac     ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr);
372da43764aSToby Isaac     if (uDof) {
373da43764aSToby Isaac       ierr = PetscSectionGetOffset(unionConeSection,uOff,&cOff);CHKERRQ(ierr);
374da43764aSToby Isaac       for (c = 0; c < dof; c++) {
375da43764aSToby Isaac         PetscInt e, eOff, eDof;
376da43764aSToby Isaac 
377da43764aSToby Isaac         e    = cone[c];
378da43764aSToby Isaac         ierr = PetscSectionGetDof(unionSection, e - pRefStart + (pEnd - pStart),&eDof);CHKERRQ(ierr);
379da43764aSToby Isaac         if (eDof) {
380da43764aSToby Isaac           ierr = PetscSectionGetOffset(unionSection, e - pRefStart + (pEnd - pStart), &eOff);CHKERRQ(ierr);
381da43764aSToby Isaac         }
382da43764aSToby Isaac         else {
383da43764aSToby Isaac           ierr = DMLabelGetValue(identityRef, e, &e);CHKERRQ(ierr);
384da43764aSToby Isaac           ierr = PetscSectionGetOffset(unionSection, e - pStart, &eOff);CHKERRQ(ierr);
385da43764aSToby Isaac         }
386da43764aSToby Isaac         unionCones[cOff + c]        = eOff;
387da43764aSToby Isaac         unionOrientations[cOff + c] = orientation[c];
388da43764aSToby Isaac       }
389da43764aSToby Isaac     }
390da43764aSToby Isaac   }
391da43764aSToby Isaac   /* get the coordinates */
392da43764aSToby Isaac   {
393da43764aSToby Isaac     PetscInt vStart, vEnd, vRefStart, vRefEnd, v, vDof, vOff;
394da43764aSToby Isaac     PetscSection KcoordsSec, KrefCoordsSec;
395da43764aSToby Isaac     Vec      KcoordsVec, KrefCoordsVec;
396da43764aSToby Isaac     PetscScalar *Kcoords;
397da43764aSToby Isaac 
398da43764aSToby Isaac     DMGetCoordinateSection(K, &KcoordsSec);CHKERRQ(ierr);
399da43764aSToby Isaac     DMGetCoordinatesLocal(K, &KcoordsVec);CHKERRQ(ierr);
400da43764aSToby Isaac     DMGetCoordinateSection(Kref, &KrefCoordsSec);CHKERRQ(ierr);
401da43764aSToby Isaac     DMGetCoordinatesLocal(Kref, &KrefCoordsVec);CHKERRQ(ierr);
402da43764aSToby Isaac 
403da43764aSToby Isaac     numVerts = numDimPoints[0];
404da43764aSToby Isaac     ierr     = PetscMalloc1(numVerts * dim,&unionCoords);CHKERRQ(ierr);
405da43764aSToby Isaac     ierr     = DMPlexGetDepthStratum(K,0,&vStart,&vEnd);CHKERRQ(ierr);
406da43764aSToby Isaac 
407da43764aSToby Isaac     offset = 0;
408da43764aSToby Isaac     for (v = vStart; v < vEnd; v++) {
409da43764aSToby Isaac       ierr = PetscSectionGetOffset(unionSection,v - pStart,&vOff);CHKERRQ(ierr);
410da43764aSToby Isaac       ierr = VecGetValuesSection(KcoordsVec, KcoordsSec, v, &Kcoords);CHKERRQ(ierr);
411da43764aSToby Isaac       for (d = 0; d < dim; d++) {
412da43764aSToby Isaac         unionCoords[offset * dim + d] = Kcoords[d];
413da43764aSToby Isaac       }
414da43764aSToby Isaac       offset++;
415da43764aSToby Isaac     }
416da43764aSToby Isaac     ierr = DMPlexGetDepthStratum(Kref,0,&vRefStart,&vRefEnd);CHKERRQ(ierr);
417da43764aSToby Isaac     for (v = vRefStart; v < vRefEnd; v++) {
418da43764aSToby Isaac       ierr = PetscSectionGetDof(unionSection,v - pRefStart + (pEnd - pStart),&vDof);CHKERRQ(ierr);
419da43764aSToby Isaac       ierr = PetscSectionGetOffset(unionSection,v - pRefStart + (pEnd - pStart),&vOff);CHKERRQ(ierr);
420da43764aSToby Isaac       ierr = VecGetValuesSection(KrefCoordsVec, KrefCoordsSec, v, &Kcoords);CHKERRQ(ierr);
421da43764aSToby Isaac       if (vDof) {
422da43764aSToby Isaac         for (d = 0; d < dim; d++) {
423da43764aSToby Isaac           unionCoords[offset * dim + d] = Kcoords[d];
424da43764aSToby Isaac         }
425da43764aSToby Isaac         offset++;
426da43764aSToby Isaac       }
427da43764aSToby Isaac     }
428da43764aSToby Isaac   }
429da43764aSToby Isaac   ierr = DMCreate(comm,ref);CHKERRQ(ierr);
430da43764aSToby Isaac   ierr = DMSetType(*ref,DMPLEX);CHKERRQ(ierr);
431da43764aSToby Isaac   ierr = DMPlexSetDimension(*ref,dim);CHKERRQ(ierr);
432da43764aSToby Isaac   ierr = DMPlexCreateFromDAG(*ref,dim,numDimPoints,coneSizes,unionCones,unionOrientations,unionCoords);CHKERRQ(ierr);
43310f7e118SToby Isaac   /* set the tree */
43410f7e118SToby Isaac   ierr = PetscSectionCreate(comm,&parentSection);CHKERRQ(ierr);
43510f7e118SToby Isaac   ierr = PetscSectionSetChart(parentSection,0,numUnionPoints);CHKERRQ(ierr);
43610f7e118SToby Isaac   for (p = pRefStart; p < pRefEnd; p++) {
43710f7e118SToby Isaac     PetscInt uDof, uOff;
43810f7e118SToby Isaac 
43910f7e118SToby Isaac     ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr);
44010f7e118SToby Isaac     ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr);
44110f7e118SToby Isaac     if (uDof) {
44210f7e118SToby Isaac       PetscSectionSetDof(parentSection,uOff,1);CHKERRQ(ierr);
44310f7e118SToby Isaac     }
44410f7e118SToby Isaac   }
44510f7e118SToby Isaac   ierr = PetscSectionSetUp(parentSection);CHKERRQ(ierr);
44610f7e118SToby Isaac   ierr = PetscSectionGetStorageSize(parentSection,&parentSize);CHKERRQ(ierr);
44710f7e118SToby Isaac   ierr = PetscMalloc2(parentSize,&parents,parentSize,&childIDs);CHKERRQ(ierr);
44810f7e118SToby Isaac   for (p = pRefStart; p < pRefEnd; p++) {
44910f7e118SToby Isaac     PetscInt uDof, uOff;
45010f7e118SToby Isaac 
45110f7e118SToby Isaac     ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr);
45210f7e118SToby Isaac     ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr);
45310f7e118SToby Isaac     if (uDof) {
45410f7e118SToby Isaac       PetscInt pOff, parent, parentU;
45510f7e118SToby Isaac       PetscSectionGetOffset(parentSection,uOff,&pOff);CHKERRQ(ierr);
45610f7e118SToby Isaac       DMLabelGetValue(identityRef,p,&parent);CHKERRQ(ierr);
45710f7e118SToby Isaac       ierr = PetscSectionGetOffset(unionSection, parent - pStart,&parentU);CHKERRQ(ierr);
45810f7e118SToby Isaac       parents[pOff] = parentU;
45910f7e118SToby Isaac       childIDs[pOff] = uOff;
46010f7e118SToby Isaac     }
46110f7e118SToby Isaac   }
462f9f063d4SToby Isaac   ierr = DMPlexSetTree_Internal(*ref,parentSection,parents,childIDs,PETSC_TRUE);CHKERRQ(ierr);
463dcbd3bf7SToby Isaac   mesh = (DM_Plex *) (*ref)->data;
464dcbd3bf7SToby Isaac   mesh->getchildsymmetry = DMPlexReferenceTreeGetChildSymmetry_Default;
46510f7e118SToby Isaac   ierr = PetscSectionDestroy(&parentSection);CHKERRQ(ierr);
46610f7e118SToby Isaac   ierr = PetscFree2(parents,childIDs);CHKERRQ(ierr);
46710f7e118SToby Isaac 
468da43764aSToby Isaac   /* clean up */
469da43764aSToby Isaac   ierr = PetscSectionDestroy(&unionSection);CHKERRQ(ierr);
470da43764aSToby Isaac   ierr = PetscSectionDestroy(&unionConeSection);CHKERRQ(ierr);
471da43764aSToby Isaac   ierr = ISDestroy(&perm);CHKERRQ(ierr);
472da43764aSToby Isaac   ierr = PetscFree(unionCoords);CHKERRQ(ierr);
473da43764aSToby Isaac   ierr = PetscFree2(unionCones,unionOrientations);CHKERRQ(ierr);
474da43764aSToby Isaac   ierr = PetscFree2(coneSizes,numDimPoints);CHKERRQ(ierr);
475da43764aSToby Isaac   ierr = DMDestroy(&K);CHKERRQ(ierr);
476da43764aSToby Isaac   ierr = DMDestroy(&Kref);CHKERRQ(ierr);
477da43764aSToby Isaac   PetscFunctionReturn(0);
478da43764aSToby Isaac }
479da43764aSToby Isaac 
480f9f063d4SToby Isaac 
481d961a43aSToby Isaac #undef __FUNCT__
482878b19aaSToby Isaac #define __FUNCT__ "DMPlexTreeSymmetrize"
483878b19aaSToby Isaac static PetscErrorCode DMPlexTreeSymmetrize(DM dm)
484878b19aaSToby Isaac {
485878b19aaSToby Isaac   DM_Plex        *mesh = (DM_Plex *)dm->data;
486878b19aaSToby Isaac   PetscSection   childSec, pSec;
487878b19aaSToby Isaac   PetscInt       p, pSize, cSize, parMax = PETSC_MIN_INT, parMin = PETSC_MAX_INT;
488878b19aaSToby Isaac   PetscInt       *offsets, *children, pStart, pEnd;
489878b19aaSToby Isaac   PetscErrorCode ierr;
490878b19aaSToby Isaac 
491878b19aaSToby Isaac   PetscFunctionBegin;
492878b19aaSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
493878b19aaSToby Isaac   ierr = PetscSectionDestroy(&mesh->childSection);CHKERRQ(ierr);
494878b19aaSToby Isaac   ierr = PetscFree(mesh->children);CHKERRQ(ierr);
495878b19aaSToby Isaac   pSec = mesh->parentSection;
496878b19aaSToby Isaac   if (!pSec) PetscFunctionReturn(0);
497878b19aaSToby Isaac   ierr = PetscSectionGetStorageSize(pSec,&pSize);CHKERRQ(ierr);
498878b19aaSToby Isaac   for (p = 0; p < pSize; p++) {
499878b19aaSToby Isaac     PetscInt par = mesh->parents[p];
500878b19aaSToby Isaac 
501878b19aaSToby Isaac     parMax = PetscMax(parMax,par+1);
502878b19aaSToby Isaac     parMin = PetscMin(parMin,par);
503878b19aaSToby Isaac   }
504878b19aaSToby Isaac   if (parMin > parMax) {
505878b19aaSToby Isaac     parMin = -1;
506878b19aaSToby Isaac     parMax = -1;
507878b19aaSToby Isaac   }
508878b19aaSToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)pSec),&childSec);CHKERRQ(ierr);
509878b19aaSToby Isaac   ierr = PetscSectionSetChart(childSec,parMin,parMax);CHKERRQ(ierr);
510878b19aaSToby Isaac   for (p = 0; p < pSize; p++) {
511878b19aaSToby Isaac     PetscInt par = mesh->parents[p];
512878b19aaSToby Isaac 
513878b19aaSToby Isaac     ierr = PetscSectionAddDof(childSec,par,1);CHKERRQ(ierr);
514878b19aaSToby Isaac   }
515878b19aaSToby Isaac   ierr = PetscSectionSetUp(childSec);CHKERRQ(ierr);
516878b19aaSToby Isaac   ierr = PetscSectionGetStorageSize(childSec,&cSize);CHKERRQ(ierr);
517878b19aaSToby Isaac   ierr = PetscMalloc1(cSize,&children);CHKERRQ(ierr);
518878b19aaSToby Isaac   ierr = PetscCalloc1(parMax-parMin,&offsets);CHKERRQ(ierr);
519878b19aaSToby Isaac   ierr = PetscSectionGetChart(pSec,&pStart,&pEnd);CHKERRQ(ierr);
520878b19aaSToby Isaac   for (p = pStart; p < pEnd; p++) {
521878b19aaSToby Isaac     PetscInt dof, off, i;
522878b19aaSToby Isaac 
523878b19aaSToby Isaac     ierr = PetscSectionGetDof(pSec,p,&dof);CHKERRQ(ierr);
524878b19aaSToby Isaac     ierr = PetscSectionGetOffset(pSec,p,&off);CHKERRQ(ierr);
525878b19aaSToby Isaac     for (i = 0; i < dof; i++) {
526878b19aaSToby Isaac       PetscInt par = mesh->parents[off + i], cOff;
527878b19aaSToby Isaac 
528878b19aaSToby Isaac       ierr = PetscSectionGetOffset(childSec,par,&cOff);CHKERRQ(ierr);
529878b19aaSToby Isaac       children[cOff + offsets[par-parMin]++] = p;
530878b19aaSToby Isaac     }
531878b19aaSToby Isaac   }
532878b19aaSToby Isaac   mesh->childSection = childSec;
533878b19aaSToby Isaac   mesh->children = children;
534878b19aaSToby Isaac   ierr = PetscFree(offsets);CHKERRQ(ierr);
535878b19aaSToby Isaac   PetscFunctionReturn(0);
536878b19aaSToby Isaac }
537878b19aaSToby Isaac 
538878b19aaSToby Isaac #undef __FUNCT__
5396dd5a8c8SToby Isaac #define __FUNCT__ "AnchorsFlatten"
5406dd5a8c8SToby Isaac static PetscErrorCode AnchorsFlatten (PetscSection section, IS is, PetscSection *sectionNew, IS *isNew)
5416dd5a8c8SToby Isaac {
5426dd5a8c8SToby Isaac   PetscInt       pStart, pEnd, size, sizeNew, i, p, *valsNew = NULL;
5436dd5a8c8SToby Isaac   const PetscInt *vals;
5446dd5a8c8SToby Isaac   PetscSection   secNew;
5456dd5a8c8SToby Isaac   PetscBool      anyNew, globalAnyNew;
5466dd5a8c8SToby Isaac   PetscBool      compress;
5476dd5a8c8SToby Isaac   PetscErrorCode ierr;
5486dd5a8c8SToby Isaac 
5496dd5a8c8SToby Isaac   PetscFunctionBegin;
5506dd5a8c8SToby Isaac   ierr = PetscSectionGetChart(section,&pStart,&pEnd);CHKERRQ(ierr);
5516dd5a8c8SToby Isaac   ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr);
5526dd5a8c8SToby Isaac   ierr = ISGetIndices(is,&vals);CHKERRQ(ierr);
5536dd5a8c8SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)section),&secNew);CHKERRQ(ierr);
5546dd5a8c8SToby Isaac   ierr = PetscSectionSetChart(secNew,pStart,pEnd);CHKERRQ(ierr);
5556dd5a8c8SToby Isaac   for (i = 0; i < size; i++) {
5566dd5a8c8SToby Isaac     PetscInt dof;
5576dd5a8c8SToby Isaac 
5586dd5a8c8SToby Isaac     p = vals[i];
5596dd5a8c8SToby Isaac     if (p < pStart || p >= pEnd) continue;
5606dd5a8c8SToby Isaac     ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
5616dd5a8c8SToby Isaac     if (dof) break;
5626dd5a8c8SToby Isaac   }
5636dd5a8c8SToby Isaac   if (i == size) {
5646dd5a8c8SToby Isaac     ierr     = PetscSectionSetUp(secNew);CHKERRQ(ierr);
5656dd5a8c8SToby Isaac     anyNew   = PETSC_FALSE;
5666dd5a8c8SToby Isaac     compress = PETSC_FALSE;
5676dd5a8c8SToby Isaac     sizeNew  = 0;
5686dd5a8c8SToby Isaac   }
5696dd5a8c8SToby Isaac   else {
5706dd5a8c8SToby Isaac     anyNew = PETSC_TRUE;
5716dd5a8c8SToby Isaac     for (p = pStart; p < pEnd; p++) {
5726dd5a8c8SToby Isaac       PetscInt dof, off;
5736dd5a8c8SToby Isaac 
5746dd5a8c8SToby Isaac       ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
5756dd5a8c8SToby Isaac       ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
5766dd5a8c8SToby Isaac       for (i = 0; i < dof; i++) {
5776dd5a8c8SToby Isaac         PetscInt q = vals[off + i], qDof = 0;
5786dd5a8c8SToby Isaac 
5796dd5a8c8SToby Isaac         if (q >= pStart && q < pEnd) {
5806dd5a8c8SToby Isaac           ierr = PetscSectionGetDof(section, q, &qDof);CHKERRQ(ierr);
5816dd5a8c8SToby Isaac         }
5826dd5a8c8SToby Isaac         if (qDof) {
5836dd5a8c8SToby Isaac           ierr = PetscSectionAddDof(secNew, p, qDof);CHKERRQ(ierr);
5846dd5a8c8SToby Isaac         }
5856dd5a8c8SToby Isaac         else {
5866dd5a8c8SToby Isaac           ierr = PetscSectionAddDof(secNew, p, 1);CHKERRQ(ierr);
5876dd5a8c8SToby Isaac         }
5886dd5a8c8SToby Isaac       }
5896dd5a8c8SToby Isaac     }
5906dd5a8c8SToby Isaac     ierr = PetscSectionSetUp(secNew);CHKERRQ(ierr);
5916dd5a8c8SToby Isaac     ierr = PetscSectionGetStorageSize(secNew,&sizeNew);CHKERRQ(ierr);
5926dd5a8c8SToby Isaac     ierr = PetscMalloc1(sizeNew,&valsNew);CHKERRQ(ierr);
5936dd5a8c8SToby Isaac     compress = PETSC_FALSE;
5946dd5a8c8SToby Isaac     for (p = pStart; p < pEnd; p++) {
5956dd5a8c8SToby Isaac       PetscInt dof, off, count, offNew, dofNew;
5966dd5a8c8SToby Isaac 
5976dd5a8c8SToby Isaac       ierr  = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr);
5986dd5a8c8SToby Isaac       ierr  = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr);
5996dd5a8c8SToby Isaac       ierr  = PetscSectionGetDof(secNew, p, &dofNew);CHKERRQ(ierr);
6006dd5a8c8SToby Isaac       ierr  = PetscSectionGetOffset(secNew, p, &offNew);CHKERRQ(ierr);
6016dd5a8c8SToby Isaac       count = 0;
6026dd5a8c8SToby Isaac       for (i = 0; i < dof; i++) {
6036dd5a8c8SToby Isaac         PetscInt q = vals[off + i], qDof = 0, qOff = 0, j;
6046dd5a8c8SToby Isaac 
6056dd5a8c8SToby Isaac         if (q >= pStart && q < pEnd) {
6066dd5a8c8SToby Isaac           ierr = PetscSectionGetDof(section, q, &qDof);CHKERRQ(ierr);
6076dd5a8c8SToby Isaac           ierr = PetscSectionGetOffset(section, q, &qOff);CHKERRQ(ierr);
6086dd5a8c8SToby Isaac         }
6096dd5a8c8SToby Isaac         if (qDof) {
6106dd5a8c8SToby Isaac           PetscInt oldCount = count;
6116dd5a8c8SToby Isaac 
6126dd5a8c8SToby Isaac           for (j = 0; j < qDof; j++) {
6136dd5a8c8SToby Isaac             PetscInt k, r = vals[qOff + j];
6146dd5a8c8SToby Isaac 
6156dd5a8c8SToby Isaac             for (k = 0; k < oldCount; k++) {
6166dd5a8c8SToby Isaac               if (valsNew[offNew + k] == r) {
6176dd5a8c8SToby Isaac                 break;
6186dd5a8c8SToby Isaac               }
6196dd5a8c8SToby Isaac             }
6206dd5a8c8SToby Isaac             if (k == oldCount) {
6216dd5a8c8SToby Isaac               valsNew[offNew + count++] = r;
6226dd5a8c8SToby Isaac             }
6236dd5a8c8SToby Isaac           }
6246dd5a8c8SToby Isaac         }
6256dd5a8c8SToby Isaac         else {
6266dd5a8c8SToby Isaac           PetscInt k, oldCount = count;
6276dd5a8c8SToby Isaac 
6286dd5a8c8SToby Isaac           for (k = 0; k < oldCount; k++) {
6296dd5a8c8SToby Isaac             if (valsNew[offNew + k] == q) {
6306dd5a8c8SToby Isaac               break;
6316dd5a8c8SToby Isaac             }
6326dd5a8c8SToby Isaac           }
6336dd5a8c8SToby Isaac           if (k == oldCount) {
6346dd5a8c8SToby Isaac             valsNew[offNew + count++] = q;
6356dd5a8c8SToby Isaac           }
6366dd5a8c8SToby Isaac         }
6376dd5a8c8SToby Isaac       }
6386dd5a8c8SToby Isaac       if (count < dofNew) {
6396dd5a8c8SToby Isaac         ierr = PetscSectionSetDof(secNew, p, count);CHKERRQ(ierr);
6406dd5a8c8SToby Isaac         compress = PETSC_TRUE;
6416dd5a8c8SToby Isaac       }
6426dd5a8c8SToby Isaac     }
6436dd5a8c8SToby Isaac   }
6446dd5a8c8SToby Isaac   ierr = ISRestoreIndices(is,&vals);CHKERRQ(ierr);
6456dd5a8c8SToby Isaac   ierr = MPI_Allreduce(&anyNew,&globalAnyNew,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)secNew));CHKERRQ(ierr);
6466dd5a8c8SToby Isaac   if (!globalAnyNew) {
6476dd5a8c8SToby Isaac     ierr = PetscSectionDestroy(&secNew);CHKERRQ(ierr);
6486dd5a8c8SToby Isaac     *sectionNew = NULL;
6496dd5a8c8SToby Isaac     *isNew = NULL;
6506dd5a8c8SToby Isaac   }
6516dd5a8c8SToby Isaac   else {
6526dd5a8c8SToby Isaac     PetscBool globalCompress;
6536dd5a8c8SToby Isaac 
6546dd5a8c8SToby Isaac     ierr = MPI_Allreduce(&compress,&globalCompress,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)secNew));CHKERRQ(ierr);
6556dd5a8c8SToby Isaac     if (compress) {
6566dd5a8c8SToby Isaac       PetscSection secComp;
6576dd5a8c8SToby Isaac       PetscInt *valsComp = NULL;
6586dd5a8c8SToby Isaac 
6596dd5a8c8SToby Isaac       ierr = PetscSectionCreate(PetscObjectComm((PetscObject)section),&secComp);CHKERRQ(ierr);
6606dd5a8c8SToby Isaac       ierr = PetscSectionSetChart(secComp,pStart,pEnd);CHKERRQ(ierr);
6616dd5a8c8SToby Isaac       for (p = pStart; p < pEnd; p++) {
6626dd5a8c8SToby Isaac         PetscInt dof;
6636dd5a8c8SToby Isaac 
6646dd5a8c8SToby Isaac         ierr = PetscSectionGetDof(secNew, p, &dof);CHKERRQ(ierr);
6656dd5a8c8SToby Isaac         ierr = PetscSectionSetDof(secComp, p, dof);CHKERRQ(ierr);
6666dd5a8c8SToby Isaac       }
6676dd5a8c8SToby Isaac       ierr = PetscSectionSetUp(secComp);CHKERRQ(ierr);
6686dd5a8c8SToby Isaac       ierr = PetscSectionGetStorageSize(secComp,&sizeNew);CHKERRQ(ierr);
6696dd5a8c8SToby Isaac       ierr = PetscMalloc1(sizeNew,&valsComp);CHKERRQ(ierr);
6706dd5a8c8SToby Isaac       for (p = pStart; p < pEnd; p++) {
6716dd5a8c8SToby Isaac         PetscInt dof, off, offNew, j;
6726dd5a8c8SToby Isaac 
6736dd5a8c8SToby Isaac         ierr = PetscSectionGetDof(secNew, p, &dof);CHKERRQ(ierr);
6746dd5a8c8SToby Isaac         ierr = PetscSectionGetOffset(secNew, p, &off);CHKERRQ(ierr);
6756dd5a8c8SToby Isaac         ierr = PetscSectionGetOffset(secComp, p, &offNew);CHKERRQ(ierr);
6766dd5a8c8SToby Isaac         for (j = 0; j < dof; j++) {
6776dd5a8c8SToby Isaac           valsComp[offNew + j] = valsNew[off + j];
6786dd5a8c8SToby Isaac         }
6796dd5a8c8SToby Isaac       }
6806dd5a8c8SToby Isaac       ierr    = PetscSectionDestroy(&secNew);CHKERRQ(ierr);
6816dd5a8c8SToby Isaac       secNew  = secComp;
6826dd5a8c8SToby Isaac       ierr    = PetscFree(valsNew);CHKERRQ(ierr);
6836dd5a8c8SToby Isaac       valsNew = valsComp;
6846dd5a8c8SToby Isaac     }
6856dd5a8c8SToby Isaac     ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),sizeNew,valsNew,PETSC_OWN_POINTER,isNew);CHKERRQ(ierr);
6866dd5a8c8SToby Isaac   }
6876dd5a8c8SToby Isaac   PetscFunctionReturn(0);
6886dd5a8c8SToby Isaac }
6896dd5a8c8SToby Isaac 
6906dd5a8c8SToby Isaac #undef __FUNCT__
69166af876cSToby Isaac #define __FUNCT__ "DMPlexComputeConstraints_Tree"
69266af876cSToby Isaac static PetscErrorCode DMPlexComputeConstraints_Tree(DM dm)
69366af876cSToby Isaac {
69466af876cSToby Isaac   PetscInt       p, pStart, pEnd, *anchors, size;
69566af876cSToby Isaac   PetscInt       aMin = PETSC_MAX_INT, aMax = PETSC_MIN_INT;
69666af876cSToby Isaac   PetscSection   aSec;
697f9f063d4SToby Isaac   DMLabel        canonLabel;
69866af876cSToby Isaac   IS             aIS;
69966af876cSToby Isaac   PetscErrorCode ierr;
70066af876cSToby Isaac 
70166af876cSToby Isaac   PetscFunctionBegin;
70266af876cSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
70366af876cSToby Isaac   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
704f9f063d4SToby Isaac   ierr = DMPlexGetLabel(dm,"canonical",&canonLabel);CHKERRQ(ierr);
70566af876cSToby Isaac   for (p = pStart; p < pEnd; p++) {
70666af876cSToby Isaac     PetscInt parent;
70766af876cSToby Isaac 
708f9f063d4SToby Isaac     if (canonLabel) {
709f9f063d4SToby Isaac       PetscInt canon;
710f9f063d4SToby Isaac 
711f9f063d4SToby Isaac       ierr = DMLabelGetValue(canonLabel,p,&canon);CHKERRQ(ierr);
712f9f063d4SToby Isaac       if (p != canon) continue;
713f9f063d4SToby Isaac     }
71466af876cSToby Isaac     ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
71566af876cSToby Isaac     if (parent != p) {
71666af876cSToby Isaac       aMin = PetscMin(aMin,p);
71766af876cSToby Isaac       aMax = PetscMax(aMax,p+1);
71866af876cSToby Isaac     }
71966af876cSToby Isaac   }
72066af876cSToby Isaac   if (aMin > aMax) {
72166af876cSToby Isaac     aMin = -1;
72266af876cSToby Isaac     aMax = -1;
72366af876cSToby Isaac   }
72466af876cSToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&aSec);CHKERRQ(ierr);
72566af876cSToby Isaac   ierr = PetscSectionSetChart(aSec,aMin,aMax);CHKERRQ(ierr);
72666af876cSToby Isaac   for (p = aMin; p < aMax; p++) {
72766af876cSToby Isaac     PetscInt parent, ancestor = p;
72866af876cSToby Isaac 
729f9f063d4SToby Isaac     if (canonLabel) {
730f9f063d4SToby Isaac       PetscInt canon;
731f9f063d4SToby Isaac 
732f9f063d4SToby Isaac       ierr = DMLabelGetValue(canonLabel,p,&canon);CHKERRQ(ierr);
733f9f063d4SToby Isaac       if (p != canon) continue;
734f9f063d4SToby Isaac     }
73566af876cSToby Isaac     ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
73666af876cSToby Isaac     while (parent != ancestor) {
73766af876cSToby Isaac       ancestor = parent;
73866af876cSToby Isaac       ierr     = DMPlexGetTreeParent(dm,ancestor,&parent,NULL);CHKERRQ(ierr);
73966af876cSToby Isaac     }
74066af876cSToby Isaac     if (ancestor != p) {
74166af876cSToby Isaac       PetscInt closureSize, *closure = NULL;
74266af876cSToby Isaac 
74366af876cSToby Isaac       ierr = DMPlexGetTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
74466af876cSToby Isaac       ierr = PetscSectionSetDof(aSec,p,closureSize);CHKERRQ(ierr);
74566af876cSToby Isaac       ierr = DMPlexRestoreTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
74666af876cSToby Isaac     }
74766af876cSToby Isaac   }
74866af876cSToby Isaac   ierr = PetscSectionSetUp(aSec);CHKERRQ(ierr);
74966af876cSToby Isaac   ierr = PetscSectionGetStorageSize(aSec,&size);CHKERRQ(ierr);
75066af876cSToby Isaac   ierr = PetscMalloc1(size,&anchors);CHKERRQ(ierr);
75166af876cSToby Isaac   for (p = aMin; p < aMax; p++) {
75266af876cSToby Isaac     PetscInt parent, ancestor = p;
75366af876cSToby Isaac 
754f9f063d4SToby Isaac     if (canonLabel) {
755f9f063d4SToby Isaac       PetscInt canon;
756f9f063d4SToby Isaac 
757f9f063d4SToby Isaac       ierr = DMLabelGetValue(canonLabel,p,&canon);CHKERRQ(ierr);
758f9f063d4SToby Isaac       if (p != canon) continue;
759f9f063d4SToby Isaac     }
76066af876cSToby Isaac     ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr);
76166af876cSToby Isaac     while (parent != ancestor) {
76266af876cSToby Isaac       ancestor = parent;
76366af876cSToby Isaac       ierr     = DMPlexGetTreeParent(dm,ancestor,&parent,NULL);CHKERRQ(ierr);
76466af876cSToby Isaac     }
76566af876cSToby Isaac     if (ancestor != p) {
76666af876cSToby Isaac       PetscInt j, closureSize, *closure = NULL, aOff;
76766af876cSToby Isaac 
76866af876cSToby Isaac       ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr);
76966af876cSToby Isaac 
77066af876cSToby Isaac       ierr = DMPlexGetTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
77166af876cSToby Isaac       for (j = 0; j < closureSize; j++) {
77266af876cSToby Isaac         anchors[aOff + j] = closure[2*j];
77366af876cSToby Isaac       }
77466af876cSToby Isaac       ierr = DMPlexRestoreTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
77566af876cSToby Isaac     }
77666af876cSToby Isaac   }
77766af876cSToby Isaac   ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm),size,anchors,PETSC_OWN_POINTER,&aIS);CHKERRQ(ierr);
7786dd5a8c8SToby Isaac   {
7796dd5a8c8SToby Isaac     PetscSection aSecNew = aSec;
7806dd5a8c8SToby Isaac     IS           aISNew  = aIS;
7816dd5a8c8SToby Isaac 
7826dd5a8c8SToby Isaac     ierr = PetscObjectReference((PetscObject)aSec);CHKERRQ(ierr);
7836dd5a8c8SToby Isaac     ierr = PetscObjectReference((PetscObject)aIS);CHKERRQ(ierr);
7846dd5a8c8SToby Isaac     while (aSecNew) {
7856dd5a8c8SToby Isaac       ierr    = PetscSectionDestroy(&aSec);CHKERRQ(ierr);
7866dd5a8c8SToby Isaac       ierr    = ISDestroy(&aIS);CHKERRQ(ierr);
7876dd5a8c8SToby Isaac       aSec    = aSecNew;
7886dd5a8c8SToby Isaac       aIS     = aISNew;
7896dd5a8c8SToby Isaac       aSecNew = NULL;
7906dd5a8c8SToby Isaac       aISNew  = NULL;
7916dd5a8c8SToby Isaac       ierr    = AnchorsFlatten(aSec,aIS,&aSecNew,&aISNew);CHKERRQ(ierr);
7926dd5a8c8SToby Isaac     }
7936dd5a8c8SToby Isaac   }
79466af876cSToby Isaac   ierr = DMPlexSetConstraints(dm,aSec,aIS);CHKERRQ(ierr);
79566af876cSToby Isaac   ierr = PetscSectionDestroy(&aSec);CHKERRQ(ierr);
79666af876cSToby Isaac   ierr = ISDestroy(&aIS);CHKERRQ(ierr);
79766af876cSToby Isaac   PetscFunctionReturn(0);
79866af876cSToby Isaac }
79966af876cSToby Isaac 
80066af876cSToby Isaac #undef __FUNCT__
801f9f063d4SToby Isaac #define __FUNCT__ "DMPlexSetTree_Internal"
802f9f063d4SToby Isaac static PetscErrorCode DMPlexSetTree_Internal(DM dm, PetscSection parentSection, PetscInt *parents, PetscInt *childIDs, PetscBool computeCanonical)
803f9f063d4SToby Isaac {
804f9f063d4SToby Isaac   DM_Plex       *mesh = (DM_Plex *)dm->data;
805f9f063d4SToby Isaac   DM             refTree;
806f9f063d4SToby Isaac   PetscInt       size;
807f9f063d4SToby Isaac   PetscErrorCode ierr;
808f9f063d4SToby Isaac 
809f9f063d4SToby Isaac   PetscFunctionBegin;
810f9f063d4SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
811f9f063d4SToby Isaac   PetscValidHeaderSpecific(parentSection, PETSC_SECTION_CLASSID, 2);
812f9f063d4SToby Isaac   ierr = PetscObjectReference((PetscObject)parentSection);CHKERRQ(ierr);
813f9f063d4SToby Isaac   ierr = PetscSectionDestroy(&mesh->parentSection);CHKERRQ(ierr);
814f9f063d4SToby Isaac   mesh->parentSection = parentSection;
815f9f063d4SToby Isaac   ierr = PetscSectionGetStorageSize(parentSection,&size);CHKERRQ(ierr);
816f9f063d4SToby Isaac   if (parents != mesh->parents) {
817f9f063d4SToby Isaac     ierr = PetscFree(mesh->parents);CHKERRQ(ierr);
818f9f063d4SToby Isaac     ierr = PetscMalloc1(size,&mesh->parents);CHKERRQ(ierr);
819f9f063d4SToby Isaac     ierr = PetscMemcpy(mesh->parents, parents, size * sizeof(*parents));CHKERRQ(ierr);
820f9f063d4SToby Isaac   }
821f9f063d4SToby Isaac   if (childIDs != mesh->childIDs) {
822f9f063d4SToby Isaac     ierr = PetscFree(mesh->childIDs);CHKERRQ(ierr);
823f9f063d4SToby Isaac     ierr = PetscMalloc1(size,&mesh->childIDs);CHKERRQ(ierr);
824f9f063d4SToby Isaac     ierr = PetscMemcpy(mesh->childIDs, childIDs, size * sizeof(*childIDs));CHKERRQ(ierr);
825f9f063d4SToby Isaac   }
826f9f063d4SToby Isaac   ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr);
827f9f063d4SToby Isaac   if (refTree) {
828f9f063d4SToby Isaac     DMLabel canonLabel;
829f9f063d4SToby Isaac 
830f9f063d4SToby Isaac     ierr = DMPlexGetLabel(refTree,"canonical",&canonLabel);CHKERRQ(ierr);
831f9f063d4SToby Isaac     if (canonLabel) {
832f9f063d4SToby Isaac       PetscInt i;
833f9f063d4SToby Isaac 
834f9f063d4SToby Isaac       for (i = 0; i < size; i++) {
835f9f063d4SToby Isaac         PetscInt canon;
836f9f063d4SToby Isaac         ierr = DMLabelGetValue(canonLabel, mesh->childIDs[i], &canon);CHKERRQ(ierr);
837f9f063d4SToby Isaac         if (canon >= 0) {
838f9f063d4SToby Isaac           mesh->childIDs[i] = canon;
839f9f063d4SToby Isaac         }
840f9f063d4SToby Isaac       }
841f9f063d4SToby Isaac     }
842f9f063d4SToby Isaac   }
843f9f063d4SToby Isaac   ierr = DMPlexTreeSymmetrize(dm);CHKERRQ(ierr);
844f9f063d4SToby Isaac   if (computeCanonical) {
845f9f063d4SToby Isaac     PetscInt d, dim;
846f9f063d4SToby Isaac 
847f9f063d4SToby Isaac     /* add the canonical label */
848f9f063d4SToby Isaac     ierr = DMPlexGetDimension(dm,&dim);CHKERRQ(ierr);
849f9f063d4SToby Isaac     ierr = DMPlexCreateLabel(dm,"canonical");CHKERRQ(ierr);
850f9f063d4SToby Isaac     for (d = 0; d <= dim; d++) {
851f9f063d4SToby Isaac       PetscInt p, dStart, dEnd, canon = -1, cNumChildren;
852f9f063d4SToby Isaac       const PetscInt *cChildren;
853f9f063d4SToby Isaac 
854f9f063d4SToby Isaac       ierr = DMPlexGetDepthStratum(dm,d,&dStart,&dEnd);CHKERRQ(ierr);
855f9f063d4SToby Isaac       for (p = dStart; p < dEnd; p++) {
856f9f063d4SToby Isaac         ierr = DMPlexGetTreeChildren(dm,p,&cNumChildren,&cChildren);CHKERRQ(ierr);
857f9f063d4SToby Isaac         if (cNumChildren) {
858f9f063d4SToby Isaac           canon = p;
859f9f063d4SToby Isaac           break;
860f9f063d4SToby Isaac         }
861f9f063d4SToby Isaac       }
862f9f063d4SToby Isaac       if (canon == -1) continue;
863f9f063d4SToby Isaac       for (p = dStart; p < dEnd; p++) {
864f9f063d4SToby Isaac         PetscInt numChildren, i;
865f9f063d4SToby Isaac         const PetscInt *children;
866f9f063d4SToby Isaac 
867f9f063d4SToby Isaac         ierr = DMPlexGetTreeChildren(dm,p,&numChildren,&children);CHKERRQ(ierr);
868f9f063d4SToby Isaac         if (numChildren) {
869f9f063d4SToby Isaac           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);
870f9f063d4SToby Isaac           ierr = DMPlexSetLabelValue(dm,"canonical",p,canon);CHKERRQ(ierr);
871f9f063d4SToby Isaac           for (i = 0; i < numChildren; i++) {
872f9f063d4SToby Isaac             ierr = DMPlexSetLabelValue(dm,"canonical",children[i],cChildren[i]);CHKERRQ(ierr);
873f9f063d4SToby Isaac           }
874f9f063d4SToby Isaac         }
875f9f063d4SToby Isaac       }
876f9f063d4SToby Isaac     }
877f9f063d4SToby Isaac   }
878f9f063d4SToby Isaac   ierr = DMPlexComputeConstraints_Tree(dm);CHKERRQ(ierr);
879f9f063d4SToby Isaac   PetscFunctionReturn(0);
880f9f063d4SToby Isaac }
881f9f063d4SToby Isaac 
882f9f063d4SToby Isaac #undef __FUNCT__
8830b7167a0SToby Isaac #define __FUNCT__ "DMPlexSetTree"
8840b7167a0SToby Isaac /*@
8850b7167a0SToby Isaac   DMPlexSetTree - set the tree that describes the hierarchy of non-conforming mesh points.  This routine also creates
8860b7167a0SToby Isaac   the point-to-point constraints determined by the tree: a point is constained to the points in the closure of its
8870b7167a0SToby Isaac   tree root.
8880b7167a0SToby Isaac 
8890b7167a0SToby Isaac   Collective on dm
8900b7167a0SToby Isaac 
8910b7167a0SToby Isaac   Input Parameters:
8920b7167a0SToby Isaac + dm - the DMPlex object
8930b7167a0SToby Isaac . parentSection - a section describing the tree: a point has a parent if it has 1 dof in the section; the section
8940b7167a0SToby Isaac                   offset indexes the parent and childID list; the reference count of parentSection is incremented
8950b7167a0SToby Isaac . parents - a list of the point parents; copied, can be destroyed
8960b7167a0SToby Isaac - childIDs - identifies the relationship of the child point to the parent point; if there is a reference tree, then
8970b7167a0SToby Isaac              the child corresponds to the point in the reference tree with index childIDs; copied, can be destroyed
8980b7167a0SToby Isaac 
8990b7167a0SToby Isaac   Level: intermediate
9000b7167a0SToby Isaac 
901b2f41788SToby Isaac .seealso: DMPlexGetTree(), DMPlexSetReferenceTree(), DMPlexSetConstraints(), DMPlexGetTreeParent(), DMPlexGetTreeChildren()
9020b7167a0SToby Isaac @*/
903b2f41788SToby Isaac PetscErrorCode DMPlexSetTree(DM dm, PetscSection parentSection, PetscInt parents[], PetscInt childIDs[])
9040b7167a0SToby Isaac {
9050b7167a0SToby Isaac   PetscErrorCode ierr;
9060b7167a0SToby Isaac 
9070b7167a0SToby Isaac   PetscFunctionBegin;
908f9f063d4SToby Isaac   ierr = DMPlexSetTree_Internal(dm,parentSection,parents,childIDs,PETSC_FALSE);CHKERRQ(ierr);
9090b7167a0SToby Isaac   PetscFunctionReturn(0);
9100b7167a0SToby Isaac }
9110b7167a0SToby Isaac 
9120b7167a0SToby Isaac #undef __FUNCT__
913b2f41788SToby Isaac #define __FUNCT__ "DMPlexGetTree"
914b2f41788SToby Isaac /*@
915b2f41788SToby Isaac   DMPlexGetTree - get the tree that describes the hierarchy of non-conforming mesh points.
916b2f41788SToby Isaac   Collective on dm
917b2f41788SToby Isaac 
918b2f41788SToby Isaac   Input Parameters:
919b2f41788SToby Isaac . dm - the DMPlex object
920b2f41788SToby Isaac 
921b2f41788SToby Isaac   Output Parameters:
922b2f41788SToby Isaac + parentSection - a section describing the tree: a point has a parent if it has 1 dof in the section; the section
923b2f41788SToby Isaac                   offset indexes the parent and childID list
924b2f41788SToby Isaac . parents - a list of the point parents
925b2f41788SToby Isaac . childIDs - identifies the relationship of the child point to the parent point; if there is a reference tree, then
926b2f41788SToby Isaac              the child corresponds to the point in the reference tree with index childID
927b2f41788SToby Isaac . childSection - the inverse of the parent section
928b2f41788SToby Isaac - children - a list of the point children
929b2f41788SToby Isaac 
930b2f41788SToby Isaac   Level: intermediate
931b2f41788SToby Isaac 
932b2f41788SToby Isaac .seealso: DMPlexSetTree(), DMPlexSetReferenceTree(), DMPlexSetConstraints(), DMPlexGetTreeParent(), DMPlexGetTreeChildren()
933b2f41788SToby Isaac @*/
934b2f41788SToby Isaac PetscErrorCode DMPlexGetTree(DM dm, PetscSection *parentSection, PetscInt *parents[], PetscInt *childIDs[], PetscSection *childSection, PetscInt *children[])
935b2f41788SToby Isaac {
936b2f41788SToby Isaac   DM_Plex        *mesh = (DM_Plex *)dm->data;
937b2f41788SToby Isaac 
938b2f41788SToby Isaac   PetscFunctionBegin;
939b2f41788SToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
940b2f41788SToby Isaac   if (parentSection) *parentSection = mesh->parentSection;
941b2f41788SToby Isaac   if (parents)       *parents       = mesh->parents;
942b2f41788SToby Isaac   if (childIDs)      *childIDs      = mesh->childIDs;
943b2f41788SToby Isaac   if (childSection)  *childSection  = mesh->childSection;
944b2f41788SToby Isaac   if (children)      *children      = mesh->children;
945b2f41788SToby Isaac   PetscFunctionReturn(0);
946b2f41788SToby Isaac }
947b2f41788SToby Isaac 
948b2f41788SToby Isaac #undef __FUNCT__
949d961a43aSToby Isaac #define __FUNCT__ "DMPlexGetTreeParent"
950d961a43aSToby Isaac /*@
951d961a43aSToby Isaac   DMPlexGetTreeParent - get the parent of a point in the tree describing the point hierarchy (not the Sieve DAG)
952d961a43aSToby Isaac 
953d961a43aSToby Isaac   Input Parameters:
954d961a43aSToby Isaac + dm - the DMPlex object
955d961a43aSToby Isaac - point - the query point
956d961a43aSToby Isaac 
957d961a43aSToby Isaac   Output Parameters:
958d961a43aSToby Isaac + parent - if not NULL, set to the parent of the point, or the point itself if the point does not have a parent
959d961a43aSToby Isaac - childID - if not NULL, set to the child ID of the point with respect to its parent, or 0 if the point
960d961a43aSToby Isaac             does not have a parent
961d961a43aSToby Isaac 
962d961a43aSToby Isaac   Level: intermediate
963d961a43aSToby Isaac 
964d961a43aSToby Isaac .seealso: DMPlexSetTree(), DMPlexGetTree(), DMPlexGetTreeChildren()
965d961a43aSToby Isaac @*/
966d961a43aSToby Isaac PetscErrorCode DMPlexGetTreeParent(DM dm, PetscInt point, PetscInt *parent, PetscInt *childID)
967d961a43aSToby Isaac {
968d961a43aSToby Isaac   DM_Plex       *mesh = (DM_Plex *)dm->data;
969d961a43aSToby Isaac   PetscSection   pSec;
970d961a43aSToby Isaac   PetscErrorCode ierr;
971d961a43aSToby Isaac 
972d961a43aSToby Isaac   PetscFunctionBegin;
973d961a43aSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
974d961a43aSToby Isaac   pSec = mesh->parentSection;
975d961a43aSToby Isaac   if (pSec && point >= pSec->pStart && point < pSec->pEnd) {
976d961a43aSToby Isaac     PetscInt dof;
977d961a43aSToby Isaac 
978d961a43aSToby Isaac     ierr = PetscSectionGetDof (pSec, point, &dof);CHKERRQ(ierr);
979d961a43aSToby Isaac     if (dof) {
980d961a43aSToby Isaac       PetscInt off;
981d961a43aSToby Isaac 
982d961a43aSToby Isaac       ierr = PetscSectionGetOffset (pSec, point, &off);CHKERRQ(ierr);
983d961a43aSToby Isaac       if (parent)  *parent = mesh->parents[off];
984d961a43aSToby Isaac       if (childID) *childID = mesh->childIDs[off];
985d961a43aSToby Isaac       PetscFunctionReturn(0);
986d961a43aSToby Isaac     }
987d961a43aSToby Isaac   }
988d961a43aSToby Isaac   if (parent) {
989d961a43aSToby Isaac     *parent = point;
990d961a43aSToby Isaac   }
991d961a43aSToby Isaac   if (childID) {
992d961a43aSToby Isaac     *childID = 0;
993d961a43aSToby Isaac   }
994d961a43aSToby Isaac   PetscFunctionReturn(0);
995d961a43aSToby Isaac }
996d961a43aSToby Isaac 
997d961a43aSToby Isaac #undef __FUNCT__
998d961a43aSToby Isaac #define __FUNCT__ "DMPlexGetTreeChildren"
999d961a43aSToby Isaac /*@C
1000d961a43aSToby Isaac   DMPlexGetTreeChildren - get the children of a point in the tree describing the point hierarchy (not the Sieve DAG)
1001d961a43aSToby Isaac 
1002d961a43aSToby Isaac   Input Parameters:
1003d961a43aSToby Isaac + dm - the DMPlex object
1004d961a43aSToby Isaac - point - the query point
1005d961a43aSToby Isaac 
1006d961a43aSToby Isaac   Output Parameters:
1007d961a43aSToby Isaac + numChildren - if not NULL, set to the number of children
1008d961a43aSToby Isaac - children - if not NULL, set to a list children, or set to NULL if the point has no children
1009d961a43aSToby Isaac 
1010d961a43aSToby Isaac   Level: intermediate
1011d961a43aSToby Isaac 
1012d961a43aSToby Isaac   Fortran Notes:
1013d961a43aSToby Isaac   Since it returns an array, this routine is only available in Fortran 90, and you must
1014d961a43aSToby Isaac   include petsc.h90 in your code.
1015d961a43aSToby Isaac 
1016d961a43aSToby Isaac .seealso: DMPlexSetTree(), DMPlexGetTree(), DMPlexGetTreeParent()
1017d961a43aSToby Isaac @*/
1018d961a43aSToby Isaac PetscErrorCode DMPlexGetTreeChildren(DM dm, PetscInt point, PetscInt *numChildren, const PetscInt *children[])
1019d961a43aSToby Isaac {
1020d961a43aSToby Isaac   DM_Plex       *mesh = (DM_Plex *)dm->data;
1021d961a43aSToby Isaac   PetscSection   childSec;
1022d961a43aSToby Isaac   PetscInt       dof = 0;
1023d961a43aSToby Isaac   PetscErrorCode ierr;
1024d961a43aSToby Isaac 
1025d961a43aSToby Isaac   PetscFunctionBegin;
1026d961a43aSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
1027d961a43aSToby Isaac   childSec = mesh->childSection;
1028d961a43aSToby Isaac   if (childSec && point >= childSec->pStart && point < childSec->pEnd) {
1029d961a43aSToby Isaac     ierr = PetscSectionGetDof (childSec, point, &dof);CHKERRQ(ierr);
1030d961a43aSToby Isaac   }
1031d961a43aSToby Isaac   if (numChildren) *numChildren = dof;
1032d961a43aSToby Isaac   if (children) {
1033d961a43aSToby Isaac     if (dof) {
1034d961a43aSToby Isaac       PetscInt off;
1035d961a43aSToby Isaac 
1036d961a43aSToby Isaac       ierr = PetscSectionGetOffset (childSec, point, &off);CHKERRQ(ierr);
1037d961a43aSToby Isaac       *children = &mesh->children[off];
1038d961a43aSToby Isaac     }
1039d961a43aSToby Isaac     else {
1040d961a43aSToby Isaac       *children = NULL;
1041d961a43aSToby Isaac     }
1042d961a43aSToby Isaac   }
1043d961a43aSToby Isaac   PetscFunctionReturn(0);
1044d961a43aSToby Isaac }
10450c37af3bSToby Isaac 
10460c37af3bSToby Isaac #undef __FUNCT__
10470c37af3bSToby Isaac #define __FUNCT__ "DMPlexComputeConstraintMatrix_ReferenceTree"
1048*ccffa8f8SToby Isaac static PetscErrorCode DMPlexComputeConstraintMatrix_ReferenceTree(DM dm)
10490c37af3bSToby Isaac {
10500c37af3bSToby Isaac   PetscDS        ds;
10510c37af3bSToby Isaac   PetscInt       spdim;
10520c37af3bSToby Isaac   PetscInt       numFields, f, c, cStart, cEnd, pStart, pEnd, conStart, conEnd;
10530c37af3bSToby Isaac   const PetscInt *anchors;
10540c37af3bSToby Isaac   PetscSection   section, cSec, aSec;
10550c37af3bSToby Isaac   Mat            cMat;
10560c37af3bSToby Isaac   PetscReal      *v0, *v0parent, *vtmp, *J, *Jparent, *invJparent, detJ, detJparent;
10570c37af3bSToby Isaac   IS             aIS;
10580c37af3bSToby Isaac   PetscErrorCode ierr;
10590c37af3bSToby Isaac 
10600c37af3bSToby Isaac   PetscFunctionBegin;
10610c37af3bSToby Isaac   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
10620c37af3bSToby Isaac   ierr = DMGetDS(dm,&ds);CHKERRQ(ierr);
10630c37af3bSToby Isaac   ierr = PetscDSGetNumFields(ds,&numFields);CHKERRQ(ierr);
10640c37af3bSToby Isaac   ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr);
10650c37af3bSToby Isaac   ierr = DMGetDefaultSection(dm,&section);CHKERRQ(ierr);
10660c37af3bSToby Isaac   ierr = DMPlexGetConstraintSection(dm,&cSec);CHKERRQ(ierr);
10670c37af3bSToby Isaac   ierr = DMPlexGetConstraints(dm,&aSec,&aIS);CHKERRQ(ierr);
10680c37af3bSToby Isaac   ierr = ISGetIndices(aIS,&anchors);CHKERRQ(ierr);
10690c37af3bSToby Isaac   ierr = PetscSectionGetChart(cSec,&conStart,&conEnd);CHKERRQ(ierr);
10700c37af3bSToby Isaac   ierr = DMPlexGetConstraintMatrix(dm,&cMat);CHKERRQ(ierr);
10710c37af3bSToby Isaac   ierr = DMPlexGetDimension(dm,&spdim);CHKERRQ(ierr);
10720c37af3bSToby Isaac   ierr = PetscMalloc6(spdim,&v0,spdim,&v0parent,spdim,&vtmp,spdim*spdim,&J,spdim*spdim,&Jparent,spdim*spdim,&invJparent);CHKERRQ(ierr);
10730c37af3bSToby Isaac 
10740c37af3bSToby Isaac   for (f = 0; f < numFields; f++) {
10750c37af3bSToby Isaac     PetscFE fe;
10760c37af3bSToby Isaac     PetscDualSpace space;
10770c37af3bSToby Isaac     PetscInt i, j, k, nPoints, offset;
10780c37af3bSToby Isaac     PetscInt fSize, fComp;
10790c37af3bSToby Isaac     PetscScalar *B = NULL;
10800c37af3bSToby Isaac     PetscReal *weights, *pointsRef, *pointsReal;
10810c37af3bSToby Isaac     Mat Amat, Bmat, Xmat;
10820c37af3bSToby Isaac 
10830c37af3bSToby Isaac     ierr = PetscDSGetDiscretization(ds,f,(PetscObject *)(&fe));CHKERRQ(ierr);
10840c37af3bSToby Isaac     ierr = PetscFEGetDualSpace(fe,&space);CHKERRQ(ierr);
10850c37af3bSToby Isaac     ierr = PetscDualSpaceGetDimension(space,&fSize);CHKERRQ(ierr);
10860c37af3bSToby Isaac     ierr = PetscFEGetNumComponents(fe,&fComp);CHKERRQ(ierr);
10870c37af3bSToby Isaac     ierr = MatCreate(PETSC_COMM_SELF,&Amat);CHKERRQ(ierr);
10880c37af3bSToby Isaac     ierr = MatSetSizes(Amat,fSize,fSize,fSize,fSize);CHKERRQ(ierr);
10890c37af3bSToby Isaac     ierr = MatSetType(Amat,MATSEQDENSE);CHKERRQ(ierr);
10900c37af3bSToby Isaac     ierr = MatSetUp(Amat);CHKERRQ(ierr);
10910c37af3bSToby Isaac     ierr = MatDuplicate(Amat,MAT_DO_NOT_COPY_VALUES,&Bmat);CHKERRQ(ierr);
10920c37af3bSToby Isaac     ierr = MatDuplicate(Amat,MAT_DO_NOT_COPY_VALUES,&Xmat);CHKERRQ(ierr);
10930c37af3bSToby Isaac     nPoints = 0;
10940c37af3bSToby Isaac     for (i = 0; i < fSize; i++) {
10950c37af3bSToby Isaac       PetscInt        qPoints;
10960c37af3bSToby Isaac       PetscQuadrature quad;
10970c37af3bSToby Isaac 
10980c37af3bSToby Isaac       ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr);
10990c37af3bSToby Isaac       ierr = PetscQuadratureGetData(quad,NULL,&qPoints,NULL,NULL);CHKERRQ(ierr);
11000c37af3bSToby Isaac       nPoints += qPoints;
11010c37af3bSToby Isaac     }
11020c37af3bSToby Isaac     ierr = PetscMalloc3(nPoints,&weights,spdim*nPoints,&pointsRef,spdim*nPoints,&pointsReal);CHKERRQ(ierr);
11030c37af3bSToby Isaac     offset = 0;
11040c37af3bSToby Isaac     for (i = 0; i < fSize; i++) {
11050c37af3bSToby Isaac       PetscInt        qPoints;
11060c37af3bSToby Isaac       const PetscReal    *p, *w;
11070c37af3bSToby Isaac       PetscQuadrature quad;
11080c37af3bSToby Isaac 
11090c37af3bSToby Isaac       ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr);
11100c37af3bSToby Isaac       ierr = PetscQuadratureGetData(quad,NULL,&qPoints,&p,&w);CHKERRQ(ierr);
11110c37af3bSToby Isaac       ierr = PetscMemcpy(weights+offset,w,qPoints*sizeof(*w));CHKERRQ(ierr);
11120c37af3bSToby Isaac       ierr = PetscMemcpy(pointsRef+spdim*offset,p,spdim*qPoints*sizeof(*p));CHKERRQ(ierr);
11130c37af3bSToby Isaac       offset += qPoints;
11140c37af3bSToby Isaac     }
11150c37af3bSToby Isaac     ierr = PetscFEGetTabulation(fe,nPoints,pointsRef,&B,NULL,NULL);CHKERRQ(ierr);
11160c37af3bSToby Isaac     offset = 0;
11170c37af3bSToby Isaac     for (i = 0; i < fSize; i++) {
11180c37af3bSToby Isaac       PetscInt        qPoints;
11190c37af3bSToby Isaac       PetscQuadrature quad;
11200c37af3bSToby Isaac 
11210c37af3bSToby Isaac       ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr);
11220c37af3bSToby Isaac       ierr = PetscQuadratureGetData(quad,NULL,&qPoints,NULL,NULL);CHKERRQ(ierr);
11230c37af3bSToby Isaac       for (j = 0; j < fSize; j++) {
11240c37af3bSToby Isaac         PetscScalar val = 0.;
11250c37af3bSToby Isaac 
11260c37af3bSToby Isaac         for (k = 0; k < qPoints; k++) {
11270c37af3bSToby Isaac           val += B[((offset + k) * fSize + j) * fComp] * weights[k];
11280c37af3bSToby Isaac         }
11290c37af3bSToby Isaac         ierr = MatSetValue(Amat,i,j,val,INSERT_VALUES);CHKERRQ(ierr);
11300c37af3bSToby Isaac       }
11310c37af3bSToby Isaac       offset += qPoints;
11320c37af3bSToby Isaac     }
11330c37af3bSToby Isaac     ierr = PetscFERestoreTabulation(fe,nPoints,pointsRef,&B,NULL,NULL);CHKERRQ(ierr);
11340c37af3bSToby Isaac     ierr = MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11350c37af3bSToby Isaac     ierr = MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11360c37af3bSToby Isaac     ierr = MatLUFactor(Amat,NULL,NULL,NULL);CHKERRQ(ierr);
11370c37af3bSToby Isaac     for (c = cStart; c < cEnd; c++) {
11380c37af3bSToby Isaac       PetscInt parent;
11390c37af3bSToby Isaac       PetscInt closureSize, closureSizeP, *closure = NULL, *closureP = NULL;
11400c37af3bSToby Isaac       PetscInt *childOffsets, *parentOffsets;
11410c37af3bSToby Isaac 
11420c37af3bSToby Isaac       ierr = DMPlexGetTreeParent(dm,c,&parent,NULL);CHKERRQ(ierr);
11430c37af3bSToby Isaac       if (parent == c) continue;
11440c37af3bSToby Isaac       ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
11450c37af3bSToby Isaac       for (i = 0; i < closureSize; i++) {
11460c37af3bSToby Isaac         PetscInt p = closure[2*i];
11470c37af3bSToby Isaac         PetscInt conDof;
11480c37af3bSToby Isaac 
11490c37af3bSToby Isaac         if (p < conStart || p >= conEnd) continue;
11500c37af3bSToby Isaac         if (numFields > 1) {
11510c37af3bSToby Isaac           ierr = PetscSectionGetFieldDof(cSec,p,f,&conDof);CHKERRQ(ierr);
11520c37af3bSToby Isaac         }
11530c37af3bSToby Isaac         else {
11540c37af3bSToby Isaac           ierr = PetscSectionGetDof(cSec,p,&conDof);CHKERRQ(ierr);
11550c37af3bSToby Isaac         }
11560c37af3bSToby Isaac         if (conDof) break;
11570c37af3bSToby Isaac       }
11580c37af3bSToby Isaac       if (i == closureSize) {
11590c37af3bSToby Isaac         ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
11600c37af3bSToby Isaac         continue;
11610c37af3bSToby Isaac       }
11620c37af3bSToby Isaac 
11630c37af3bSToby Isaac       ierr = DMPlexComputeCellGeometry(dm, c, v0, J, NULL, &detJ);CHKERRQ(ierr);
11640c37af3bSToby Isaac       ierr = DMPlexComputeCellGeometry(dm, parent, v0parent, Jparent, invJparent, &detJparent);CHKERRQ(ierr);
11650c37af3bSToby Isaac       for (i = 0; i < nPoints; i++) {
11660c37af3bSToby Isaac         CoordinatesRefToReal(spdim, spdim, v0, J, &pointsRef[i*spdim],vtmp);CHKERRQ(ierr);
11670c37af3bSToby Isaac         CoordinatesRealToRef(spdim, spdim, v0parent, invJparent, vtmp, &pointsReal[i*spdim]);CHKERRQ(ierr);
11680c37af3bSToby Isaac       }
11690c37af3bSToby Isaac       ierr = PetscFEGetTabulation(fe,nPoints,pointsReal,&B,NULL,NULL);CHKERRQ(ierr);
11700c37af3bSToby Isaac       offset = 0;
11710c37af3bSToby Isaac       for (i = 0; i < fSize; i++) {
11720c37af3bSToby Isaac         PetscInt        qPoints;
11730c37af3bSToby Isaac         PetscQuadrature quad;
11740c37af3bSToby Isaac 
11750c37af3bSToby Isaac         ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr);
11760c37af3bSToby Isaac         ierr = PetscQuadratureGetData(quad,NULL,&qPoints,NULL,NULL);CHKERRQ(ierr);
11770c37af3bSToby Isaac         for (j = 0; j < fSize; j++) {
11780c37af3bSToby Isaac           PetscScalar val = 0.;
11790c37af3bSToby Isaac 
11800c37af3bSToby Isaac           for (k = 0; k < qPoints; k++) {
11810c37af3bSToby Isaac             val += B[((offset + k) * fSize + j) * fComp] * weights[k];
11820c37af3bSToby Isaac           }
11830c37af3bSToby Isaac           MatSetValue(Bmat,i,j,val,INSERT_VALUES);CHKERRQ(ierr);
11840c37af3bSToby Isaac         }
11850c37af3bSToby Isaac         offset += qPoints;
11860c37af3bSToby Isaac       }
11870c37af3bSToby Isaac       ierr = PetscFERestoreTabulation(fe,nPoints,pointsReal,&B,NULL,NULL);CHKERRQ(ierr);
11880c37af3bSToby Isaac       ierr = MatAssemblyBegin(Bmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11890c37af3bSToby Isaac       ierr = MatAssemblyEnd(Bmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
11900c37af3bSToby Isaac       ierr = MatMatSolve(Amat,Bmat,Xmat);CHKERRQ(ierr);
11910c37af3bSToby Isaac       ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSizeP,&closureP);CHKERRQ(ierr);
11920c37af3bSToby Isaac       ierr = PetscMalloc2(closureSize+1,&childOffsets,closureSizeP+1,&parentOffsets);CHKERRQ(ierr);
11930c37af3bSToby Isaac       childOffsets[0] = 0;
11940c37af3bSToby Isaac       for (i = 0; i < closureSize; i++) {
11950c37af3bSToby Isaac         PetscInt p = closure[2*i];
11960c37af3bSToby Isaac         PetscInt dof;
11970c37af3bSToby Isaac 
11980c37af3bSToby Isaac         if (numFields > 1) {
11990c37af3bSToby Isaac           ierr = PetscSectionGetFieldDof(section,p,f,&dof);CHKERRQ(ierr);
12000c37af3bSToby Isaac         }
12010c37af3bSToby Isaac         else {
12020c37af3bSToby Isaac           ierr = PetscSectionGetDof(section,p,&dof);CHKERRQ(ierr);
12030c37af3bSToby Isaac         }
12040c37af3bSToby Isaac         childOffsets[i+1]=childOffsets[i]+dof / fComp;
12050c37af3bSToby Isaac       }
12060c37af3bSToby Isaac       parentOffsets[0] = 0;
12070c37af3bSToby Isaac       for (i = 0; i < closureSizeP; i++) {
12080c37af3bSToby Isaac         PetscInt p = closureP[2*i];
12090c37af3bSToby Isaac         PetscInt dof;
12100c37af3bSToby Isaac 
12110c37af3bSToby Isaac         if (numFields > 1) {
12120c37af3bSToby Isaac           ierr = PetscSectionGetFieldDof(section,p,f,&dof);CHKERRQ(ierr);
12130c37af3bSToby Isaac         }
12140c37af3bSToby Isaac         else {
12150c37af3bSToby Isaac           ierr = PetscSectionGetDof(section,p,&dof);CHKERRQ(ierr);
12160c37af3bSToby Isaac         }
12170c37af3bSToby Isaac         parentOffsets[i+1]=parentOffsets[i]+dof / fComp;
12180c37af3bSToby Isaac       }
12190c37af3bSToby Isaac       for (i = 0; i < closureSize; i++) {
12200c37af3bSToby Isaac         PetscInt conDof, conOff, aDof, aOff;
12210c37af3bSToby Isaac         PetscInt p = closure[2*i];
12220c37af3bSToby Isaac         PetscInt o = closure[2*i+1];
12230c37af3bSToby Isaac 
12240c37af3bSToby Isaac         if (p < conStart || p >= conEnd) continue;
12250c37af3bSToby Isaac         if (numFields > 1) {
12260c37af3bSToby Isaac           ierr = PetscSectionGetFieldDof(cSec,p,f,&conDof);CHKERRQ(ierr);
12270c37af3bSToby Isaac           ierr = PetscSectionGetFieldOffset(cSec,p,f,&conOff);CHKERRQ(ierr);
12280c37af3bSToby Isaac         }
12290c37af3bSToby Isaac         else {
12300c37af3bSToby Isaac           ierr = PetscSectionGetDof(cSec,p,&conDof);CHKERRQ(ierr);
12310c37af3bSToby Isaac           ierr = PetscSectionGetOffset(cSec,p,&conOff);CHKERRQ(ierr);
12320c37af3bSToby Isaac         }
12330c37af3bSToby Isaac         if (!conDof) continue;
12340c37af3bSToby Isaac         ierr = PetscSectionGetDof(aSec,p,&aDof);CHKERRQ(ierr);
12350c37af3bSToby Isaac         ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr);
12360c37af3bSToby Isaac         for (k = 0; k < aDof; k++) {
12370c37af3bSToby Isaac           PetscInt a = anchors[aOff + k];
12380c37af3bSToby Isaac           PetscInt aSecDof, aSecOff;
12390c37af3bSToby Isaac 
12400c37af3bSToby Isaac           if (numFields > 1) {
12410c37af3bSToby Isaac             ierr = PetscSectionGetFieldDof(section,a,f,&aSecDof);CHKERRQ(ierr);
12420c37af3bSToby Isaac             ierr = PetscSectionGetFieldOffset(section,a,f,&aSecOff);CHKERRQ(ierr);
12430c37af3bSToby Isaac           }
12440c37af3bSToby Isaac           else {
12450c37af3bSToby Isaac             ierr = PetscSectionGetDof(section,a,&aSecDof);CHKERRQ(ierr);
12460c37af3bSToby Isaac             ierr = PetscSectionGetOffset(section,a,&aSecOff);CHKERRQ(ierr);
12470c37af3bSToby Isaac           }
12480c37af3bSToby Isaac           if (!aSecDof) continue;
12490c37af3bSToby Isaac 
12500c37af3bSToby Isaac           for (j = 0; j < closureSizeP; j++) {
12510c37af3bSToby Isaac             PetscInt q = closureP[2*j];
12520c37af3bSToby Isaac             PetscInt oq = closureP[2*j+1];
12530c37af3bSToby Isaac 
12540c37af3bSToby Isaac             if (q == a) {
12550c37af3bSToby Isaac               PetscInt r, s, t;
12560c37af3bSToby Isaac 
12570c37af3bSToby Isaac               for (r = childOffsets[i]; r < childOffsets[i+1]; r++) {
12580c37af3bSToby Isaac                 for (s = parentOffsets[j]; s < parentOffsets[j+1]; s++) {
12590c37af3bSToby Isaac                   PetscScalar val;
12600c37af3bSToby Isaac                   PetscInt insertCol, insertRow;
12610c37af3bSToby Isaac 
12620c37af3bSToby Isaac                   ierr = MatGetValue(Xmat,r,s,&val);CHKERRQ(ierr);
12630c37af3bSToby Isaac                   if (o >= 0) {
12640c37af3bSToby Isaac                     insertRow = conOff + fComp * (r - childOffsets[i]);
12650c37af3bSToby Isaac                   }
12660c37af3bSToby Isaac                   else {
12670c37af3bSToby Isaac                     insertRow = conOff + fComp * (childOffsets[i + 1] - 1 - r);
12680c37af3bSToby Isaac                   }
12690c37af3bSToby Isaac                   if (oq >= 0) {
12700c37af3bSToby Isaac                     insertCol = aSecOff + fComp * (s - parentOffsets[j]);
12710c37af3bSToby Isaac                   }
12720c37af3bSToby Isaac                   else {
12730c37af3bSToby Isaac                     insertCol = aSecOff + fComp * (parentOffsets[j + 1] - 1 - s);
12740c37af3bSToby Isaac                   }
12750c37af3bSToby Isaac                   for (t = 0; t < fComp; t++) {
12760c37af3bSToby Isaac                     ierr = MatSetValue(cMat,insertRow + t,insertCol + t,val,INSERT_VALUES);CHKERRQ(ierr);
12770c37af3bSToby Isaac                   }
12780c37af3bSToby Isaac                 }
12790c37af3bSToby Isaac               }
12800c37af3bSToby Isaac             }
12810c37af3bSToby Isaac           }
12820c37af3bSToby Isaac         }
12830c37af3bSToby Isaac       }
12840c37af3bSToby Isaac       ierr = PetscFree2(childOffsets,parentOffsets);CHKERRQ(ierr);
12850c37af3bSToby Isaac       ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
12860c37af3bSToby Isaac       ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSizeP,&closureP);CHKERRQ(ierr);
12870c37af3bSToby Isaac     }
12880c37af3bSToby Isaac     ierr = MatDestroy(&Amat);CHKERRQ(ierr);
12890c37af3bSToby Isaac     ierr = MatDestroy(&Bmat);CHKERRQ(ierr);
12900c37af3bSToby Isaac     ierr = MatDestroy(&Xmat);CHKERRQ(ierr);
12910c37af3bSToby Isaac     ierr = PetscFree3(weights,pointsRef,pointsReal);CHKERRQ(ierr);
12920c37af3bSToby Isaac   }
12930c37af3bSToby Isaac   ierr = MatAssemblyBegin(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12940c37af3bSToby Isaac   ierr = MatAssemblyEnd(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12950c37af3bSToby Isaac   ierr = PetscFree6(v0,v0parent,vtmp,J,Jparent,invJparent);CHKERRQ(ierr);
12960c37af3bSToby Isaac   ierr = ISRestoreIndices(aIS,&anchors);CHKERRQ(ierr);
12970c37af3bSToby Isaac 
12980c37af3bSToby Isaac   PetscFunctionReturn(0);
12990c37af3bSToby Isaac }
130095a0b26dSToby Isaac 
130195a0b26dSToby Isaac #undef __FUNCT__
130295a0b26dSToby Isaac #define __FUNCT__ "DMPlexComputeConstraintMatrix_Tree"
130395a0b26dSToby Isaac PetscErrorCode DMPlexComputeConstraintMatrix_Tree(DM dm)
130495a0b26dSToby Isaac {
130595a0b26dSToby Isaac   DM             refTree;
130695a0b26dSToby Isaac   PetscDS        ds;
130795a0b26dSToby Isaac   Mat            refCmat, cMat;
130895a0b26dSToby Isaac   PetscInt       numFields, f, pRefStart, pRefEnd, p, *rows, *cols, maxDof, maxAnDof, *perm, *iperm, pStart, pEnd, conStart, conEnd, **refPointFieldN;
130995a0b26dSToby Isaac   PetscScalar ***refPointFieldMats, *pointWork;
131095a0b26dSToby Isaac   PetscSection   refConSec, conSec, refAnSec, anSec, section, refSection;
131195a0b26dSToby Isaac   IS             refAnIS, anIS;
131295a0b26dSToby Isaac   const PetscInt *refAnchors, *anchors;
131395a0b26dSToby Isaac   PetscErrorCode ierr;
131495a0b26dSToby Isaac 
131595a0b26dSToby Isaac   PetscFunctionBegin;
131695a0b26dSToby Isaac   PetscValidHeaderSpecific(dm, DM_CLASSID, 1);
131795a0b26dSToby Isaac   ierr = DMGetDS(dm,&ds);CHKERRQ(ierr);
131895a0b26dSToby Isaac   ierr = PetscDSGetNumFields(ds,&numFields);CHKERRQ(ierr);
131995a0b26dSToby Isaac   ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr);
132095a0b26dSToby Isaac   ierr = DMSetDS(refTree,ds);CHKERRQ(ierr);
132195a0b26dSToby Isaac   ierr = DMPlexComputeConstraintMatrix_ReferenceTree(refTree);CHKERRQ(ierr);
132295a0b26dSToby Isaac   ierr = DMPlexGetConstraintMatrix(refTree,&refCmat);CHKERRQ(ierr);
132395a0b26dSToby Isaac   ierr = DMPlexGetChart(refTree,&pRefStart,&pRefEnd);CHKERRQ(ierr);
132495a0b26dSToby Isaac   ierr = DMPlexGetConstraints(refTree,&refAnSec,&refAnIS);CHKERRQ(ierr);
132595a0b26dSToby Isaac   ierr = DMPlexGetConstraints(dm,&anSec,&anIS);CHKERRQ(ierr);
132695a0b26dSToby Isaac   ierr = ISGetIndices(refAnIS,&refAnchors);CHKERRQ(ierr);
132795a0b26dSToby Isaac   ierr = ISGetIndices(anIS,&anchors);CHKERRQ(ierr);
132895a0b26dSToby Isaac   ierr = DMPlexGetConstraintSection(refTree,&refConSec);CHKERRQ(ierr);
132995a0b26dSToby Isaac   ierr = DMPlexGetConstraintSection(dm,&conSec);CHKERRQ(ierr);
133095a0b26dSToby Isaac   ierr = DMGetDefaultSection(refTree,&refSection);CHKERRQ(ierr);
133195a0b26dSToby Isaac   ierr = DMGetDefaultSection(dm,&section);CHKERRQ(ierr);
133295a0b26dSToby Isaac   ierr = PetscSectionGetChart(refConSec,&pRefStart,&pRefEnd);CHKERRQ(ierr);
133395a0b26dSToby Isaac   ierr = PetscSectionGetChart(conSec,&conStart,&conEnd);CHKERRQ(ierr);
133495a0b26dSToby Isaac   ierr = PetscMalloc1(pRefEnd-pRefStart,&refPointFieldMats);CHKERRQ(ierr);
133595a0b26dSToby Isaac   ierr = PetscMalloc1(pRefEnd-pRefStart,&refPointFieldN);CHKERRQ(ierr);
133695a0b26dSToby Isaac   ierr = PetscSectionGetMaxDof(refConSec,&maxDof);CHKERRQ(ierr);
133795a0b26dSToby Isaac   ierr = PetscSectionGetMaxDof(refAnSec,&maxAnDof);CHKERRQ(ierr);
133895a0b26dSToby Isaac   ierr = PetscMalloc1(maxDof,&rows);CHKERRQ(ierr);
133995a0b26dSToby Isaac   ierr = PetscMalloc1(maxDof*maxAnDof,&cols);CHKERRQ(ierr);
134095a0b26dSToby Isaac   ierr = PetscMalloc1(maxDof*maxDof*maxAnDof,&pointWork);CHKERRQ(ierr);
134195a0b26dSToby Isaac 
134295a0b26dSToby Isaac   /* step 1: get submats for every constrained point in the reference tree */
134395a0b26dSToby Isaac   for (p = pRefStart; p < pRefEnd; p++) {
134495a0b26dSToby Isaac     PetscInt parent, closureSize, *closure = NULL, pDof;
134595a0b26dSToby Isaac 
134695a0b26dSToby Isaac     ierr = DMPlexGetTreeParent(refTree,p,&parent,NULL);CHKERRQ(ierr);
134795a0b26dSToby Isaac     ierr = PetscSectionGetDof(refConSec,p,&pDof);CHKERRQ(ierr);
134895a0b26dSToby Isaac     if (!pDof || parent == p) continue;
134995a0b26dSToby Isaac 
135095a0b26dSToby Isaac     ierr = PetscMalloc1(numFields,&refPointFieldMats[p-pRefStart]);CHKERRQ(ierr);
135195a0b26dSToby Isaac     ierr = PetscCalloc1(numFields,&refPointFieldN[p-pRefStart]);CHKERRQ(ierr);
135295a0b26dSToby Isaac     ierr = DMPlexGetTransitiveClosure(refTree,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
135395a0b26dSToby Isaac     for (f = 0; f < numFields; f++) {
135495a0b26dSToby Isaac       PetscInt cDof, cOff, numCols, r, i, fComp;
135595a0b26dSToby Isaac       PetscFE fe;
135695a0b26dSToby Isaac 
135795a0b26dSToby Isaac       ierr = PetscDSGetDiscretization(ds,f,(PetscObject *) &fe);CHKERRQ(ierr);
135895a0b26dSToby Isaac       ierr = PetscFEGetNumComponents(fe,&fComp);CHKERRQ(ierr);
135995a0b26dSToby Isaac 
136095a0b26dSToby Isaac       if (numFields > 1) {
136195a0b26dSToby Isaac         ierr = PetscSectionGetFieldDof(refConSec,p,f,&cDof);CHKERRQ(ierr);
136295a0b26dSToby Isaac         ierr = PetscSectionGetFieldOffset(refConSec,p,f,&cOff);CHKERRQ(ierr);
136395a0b26dSToby Isaac       }
136495a0b26dSToby Isaac       else {
136595a0b26dSToby Isaac         ierr = PetscSectionGetDof(refConSec,p,&cDof);CHKERRQ(ierr);
136695a0b26dSToby Isaac         ierr = PetscSectionGetOffset(refConSec,p,&cOff);CHKERRQ(ierr);
136795a0b26dSToby Isaac       }
136895a0b26dSToby Isaac 
136995a0b26dSToby Isaac       if (!cDof) continue;
137095a0b26dSToby Isaac       for (r = 0; r < cDof; r++) {
137195a0b26dSToby Isaac         rows[r] = cOff + r;
137295a0b26dSToby Isaac       }
137395a0b26dSToby Isaac       numCols = 0;
137495a0b26dSToby Isaac       for (i = 0; i < closureSize; i++) {
137595a0b26dSToby Isaac         PetscInt q = closure[2*i];
137695a0b26dSToby Isaac         PetscInt o = closure[2*i+1];
137795a0b26dSToby Isaac         PetscInt aDof, aOff, j;
137895a0b26dSToby Isaac 
137995a0b26dSToby Isaac         if (numFields > 1) {
138095a0b26dSToby Isaac           ierr = PetscSectionGetFieldDof(refSection,q,f,&aDof);CHKERRQ(ierr);
138195a0b26dSToby Isaac           ierr = PetscSectionGetFieldOffset(refSection,q,f,&aOff);CHKERRQ(ierr);
138295a0b26dSToby Isaac         }
138395a0b26dSToby Isaac         else {
138495a0b26dSToby Isaac           ierr = PetscSectionGetDof(refSection,q,&aDof);CHKERRQ(ierr);
138595a0b26dSToby Isaac           ierr = PetscSectionGetOffset(refSection,q,&aOff);CHKERRQ(ierr);
138695a0b26dSToby Isaac         }
138795a0b26dSToby Isaac 
138895a0b26dSToby Isaac         for (j = 0; j < aDof; j++) {
138995a0b26dSToby Isaac           PetscInt node = (o >= 0) ? (j / fComp) : ((aDof - 1 - j) / fComp);
139095a0b26dSToby Isaac           PetscInt comp = (j % fComp);
139195a0b26dSToby Isaac 
139295a0b26dSToby Isaac           cols[numCols++] = aOff + node * fComp + comp;
139395a0b26dSToby Isaac         }
139495a0b26dSToby Isaac       }
139595a0b26dSToby Isaac       refPointFieldN[p-pRefStart][f] = numCols;
139695a0b26dSToby Isaac       ierr = PetscMalloc1(cDof*numCols,&refPointFieldMats[p-pRefStart][f]);CHKERRQ(ierr);
139795a0b26dSToby Isaac       ierr = MatGetValues(refCmat,cDof,rows,numCols,cols,refPointFieldMats[p-pRefStart][f]);CHKERRQ(ierr);
139895a0b26dSToby Isaac     }
139995a0b26dSToby Isaac     ierr = DMPlexRestoreTransitiveClosure(refTree,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
140095a0b26dSToby Isaac   }
140195a0b26dSToby Isaac 
140295a0b26dSToby Isaac   /* step 2: compute the preorder */
140395a0b26dSToby Isaac   ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr);
140495a0b26dSToby Isaac   ierr = PetscMalloc2(pEnd-pStart,&perm,pEnd-pStart,&iperm);CHKERRQ(ierr);
140595a0b26dSToby Isaac   for (p = pStart; p < pEnd; p++) {
140695a0b26dSToby Isaac     perm[p - pStart] = p;
140795a0b26dSToby Isaac     iperm[p - pStart] = p-pStart;
140895a0b26dSToby Isaac   }
140995a0b26dSToby Isaac   for (p = 0; p < pEnd - pStart;) {
141095a0b26dSToby Isaac     PetscInt point = perm[p];
141195a0b26dSToby Isaac     PetscInt parent;
141295a0b26dSToby Isaac 
141395a0b26dSToby Isaac     ierr = DMPlexGetTreeParent(dm,point,&parent,NULL);CHKERRQ(ierr);
141495a0b26dSToby Isaac     if (parent == point) {
141595a0b26dSToby Isaac       p++;
141695a0b26dSToby Isaac     }
141795a0b26dSToby Isaac     else {
141895a0b26dSToby Isaac       PetscInt size, closureSize, *closure = NULL, i;
141995a0b26dSToby Isaac 
142095a0b26dSToby Isaac       ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
142195a0b26dSToby Isaac       for (i = 0; i < closureSize; i++) {
142295a0b26dSToby Isaac         PetscInt q = closure[2*i];
142395a0b26dSToby Isaac         if (iperm[q-pStart] > iperm[point-pStart]) {
142495a0b26dSToby Isaac           /* swap */
142595a0b26dSToby Isaac           perm[p]               = q;
142695a0b26dSToby Isaac           perm[iperm[q-pStart]] = point;
142795a0b26dSToby Isaac           iperm[point-pStart]   = iperm[q-pStart];
142895a0b26dSToby Isaac           iperm[q-pStart]       = p;
142995a0b26dSToby Isaac           break;
143095a0b26dSToby Isaac         }
143195a0b26dSToby Isaac       }
143295a0b26dSToby Isaac       size = closureSize;
143395a0b26dSToby Isaac       ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
143495a0b26dSToby Isaac       if (i == size) {
143595a0b26dSToby Isaac         p++;
143695a0b26dSToby Isaac       }
143795a0b26dSToby Isaac     }
143895a0b26dSToby Isaac   }
143995a0b26dSToby Isaac 
144095a0b26dSToby Isaac   /* step 3: fill the constraint matrix */
144195a0b26dSToby Isaac   /* we are going to use a preorder progressive fill strategy.  Mat doesn't
144295a0b26dSToby Isaac    * allow progressive fill without assembly, so we are going to set up the
144395a0b26dSToby Isaac    * values outside of the Mat first.
144495a0b26dSToby Isaac    */
144595a0b26dSToby Isaac   {
144695a0b26dSToby Isaac     PetscInt nRows, row, nnz;
144795a0b26dSToby Isaac     PetscBool done;
144895a0b26dSToby Isaac     const PetscInt *ia, *ja;
144995a0b26dSToby Isaac     PetscScalar *vals;
145095a0b26dSToby Isaac 
145195a0b26dSToby Isaac     ierr = DMPlexGetConstraintMatrix(dm,&cMat);CHKERRQ(ierr);
145295a0b26dSToby Isaac     ierr = MatGetRowIJ(cMat,0,PETSC_FALSE,PETSC_FALSE,&nRows,&ia,&ja,&done);CHKERRQ(ierr);
145395a0b26dSToby Isaac     if (!done) SETERRQ(PetscObjectComm((PetscObject)cMat),PETSC_ERR_PLIB,"Could not get RowIJ of constraint matrix");
145495a0b26dSToby Isaac     nnz  = ia[nRows];
145595a0b26dSToby Isaac     /* malloc and then zero rows right before we fill them: this way valgrind
145695a0b26dSToby Isaac      * can tell if we are doing progressive fill in the wrong order */
145795a0b26dSToby Isaac     ierr = PetscMalloc1(nnz,&vals);CHKERRQ(ierr);
145895a0b26dSToby Isaac     for (p = 0; p < pEnd - pStart; p++) {
145995a0b26dSToby Isaac       PetscInt parent, childid, closureSize, *closure = NULL;
146095a0b26dSToby Isaac       PetscInt point = perm[p], pointDof;
146195a0b26dSToby Isaac 
146295a0b26dSToby Isaac       ierr = DMPlexGetTreeParent(dm,point,&parent,&childid);CHKERRQ(ierr);
146395a0b26dSToby Isaac       if ((point < conStart) || (point >= conEnd) || (parent == point)) continue;
146495a0b26dSToby Isaac       ierr = PetscSectionGetDof(conSec,point,&pointDof);CHKERRQ(ierr);
146595a0b26dSToby Isaac       if (!pointDof) continue;
146695a0b26dSToby Isaac       ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
146795a0b26dSToby Isaac       for (f = 0; f < numFields; f++) {
146895a0b26dSToby Isaac         PetscInt cDof, cOff, numCols, numFillCols, i, r, fComp, matOffset, offset;
146995a0b26dSToby Isaac         PetscScalar *pointMat;
147095a0b26dSToby Isaac         PetscFE fe;
147195a0b26dSToby Isaac 
147295a0b26dSToby Isaac         ierr = PetscDSGetDiscretization(ds,f,(PetscObject *) &fe);CHKERRQ(ierr);
147395a0b26dSToby Isaac         ierr = PetscFEGetNumComponents(fe,&fComp);CHKERRQ(ierr);
147495a0b26dSToby Isaac 
147595a0b26dSToby Isaac         if (numFields > 1) {
147695a0b26dSToby Isaac           ierr = PetscSectionGetFieldDof(conSec,point,f,&cDof);CHKERRQ(ierr);
147795a0b26dSToby Isaac           ierr = PetscSectionGetFieldOffset(conSec,point,f,&cOff);CHKERRQ(ierr);
147895a0b26dSToby Isaac         }
147995a0b26dSToby Isaac         else {
148095a0b26dSToby Isaac           ierr = PetscSectionGetDof(conSec,point,&cDof);CHKERRQ(ierr);
148195a0b26dSToby Isaac           ierr = PetscSectionGetOffset(conSec,point,&cOff);CHKERRQ(ierr);
148295a0b26dSToby Isaac         }
148395a0b26dSToby Isaac         if (!cDof) continue;
148495a0b26dSToby Isaac 
148595a0b26dSToby Isaac         /* make sure that every row for this point is the same size */
148695a0b26dSToby Isaac #if defined(PETSC_USE_DEBUG)
148795a0b26dSToby Isaac         for (r = 0; r < cDof; r++) {
148895a0b26dSToby Isaac           if (cDof > 1 && r) {
148995a0b26dSToby Isaac             if ((ia[cOff+r+1]-ia[cOff+r]) != (ia[cOff+r]-ia[cOff+r-1])) {
149095a0b26dSToby Isaac               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]));
149195a0b26dSToby Isaac             }
149295a0b26dSToby Isaac           }
149395a0b26dSToby Isaac         }
149495a0b26dSToby Isaac #endif
149595a0b26dSToby Isaac         /* zero rows */
149695a0b26dSToby Isaac         for (i = ia[cOff] ; i< ia[cOff+cDof];i++) {
149795a0b26dSToby Isaac           vals[i] = 0.;
149895a0b26dSToby Isaac         }
149995a0b26dSToby Isaac         matOffset = ia[cOff];
150095a0b26dSToby Isaac         numFillCols = ia[cOff+1] - matOffset;
150195a0b26dSToby Isaac         pointMat = refPointFieldMats[childid-pRefStart][f];
150295a0b26dSToby Isaac         numCols = refPointFieldN[childid-pRefStart][f];
150395a0b26dSToby Isaac         offset = 0;
150495a0b26dSToby Isaac         for (i = 0; i < closureSize; i++) {
150595a0b26dSToby Isaac           PetscInt q = closure[2*i];
150695a0b26dSToby Isaac           PetscInt o = closure[2*i+1];
150795a0b26dSToby Isaac           PetscInt aDof, aOff, j, k, qConDof, qConOff;
150895a0b26dSToby Isaac 
150995a0b26dSToby Isaac           qConDof = qConOff = 0;
151095a0b26dSToby Isaac           if (numFields > 1) {
151195a0b26dSToby Isaac             ierr = PetscSectionGetFieldDof(section,q,f,&aDof);CHKERRQ(ierr);
151295a0b26dSToby Isaac             ierr = PetscSectionGetFieldOffset(section,q,f,&aOff);CHKERRQ(ierr);
151395a0b26dSToby Isaac             if (q >= conStart && q < conEnd) {
151495a0b26dSToby Isaac               ierr = PetscSectionGetFieldDof(conSec,q,f,&qConDof);CHKERRQ(ierr);
151595a0b26dSToby Isaac               ierr = PetscSectionGetFieldOffset(conSec,q,f,&qConOff);CHKERRQ(ierr);
151695a0b26dSToby Isaac             }
151795a0b26dSToby Isaac           }
151895a0b26dSToby Isaac           else {
151995a0b26dSToby Isaac             ierr = PetscSectionGetDof(section,q,&aDof);CHKERRQ(ierr);
152095a0b26dSToby Isaac             ierr = PetscSectionGetOffset(section,q,&aOff);CHKERRQ(ierr);
152195a0b26dSToby Isaac             if (q >= conStart && q < conEnd) {
152295a0b26dSToby Isaac               ierr = PetscSectionGetDof(conSec,q,&qConDof);CHKERRQ(ierr);
152395a0b26dSToby Isaac               ierr = PetscSectionGetOffset(conSec,q,&qConOff);CHKERRQ(ierr);
152495a0b26dSToby Isaac             }
152595a0b26dSToby Isaac           }
152695a0b26dSToby Isaac           if (!aDof) continue;
152795a0b26dSToby Isaac           if (qConDof) {
152895a0b26dSToby Isaac             /* this point has anchors: its rows of the matrix should already
152995a0b26dSToby Isaac              * be filled, thanks to preordering */
153095a0b26dSToby Isaac             /* first multiply into pointWork, then set in matrix */
153195a0b26dSToby Isaac             PetscInt aMatOffset = ia[qConOff];
153295a0b26dSToby Isaac             PetscInt aNumFillCols = ia[qConOff + 1] - aMatOffset;
153395a0b26dSToby Isaac             for (r = 0; r < cDof; r++) {
153495a0b26dSToby Isaac               for (j = 0; j < aNumFillCols; j++) {
153595a0b26dSToby Isaac                 PetscScalar inVal = 0;
153695a0b26dSToby Isaac                 for (k = 0; k < aDof; k++) {
153795a0b26dSToby Isaac                   PetscInt node = (o >= 0) ? (k / fComp) : ((aDof - 1 - k) / fComp);
153895a0b26dSToby Isaac                   PetscInt comp = (k % fComp);
153995a0b26dSToby Isaac                   PetscInt col  = node * fComp + comp;
154095a0b26dSToby Isaac 
154195a0b26dSToby Isaac                   inVal += pointMat[r * numCols + offset + col] * vals[aMatOffset + aNumFillCols * k + j];
154295a0b26dSToby Isaac                 }
154395a0b26dSToby Isaac                 pointWork[r * aNumFillCols + j] = inVal;
154495a0b26dSToby Isaac               }
154595a0b26dSToby Isaac             }
154695a0b26dSToby Isaac             /* assume that the columns are sorted, spend less time searching */
154795a0b26dSToby Isaac             for (j = 0, k = 0; j < aNumFillCols; j++) {
154895a0b26dSToby Isaac               PetscInt col = ja[aMatOffset + j];
154995a0b26dSToby Isaac               for (;k < numFillCols; k++) {
155095a0b26dSToby Isaac                 if (ja[matOffset + k] == col) {
155195a0b26dSToby Isaac                   break;
155295a0b26dSToby Isaac                 }
155395a0b26dSToby Isaac               }
155495a0b26dSToby Isaac               if (k == numFillCols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"No nonzero space for (%d, %d)", cOff, col);
155595a0b26dSToby Isaac               for (r = 0; r < cDof; r++) {
155695a0b26dSToby Isaac                 vals[matOffset + numFillCols * r + k] = pointWork[r * aNumFillCols + j];
155795a0b26dSToby Isaac               }
155895a0b26dSToby Isaac             }
155995a0b26dSToby Isaac           }
156095a0b26dSToby Isaac           else {
156195a0b26dSToby Isaac             /* find where to put this portion of pointMat into the matrix */
156295a0b26dSToby Isaac             for (k = 0; k < numFillCols; k++) {
156395a0b26dSToby Isaac               if (ja[matOffset + k] == aOff) {
156495a0b26dSToby Isaac                 break;
156595a0b26dSToby Isaac               }
156695a0b26dSToby Isaac             }
156795a0b26dSToby Isaac             if (k == numFillCols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"No nonzero space for (%d, %d)", cOff, aOff);
156895a0b26dSToby Isaac             for (j = 0; j < aDof; j++) {
156995a0b26dSToby Isaac               PetscInt node = (o >= 0) ? (j / fComp) : ((aDof - 1 - j) / fComp);
157095a0b26dSToby Isaac               PetscInt comp = (j % fComp);
157195a0b26dSToby Isaac               PetscInt col  = node * fComp + comp;
157295a0b26dSToby Isaac               for (r = 0; r < cDof; r++) {
157395a0b26dSToby Isaac                 vals[matOffset + numFillCols * r + k + col] += pointMat[r * numCols + offset + col];
157495a0b26dSToby Isaac               }
157595a0b26dSToby Isaac             }
157695a0b26dSToby Isaac           }
157795a0b26dSToby Isaac           offset += aDof;
157895a0b26dSToby Isaac         }
157995a0b26dSToby Isaac       }
158095a0b26dSToby Isaac       ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr);
158195a0b26dSToby Isaac     }
158295a0b26dSToby Isaac     for (row = 0; row < nRows; row++) {
158395a0b26dSToby Isaac       ierr = MatSetValues(cMat,1,&row,ia[row+1]-ia[row],&ja[ia[row]],&vals[ia[row]],INSERT_VALUES);CHKERRQ(ierr);
158495a0b26dSToby Isaac     }
158595a0b26dSToby Isaac     ierr = MatRestoreRowIJ(cMat,0,PETSC_FALSE,PETSC_FALSE,&nRows,&ia,&ja,&done);CHKERRQ(ierr);
158695a0b26dSToby Isaac     if (!done) SETERRQ(PetscObjectComm((PetscObject)cMat),PETSC_ERR_PLIB,"Could not restore RowIJ of constraint matrix");
158795a0b26dSToby Isaac     ierr = MatAssemblyBegin(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
158895a0b26dSToby Isaac     ierr = MatAssemblyEnd(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
158995a0b26dSToby Isaac     ierr = PetscFree(vals);CHKERRQ(ierr);
159095a0b26dSToby Isaac   }
159195a0b26dSToby Isaac 
159295a0b26dSToby Isaac   /* clean up */
159395a0b26dSToby Isaac   ierr = ISRestoreIndices(refAnIS,&refAnchors);CHKERRQ(ierr);
159495a0b26dSToby Isaac   ierr = ISRestoreIndices(anIS,&anchors);CHKERRQ(ierr);
159595a0b26dSToby Isaac   ierr = PetscFree2(perm,iperm);CHKERRQ(ierr);
159695a0b26dSToby Isaac   ierr = PetscFree(rows);CHKERRQ(ierr);
159795a0b26dSToby Isaac   ierr = PetscFree(cols);CHKERRQ(ierr);
159895a0b26dSToby Isaac   ierr = PetscFree(pointWork);CHKERRQ(ierr);
159995a0b26dSToby Isaac   for (p = pRefStart; p < pRefEnd; p++) {
160095a0b26dSToby Isaac     PetscInt parent, pDof;
160195a0b26dSToby Isaac 
160295a0b26dSToby Isaac     ierr = DMPlexGetTreeParent(refTree,p,&parent,NULL);CHKERRQ(ierr);
160395a0b26dSToby Isaac     ierr = PetscSectionGetDof(refConSec,p,&pDof);CHKERRQ(ierr);
160495a0b26dSToby Isaac     if (!pDof || parent == p) continue;
160595a0b26dSToby Isaac 
160695a0b26dSToby Isaac     for (f = 0; f < numFields; f++) {
160795a0b26dSToby Isaac       PetscInt cDof;
160895a0b26dSToby Isaac 
160995a0b26dSToby Isaac       if (numFields > 1) {
161095a0b26dSToby Isaac         ierr = PetscSectionGetFieldDof(refConSec,p,f,&cDof);CHKERRQ(ierr);
161195a0b26dSToby Isaac       }
161295a0b26dSToby Isaac       else {
161395a0b26dSToby Isaac         ierr = PetscSectionGetDof(refConSec,p,&cDof);CHKERRQ(ierr);
161495a0b26dSToby Isaac       }
161595a0b26dSToby Isaac 
161695a0b26dSToby Isaac       if (!cDof) continue;
161795a0b26dSToby Isaac       ierr = PetscFree(refPointFieldMats[p - pRefStart][f]);CHKERRQ(ierr);
161895a0b26dSToby Isaac     }
161995a0b26dSToby Isaac     ierr = PetscFree(refPointFieldMats[p - pRefStart]);CHKERRQ(ierr);
162095a0b26dSToby Isaac     ierr = PetscFree(refPointFieldN[p - pRefStart]);CHKERRQ(ierr);
162195a0b26dSToby Isaac   }
162295a0b26dSToby Isaac   ierr = PetscFree(refPointFieldMats);CHKERRQ(ierr);
162395a0b26dSToby Isaac   ierr = PetscFree(refPointFieldN);CHKERRQ(ierr);
162495a0b26dSToby Isaac   PetscFunctionReturn(0);
162595a0b26dSToby Isaac }
162695a0b26dSToby Isaac 
16276f5f1567SToby Isaac #undef __FUNCT__
16286f5f1567SToby Isaac #define __FUNCT__ "DMPlexTreeRefineCell"
16296f5f1567SToby Isaac /* refine a single cell on rank 0: this is not intended to provide good local refinement, only to create an example of
16306f5f1567SToby Isaac  * a non-conforming mesh.  Local refinement comes later */
16316f5f1567SToby Isaac PetscErrorCode DMPlexTreeRefineCell (DM dm, PetscInt cell, DM *ncdm)
16326f5f1567SToby Isaac {
16336f5f1567SToby Isaac   DM K;
16346f5f1567SToby Isaac   PetscInt rank;
16356f5f1567SToby Isaac   PetscInt dim, *pNewStart, *pNewEnd, *pNewCount, *pOldStart, *pOldEnd, offset, d, pStart, pEnd;
16366f5f1567SToby Isaac   PetscInt numNewCones, *newConeSizes, *newCones, *newOrientations;
16376f5f1567SToby Isaac   PetscInt *Kembedding;
16386f5f1567SToby Isaac   PetscInt *cellClosure=NULL, nc;
16396f5f1567SToby Isaac   PetscScalar *newVertexCoords;
16406f5f1567SToby Isaac   PetscInt numPointsWithParents, *parents, *childIDs, *perm, *iperm, *preOrient, pOffset;
16416f5f1567SToby Isaac   PetscSection parentSection;
16426f5f1567SToby Isaac   PetscErrorCode ierr;
16436f5f1567SToby Isaac 
16446f5f1567SToby Isaac   PetscFunctionBegin;
16456f5f1567SToby Isaac   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr);
16466f5f1567SToby Isaac   ierr = DMPlexGetDimension(dm,&dim);CHKERRQ(ierr);
16476f5f1567SToby Isaac   ierr = DMPlexCreate(PetscObjectComm((PetscObject)dm), ncdm);CHKERRQ(ierr);
16486f5f1567SToby Isaac   ierr = DMPlexSetDimension(*ncdm,dim);CHKERRQ(ierr);
16496f5f1567SToby Isaac 
16506f5f1567SToby Isaac   ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr);
16516f5f1567SToby Isaac   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&parentSection);CHKERRQ(ierr);
16526f5f1567SToby Isaac   ierr = DMPlexGetReferenceTree(dm,&K);CHKERRQ(ierr);
16536f5f1567SToby Isaac   if (!rank) {
16546f5f1567SToby Isaac     /* compute the new charts */
16556f5f1567SToby Isaac     ierr = PetscMalloc5(dim+1,&pNewCount,dim+1,&pNewStart,dim+1,&pNewEnd,dim+1,&pOldStart,dim+1,&pOldEnd);CHKERRQ(ierr);
16566f5f1567SToby Isaac     offset = 0;
16576f5f1567SToby Isaac     for (d = 0; d <= dim; d++) {
16586f5f1567SToby Isaac       PetscInt pOldCount, kStart, kEnd, k;
16596f5f1567SToby Isaac 
16606f5f1567SToby Isaac       pNewStart[d] = offset;
16616f5f1567SToby Isaac       ierr = DMPlexGetHeightStratum(dm,d,&pOldStart[d],&pOldEnd[d]);CHKERRQ(ierr);
16626f5f1567SToby Isaac       ierr = DMPlexGetHeightStratum(K,d,&kStart,&kEnd);CHKERRQ(ierr);
16636f5f1567SToby Isaac       pOldCount = pOldEnd[d] - pOldStart[d];
16646f5f1567SToby Isaac       /* adding the new points */
16656f5f1567SToby Isaac       pNewCount[d] = pOldCount + kEnd - kStart;
16666f5f1567SToby Isaac       if (!d) {
16676f5f1567SToby Isaac         /* removing the cell */
16686f5f1567SToby Isaac         pNewCount[d]--;
16696f5f1567SToby Isaac       }
16706f5f1567SToby Isaac       for (k = kStart; k < kEnd; k++) {
16716f5f1567SToby Isaac         PetscInt parent;
16726f5f1567SToby Isaac         ierr = DMPlexGetTreeParent(K,k,&parent,NULL);CHKERRQ(ierr);
16736f5f1567SToby Isaac         if (parent == k) {
16746f5f1567SToby Isaac           /* avoid double counting points that won't actually be new */
16756f5f1567SToby Isaac           pNewCount[d]--;
16766f5f1567SToby Isaac         }
16776f5f1567SToby Isaac       }
16786f5f1567SToby Isaac       pNewEnd[d] = pNewStart[d] + pNewCount[d];
16796f5f1567SToby Isaac       offset = pNewEnd[d];
16806f5f1567SToby Isaac 
16816f5f1567SToby Isaac     }
16826f5f1567SToby Isaac     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]);
16836f5f1567SToby Isaac     /* get the current closure of the cell that we are removing */
16846f5f1567SToby Isaac     ierr = DMPlexGetTransitiveClosure(dm,cell,PETSC_TRUE,&nc,&cellClosure);CHKERRQ(ierr);
16856f5f1567SToby Isaac 
16866f5f1567SToby Isaac     ierr = PetscMalloc1(pNewEnd[dim],&newConeSizes);CHKERRQ(ierr);
16876f5f1567SToby Isaac     {
16886f5f1567SToby Isaac       PetscInt kStart, kEnd, k, closureSizeK, *closureK = NULL, j;
16896f5f1567SToby Isaac 
16906f5f1567SToby Isaac       ierr = DMPlexGetChart(K,&kStart,&kEnd);CHKERRQ(ierr);
16916f5f1567SToby Isaac       ierr = PetscMalloc4(kEnd-kStart,&Kembedding,kEnd-kStart,&perm,kEnd-kStart,&iperm,kEnd-kStart,&preOrient);CHKERRQ(ierr);
16926f5f1567SToby Isaac 
16936f5f1567SToby Isaac       for (k = kStart; k < kEnd; k++) {
16946f5f1567SToby Isaac         perm[k - kStart] = k;
16956f5f1567SToby Isaac         iperm [k - kStart] = k - kStart;
16966f5f1567SToby Isaac         preOrient[k - kStart] = 0;
16976f5f1567SToby Isaac       }
16986f5f1567SToby Isaac 
16996f5f1567SToby Isaac       ierr = DMPlexGetTransitiveClosure(K,0,PETSC_TRUE,&closureSizeK,&closureK);CHKERRQ(ierr);
17006f5f1567SToby Isaac       for (j = 1; j < closureSizeK; j++) {
17016f5f1567SToby Isaac         PetscInt parentOrientA = closureK[2*j+1];
17026f5f1567SToby Isaac         PetscInt parentOrientB = cellClosure[2*j+1];
17036f5f1567SToby Isaac         PetscInt p, q;
17046f5f1567SToby Isaac 
17056f5f1567SToby Isaac         p = closureK[2*j];
17066f5f1567SToby Isaac         q = cellClosure[2*j];
17076f5f1567SToby Isaac         for (d = 0; d <= dim; d++) {
17086f5f1567SToby Isaac           if (q >= pOldStart[d] && q < pOldEnd[d]) {
17096f5f1567SToby Isaac             Kembedding[p] = (q - pOldStart[d]) + pNewStart[d];
17106f5f1567SToby Isaac           }
17116f5f1567SToby Isaac         }
17126f5f1567SToby Isaac         if (parentOrientA != parentOrientB) {
17136f5f1567SToby Isaac           PetscInt numChildren, i;
17146f5f1567SToby Isaac           const PetscInt *children;
17156f5f1567SToby Isaac 
17166f5f1567SToby Isaac           ierr = DMPlexGetTreeChildren(K,p,&numChildren,&children);CHKERRQ(ierr);
17176f5f1567SToby Isaac           for (i = 0; i < numChildren; i++) {
17186f5f1567SToby Isaac             PetscInt kPerm, oPerm;
17196f5f1567SToby Isaac 
17206f5f1567SToby Isaac             k    = children[i];
17216f5f1567SToby Isaac             ierr = DMPlexReferenceTreeGetChildSymmetry(K,p,parentOrientA,0,k,parentOrientB,&oPerm,&kPerm);CHKERRQ(ierr);
17226f5f1567SToby Isaac             /* perm = what refTree position I'm in */
17236f5f1567SToby Isaac             perm[kPerm-kStart]      = k;
17246f5f1567SToby Isaac             /* iperm = who is at this position */
17256f5f1567SToby Isaac             iperm[k-kStart]         = kPerm-kStart;
17266f5f1567SToby Isaac             preOrient[kPerm-kStart] = oPerm;
17276f5f1567SToby Isaac           }
17286f5f1567SToby Isaac         }
17296f5f1567SToby Isaac       }
17306f5f1567SToby Isaac       ierr = DMPlexRestoreTransitiveClosure(K,0,PETSC_TRUE,&closureSizeK,&closureK);CHKERRQ(ierr);
17316f5f1567SToby Isaac     }
17326f5f1567SToby Isaac     ierr = PetscSectionSetChart(parentSection,0,pNewEnd[dim]);CHKERRQ(ierr);
17336f5f1567SToby Isaac     offset = 0;
17346f5f1567SToby Isaac     numNewCones = 0;
17356f5f1567SToby Isaac     for (d = 0; d <= dim; d++) {
17366f5f1567SToby Isaac       PetscInt kStart, kEnd, k;
17376f5f1567SToby Isaac       PetscInt p;
17386f5f1567SToby Isaac       PetscInt size;
17396f5f1567SToby Isaac 
17406f5f1567SToby Isaac       for (p = pOldStart[d]; p < pOldEnd[d]; p++) {
17416f5f1567SToby Isaac         /* skip cell 0 */
17426f5f1567SToby Isaac         if (p == cell) continue;
17436f5f1567SToby Isaac         /* old cones to new cones */
17446f5f1567SToby Isaac         ierr = DMPlexGetConeSize(dm,p,&size);CHKERRQ(ierr);
17456f5f1567SToby Isaac         newConeSizes[offset++] = size;
17466f5f1567SToby Isaac         numNewCones += size;
17476f5f1567SToby Isaac       }
17486f5f1567SToby Isaac 
17496f5f1567SToby Isaac       ierr = DMPlexGetHeightStratum(K,d,&kStart,&kEnd);CHKERRQ(ierr);
17506f5f1567SToby Isaac       for (k = kStart; k < kEnd; k++) {
17516f5f1567SToby Isaac         PetscInt kParent;
17526f5f1567SToby Isaac 
17536f5f1567SToby Isaac         ierr = DMPlexGetTreeParent(K,k,&kParent,NULL);CHKERRQ(ierr);
17546f5f1567SToby Isaac         if (kParent != k) {
17556f5f1567SToby Isaac           Kembedding[k] = offset;
17566f5f1567SToby Isaac           ierr = DMPlexGetConeSize(K,k,&size);CHKERRQ(ierr);
17576f5f1567SToby Isaac           newConeSizes[offset++] = size;
17586f5f1567SToby Isaac           numNewCones += size;
17596f5f1567SToby Isaac           if (kParent != 0) {
17606f5f1567SToby Isaac             ierr = PetscSectionSetDof(parentSection,Kembedding[k],1);CHKERRQ(ierr);
17616f5f1567SToby Isaac           }
17626f5f1567SToby Isaac         }
17636f5f1567SToby Isaac       }
17646f5f1567SToby Isaac     }
17656f5f1567SToby Isaac 
17666f5f1567SToby Isaac     ierr = PetscSectionSetUp(parentSection);CHKERRQ(ierr);
17676f5f1567SToby Isaac     ierr = PetscSectionGetStorageSize(parentSection,&numPointsWithParents);CHKERRQ(ierr);
17686f5f1567SToby Isaac     ierr = PetscMalloc2(numNewCones,&newCones,numNewCones,&newOrientations);CHKERRQ(ierr);
17696f5f1567SToby Isaac     ierr = PetscMalloc2(numPointsWithParents,&parents,numPointsWithParents,&childIDs);CHKERRQ(ierr);
17706f5f1567SToby Isaac 
17716f5f1567SToby Isaac     /* fill new cones */
17726f5f1567SToby Isaac     offset = 0;
17736f5f1567SToby Isaac     for (d = 0; d <= dim; d++) {
17746f5f1567SToby Isaac       PetscInt kStart, kEnd, k, l;
17756f5f1567SToby Isaac       PetscInt p;
17766f5f1567SToby Isaac       PetscInt size;
17776f5f1567SToby Isaac       const PetscInt *cone, *orientation;
17786f5f1567SToby Isaac 
17796f5f1567SToby Isaac       for (p = pOldStart[d]; p < pOldEnd[d]; p++) {
17806f5f1567SToby Isaac         /* skip cell 0 */
17816f5f1567SToby Isaac         if (p == cell) continue;
17826f5f1567SToby Isaac         /* old cones to new cones */
17836f5f1567SToby Isaac         ierr = DMPlexGetConeSize(dm,p,&size);CHKERRQ(ierr);
17846f5f1567SToby Isaac         ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr);
17856f5f1567SToby Isaac         ierr = DMPlexGetConeOrientation(dm,p,&orientation);CHKERRQ(ierr);
17866f5f1567SToby Isaac         for (l = 0; l < size; l++) {
17876f5f1567SToby Isaac           newCones[offset]          = (cone[l] - pOldStart[d + 1]) + pNewStart[d + 1];
17886f5f1567SToby Isaac           newOrientations[offset++] = orientation[l];
17896f5f1567SToby Isaac         }
17906f5f1567SToby Isaac       }
17916f5f1567SToby Isaac 
17926f5f1567SToby Isaac       ierr = DMPlexGetHeightStratum(K,d,&kStart,&kEnd);CHKERRQ(ierr);
17936f5f1567SToby Isaac       for (k = kStart; k < kEnd; k++) {
17946f5f1567SToby Isaac         PetscInt kPerm = perm[k], kParent;
17956f5f1567SToby Isaac         PetscInt preO  = preOrient[k];
17966f5f1567SToby Isaac 
17976f5f1567SToby Isaac         ierr = DMPlexGetTreeParent(K,k,&kParent,NULL);CHKERRQ(ierr);
17986f5f1567SToby Isaac         if (kParent != k) {
17996f5f1567SToby Isaac           /* embed new cones */
18006f5f1567SToby Isaac           ierr = DMPlexGetConeSize(K,k,&size);CHKERRQ(ierr);
18016f5f1567SToby Isaac           ierr = DMPlexGetCone(K,kPerm,&cone);CHKERRQ(ierr);
18026f5f1567SToby Isaac           ierr = DMPlexGetConeOrientation(K,kPerm,&orientation);CHKERRQ(ierr);
18036f5f1567SToby Isaac           for (l = 0; l < size; l++) {
18046f5f1567SToby Isaac             PetscInt q, m = (preO >= 0) ? ((preO + l) % size) : ((size -(preO + 1) - l) % size);
18056f5f1567SToby Isaac             PetscInt newO, lSize, oTrue;
18066f5f1567SToby Isaac 
18076f5f1567SToby Isaac             q                         = iperm[cone[m]];
18086f5f1567SToby Isaac             newCones[offset]          = Kembedding[q];
18096f5f1567SToby Isaac             ierr                      = DMPlexGetConeSize(K,q,&lSize);CHKERRQ(ierr);
18106f5f1567SToby Isaac             oTrue                     = orientation[m];
18116f5f1567SToby Isaac             oTrue                     = ((!lSize) || (preOrient[k] >= 0)) ? oTrue : -(oTrue + 2);
18126f5f1567SToby Isaac             newO                      = DihedralCompose(lSize,oTrue,preOrient[q]);
18136f5f1567SToby Isaac             newOrientations[offset++] = newO;
18146f5f1567SToby Isaac           }
18156f5f1567SToby Isaac           if (kParent != 0) {
18166f5f1567SToby Isaac             PetscInt newPoint = Kembedding[kParent];
18176f5f1567SToby Isaac             ierr              = PetscSectionGetOffset(parentSection,Kembedding[k],&pOffset);CHKERRQ(ierr);
18186f5f1567SToby Isaac             parents[pOffset]  = newPoint;
18196f5f1567SToby Isaac             childIDs[pOffset] = k;
18206f5f1567SToby Isaac           }
18216f5f1567SToby Isaac         }
18226f5f1567SToby Isaac       }
18236f5f1567SToby Isaac     }
18246f5f1567SToby Isaac 
18256f5f1567SToby Isaac     ierr = PetscMalloc1(dim*(pNewEnd[dim]-pNewStart[dim]),&newVertexCoords);CHKERRQ(ierr);
18266f5f1567SToby Isaac 
18276f5f1567SToby Isaac     /* fill coordinates */
18286f5f1567SToby Isaac     offset = 0;
18296f5f1567SToby Isaac     {
18306f5f1567SToby Isaac       PetscInt k, kStart, kEnd, l;
18316f5f1567SToby Isaac       PetscSection vSection;
18326f5f1567SToby Isaac       PetscInt v;
18336f5f1567SToby Isaac       Vec coords;
18346f5f1567SToby Isaac       PetscScalar *coordvals;
18356f5f1567SToby Isaac       PetscInt dof, off;
18366f5f1567SToby Isaac       PetscScalar v0[3], J[9], detJ;
18376f5f1567SToby Isaac 
18386f5f1567SToby Isaac #if defined(PETSC_USE_DEBUG)
18396f5f1567SToby Isaac       ierr = DMPlexGetHeightStratum(K,0,&kStart,&kEnd);CHKERRQ(ierr);
18406f5f1567SToby Isaac       for (k = kStart; k < kEnd; k++) {
18416f5f1567SToby Isaac         ierr = DMPlexComputeCellGeometry(K, k, v0, J, NULL, &detJ);CHKERRQ(ierr);
18426f5f1567SToby Isaac         if (detJ <= 0.) SETERRQ1 (PETSC_COMM_SELF,PETSC_ERR_PLIB,"reference tree cell %d has bad determinant",k);
18436f5f1567SToby Isaac       }
18446f5f1567SToby Isaac #endif
18456f5f1567SToby Isaac       ierr = DMPlexComputeCellGeometry(dm, cell, v0, J, NULL, &detJ);CHKERRQ(ierr);
18466f5f1567SToby Isaac       ierr = DMGetCoordinateSection(dm,&vSection);CHKERRQ(ierr);
18476f5f1567SToby Isaac       ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr);
18486f5f1567SToby Isaac       ierr = VecGetArray(coords,&coordvals);CHKERRQ(ierr);
18496f5f1567SToby Isaac       for (v = pOldStart[dim]; v < pOldEnd[dim]; v++) {
18506f5f1567SToby Isaac 
18516f5f1567SToby Isaac         ierr = PetscSectionGetDof(vSection,v,&dof);CHKERRQ(ierr);
18526f5f1567SToby Isaac         ierr = PetscSectionGetOffset(vSection,v,&off);CHKERRQ(ierr);
18536f5f1567SToby Isaac         for (l = 0; l < dof; l++) {
18546f5f1567SToby Isaac           newVertexCoords[offset++] = coordvals[off + l];
18556f5f1567SToby Isaac         }
18566f5f1567SToby Isaac       }
18576f5f1567SToby Isaac       ierr = VecRestoreArray(coords,&coordvals);CHKERRQ(ierr);
18586f5f1567SToby Isaac 
18596f5f1567SToby Isaac       ierr = DMGetCoordinateSection(K,&vSection);CHKERRQ(ierr);
18606f5f1567SToby Isaac       ierr = DMGetCoordinatesLocal(K,&coords);CHKERRQ(ierr);
18616f5f1567SToby Isaac       ierr = VecGetArray(coords,&coordvals);CHKERRQ(ierr);
18626f5f1567SToby Isaac       ierr = DMPlexGetDepthStratum(K,0,&kStart,&kEnd);CHKERRQ(ierr);
18636f5f1567SToby Isaac       for (v = kStart; v < kEnd; v++) {
18646f5f1567SToby Isaac         PetscInt vPerm = perm[v];
18656f5f1567SToby Isaac         PetscInt kParent;
18666f5f1567SToby Isaac 
18676f5f1567SToby Isaac         ierr = DMPlexGetTreeParent(K,v,&kParent,NULL);CHKERRQ(ierr);
18686f5f1567SToby Isaac         if (kParent != v) {
18696f5f1567SToby Isaac           /* this is a new vertex */
18706f5f1567SToby Isaac           ierr = PetscSectionGetOffset(vSection,vPerm,&off);CHKERRQ(ierr);
18716f5f1567SToby Isaac           CoordinatesRefToReal(dim, dim, v0, J, &coordvals[off],&newVertexCoords[offset]);CHKERRQ(ierr);
18726f5f1567SToby Isaac           offset += dim;
18736f5f1567SToby Isaac         }
18746f5f1567SToby Isaac       }
18756f5f1567SToby Isaac       ierr = VecRestoreArray(coords,&coordvals);CHKERRQ(ierr);
18766f5f1567SToby Isaac     }
18776f5f1567SToby Isaac 
18786f5f1567SToby Isaac     /* need to reverse the order of pNewCount: vertices first, cells last */
18796f5f1567SToby Isaac     for (d = 0; d < (dim + 1) / 2; d++) {
18806f5f1567SToby Isaac       PetscInt tmp;
18816f5f1567SToby Isaac 
18826f5f1567SToby Isaac       tmp = pNewCount[d];
18836f5f1567SToby Isaac       pNewCount[d] = pNewCount[dim - d];
18846f5f1567SToby Isaac       pNewCount[dim - d] = tmp;
18856f5f1567SToby Isaac     }
18866f5f1567SToby Isaac 
18876f5f1567SToby Isaac     ierr = DMPlexCreateFromDAG(*ncdm,dim,pNewCount,newConeSizes,newCones,newOrientations,newVertexCoords);CHKERRQ(ierr);
18886f5f1567SToby Isaac     ierr = DMPlexSetReferenceTree(*ncdm,K);CHKERRQ(ierr);
18896f5f1567SToby Isaac     ierr = DMPlexSetTree(*ncdm,parentSection,parents,childIDs);CHKERRQ(ierr);
18906f5f1567SToby Isaac 
18916f5f1567SToby Isaac     /* clean up */
18926f5f1567SToby Isaac     ierr = DMPlexRestoreTransitiveClosure(dm,cell,PETSC_TRUE,&nc,&cellClosure);CHKERRQ(ierr);
18936f5f1567SToby Isaac     ierr = PetscFree5(pNewCount,pNewStart,pNewEnd,pOldStart,pOldEnd);CHKERRQ(ierr);
18946f5f1567SToby Isaac     ierr = PetscFree(newConeSizes);CHKERRQ(ierr);
18956f5f1567SToby Isaac     ierr = PetscFree2(newCones,newOrientations);CHKERRQ(ierr);
18966f5f1567SToby Isaac     ierr = PetscFree(newVertexCoords);CHKERRQ(ierr);
18976f5f1567SToby Isaac     ierr = PetscFree2(parents,childIDs);CHKERRQ(ierr);
18986f5f1567SToby Isaac     ierr = PetscFree4(Kembedding,perm,iperm,preOrient);CHKERRQ(ierr);
18996f5f1567SToby Isaac   }
19006f5f1567SToby Isaac   else {
19016f5f1567SToby Isaac     PetscInt    p, counts[4];
19026f5f1567SToby Isaac     PetscInt    *coneSizes, *cones, *orientations;
19036f5f1567SToby Isaac     Vec         coordVec;
19046f5f1567SToby Isaac     PetscScalar *coords;
19056f5f1567SToby Isaac 
19066f5f1567SToby Isaac     for (d = 0; d <= dim; d++) {
19076f5f1567SToby Isaac       PetscInt dStart, dEnd;
19086f5f1567SToby Isaac 
19096f5f1567SToby Isaac       ierr = DMPlexGetDepthStratum(dm,d,&dStart,&dEnd);CHKERRQ(ierr);
19106f5f1567SToby Isaac       counts[d] = dEnd - dStart;
19116f5f1567SToby Isaac     }
19126f5f1567SToby Isaac     ierr = PetscMalloc1(pEnd-pStart,&coneSizes);CHKERRQ(ierr);
19136f5f1567SToby Isaac     for (p = pStart; p < pEnd; p++) {
19146f5f1567SToby Isaac       ierr = DMPlexGetConeSize(dm,p,&coneSizes[p-pStart]);CHKERRQ(ierr);
19156f5f1567SToby Isaac     }
19166f5f1567SToby Isaac     ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr);
19176f5f1567SToby Isaac     ierr = DMPlexGetConeOrientations(dm, &orientations);CHKERRQ(ierr);
19186f5f1567SToby Isaac     ierr = DMGetCoordinatesLocal(dm,&coordVec);CHKERRQ(ierr);
19196f5f1567SToby Isaac     ierr = VecGetArray(coordVec,&coords);CHKERRQ(ierr);
19206f5f1567SToby Isaac 
19216f5f1567SToby Isaac     ierr = PetscSectionSetChart(parentSection,pStart,pEnd);CHKERRQ(ierr);
19226f5f1567SToby Isaac     ierr = PetscSectionSetUp(parentSection);CHKERRQ(ierr);
19236f5f1567SToby Isaac     ierr = DMPlexCreateFromDAG(*ncdm,dim,counts,coneSizes,cones,orientations,NULL);CHKERRQ(ierr);
19246f5f1567SToby Isaac     ierr = DMPlexSetReferenceTree(*ncdm,K);CHKERRQ(ierr);
19256f5f1567SToby Isaac     ierr = DMPlexSetTree(*ncdm,parentSection,NULL,NULL);CHKERRQ(ierr);
19266f5f1567SToby Isaac     ierr = VecRestoreArray(coordVec,&coords);CHKERRQ(ierr);
19276f5f1567SToby Isaac   }
19286f5f1567SToby Isaac   ierr = PetscSectionDestroy(&parentSection);CHKERRQ(ierr);
19296f5f1567SToby Isaac 
19306f5f1567SToby Isaac   PetscFunctionReturn(0);
19316f5f1567SToby Isaac }
1932