1af0996ceSBarry Smith #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2d6a7ad0dSToby Isaac #include <../src/sys/utils/hash.h> 3af0996ceSBarry Smith #include <petsc/private/isimpl.h> 4af0996ceSBarry Smith #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); 3247a1df27SMatthew G. Knepley if (ref) {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 221776742edSToby Isaac static PetscErrorCode DMPlexSetTree_Internal(DM,PetscSection,PetscInt*,PetscInt*,PetscBool,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; 256e228b242SToby Isaac #if 1 257e228b242SToby Isaac comm = PETSC_COMM_SELF; 258e228b242SToby Isaac #endif 259da43764aSToby Isaac /* create a reference element */ 260da43764aSToby Isaac ierr = DMPlexCreateReferenceCell(comm, dim, simplex, &K);CHKERRQ(ierr); 261da43764aSToby Isaac ierr = DMPlexCreateLabel(K, "identity");CHKERRQ(ierr); 262da43764aSToby Isaac ierr = DMPlexGetLabel(K, "identity", &identity);CHKERRQ(ierr); 263da43764aSToby Isaac ierr = DMPlexGetChart(K, &pStart, &pEnd);CHKERRQ(ierr); 264da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 265da43764aSToby Isaac ierr = DMLabelSetValue(identity, p, p);CHKERRQ(ierr); 266da43764aSToby Isaac } 267da43764aSToby Isaac /* refine it */ 268da43764aSToby Isaac ierr = DMRefine(K,comm,&Kref);CHKERRQ(ierr); 269da43764aSToby Isaac 270da43764aSToby Isaac /* the reference tree is the union of these two, without duplicating 271da43764aSToby Isaac * points that appear in both */ 272da43764aSToby Isaac ierr = DMPlexGetLabel(Kref, "identity", &identityRef);CHKERRQ(ierr); 273da43764aSToby Isaac ierr = DMPlexGetChart(Kref, &pRefStart, &pRefEnd);CHKERRQ(ierr); 274da43764aSToby Isaac ierr = PetscSectionCreate(comm, &unionSection);CHKERRQ(ierr); 275da43764aSToby Isaac ierr = PetscSectionSetChart(unionSection, 0, (pEnd - pStart) + (pRefEnd - pRefStart));CHKERRQ(ierr); 276da43764aSToby Isaac /* count points that will go in the union */ 277da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 278da43764aSToby Isaac ierr = PetscSectionSetDof(unionSection, p - pStart, 1);CHKERRQ(ierr); 279da43764aSToby Isaac } 280da43764aSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 281da43764aSToby Isaac PetscInt q, qSize; 282da43764aSToby Isaac ierr = DMLabelGetValue(identityRef, p, &q);CHKERRQ(ierr); 283da43764aSToby Isaac ierr = DMLabelGetStratumSize(identityRef, q, &qSize);CHKERRQ(ierr); 284da43764aSToby Isaac if (qSize > 1) { 285da43764aSToby Isaac ierr = PetscSectionSetDof(unionSection, p - pRefStart + (pEnd - pStart), 1);CHKERRQ(ierr); 286da43764aSToby Isaac } 287da43764aSToby Isaac } 288854ce69bSBarry Smith ierr = PetscMalloc1(pEnd - pStart + pRefEnd - pRefStart,&permvals);CHKERRQ(ierr); 289da43764aSToby Isaac offset = 0; 290da43764aSToby Isaac /* stratify points in the union by topological dimension */ 291da43764aSToby Isaac for (d = 0; d <= dim; d++) { 292da43764aSToby Isaac PetscInt cStart, cEnd, c; 293da43764aSToby Isaac 294da43764aSToby Isaac ierr = DMPlexGetHeightStratum(K, d, &cStart, &cEnd);CHKERRQ(ierr); 295da43764aSToby Isaac for (c = cStart; c < cEnd; c++) { 296da43764aSToby Isaac permvals[offset++] = c; 297da43764aSToby Isaac } 298da43764aSToby Isaac 299da43764aSToby Isaac ierr = DMPlexGetHeightStratum(Kref, d, &cStart, &cEnd);CHKERRQ(ierr); 300da43764aSToby Isaac for (c = cStart; c < cEnd; c++) { 301da43764aSToby Isaac permvals[offset++] = c + (pEnd - pStart); 302da43764aSToby Isaac } 303da43764aSToby Isaac } 304da43764aSToby Isaac ierr = ISCreateGeneral(comm, (pEnd - pStart) + (pRefEnd - pRefStart), permvals, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr); 305da43764aSToby Isaac ierr = PetscSectionSetPermutation(unionSection,perm);CHKERRQ(ierr); 306da43764aSToby Isaac ierr = PetscSectionSetUp(unionSection);CHKERRQ(ierr); 307da43764aSToby Isaac ierr = PetscSectionGetStorageSize(unionSection,&numUnionPoints);CHKERRQ(ierr); 308da43764aSToby Isaac ierr = PetscMalloc2(numUnionPoints,&coneSizes,dim+1,&numDimPoints);CHKERRQ(ierr); 309da43764aSToby Isaac /* count dimension points */ 310da43764aSToby Isaac for (d = 0; d <= dim; d++) { 311da43764aSToby Isaac PetscInt cStart, cOff, cOff2; 312da43764aSToby Isaac ierr = DMPlexGetHeightStratum(K,d,&cStart,NULL);CHKERRQ(ierr); 313da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,cStart-pStart,&cOff);CHKERRQ(ierr); 314da43764aSToby Isaac if (d < dim) { 315da43764aSToby Isaac ierr = DMPlexGetHeightStratum(K,d+1,&cStart,NULL);CHKERRQ(ierr); 316da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,cStart-pStart,&cOff2);CHKERRQ(ierr); 317da43764aSToby Isaac } 318da43764aSToby Isaac else { 319da43764aSToby Isaac cOff2 = numUnionPoints; 320da43764aSToby Isaac } 321da43764aSToby Isaac numDimPoints[dim - d] = cOff2 - cOff; 322da43764aSToby Isaac } 323da43764aSToby Isaac ierr = PetscSectionCreate(comm, &unionConeSection);CHKERRQ(ierr); 324da43764aSToby Isaac ierr = PetscSectionSetChart(unionConeSection, 0, numUnionPoints);CHKERRQ(ierr); 325da43764aSToby Isaac /* count the cones in the union */ 326da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 327da43764aSToby Isaac PetscInt dof, uOff; 328da43764aSToby Isaac 329da43764aSToby Isaac ierr = DMPlexGetConeSize(K, p, &dof);CHKERRQ(ierr); 330da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pStart,&uOff);CHKERRQ(ierr); 331da43764aSToby Isaac ierr = PetscSectionSetDof(unionConeSection, uOff, dof);CHKERRQ(ierr); 332da43764aSToby Isaac coneSizes[uOff] = dof; 333da43764aSToby Isaac } 334da43764aSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 335da43764aSToby Isaac PetscInt dof, uDof, uOff; 336da43764aSToby Isaac 337da43764aSToby Isaac ierr = DMPlexGetConeSize(Kref, p, &dof);CHKERRQ(ierr); 338da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 339da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 340da43764aSToby Isaac if (uDof) { 341da43764aSToby Isaac ierr = PetscSectionSetDof(unionConeSection, uOff, dof);CHKERRQ(ierr); 342da43764aSToby Isaac coneSizes[uOff] = dof; 343da43764aSToby Isaac } 344da43764aSToby Isaac } 345da43764aSToby Isaac ierr = PetscSectionSetUp(unionConeSection);CHKERRQ(ierr); 346da43764aSToby Isaac ierr = PetscSectionGetStorageSize(unionConeSection,&numCones);CHKERRQ(ierr); 347da43764aSToby Isaac ierr = PetscMalloc2(numCones,&unionCones,numCones,&unionOrientations);CHKERRQ(ierr); 348da43764aSToby Isaac /* write the cones in the union */ 349da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 350da43764aSToby Isaac PetscInt dof, uOff, c, cOff; 351da43764aSToby Isaac const PetscInt *cone, *orientation; 352da43764aSToby Isaac 353da43764aSToby Isaac ierr = DMPlexGetConeSize(K, p, &dof);CHKERRQ(ierr); 354da43764aSToby Isaac ierr = DMPlexGetCone(K, p, &cone);CHKERRQ(ierr); 355da43764aSToby Isaac ierr = DMPlexGetConeOrientation(K, p, &orientation);CHKERRQ(ierr); 356da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pStart,&uOff);CHKERRQ(ierr); 357da43764aSToby Isaac ierr = PetscSectionGetOffset(unionConeSection,uOff,&cOff);CHKERRQ(ierr); 358da43764aSToby Isaac for (c = 0; c < dof; c++) { 359da43764aSToby Isaac PetscInt e, eOff; 360da43764aSToby Isaac e = cone[c]; 361da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, e - pStart, &eOff);CHKERRQ(ierr); 362da43764aSToby Isaac unionCones[cOff + c] = eOff; 363da43764aSToby Isaac unionOrientations[cOff + c] = orientation[c]; 364da43764aSToby Isaac } 365da43764aSToby Isaac } 366da43764aSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 367da43764aSToby Isaac PetscInt dof, uDof, uOff, c, cOff; 368da43764aSToby Isaac const PetscInt *cone, *orientation; 369da43764aSToby Isaac 370da43764aSToby Isaac ierr = DMPlexGetConeSize(Kref, p, &dof);CHKERRQ(ierr); 371da43764aSToby Isaac ierr = DMPlexGetCone(Kref, p, &cone);CHKERRQ(ierr); 372da43764aSToby Isaac ierr = DMPlexGetConeOrientation(Kref, p, &orientation);CHKERRQ(ierr); 373da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 374da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 375da43764aSToby Isaac if (uDof) { 376da43764aSToby Isaac ierr = PetscSectionGetOffset(unionConeSection,uOff,&cOff);CHKERRQ(ierr); 377da43764aSToby Isaac for (c = 0; c < dof; c++) { 378da43764aSToby Isaac PetscInt e, eOff, eDof; 379da43764aSToby Isaac 380da43764aSToby Isaac e = cone[c]; 381da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection, e - pRefStart + (pEnd - pStart),&eDof);CHKERRQ(ierr); 382da43764aSToby Isaac if (eDof) { 383da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, e - pRefStart + (pEnd - pStart), &eOff);CHKERRQ(ierr); 384da43764aSToby Isaac } 385da43764aSToby Isaac else { 386da43764aSToby Isaac ierr = DMLabelGetValue(identityRef, e, &e);CHKERRQ(ierr); 387da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, e - pStart, &eOff);CHKERRQ(ierr); 388da43764aSToby Isaac } 389da43764aSToby Isaac unionCones[cOff + c] = eOff; 390da43764aSToby Isaac unionOrientations[cOff + c] = orientation[c]; 391da43764aSToby Isaac } 392da43764aSToby Isaac } 393da43764aSToby Isaac } 394da43764aSToby Isaac /* get the coordinates */ 395da43764aSToby Isaac { 396da43764aSToby Isaac PetscInt vStart, vEnd, vRefStart, vRefEnd, v, vDof, vOff; 397da43764aSToby Isaac PetscSection KcoordsSec, KrefCoordsSec; 398da43764aSToby Isaac Vec KcoordsVec, KrefCoordsVec; 399da43764aSToby Isaac PetscScalar *Kcoords; 400da43764aSToby Isaac 401da43764aSToby Isaac DMGetCoordinateSection(K, &KcoordsSec);CHKERRQ(ierr); 402da43764aSToby Isaac DMGetCoordinatesLocal(K, &KcoordsVec);CHKERRQ(ierr); 403da43764aSToby Isaac DMGetCoordinateSection(Kref, &KrefCoordsSec);CHKERRQ(ierr); 404da43764aSToby Isaac DMGetCoordinatesLocal(Kref, &KrefCoordsVec);CHKERRQ(ierr); 405da43764aSToby Isaac 406da43764aSToby Isaac numVerts = numDimPoints[0]; 407da43764aSToby Isaac ierr = PetscMalloc1(numVerts * dim,&unionCoords);CHKERRQ(ierr); 408da43764aSToby Isaac ierr = DMPlexGetDepthStratum(K,0,&vStart,&vEnd);CHKERRQ(ierr); 409da43764aSToby Isaac 410da43764aSToby Isaac offset = 0; 411da43764aSToby Isaac for (v = vStart; v < vEnd; v++) { 412da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,v - pStart,&vOff);CHKERRQ(ierr); 413da43764aSToby Isaac ierr = VecGetValuesSection(KcoordsVec, KcoordsSec, v, &Kcoords);CHKERRQ(ierr); 414da43764aSToby Isaac for (d = 0; d < dim; d++) { 415da43764aSToby Isaac unionCoords[offset * dim + d] = Kcoords[d]; 416da43764aSToby Isaac } 417da43764aSToby Isaac offset++; 418da43764aSToby Isaac } 419da43764aSToby Isaac ierr = DMPlexGetDepthStratum(Kref,0,&vRefStart,&vRefEnd);CHKERRQ(ierr); 420da43764aSToby Isaac for (v = vRefStart; v < vRefEnd; v++) { 421da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection,v - pRefStart + (pEnd - pStart),&vDof);CHKERRQ(ierr); 422da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,v - pRefStart + (pEnd - pStart),&vOff);CHKERRQ(ierr); 423da43764aSToby Isaac ierr = VecGetValuesSection(KrefCoordsVec, KrefCoordsSec, v, &Kcoords);CHKERRQ(ierr); 424da43764aSToby Isaac if (vDof) { 425da43764aSToby Isaac for (d = 0; d < dim; d++) { 426da43764aSToby Isaac unionCoords[offset * dim + d] = Kcoords[d]; 427da43764aSToby Isaac } 428da43764aSToby Isaac offset++; 429da43764aSToby Isaac } 430da43764aSToby Isaac } 431da43764aSToby Isaac } 432da43764aSToby Isaac ierr = DMCreate(comm,ref);CHKERRQ(ierr); 433da43764aSToby Isaac ierr = DMSetType(*ref,DMPLEX);CHKERRQ(ierr); 43428f4b327SMatthew G. Knepley ierr = DMSetDimension(*ref,dim);CHKERRQ(ierr); 435da43764aSToby Isaac ierr = DMPlexCreateFromDAG(*ref,dim,numDimPoints,coneSizes,unionCones,unionOrientations,unionCoords);CHKERRQ(ierr); 43610f7e118SToby Isaac /* set the tree */ 43710f7e118SToby Isaac ierr = PetscSectionCreate(comm,&parentSection);CHKERRQ(ierr); 43810f7e118SToby Isaac ierr = PetscSectionSetChart(parentSection,0,numUnionPoints);CHKERRQ(ierr); 43910f7e118SToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 44010f7e118SToby Isaac PetscInt uDof, uOff; 44110f7e118SToby Isaac 44210f7e118SToby Isaac ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 44310f7e118SToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 44410f7e118SToby Isaac if (uDof) { 44510f7e118SToby Isaac PetscSectionSetDof(parentSection,uOff,1);CHKERRQ(ierr); 44610f7e118SToby Isaac } 44710f7e118SToby Isaac } 44810f7e118SToby Isaac ierr = PetscSectionSetUp(parentSection);CHKERRQ(ierr); 44910f7e118SToby Isaac ierr = PetscSectionGetStorageSize(parentSection,&parentSize);CHKERRQ(ierr); 45010f7e118SToby Isaac ierr = PetscMalloc2(parentSize,&parents,parentSize,&childIDs);CHKERRQ(ierr); 45110f7e118SToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 45210f7e118SToby Isaac PetscInt uDof, uOff; 45310f7e118SToby Isaac 45410f7e118SToby Isaac ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 45510f7e118SToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 45610f7e118SToby Isaac if (uDof) { 45710f7e118SToby Isaac PetscInt pOff, parent, parentU; 45810f7e118SToby Isaac PetscSectionGetOffset(parentSection,uOff,&pOff);CHKERRQ(ierr); 45910f7e118SToby Isaac DMLabelGetValue(identityRef,p,&parent);CHKERRQ(ierr); 46010f7e118SToby Isaac ierr = PetscSectionGetOffset(unionSection, parent - pStart,&parentU);CHKERRQ(ierr); 46110f7e118SToby Isaac parents[pOff] = parentU; 46210f7e118SToby Isaac childIDs[pOff] = uOff; 46310f7e118SToby Isaac } 46410f7e118SToby Isaac } 465776742edSToby Isaac ierr = DMPlexSetTree_Internal(*ref,parentSection,parents,childIDs,PETSC_TRUE,PETSC_FALSE);CHKERRQ(ierr); 466dcbd3bf7SToby Isaac mesh = (DM_Plex *) (*ref)->data; 467dcbd3bf7SToby Isaac mesh->getchildsymmetry = DMPlexReferenceTreeGetChildSymmetry_Default; 46810f7e118SToby Isaac ierr = PetscSectionDestroy(&parentSection);CHKERRQ(ierr); 46910f7e118SToby Isaac ierr = PetscFree2(parents,childIDs);CHKERRQ(ierr); 47010f7e118SToby Isaac 471da43764aSToby Isaac /* clean up */ 472da43764aSToby Isaac ierr = PetscSectionDestroy(&unionSection);CHKERRQ(ierr); 473da43764aSToby Isaac ierr = PetscSectionDestroy(&unionConeSection);CHKERRQ(ierr); 474da43764aSToby Isaac ierr = ISDestroy(&perm);CHKERRQ(ierr); 475da43764aSToby Isaac ierr = PetscFree(unionCoords);CHKERRQ(ierr); 476da43764aSToby Isaac ierr = PetscFree2(unionCones,unionOrientations);CHKERRQ(ierr); 477da43764aSToby Isaac ierr = PetscFree2(coneSizes,numDimPoints);CHKERRQ(ierr); 478da43764aSToby Isaac ierr = DMDestroy(&K);CHKERRQ(ierr); 479da43764aSToby Isaac ierr = DMDestroy(&Kref);CHKERRQ(ierr); 480da43764aSToby Isaac PetscFunctionReturn(0); 481da43764aSToby Isaac } 482da43764aSToby Isaac 483d961a43aSToby Isaac #undef __FUNCT__ 484878b19aaSToby Isaac #define __FUNCT__ "DMPlexTreeSymmetrize" 485878b19aaSToby Isaac static PetscErrorCode DMPlexTreeSymmetrize(DM dm) 486878b19aaSToby Isaac { 487878b19aaSToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 488878b19aaSToby Isaac PetscSection childSec, pSec; 489878b19aaSToby Isaac PetscInt p, pSize, cSize, parMax = PETSC_MIN_INT, parMin = PETSC_MAX_INT; 490878b19aaSToby Isaac PetscInt *offsets, *children, pStart, pEnd; 491878b19aaSToby Isaac PetscErrorCode ierr; 492878b19aaSToby Isaac 493878b19aaSToby Isaac PetscFunctionBegin; 494878b19aaSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 495878b19aaSToby Isaac ierr = PetscSectionDestroy(&mesh->childSection);CHKERRQ(ierr); 496878b19aaSToby Isaac ierr = PetscFree(mesh->children);CHKERRQ(ierr); 497878b19aaSToby Isaac pSec = mesh->parentSection; 498878b19aaSToby Isaac if (!pSec) PetscFunctionReturn(0); 499878b19aaSToby Isaac ierr = PetscSectionGetStorageSize(pSec,&pSize);CHKERRQ(ierr); 500878b19aaSToby Isaac for (p = 0; p < pSize; p++) { 501878b19aaSToby Isaac PetscInt par = mesh->parents[p]; 502878b19aaSToby Isaac 503878b19aaSToby Isaac parMax = PetscMax(parMax,par+1); 504878b19aaSToby Isaac parMin = PetscMin(parMin,par); 505878b19aaSToby Isaac } 506878b19aaSToby Isaac if (parMin > parMax) { 507878b19aaSToby Isaac parMin = -1; 508878b19aaSToby Isaac parMax = -1; 509878b19aaSToby Isaac } 510878b19aaSToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)pSec),&childSec);CHKERRQ(ierr); 511878b19aaSToby Isaac ierr = PetscSectionSetChart(childSec,parMin,parMax);CHKERRQ(ierr); 512878b19aaSToby Isaac for (p = 0; p < pSize; p++) { 513878b19aaSToby Isaac PetscInt par = mesh->parents[p]; 514878b19aaSToby Isaac 515878b19aaSToby Isaac ierr = PetscSectionAddDof(childSec,par,1);CHKERRQ(ierr); 516878b19aaSToby Isaac } 517878b19aaSToby Isaac ierr = PetscSectionSetUp(childSec);CHKERRQ(ierr); 518878b19aaSToby Isaac ierr = PetscSectionGetStorageSize(childSec,&cSize);CHKERRQ(ierr); 519878b19aaSToby Isaac ierr = PetscMalloc1(cSize,&children);CHKERRQ(ierr); 520878b19aaSToby Isaac ierr = PetscCalloc1(parMax-parMin,&offsets);CHKERRQ(ierr); 521878b19aaSToby Isaac ierr = PetscSectionGetChart(pSec,&pStart,&pEnd);CHKERRQ(ierr); 522878b19aaSToby Isaac for (p = pStart; p < pEnd; p++) { 523878b19aaSToby Isaac PetscInt dof, off, i; 524878b19aaSToby Isaac 525878b19aaSToby Isaac ierr = PetscSectionGetDof(pSec,p,&dof);CHKERRQ(ierr); 526878b19aaSToby Isaac ierr = PetscSectionGetOffset(pSec,p,&off);CHKERRQ(ierr); 527878b19aaSToby Isaac for (i = 0; i < dof; i++) { 528878b19aaSToby Isaac PetscInt par = mesh->parents[off + i], cOff; 529878b19aaSToby Isaac 530878b19aaSToby Isaac ierr = PetscSectionGetOffset(childSec,par,&cOff);CHKERRQ(ierr); 531878b19aaSToby Isaac children[cOff + offsets[par-parMin]++] = p; 532878b19aaSToby Isaac } 533878b19aaSToby Isaac } 534878b19aaSToby Isaac mesh->childSection = childSec; 535878b19aaSToby Isaac mesh->children = children; 536878b19aaSToby Isaac ierr = PetscFree(offsets);CHKERRQ(ierr); 537878b19aaSToby Isaac PetscFunctionReturn(0); 538878b19aaSToby Isaac } 539878b19aaSToby Isaac 540878b19aaSToby Isaac #undef __FUNCT__ 5416dd5a8c8SToby Isaac #define __FUNCT__ "AnchorsFlatten" 5426dd5a8c8SToby Isaac static PetscErrorCode AnchorsFlatten (PetscSection section, IS is, PetscSection *sectionNew, IS *isNew) 5436dd5a8c8SToby Isaac { 5446dd5a8c8SToby Isaac PetscInt pStart, pEnd, size, sizeNew, i, p, *valsNew = NULL; 5456dd5a8c8SToby Isaac const PetscInt *vals; 5466dd5a8c8SToby Isaac PetscSection secNew; 5476dd5a8c8SToby Isaac PetscBool anyNew, globalAnyNew; 5486dd5a8c8SToby Isaac PetscBool compress; 5496dd5a8c8SToby Isaac PetscErrorCode ierr; 5506dd5a8c8SToby Isaac 5516dd5a8c8SToby Isaac PetscFunctionBegin; 5526dd5a8c8SToby Isaac ierr = PetscSectionGetChart(section,&pStart,&pEnd);CHKERRQ(ierr); 5536dd5a8c8SToby Isaac ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr); 5546dd5a8c8SToby Isaac ierr = ISGetIndices(is,&vals);CHKERRQ(ierr); 5556dd5a8c8SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)section),&secNew);CHKERRQ(ierr); 5566dd5a8c8SToby Isaac ierr = PetscSectionSetChart(secNew,pStart,pEnd);CHKERRQ(ierr); 5576dd5a8c8SToby Isaac for (i = 0; i < size; i++) { 5586dd5a8c8SToby Isaac PetscInt dof; 5596dd5a8c8SToby Isaac 5606dd5a8c8SToby Isaac p = vals[i]; 5616dd5a8c8SToby Isaac if (p < pStart || p >= pEnd) continue; 5626dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 5636dd5a8c8SToby Isaac if (dof) break; 5646dd5a8c8SToby Isaac } 5656dd5a8c8SToby Isaac if (i == size) { 5666dd5a8c8SToby Isaac ierr = PetscSectionSetUp(secNew);CHKERRQ(ierr); 5676dd5a8c8SToby Isaac anyNew = PETSC_FALSE; 5686dd5a8c8SToby Isaac compress = PETSC_FALSE; 5696dd5a8c8SToby Isaac sizeNew = 0; 5706dd5a8c8SToby Isaac } 5716dd5a8c8SToby Isaac else { 5726dd5a8c8SToby Isaac anyNew = PETSC_TRUE; 5736dd5a8c8SToby Isaac for (p = pStart; p < pEnd; p++) { 5746dd5a8c8SToby Isaac PetscInt dof, off; 5756dd5a8c8SToby Isaac 5766dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 5776dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 5786dd5a8c8SToby Isaac for (i = 0; i < dof; i++) { 5796dd5a8c8SToby Isaac PetscInt q = vals[off + i], qDof = 0; 5806dd5a8c8SToby Isaac 5816dd5a8c8SToby Isaac if (q >= pStart && q < pEnd) { 5826dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, q, &qDof);CHKERRQ(ierr); 5836dd5a8c8SToby Isaac } 5846dd5a8c8SToby Isaac if (qDof) { 5856dd5a8c8SToby Isaac ierr = PetscSectionAddDof(secNew, p, qDof);CHKERRQ(ierr); 5866dd5a8c8SToby Isaac } 5876dd5a8c8SToby Isaac else { 5886dd5a8c8SToby Isaac ierr = PetscSectionAddDof(secNew, p, 1);CHKERRQ(ierr); 5896dd5a8c8SToby Isaac } 5906dd5a8c8SToby Isaac } 5916dd5a8c8SToby Isaac } 5926dd5a8c8SToby Isaac ierr = PetscSectionSetUp(secNew);CHKERRQ(ierr); 5936dd5a8c8SToby Isaac ierr = PetscSectionGetStorageSize(secNew,&sizeNew);CHKERRQ(ierr); 5946dd5a8c8SToby Isaac ierr = PetscMalloc1(sizeNew,&valsNew);CHKERRQ(ierr); 5956dd5a8c8SToby Isaac compress = PETSC_FALSE; 5966dd5a8c8SToby Isaac for (p = pStart; p < pEnd; p++) { 5976dd5a8c8SToby Isaac PetscInt dof, off, count, offNew, dofNew; 5986dd5a8c8SToby Isaac 5996dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 6006dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 6016dd5a8c8SToby Isaac ierr = PetscSectionGetDof(secNew, p, &dofNew);CHKERRQ(ierr); 6026dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(secNew, p, &offNew);CHKERRQ(ierr); 6036dd5a8c8SToby Isaac count = 0; 6046dd5a8c8SToby Isaac for (i = 0; i < dof; i++) { 6056dd5a8c8SToby Isaac PetscInt q = vals[off + i], qDof = 0, qOff = 0, j; 6066dd5a8c8SToby Isaac 6076dd5a8c8SToby Isaac if (q >= pStart && q < pEnd) { 6086dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, q, &qDof);CHKERRQ(ierr); 6096dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(section, q, &qOff);CHKERRQ(ierr); 6106dd5a8c8SToby Isaac } 6116dd5a8c8SToby Isaac if (qDof) { 6126dd5a8c8SToby Isaac PetscInt oldCount = count; 6136dd5a8c8SToby Isaac 6146dd5a8c8SToby Isaac for (j = 0; j < qDof; j++) { 6156dd5a8c8SToby Isaac PetscInt k, r = vals[qOff + j]; 6166dd5a8c8SToby Isaac 6176dd5a8c8SToby Isaac for (k = 0; k < oldCount; k++) { 6186dd5a8c8SToby Isaac if (valsNew[offNew + k] == r) { 6196dd5a8c8SToby Isaac break; 6206dd5a8c8SToby Isaac } 6216dd5a8c8SToby Isaac } 6226dd5a8c8SToby Isaac if (k == oldCount) { 6236dd5a8c8SToby Isaac valsNew[offNew + count++] = r; 6246dd5a8c8SToby Isaac } 6256dd5a8c8SToby Isaac } 6266dd5a8c8SToby Isaac } 6276dd5a8c8SToby Isaac else { 6286dd5a8c8SToby Isaac PetscInt k, oldCount = count; 6296dd5a8c8SToby Isaac 6306dd5a8c8SToby Isaac for (k = 0; k < oldCount; k++) { 6316dd5a8c8SToby Isaac if (valsNew[offNew + k] == q) { 6326dd5a8c8SToby Isaac break; 6336dd5a8c8SToby Isaac } 6346dd5a8c8SToby Isaac } 6356dd5a8c8SToby Isaac if (k == oldCount) { 6366dd5a8c8SToby Isaac valsNew[offNew + count++] = q; 6376dd5a8c8SToby Isaac } 6386dd5a8c8SToby Isaac } 6396dd5a8c8SToby Isaac } 6406dd5a8c8SToby Isaac if (count < dofNew) { 6416dd5a8c8SToby Isaac ierr = PetscSectionSetDof(secNew, p, count);CHKERRQ(ierr); 6426dd5a8c8SToby Isaac compress = PETSC_TRUE; 6436dd5a8c8SToby Isaac } 6446dd5a8c8SToby Isaac } 6456dd5a8c8SToby Isaac } 6466dd5a8c8SToby Isaac ierr = ISRestoreIndices(is,&vals);CHKERRQ(ierr); 6476dd5a8c8SToby Isaac ierr = MPI_Allreduce(&anyNew,&globalAnyNew,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)secNew));CHKERRQ(ierr); 6486dd5a8c8SToby Isaac if (!globalAnyNew) { 6496dd5a8c8SToby Isaac ierr = PetscSectionDestroy(&secNew);CHKERRQ(ierr); 6506dd5a8c8SToby Isaac *sectionNew = NULL; 6516dd5a8c8SToby Isaac *isNew = NULL; 6526dd5a8c8SToby Isaac } 6536dd5a8c8SToby Isaac else { 6546dd5a8c8SToby Isaac PetscBool globalCompress; 6556dd5a8c8SToby Isaac 6566dd5a8c8SToby Isaac ierr = MPI_Allreduce(&compress,&globalCompress,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)secNew));CHKERRQ(ierr); 6576dd5a8c8SToby Isaac if (compress) { 6586dd5a8c8SToby Isaac PetscSection secComp; 6596dd5a8c8SToby Isaac PetscInt *valsComp = NULL; 6606dd5a8c8SToby Isaac 6616dd5a8c8SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)section),&secComp);CHKERRQ(ierr); 6626dd5a8c8SToby Isaac ierr = PetscSectionSetChart(secComp,pStart,pEnd);CHKERRQ(ierr); 6636dd5a8c8SToby Isaac for (p = pStart; p < pEnd; p++) { 6646dd5a8c8SToby Isaac PetscInt dof; 6656dd5a8c8SToby Isaac 6666dd5a8c8SToby Isaac ierr = PetscSectionGetDof(secNew, p, &dof);CHKERRQ(ierr); 6676dd5a8c8SToby Isaac ierr = PetscSectionSetDof(secComp, p, dof);CHKERRQ(ierr); 6686dd5a8c8SToby Isaac } 6696dd5a8c8SToby Isaac ierr = PetscSectionSetUp(secComp);CHKERRQ(ierr); 6706dd5a8c8SToby Isaac ierr = PetscSectionGetStorageSize(secComp,&sizeNew);CHKERRQ(ierr); 6716dd5a8c8SToby Isaac ierr = PetscMalloc1(sizeNew,&valsComp);CHKERRQ(ierr); 6726dd5a8c8SToby Isaac for (p = pStart; p < pEnd; p++) { 6736dd5a8c8SToby Isaac PetscInt dof, off, offNew, j; 6746dd5a8c8SToby Isaac 6756dd5a8c8SToby Isaac ierr = PetscSectionGetDof(secNew, p, &dof);CHKERRQ(ierr); 6766dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(secNew, p, &off);CHKERRQ(ierr); 6776dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(secComp, p, &offNew);CHKERRQ(ierr); 6786dd5a8c8SToby Isaac for (j = 0; j < dof; j++) { 6796dd5a8c8SToby Isaac valsComp[offNew + j] = valsNew[off + j]; 6806dd5a8c8SToby Isaac } 6816dd5a8c8SToby Isaac } 6826dd5a8c8SToby Isaac ierr = PetscSectionDestroy(&secNew);CHKERRQ(ierr); 6836dd5a8c8SToby Isaac secNew = secComp; 6846dd5a8c8SToby Isaac ierr = PetscFree(valsNew);CHKERRQ(ierr); 6856dd5a8c8SToby Isaac valsNew = valsComp; 6866dd5a8c8SToby Isaac } 6876dd5a8c8SToby Isaac ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),sizeNew,valsNew,PETSC_OWN_POINTER,isNew);CHKERRQ(ierr); 6886dd5a8c8SToby Isaac } 6896dd5a8c8SToby Isaac PetscFunctionReturn(0); 6906dd5a8c8SToby Isaac } 6916dd5a8c8SToby Isaac 6926dd5a8c8SToby Isaac #undef __FUNCT__ 693f7c74593SToby Isaac #define __FUNCT__ "DMPlexCreateAnchors_Tree" 694f7c74593SToby Isaac static PetscErrorCode DMPlexCreateAnchors_Tree(DM dm) 69566af876cSToby Isaac { 69666af876cSToby Isaac PetscInt p, pStart, pEnd, *anchors, size; 69766af876cSToby Isaac PetscInt aMin = PETSC_MAX_INT, aMax = PETSC_MIN_INT; 69866af876cSToby Isaac PetscSection aSec; 699f9f063d4SToby Isaac DMLabel canonLabel; 70066af876cSToby Isaac IS aIS; 70166af876cSToby Isaac PetscErrorCode ierr; 70266af876cSToby Isaac 70366af876cSToby Isaac PetscFunctionBegin; 70466af876cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 70566af876cSToby Isaac ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 706f9f063d4SToby Isaac ierr = DMPlexGetLabel(dm,"canonical",&canonLabel);CHKERRQ(ierr); 70766af876cSToby Isaac for (p = pStart; p < pEnd; p++) { 70866af876cSToby Isaac PetscInt parent; 70966af876cSToby Isaac 710f9f063d4SToby Isaac if (canonLabel) { 711f9f063d4SToby Isaac PetscInt canon; 712f9f063d4SToby Isaac 713f9f063d4SToby Isaac ierr = DMLabelGetValue(canonLabel,p,&canon);CHKERRQ(ierr); 714f9f063d4SToby Isaac if (p != canon) continue; 715f9f063d4SToby Isaac } 71666af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 71766af876cSToby Isaac if (parent != p) { 71866af876cSToby Isaac aMin = PetscMin(aMin,p); 71966af876cSToby Isaac aMax = PetscMax(aMax,p+1); 72066af876cSToby Isaac } 72166af876cSToby Isaac } 72266af876cSToby Isaac if (aMin > aMax) { 72366af876cSToby Isaac aMin = -1; 72466af876cSToby Isaac aMax = -1; 72566af876cSToby Isaac } 726e228b242SToby Isaac ierr = PetscSectionCreate(PETSC_COMM_SELF,&aSec);CHKERRQ(ierr); 72766af876cSToby Isaac ierr = PetscSectionSetChart(aSec,aMin,aMax);CHKERRQ(ierr); 72866af876cSToby Isaac for (p = aMin; p < aMax; p++) { 72966af876cSToby Isaac PetscInt parent, ancestor = p; 73066af876cSToby Isaac 731f9f063d4SToby Isaac if (canonLabel) { 732f9f063d4SToby Isaac PetscInt canon; 733f9f063d4SToby Isaac 734f9f063d4SToby Isaac ierr = DMLabelGetValue(canonLabel,p,&canon);CHKERRQ(ierr); 735f9f063d4SToby Isaac if (p != canon) continue; 736f9f063d4SToby Isaac } 73766af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 73866af876cSToby Isaac while (parent != ancestor) { 73966af876cSToby Isaac ancestor = parent; 74066af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,ancestor,&parent,NULL);CHKERRQ(ierr); 74166af876cSToby Isaac } 74266af876cSToby Isaac if (ancestor != p) { 74366af876cSToby Isaac PetscInt closureSize, *closure = NULL; 74466af876cSToby Isaac 74566af876cSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 74666af876cSToby Isaac ierr = PetscSectionSetDof(aSec,p,closureSize);CHKERRQ(ierr); 74766af876cSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 74866af876cSToby Isaac } 74966af876cSToby Isaac } 75066af876cSToby Isaac ierr = PetscSectionSetUp(aSec);CHKERRQ(ierr); 75166af876cSToby Isaac ierr = PetscSectionGetStorageSize(aSec,&size);CHKERRQ(ierr); 75266af876cSToby Isaac ierr = PetscMalloc1(size,&anchors);CHKERRQ(ierr); 75366af876cSToby Isaac for (p = aMin; p < aMax; p++) { 75466af876cSToby Isaac PetscInt parent, ancestor = p; 75566af876cSToby Isaac 756f9f063d4SToby Isaac if (canonLabel) { 757f9f063d4SToby Isaac PetscInt canon; 758f9f063d4SToby Isaac 759f9f063d4SToby Isaac ierr = DMLabelGetValue(canonLabel,p,&canon);CHKERRQ(ierr); 760f9f063d4SToby Isaac if (p != canon) continue; 761f9f063d4SToby Isaac } 76266af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 76366af876cSToby Isaac while (parent != ancestor) { 76466af876cSToby Isaac ancestor = parent; 76566af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,ancestor,&parent,NULL);CHKERRQ(ierr); 76666af876cSToby Isaac } 76766af876cSToby Isaac if (ancestor != p) { 76866af876cSToby Isaac PetscInt j, closureSize, *closure = NULL, aOff; 76966af876cSToby Isaac 77066af876cSToby Isaac ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr); 77166af876cSToby Isaac 77266af876cSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 77366af876cSToby Isaac for (j = 0; j < closureSize; j++) { 77466af876cSToby Isaac anchors[aOff + j] = closure[2*j]; 77566af876cSToby Isaac } 77666af876cSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 77766af876cSToby Isaac } 77866af876cSToby Isaac } 779e228b242SToby Isaac ierr = ISCreateGeneral(PETSC_COMM_SELF,size,anchors,PETSC_OWN_POINTER,&aIS);CHKERRQ(ierr); 7806dd5a8c8SToby Isaac { 7816dd5a8c8SToby Isaac PetscSection aSecNew = aSec; 7826dd5a8c8SToby Isaac IS aISNew = aIS; 7836dd5a8c8SToby Isaac 7846dd5a8c8SToby Isaac ierr = PetscObjectReference((PetscObject)aSec);CHKERRQ(ierr); 7856dd5a8c8SToby Isaac ierr = PetscObjectReference((PetscObject)aIS);CHKERRQ(ierr); 7866dd5a8c8SToby Isaac while (aSecNew) { 7876dd5a8c8SToby Isaac ierr = PetscSectionDestroy(&aSec);CHKERRQ(ierr); 7886dd5a8c8SToby Isaac ierr = ISDestroy(&aIS);CHKERRQ(ierr); 7896dd5a8c8SToby Isaac aSec = aSecNew; 7906dd5a8c8SToby Isaac aIS = aISNew; 7916dd5a8c8SToby Isaac aSecNew = NULL; 7926dd5a8c8SToby Isaac aISNew = NULL; 7936dd5a8c8SToby Isaac ierr = AnchorsFlatten(aSec,aIS,&aSecNew,&aISNew);CHKERRQ(ierr); 7946dd5a8c8SToby Isaac } 7956dd5a8c8SToby Isaac } 796a17985deSToby Isaac ierr = DMPlexSetAnchors(dm,aSec,aIS);CHKERRQ(ierr); 79766af876cSToby Isaac ierr = PetscSectionDestroy(&aSec);CHKERRQ(ierr); 79866af876cSToby Isaac ierr = ISDestroy(&aIS);CHKERRQ(ierr); 79966af876cSToby Isaac PetscFunctionReturn(0); 80066af876cSToby Isaac } 80166af876cSToby Isaac 80266af876cSToby Isaac #undef __FUNCT__ 803*6461c1adSToby Isaac #define __FUNCT__ "DMPlexGetTrueSupportSize" 804*6461c1adSToby Isaac static PetscErrorCode DMPlexGetTrueSupportSize(DM dm,PetscInt p,PetscInt *dof,PetscInt *numTrueSupp) 805*6461c1adSToby Isaac { 806*6461c1adSToby Isaac PetscErrorCode ierr; 807*6461c1adSToby Isaac 808*6461c1adSToby Isaac PetscFunctionBegin; 809*6461c1adSToby Isaac if (numTrueSupp[p] == -1) { 810*6461c1adSToby Isaac PetscInt i, alldof; 811*6461c1adSToby Isaac const PetscInt *supp; 812*6461c1adSToby Isaac PetscInt count = 0; 813*6461c1adSToby Isaac 814*6461c1adSToby Isaac ierr = DMPlexGetSupportSize(dm,p,&alldof);CHKERRQ(ierr); 815*6461c1adSToby Isaac ierr = DMPlexGetSupport(dm,p,&supp);CHKERRQ(ierr); 816*6461c1adSToby Isaac for (i = 0; i < alldof; i++) { 817*6461c1adSToby Isaac PetscInt q = supp[i], numCones, j; 818*6461c1adSToby Isaac const PetscInt *cone; 819*6461c1adSToby Isaac 820*6461c1adSToby Isaac ierr = DMPlexGetConeSize(dm,q,&numCones);CHKERRQ(ierr); 821*6461c1adSToby Isaac ierr = DMPlexGetCone(dm,q,&cone);CHKERRQ(ierr); 822*6461c1adSToby Isaac for (j = 0; j < numCones; j++) { 823*6461c1adSToby Isaac if (cone[j] == p) break; 824*6461c1adSToby Isaac } 825*6461c1adSToby Isaac if (j < numCones) count++; 826*6461c1adSToby Isaac } 827*6461c1adSToby Isaac numTrueSupp[p] = count; 828*6461c1adSToby Isaac } 829*6461c1adSToby Isaac *dof = numTrueSupp[p]; 830*6461c1adSToby Isaac PetscFunctionReturn(0); 831*6461c1adSToby Isaac } 832*6461c1adSToby Isaac 833*6461c1adSToby Isaac #undef __FUNCT__ 834776742edSToby Isaac #define __FUNCT__ "DMPlexTreeExchangeSupports" 835776742edSToby Isaac static PetscErrorCode DMPlexTreeExchangeSupports(DM dm) 836776742edSToby Isaac { 837776742edSToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 838776742edSToby Isaac PetscSection newSupportSection; 839776742edSToby Isaac PetscInt newSize, *newSupports, pStart, pEnd, p, d, depth; 840*6461c1adSToby Isaac PetscInt *numTrueSupp; 841776742edSToby Isaac PetscInt *offsets; 842776742edSToby Isaac PetscErrorCode ierr; 843776742edSToby Isaac 844776742edSToby Isaac PetscFunctionBegin; 845776742edSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 846776742edSToby Isaac /* symmetrize the hierarchy */ 847776742edSToby Isaac ierr = DMPlexGetDepth(dm,&depth);CHKERRQ(ierr); 848e228b242SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)(mesh->supportSection)),&newSupportSection);CHKERRQ(ierr); 849776742edSToby Isaac ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 850776742edSToby Isaac ierr = PetscSectionSetChart(newSupportSection,pStart,pEnd);CHKERRQ(ierr); 851776742edSToby Isaac ierr = PetscCalloc1(pEnd,&offsets);CHKERRQ(ierr); 852*6461c1adSToby Isaac ierr = PetscMalloc1(pEnd,&numTrueSupp);CHKERRQ(ierr); 853*6461c1adSToby Isaac for (p = 0; p < pEnd; p++) numTrueSupp[p] = -1; 854*6461c1adSToby Isaac /* if a point is in the (true) support of q, it should be in the support of 855776742edSToby Isaac * parent(q) */ 856776742edSToby Isaac for (d = 0; d <= depth; d++) { 857776742edSToby Isaac ierr = DMPlexGetHeightStratum(dm,d,&pStart,&pEnd);CHKERRQ(ierr); 858776742edSToby Isaac for (p = pStart; p < pEnd; ++p) { 859776742edSToby Isaac PetscInt dof, q, qdof, parent; 860776742edSToby Isaac 861*6461c1adSToby Isaac ierr = DMPlexGetTrueSupportSize(dm,p,&dof,numTrueSupp);CHKERRQ(ierr); 862776742edSToby Isaac ierr = PetscSectionAddDof(newSupportSection, p, dof);CHKERRQ(ierr); 863776742edSToby Isaac q = p; 864776742edSToby Isaac ierr = DMPlexGetTreeParent(dm,q,&parent,NULL);CHKERRQ(ierr); 865776742edSToby Isaac while (parent != q && parent >= pStart && parent < pEnd) { 866776742edSToby Isaac q = parent; 867776742edSToby Isaac 868*6461c1adSToby Isaac ierr = DMPlexGetTrueSupportSize(dm,q,&qdof,numTrueSupp);CHKERRQ(ierr); 869776742edSToby Isaac ierr = PetscSectionAddDof(newSupportSection,p,qdof);CHKERRQ(ierr); 870776742edSToby Isaac ierr = PetscSectionAddDof(newSupportSection,q,dof);CHKERRQ(ierr); 871776742edSToby Isaac ierr = DMPlexGetTreeParent(dm,q,&parent,NULL);CHKERRQ(ierr); 872776742edSToby Isaac } 873776742edSToby Isaac } 874776742edSToby Isaac } 875776742edSToby Isaac ierr = PetscSectionSetUp(newSupportSection);CHKERRQ(ierr); 876776742edSToby Isaac ierr = PetscSectionGetStorageSize(newSupportSection,&newSize);CHKERRQ(ierr); 877776742edSToby Isaac ierr = PetscMalloc1(newSize,&newSupports);CHKERRQ(ierr); 878776742edSToby Isaac for (d = 0; d <= depth; d++) { 879776742edSToby Isaac ierr = DMPlexGetHeightStratum(dm,d,&pStart,&pEnd);CHKERRQ(ierr); 880776742edSToby Isaac for (p = pStart; p < pEnd; p++) { 881776742edSToby Isaac PetscInt dof, off, q, qdof, qoff, newDof, newOff, newqOff, i, parent; 882776742edSToby Isaac 883776742edSToby Isaac ierr = PetscSectionGetDof(mesh->supportSection, p, &dof);CHKERRQ(ierr); 884776742edSToby Isaac ierr = PetscSectionGetOffset(mesh->supportSection, p, &off);CHKERRQ(ierr); 885776742edSToby Isaac ierr = PetscSectionGetDof(newSupportSection, p, &newDof);CHKERRQ(ierr); 886776742edSToby Isaac ierr = PetscSectionGetOffset(newSupportSection, p, &newOff);CHKERRQ(ierr); 887776742edSToby Isaac for (i = 0; i < dof; i++) { 888*6461c1adSToby Isaac PetscInt numCones, j; 889*6461c1adSToby Isaac const PetscInt *cone; 890*6461c1adSToby Isaac PetscInt q = mesh->supports[off + i]; 891*6461c1adSToby Isaac 892*6461c1adSToby Isaac ierr = DMPlexGetConeSize(dm,q,&numCones);CHKERRQ(ierr); 893*6461c1adSToby Isaac ierr = DMPlexGetCone(dm,q,&cone);CHKERRQ(ierr); 894*6461c1adSToby Isaac for (j = 0; j < numCones; j++) { 895*6461c1adSToby Isaac if (cone[j] == p) break; 896*6461c1adSToby Isaac } 897*6461c1adSToby Isaac if (j < numCones) newSupports[newOff+offsets[p]++] = q; 898776742edSToby Isaac } 899776742edSToby Isaac mesh->maxSupportSize = PetscMax(mesh->maxSupportSize,newDof); 900776742edSToby Isaac 901776742edSToby Isaac q = p; 902776742edSToby Isaac ierr = DMPlexGetTreeParent(dm,q,&parent,NULL);CHKERRQ(ierr); 903776742edSToby Isaac while (parent != q && parent >= pStart && parent < pEnd) { 904776742edSToby Isaac q = parent; 905776742edSToby Isaac ierr = PetscSectionGetDof(mesh->supportSection, q, &qdof);CHKERRQ(ierr); 906776742edSToby Isaac ierr = PetscSectionGetOffset(mesh->supportSection, q, &qoff);CHKERRQ(ierr); 907776742edSToby Isaac ierr = PetscSectionGetOffset(newSupportSection, q, &newqOff);CHKERRQ(ierr); 908776742edSToby Isaac for (i = 0; i < qdof; i++) { 909*6461c1adSToby Isaac PetscInt numCones, j; 910*6461c1adSToby Isaac const PetscInt *cone; 911*6461c1adSToby Isaac PetscInt r = mesh->supports[qoff + i]; 912*6461c1adSToby Isaac 913*6461c1adSToby Isaac ierr = DMPlexGetConeSize(dm,r,&numCones);CHKERRQ(ierr); 914*6461c1adSToby Isaac ierr = DMPlexGetCone(dm,r,&cone);CHKERRQ(ierr); 915*6461c1adSToby Isaac for (j = 0; j < numCones; j++) { 916*6461c1adSToby Isaac if (cone[j] == q) break; 917*6461c1adSToby Isaac } 918*6461c1adSToby Isaac if (j < numCones) newSupports[newOff+offsets[p]++] = r; 919776742edSToby Isaac } 920776742edSToby Isaac for (i = 0; i < dof; i++) { 921*6461c1adSToby Isaac PetscInt numCones, j; 922*6461c1adSToby Isaac const PetscInt *cone; 923*6461c1adSToby Isaac PetscInt r = mesh->supports[off + i]; 924*6461c1adSToby Isaac 925*6461c1adSToby Isaac ierr = DMPlexGetConeSize(dm,r,&numCones);CHKERRQ(ierr); 926*6461c1adSToby Isaac ierr = DMPlexGetCone(dm,r,&cone);CHKERRQ(ierr); 927*6461c1adSToby Isaac for (j = 0; j < numCones; j++) { 928*6461c1adSToby Isaac if (cone[j] == p) break; 929*6461c1adSToby Isaac } 930*6461c1adSToby Isaac if (j < numCones) newSupports[newqOff+offsets[q]++] = r; 931776742edSToby Isaac } 932776742edSToby Isaac ierr = DMPlexGetTreeParent(dm,q,&parent,NULL);CHKERRQ(ierr); 933776742edSToby Isaac } 934776742edSToby Isaac } 935776742edSToby Isaac } 936776742edSToby Isaac ierr = PetscSectionDestroy(&mesh->supportSection);CHKERRQ(ierr); 937776742edSToby Isaac mesh->supportSection = newSupportSection; 938776742edSToby Isaac ierr = PetscFree(mesh->supports);CHKERRQ(ierr); 939776742edSToby Isaac mesh->supports = newSupports; 940776742edSToby Isaac ierr = PetscFree(offsets);CHKERRQ(ierr); 941*6461c1adSToby Isaac ierr = PetscFree(numTrueSupp);CHKERRQ(ierr); 942776742edSToby Isaac 943776742edSToby Isaac PetscFunctionReturn(0); 944776742edSToby Isaac } 945776742edSToby Isaac 946f7c74593SToby Isaac static PetscErrorCode DMPlexComputeAnchorMatrix_Tree_Direct(DM,PetscSection,PetscSection,Mat); 947f7c74593SToby Isaac static PetscErrorCode DMPlexComputeAnchorMatrix_Tree_FromReference(DM,PetscSection,PetscSection,Mat); 948f7c74593SToby Isaac 949776742edSToby Isaac #undef __FUNCT__ 950f9f063d4SToby Isaac #define __FUNCT__ "DMPlexSetTree_Internal" 951776742edSToby Isaac static PetscErrorCode DMPlexSetTree_Internal(DM dm, PetscSection parentSection, PetscInt *parents, PetscInt *childIDs, PetscBool computeCanonical, PetscBool exchangeSupports) 952f9f063d4SToby Isaac { 953f9f063d4SToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 954f9f063d4SToby Isaac DM refTree; 955f9f063d4SToby Isaac PetscInt size; 956f9f063d4SToby Isaac PetscErrorCode ierr; 957f9f063d4SToby Isaac 958f9f063d4SToby Isaac PetscFunctionBegin; 959f9f063d4SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 960f9f063d4SToby Isaac PetscValidHeaderSpecific(parentSection, PETSC_SECTION_CLASSID, 2); 961f9f063d4SToby Isaac ierr = PetscObjectReference((PetscObject)parentSection);CHKERRQ(ierr); 962f9f063d4SToby Isaac ierr = PetscSectionDestroy(&mesh->parentSection);CHKERRQ(ierr); 963f9f063d4SToby Isaac mesh->parentSection = parentSection; 964f9f063d4SToby Isaac ierr = PetscSectionGetStorageSize(parentSection,&size);CHKERRQ(ierr); 965f9f063d4SToby Isaac if (parents != mesh->parents) { 966f9f063d4SToby Isaac ierr = PetscFree(mesh->parents);CHKERRQ(ierr); 967f9f063d4SToby Isaac ierr = PetscMalloc1(size,&mesh->parents);CHKERRQ(ierr); 968f9f063d4SToby Isaac ierr = PetscMemcpy(mesh->parents, parents, size * sizeof(*parents));CHKERRQ(ierr); 969f9f063d4SToby Isaac } 970f9f063d4SToby Isaac if (childIDs != mesh->childIDs) { 971f9f063d4SToby Isaac ierr = PetscFree(mesh->childIDs);CHKERRQ(ierr); 972f9f063d4SToby Isaac ierr = PetscMalloc1(size,&mesh->childIDs);CHKERRQ(ierr); 973f9f063d4SToby Isaac ierr = PetscMemcpy(mesh->childIDs, childIDs, size * sizeof(*childIDs));CHKERRQ(ierr); 974f9f063d4SToby Isaac } 975f9f063d4SToby Isaac ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr); 976f9f063d4SToby Isaac if (refTree) { 977f9f063d4SToby Isaac DMLabel canonLabel; 978f9f063d4SToby Isaac 979f9f063d4SToby Isaac ierr = DMPlexGetLabel(refTree,"canonical",&canonLabel);CHKERRQ(ierr); 980f9f063d4SToby Isaac if (canonLabel) { 981f9f063d4SToby Isaac PetscInt i; 982f9f063d4SToby Isaac 983f9f063d4SToby Isaac for (i = 0; i < size; i++) { 984f9f063d4SToby Isaac PetscInt canon; 985f9f063d4SToby Isaac ierr = DMLabelGetValue(canonLabel, mesh->childIDs[i], &canon);CHKERRQ(ierr); 986f9f063d4SToby Isaac if (canon >= 0) { 987f9f063d4SToby Isaac mesh->childIDs[i] = canon; 988f9f063d4SToby Isaac } 989f9f063d4SToby Isaac } 990f9f063d4SToby Isaac } 991f7c74593SToby Isaac mesh->computeanchormatrix = DMPlexComputeAnchorMatrix_Tree_FromReference; 992f7c74593SToby Isaac } 993f7c74593SToby Isaac else { 994f7c74593SToby Isaac mesh->computeanchormatrix = DMPlexComputeAnchorMatrix_Tree_Direct; 995f9f063d4SToby Isaac } 996f9f063d4SToby Isaac ierr = DMPlexTreeSymmetrize(dm);CHKERRQ(ierr); 997f9f063d4SToby Isaac if (computeCanonical) { 998f9f063d4SToby Isaac PetscInt d, dim; 999f9f063d4SToby Isaac 1000f9f063d4SToby Isaac /* add the canonical label */ 100128f4b327SMatthew G. Knepley ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 1002f9f063d4SToby Isaac ierr = DMPlexCreateLabel(dm,"canonical");CHKERRQ(ierr); 1003f9f063d4SToby Isaac for (d = 0; d <= dim; d++) { 1004f9f063d4SToby Isaac PetscInt p, dStart, dEnd, canon = -1, cNumChildren; 1005f9f063d4SToby Isaac const PetscInt *cChildren; 1006f9f063d4SToby Isaac 1007f9f063d4SToby Isaac ierr = DMPlexGetDepthStratum(dm,d,&dStart,&dEnd);CHKERRQ(ierr); 1008f9f063d4SToby Isaac for (p = dStart; p < dEnd; p++) { 1009f9f063d4SToby Isaac ierr = DMPlexGetTreeChildren(dm,p,&cNumChildren,&cChildren);CHKERRQ(ierr); 1010f9f063d4SToby Isaac if (cNumChildren) { 1011f9f063d4SToby Isaac canon = p; 1012f9f063d4SToby Isaac break; 1013f9f063d4SToby Isaac } 1014f9f063d4SToby Isaac } 1015f9f063d4SToby Isaac if (canon == -1) continue; 1016f9f063d4SToby Isaac for (p = dStart; p < dEnd; p++) { 1017f9f063d4SToby Isaac PetscInt numChildren, i; 1018f9f063d4SToby Isaac const PetscInt *children; 1019f9f063d4SToby Isaac 1020f9f063d4SToby Isaac ierr = DMPlexGetTreeChildren(dm,p,&numChildren,&children);CHKERRQ(ierr); 1021f9f063d4SToby Isaac if (numChildren) { 1022f9f063d4SToby 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); 1023f9f063d4SToby Isaac ierr = DMPlexSetLabelValue(dm,"canonical",p,canon);CHKERRQ(ierr); 1024f9f063d4SToby Isaac for (i = 0; i < numChildren; i++) { 1025f9f063d4SToby Isaac ierr = DMPlexSetLabelValue(dm,"canonical",children[i],cChildren[i]);CHKERRQ(ierr); 1026f9f063d4SToby Isaac } 1027f9f063d4SToby Isaac } 1028f9f063d4SToby Isaac } 1029f9f063d4SToby Isaac } 1030f9f063d4SToby Isaac } 1031776742edSToby Isaac if (exchangeSupports) { 1032776742edSToby Isaac ierr = DMPlexTreeExchangeSupports(dm);CHKERRQ(ierr); 1033776742edSToby Isaac } 1034f7c74593SToby Isaac mesh->createanchors = DMPlexCreateAnchors_Tree; 1035f7c74593SToby Isaac /* reset anchors */ 1036f7c74593SToby Isaac ierr = DMPlexSetAnchors(dm,NULL,NULL);CHKERRQ(ierr); 1037f9f063d4SToby Isaac PetscFunctionReturn(0); 1038f9f063d4SToby Isaac } 1039f9f063d4SToby Isaac 1040f9f063d4SToby Isaac #undef __FUNCT__ 10410b7167a0SToby Isaac #define __FUNCT__ "DMPlexSetTree" 10420b7167a0SToby Isaac /*@ 10430b7167a0SToby Isaac DMPlexSetTree - set the tree that describes the hierarchy of non-conforming mesh points. This routine also creates 10440b7167a0SToby Isaac the point-to-point constraints determined by the tree: a point is constained to the points in the closure of its 10450b7167a0SToby Isaac tree root. 10460b7167a0SToby Isaac 10470b7167a0SToby Isaac Collective on dm 10480b7167a0SToby Isaac 10490b7167a0SToby Isaac Input Parameters: 10500b7167a0SToby Isaac + dm - the DMPlex object 10510b7167a0SToby Isaac . parentSection - a section describing the tree: a point has a parent if it has 1 dof in the section; the section 10520b7167a0SToby Isaac offset indexes the parent and childID list; the reference count of parentSection is incremented 10530b7167a0SToby Isaac . parents - a list of the point parents; copied, can be destroyed 10540b7167a0SToby Isaac - childIDs - identifies the relationship of the child point to the parent point; if there is a reference tree, then 10550b7167a0SToby Isaac the child corresponds to the point in the reference tree with index childIDs; copied, can be destroyed 10560b7167a0SToby Isaac 10570b7167a0SToby Isaac Level: intermediate 10580b7167a0SToby Isaac 1059a17985deSToby Isaac .seealso: DMPlexGetTree(), DMPlexSetReferenceTree(), DMPlexSetAnchors(), DMPlexGetTreeParent(), DMPlexGetTreeChildren() 10600b7167a0SToby Isaac @*/ 1061b2f41788SToby Isaac PetscErrorCode DMPlexSetTree(DM dm, PetscSection parentSection, PetscInt parents[], PetscInt childIDs[]) 10620b7167a0SToby Isaac { 10630b7167a0SToby Isaac PetscErrorCode ierr; 10640b7167a0SToby Isaac 10650b7167a0SToby Isaac PetscFunctionBegin; 1066776742edSToby Isaac ierr = DMPlexSetTree_Internal(dm,parentSection,parents,childIDs,PETSC_FALSE,PETSC_TRUE);CHKERRQ(ierr); 10670b7167a0SToby Isaac PetscFunctionReturn(0); 10680b7167a0SToby Isaac } 10690b7167a0SToby Isaac 10700b7167a0SToby Isaac #undef __FUNCT__ 1071b2f41788SToby Isaac #define __FUNCT__ "DMPlexGetTree" 1072b2f41788SToby Isaac /*@ 1073b2f41788SToby Isaac DMPlexGetTree - get the tree that describes the hierarchy of non-conforming mesh points. 1074b2f41788SToby Isaac Collective on dm 1075b2f41788SToby Isaac 1076b2f41788SToby Isaac Input Parameters: 1077b2f41788SToby Isaac . dm - the DMPlex object 1078b2f41788SToby Isaac 1079b2f41788SToby Isaac Output Parameters: 1080b2f41788SToby Isaac + parentSection - a section describing the tree: a point has a parent if it has 1 dof in the section; the section 1081b2f41788SToby Isaac offset indexes the parent and childID list 1082b2f41788SToby Isaac . parents - a list of the point parents 1083b2f41788SToby Isaac . childIDs - identifies the relationship of the child point to the parent point; if there is a reference tree, then 1084b2f41788SToby Isaac the child corresponds to the point in the reference tree with index childID 1085b2f41788SToby Isaac . childSection - the inverse of the parent section 1086b2f41788SToby Isaac - children - a list of the point children 1087b2f41788SToby Isaac 1088b2f41788SToby Isaac Level: intermediate 1089b2f41788SToby Isaac 1090a17985deSToby Isaac .seealso: DMPlexSetTree(), DMPlexSetReferenceTree(), DMPlexSetAnchors(), DMPlexGetTreeParent(), DMPlexGetTreeChildren() 1091b2f41788SToby Isaac @*/ 1092b2f41788SToby Isaac PetscErrorCode DMPlexGetTree(DM dm, PetscSection *parentSection, PetscInt *parents[], PetscInt *childIDs[], PetscSection *childSection, PetscInt *children[]) 1093b2f41788SToby Isaac { 1094b2f41788SToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 1095b2f41788SToby Isaac 1096b2f41788SToby Isaac PetscFunctionBegin; 1097b2f41788SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1098b2f41788SToby Isaac if (parentSection) *parentSection = mesh->parentSection; 1099b2f41788SToby Isaac if (parents) *parents = mesh->parents; 1100b2f41788SToby Isaac if (childIDs) *childIDs = mesh->childIDs; 1101b2f41788SToby Isaac if (childSection) *childSection = mesh->childSection; 1102b2f41788SToby Isaac if (children) *children = mesh->children; 1103b2f41788SToby Isaac PetscFunctionReturn(0); 1104b2f41788SToby Isaac } 1105b2f41788SToby Isaac 1106b2f41788SToby Isaac #undef __FUNCT__ 1107d961a43aSToby Isaac #define __FUNCT__ "DMPlexGetTreeParent" 1108d961a43aSToby Isaac /*@ 1109d961a43aSToby Isaac DMPlexGetTreeParent - get the parent of a point in the tree describing the point hierarchy (not the Sieve DAG) 1110d961a43aSToby Isaac 1111d961a43aSToby Isaac Input Parameters: 1112d961a43aSToby Isaac + dm - the DMPlex object 1113d961a43aSToby Isaac - point - the query point 1114d961a43aSToby Isaac 1115d961a43aSToby Isaac Output Parameters: 1116d961a43aSToby Isaac + parent - if not NULL, set to the parent of the point, or the point itself if the point does not have a parent 1117d961a43aSToby Isaac - childID - if not NULL, set to the child ID of the point with respect to its parent, or 0 if the point 1118d961a43aSToby Isaac does not have a parent 1119d961a43aSToby Isaac 1120d961a43aSToby Isaac Level: intermediate 1121d961a43aSToby Isaac 1122d961a43aSToby Isaac .seealso: DMPlexSetTree(), DMPlexGetTree(), DMPlexGetTreeChildren() 1123d961a43aSToby Isaac @*/ 1124d961a43aSToby Isaac PetscErrorCode DMPlexGetTreeParent(DM dm, PetscInt point, PetscInt *parent, PetscInt *childID) 1125d961a43aSToby Isaac { 1126d961a43aSToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 1127d961a43aSToby Isaac PetscSection pSec; 1128d961a43aSToby Isaac PetscErrorCode ierr; 1129d961a43aSToby Isaac 1130d961a43aSToby Isaac PetscFunctionBegin; 1131d961a43aSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1132d961a43aSToby Isaac pSec = mesh->parentSection; 1133d961a43aSToby Isaac if (pSec && point >= pSec->pStart && point < pSec->pEnd) { 1134d961a43aSToby Isaac PetscInt dof; 1135d961a43aSToby Isaac 1136d961a43aSToby Isaac ierr = PetscSectionGetDof (pSec, point, &dof);CHKERRQ(ierr); 1137d961a43aSToby Isaac if (dof) { 1138d961a43aSToby Isaac PetscInt off; 1139d961a43aSToby Isaac 1140d961a43aSToby Isaac ierr = PetscSectionGetOffset (pSec, point, &off);CHKERRQ(ierr); 1141d961a43aSToby Isaac if (parent) *parent = mesh->parents[off]; 1142d961a43aSToby Isaac if (childID) *childID = mesh->childIDs[off]; 1143d961a43aSToby Isaac PetscFunctionReturn(0); 1144d961a43aSToby Isaac } 1145d961a43aSToby Isaac } 1146d961a43aSToby Isaac if (parent) { 1147d961a43aSToby Isaac *parent = point; 1148d961a43aSToby Isaac } 1149d961a43aSToby Isaac if (childID) { 1150d961a43aSToby Isaac *childID = 0; 1151d961a43aSToby Isaac } 1152d961a43aSToby Isaac PetscFunctionReturn(0); 1153d961a43aSToby Isaac } 1154d961a43aSToby Isaac 1155d961a43aSToby Isaac #undef __FUNCT__ 1156d961a43aSToby Isaac #define __FUNCT__ "DMPlexGetTreeChildren" 1157d961a43aSToby Isaac /*@C 1158d961a43aSToby Isaac DMPlexGetTreeChildren - get the children of a point in the tree describing the point hierarchy (not the Sieve DAG) 1159d961a43aSToby Isaac 1160d961a43aSToby Isaac Input Parameters: 1161d961a43aSToby Isaac + dm - the DMPlex object 1162d961a43aSToby Isaac - point - the query point 1163d961a43aSToby Isaac 1164d961a43aSToby Isaac Output Parameters: 1165d961a43aSToby Isaac + numChildren - if not NULL, set to the number of children 1166d961a43aSToby Isaac - children - if not NULL, set to a list children, or set to NULL if the point has no children 1167d961a43aSToby Isaac 1168d961a43aSToby Isaac Level: intermediate 1169d961a43aSToby Isaac 1170d961a43aSToby Isaac Fortran Notes: 1171d961a43aSToby Isaac Since it returns an array, this routine is only available in Fortran 90, and you must 1172d961a43aSToby Isaac include petsc.h90 in your code. 1173d961a43aSToby Isaac 1174d961a43aSToby Isaac .seealso: DMPlexSetTree(), DMPlexGetTree(), DMPlexGetTreeParent() 1175d961a43aSToby Isaac @*/ 1176d961a43aSToby Isaac PetscErrorCode DMPlexGetTreeChildren(DM dm, PetscInt point, PetscInt *numChildren, const PetscInt *children[]) 1177d961a43aSToby Isaac { 1178d961a43aSToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 1179d961a43aSToby Isaac PetscSection childSec; 1180d961a43aSToby Isaac PetscInt dof = 0; 1181d961a43aSToby Isaac PetscErrorCode ierr; 1182d961a43aSToby Isaac 1183d961a43aSToby Isaac PetscFunctionBegin; 1184d961a43aSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 1185d961a43aSToby Isaac childSec = mesh->childSection; 1186d961a43aSToby Isaac if (childSec && point >= childSec->pStart && point < childSec->pEnd) { 1187d961a43aSToby Isaac ierr = PetscSectionGetDof (childSec, point, &dof);CHKERRQ(ierr); 1188d961a43aSToby Isaac } 1189d961a43aSToby Isaac if (numChildren) *numChildren = dof; 1190d961a43aSToby Isaac if (children) { 1191d961a43aSToby Isaac if (dof) { 1192d961a43aSToby Isaac PetscInt off; 1193d961a43aSToby Isaac 1194d961a43aSToby Isaac ierr = PetscSectionGetOffset (childSec, point, &off);CHKERRQ(ierr); 1195d961a43aSToby Isaac *children = &mesh->children[off]; 1196d961a43aSToby Isaac } 1197d961a43aSToby Isaac else { 1198d961a43aSToby Isaac *children = NULL; 1199d961a43aSToby Isaac } 1200d961a43aSToby Isaac } 1201d961a43aSToby Isaac PetscFunctionReturn(0); 1202d961a43aSToby Isaac } 12030c37af3bSToby Isaac 12040c37af3bSToby Isaac #undef __FUNCT__ 1205f7c74593SToby Isaac #define __FUNCT__ "DMPlexComputeAnchorMatrix_Tree_Direct" 1206f7c74593SToby Isaac static PetscErrorCode DMPlexComputeAnchorMatrix_Tree_Direct(DM dm, PetscSection section, PetscSection cSec, Mat cMat) 12070c37af3bSToby Isaac { 12080c37af3bSToby Isaac PetscDS ds; 12090c37af3bSToby Isaac PetscInt spdim; 12100c37af3bSToby Isaac PetscInt numFields, f, c, cStart, cEnd, pStart, pEnd, conStart, conEnd; 12110c37af3bSToby Isaac const PetscInt *anchors; 1212f7c74593SToby Isaac PetscSection aSec; 12130c37af3bSToby Isaac PetscReal *v0, *v0parent, *vtmp, *J, *Jparent, *invJparent, detJ, detJparent; 12140c37af3bSToby Isaac IS aIS; 12150c37af3bSToby Isaac PetscErrorCode ierr; 12160c37af3bSToby Isaac 12170c37af3bSToby Isaac PetscFunctionBegin; 12180c37af3bSToby Isaac ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 12190c37af3bSToby Isaac ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); 12200c37af3bSToby Isaac ierr = PetscDSGetNumFields(ds,&numFields);CHKERRQ(ierr); 12210c37af3bSToby Isaac ierr = DMPlexGetHeightStratum(dm,0,&cStart,&cEnd);CHKERRQ(ierr); 1222a17985deSToby Isaac ierr = DMPlexGetAnchors(dm,&aSec,&aIS);CHKERRQ(ierr); 12230c37af3bSToby Isaac ierr = ISGetIndices(aIS,&anchors);CHKERRQ(ierr); 12240c37af3bSToby Isaac ierr = PetscSectionGetChart(cSec,&conStart,&conEnd);CHKERRQ(ierr); 122528f4b327SMatthew G. Knepley ierr = DMGetDimension(dm,&spdim);CHKERRQ(ierr); 12260c37af3bSToby Isaac ierr = PetscMalloc6(spdim,&v0,spdim,&v0parent,spdim,&vtmp,spdim*spdim,&J,spdim*spdim,&Jparent,spdim*spdim,&invJparent);CHKERRQ(ierr); 12270c37af3bSToby Isaac 12280c37af3bSToby Isaac for (f = 0; f < numFields; f++) { 12290dd1b1feSToby Isaac PetscObject disc; 12300dd1b1feSToby Isaac PetscFE fe = NULL; 12310dd1b1feSToby Isaac PetscFV fv = NULL; 12320dd1b1feSToby Isaac PetscClassId id; 12330c37af3bSToby Isaac PetscDualSpace space; 12340c37af3bSToby Isaac PetscInt i, j, k, nPoints, offset; 12350c37af3bSToby Isaac PetscInt fSize, fComp; 1236c111c6b7SMatthew G. Knepley PetscReal *B = NULL; 12370c37af3bSToby Isaac PetscReal *weights, *pointsRef, *pointsReal; 12380c37af3bSToby Isaac Mat Amat, Bmat, Xmat; 12390c37af3bSToby Isaac 12400dd1b1feSToby Isaac ierr = PetscDSGetDiscretization(ds,f,&disc);CHKERRQ(ierr); 12410dd1b1feSToby Isaac ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 12420dd1b1feSToby Isaac if (id == PETSCFE_CLASSID) { 12430dd1b1feSToby Isaac fe = (PetscFE) disc; 12440c37af3bSToby Isaac ierr = PetscFEGetDualSpace(fe,&space);CHKERRQ(ierr); 12450c37af3bSToby Isaac ierr = PetscDualSpaceGetDimension(space,&fSize);CHKERRQ(ierr); 12460c37af3bSToby Isaac ierr = PetscFEGetNumComponents(fe,&fComp);CHKERRQ(ierr); 12470dd1b1feSToby Isaac } 12480dd1b1feSToby Isaac else if (id == PETSCFV_CLASSID) { 12490dd1b1feSToby Isaac fv = (PetscFV) disc; 12500dd1b1feSToby Isaac ierr = PetscFVGetDualSpace(fv,&space);CHKERRQ(ierr); 12510dd1b1feSToby Isaac ierr = PetscDualSpaceGetDimension(space,&fSize);CHKERRQ(ierr); 12520dd1b1feSToby Isaac ierr = PetscFVGetNumComponents(fv,&fComp);CHKERRQ(ierr); 12530dd1b1feSToby Isaac } 12540dd1b1feSToby Isaac else SETERRQ1(PetscObjectComm(disc),PETSC_ERR_ARG_UNKNOWN_TYPE, "PetscDS discretization id %d not recognized.", id); 12550dd1b1feSToby Isaac 12560c37af3bSToby Isaac ierr = MatCreate(PETSC_COMM_SELF,&Amat);CHKERRQ(ierr); 12570c37af3bSToby Isaac ierr = MatSetSizes(Amat,fSize,fSize,fSize,fSize);CHKERRQ(ierr); 12580c37af3bSToby Isaac ierr = MatSetType(Amat,MATSEQDENSE);CHKERRQ(ierr); 12590c37af3bSToby Isaac ierr = MatSetUp(Amat);CHKERRQ(ierr); 12600c37af3bSToby Isaac ierr = MatDuplicate(Amat,MAT_DO_NOT_COPY_VALUES,&Bmat);CHKERRQ(ierr); 12610c37af3bSToby Isaac ierr = MatDuplicate(Amat,MAT_DO_NOT_COPY_VALUES,&Xmat);CHKERRQ(ierr); 12620c37af3bSToby Isaac nPoints = 0; 12630c37af3bSToby Isaac for (i = 0; i < fSize; i++) { 12640c37af3bSToby Isaac PetscInt qPoints; 12650c37af3bSToby Isaac PetscQuadrature quad; 12660c37af3bSToby Isaac 12670c37af3bSToby Isaac ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr); 12680c37af3bSToby Isaac ierr = PetscQuadratureGetData(quad,NULL,&qPoints,NULL,NULL);CHKERRQ(ierr); 12690c37af3bSToby Isaac nPoints += qPoints; 12700c37af3bSToby Isaac } 12710c37af3bSToby Isaac ierr = PetscMalloc3(nPoints,&weights,spdim*nPoints,&pointsRef,spdim*nPoints,&pointsReal);CHKERRQ(ierr); 12720c37af3bSToby Isaac offset = 0; 12730c37af3bSToby Isaac for (i = 0; i < fSize; i++) { 12740c37af3bSToby Isaac PetscInt qPoints; 12750c37af3bSToby Isaac const PetscReal *p, *w; 12760c37af3bSToby Isaac PetscQuadrature quad; 12770c37af3bSToby Isaac 12780c37af3bSToby Isaac ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr); 12790c37af3bSToby Isaac ierr = PetscQuadratureGetData(quad,NULL,&qPoints,&p,&w);CHKERRQ(ierr); 12800c37af3bSToby Isaac ierr = PetscMemcpy(weights+offset,w,qPoints*sizeof(*w));CHKERRQ(ierr); 12810c37af3bSToby Isaac ierr = PetscMemcpy(pointsRef+spdim*offset,p,spdim*qPoints*sizeof(*p));CHKERRQ(ierr); 12820c37af3bSToby Isaac offset += qPoints; 12830c37af3bSToby Isaac } 12840dd1b1feSToby Isaac if (id == PETSCFE_CLASSID) { 12850c37af3bSToby Isaac ierr = PetscFEGetTabulation(fe,nPoints,pointsRef,&B,NULL,NULL);CHKERRQ(ierr); 12860dd1b1feSToby Isaac } 12870dd1b1feSToby Isaac else { 12880dd1b1feSToby Isaac ierr = PetscFVGetTabulation(fv,nPoints,pointsRef,&B,NULL,NULL);CHKERRQ(ierr); 12890dd1b1feSToby Isaac } 12900c37af3bSToby Isaac offset = 0; 12910c37af3bSToby Isaac for (i = 0; i < fSize; i++) { 12920c37af3bSToby Isaac PetscInt qPoints; 12930c37af3bSToby Isaac PetscQuadrature quad; 12940c37af3bSToby Isaac 12950c37af3bSToby Isaac ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr); 12960c37af3bSToby Isaac ierr = PetscQuadratureGetData(quad,NULL,&qPoints,NULL,NULL);CHKERRQ(ierr); 12970c37af3bSToby Isaac for (j = 0; j < fSize; j++) { 12980c37af3bSToby Isaac PetscScalar val = 0.; 12990c37af3bSToby Isaac 13000c37af3bSToby Isaac for (k = 0; k < qPoints; k++) { 13010c37af3bSToby Isaac val += B[((offset + k) * fSize + j) * fComp] * weights[k]; 13020c37af3bSToby Isaac } 13030c37af3bSToby Isaac ierr = MatSetValue(Amat,i,j,val,INSERT_VALUES);CHKERRQ(ierr); 13040c37af3bSToby Isaac } 13050c37af3bSToby Isaac offset += qPoints; 13060c37af3bSToby Isaac } 13070dd1b1feSToby Isaac if (id == PETSCFE_CLASSID) { 13080c37af3bSToby Isaac ierr = PetscFERestoreTabulation(fe,nPoints,pointsRef,&B,NULL,NULL);CHKERRQ(ierr); 13090dd1b1feSToby Isaac } 13100dd1b1feSToby Isaac else { 13110dd1b1feSToby Isaac ierr = PetscFVRestoreTabulation(fv,nPoints,pointsRef,&B,NULL,NULL);CHKERRQ(ierr); 13120dd1b1feSToby Isaac } 13130c37af3bSToby Isaac ierr = MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13140c37af3bSToby Isaac ierr = MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13150c37af3bSToby Isaac ierr = MatLUFactor(Amat,NULL,NULL,NULL);CHKERRQ(ierr); 13160c37af3bSToby Isaac for (c = cStart; c < cEnd; c++) { 13170c37af3bSToby Isaac PetscInt parent; 13180c37af3bSToby Isaac PetscInt closureSize, closureSizeP, *closure = NULL, *closureP = NULL; 13190c37af3bSToby Isaac PetscInt *childOffsets, *parentOffsets; 13200c37af3bSToby Isaac 13210c37af3bSToby Isaac ierr = DMPlexGetTreeParent(dm,c,&parent,NULL);CHKERRQ(ierr); 13220c37af3bSToby Isaac if (parent == c) continue; 13230c37af3bSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 13240c37af3bSToby Isaac for (i = 0; i < closureSize; i++) { 13250c37af3bSToby Isaac PetscInt p = closure[2*i]; 13260c37af3bSToby Isaac PetscInt conDof; 13270c37af3bSToby Isaac 13280c37af3bSToby Isaac if (p < conStart || p >= conEnd) continue; 13290c37af3bSToby Isaac if (numFields > 1) { 13300c37af3bSToby Isaac ierr = PetscSectionGetFieldDof(cSec,p,f,&conDof);CHKERRQ(ierr); 13310c37af3bSToby Isaac } 13320c37af3bSToby Isaac else { 13330c37af3bSToby Isaac ierr = PetscSectionGetDof(cSec,p,&conDof);CHKERRQ(ierr); 13340c37af3bSToby Isaac } 13350c37af3bSToby Isaac if (conDof) break; 13360c37af3bSToby Isaac } 13370c37af3bSToby Isaac if (i == closureSize) { 13380c37af3bSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 13390c37af3bSToby Isaac continue; 13400c37af3bSToby Isaac } 13410c37af3bSToby Isaac 134273a7f2aaSMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, c, NULL, v0, J, NULL, &detJ);CHKERRQ(ierr); 134373a7f2aaSMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, parent, NULL, v0parent, Jparent, invJparent, &detJparent);CHKERRQ(ierr); 13440c37af3bSToby Isaac for (i = 0; i < nPoints; i++) { 13450c37af3bSToby Isaac CoordinatesRefToReal(spdim, spdim, v0, J, &pointsRef[i*spdim],vtmp);CHKERRQ(ierr); 13460c37af3bSToby Isaac CoordinatesRealToRef(spdim, spdim, v0parent, invJparent, vtmp, &pointsReal[i*spdim]);CHKERRQ(ierr); 13470c37af3bSToby Isaac } 13480dd1b1feSToby Isaac if (id == PETSCFE_CLASSID) { 13490c37af3bSToby Isaac ierr = PetscFEGetTabulation(fe,nPoints,pointsReal,&B,NULL,NULL);CHKERRQ(ierr); 13500dd1b1feSToby Isaac } 13510dd1b1feSToby Isaac else { 13520dd1b1feSToby Isaac ierr = PetscFVGetTabulation(fv,nPoints,pointsReal,&B,NULL,NULL);CHKERRQ(ierr); 13530dd1b1feSToby Isaac } 13540c37af3bSToby Isaac offset = 0; 13550c37af3bSToby Isaac for (i = 0; i < fSize; i++) { 13560c37af3bSToby Isaac PetscInt qPoints; 13570c37af3bSToby Isaac PetscQuadrature quad; 13580c37af3bSToby Isaac 13590c37af3bSToby Isaac ierr = PetscDualSpaceGetFunctional(space,i,&quad);CHKERRQ(ierr); 13600c37af3bSToby Isaac ierr = PetscQuadratureGetData(quad,NULL,&qPoints,NULL,NULL);CHKERRQ(ierr); 13610c37af3bSToby Isaac for (j = 0; j < fSize; j++) { 13620c37af3bSToby Isaac PetscScalar val = 0.; 13630c37af3bSToby Isaac 13640c37af3bSToby Isaac for (k = 0; k < qPoints; k++) { 13650c37af3bSToby Isaac val += B[((offset + k) * fSize + j) * fComp] * weights[k]; 13660c37af3bSToby Isaac } 13670c37af3bSToby Isaac MatSetValue(Bmat,i,j,val,INSERT_VALUES);CHKERRQ(ierr); 13680c37af3bSToby Isaac } 13690c37af3bSToby Isaac offset += qPoints; 13700c37af3bSToby Isaac } 13710dd1b1feSToby Isaac if (id == PETSCFE_CLASSID) { 13720c37af3bSToby Isaac ierr = PetscFERestoreTabulation(fe,nPoints,pointsReal,&B,NULL,NULL);CHKERRQ(ierr); 13730dd1b1feSToby Isaac } 13740dd1b1feSToby Isaac else { 13750dd1b1feSToby Isaac ierr = PetscFVRestoreTabulation(fv,nPoints,pointsReal,&B,NULL,NULL);CHKERRQ(ierr); 13760dd1b1feSToby Isaac } 13770c37af3bSToby Isaac ierr = MatAssemblyBegin(Bmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13780c37af3bSToby Isaac ierr = MatAssemblyEnd(Bmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 13790c37af3bSToby Isaac ierr = MatMatSolve(Amat,Bmat,Xmat);CHKERRQ(ierr); 13800c37af3bSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSizeP,&closureP);CHKERRQ(ierr); 13810c37af3bSToby Isaac ierr = PetscMalloc2(closureSize+1,&childOffsets,closureSizeP+1,&parentOffsets);CHKERRQ(ierr); 13820c37af3bSToby Isaac childOffsets[0] = 0; 13830c37af3bSToby Isaac for (i = 0; i < closureSize; i++) { 13840c37af3bSToby Isaac PetscInt p = closure[2*i]; 13850c37af3bSToby Isaac PetscInt dof; 13860c37af3bSToby Isaac 13870c37af3bSToby Isaac if (numFields > 1) { 13880c37af3bSToby Isaac ierr = PetscSectionGetFieldDof(section,p,f,&dof);CHKERRQ(ierr); 13890c37af3bSToby Isaac } 13900c37af3bSToby Isaac else { 13910c37af3bSToby Isaac ierr = PetscSectionGetDof(section,p,&dof);CHKERRQ(ierr); 13920c37af3bSToby Isaac } 13930c37af3bSToby Isaac childOffsets[i+1]=childOffsets[i]+dof / fComp; 13940c37af3bSToby Isaac } 13950c37af3bSToby Isaac parentOffsets[0] = 0; 13960c37af3bSToby Isaac for (i = 0; i < closureSizeP; i++) { 13970c37af3bSToby Isaac PetscInt p = closureP[2*i]; 13980c37af3bSToby Isaac PetscInt dof; 13990c37af3bSToby Isaac 14000c37af3bSToby Isaac if (numFields > 1) { 14010c37af3bSToby Isaac ierr = PetscSectionGetFieldDof(section,p,f,&dof);CHKERRQ(ierr); 14020c37af3bSToby Isaac } 14030c37af3bSToby Isaac else { 14040c37af3bSToby Isaac ierr = PetscSectionGetDof(section,p,&dof);CHKERRQ(ierr); 14050c37af3bSToby Isaac } 14060c37af3bSToby Isaac parentOffsets[i+1]=parentOffsets[i]+dof / fComp; 14070c37af3bSToby Isaac } 14080c37af3bSToby Isaac for (i = 0; i < closureSize; i++) { 14090c37af3bSToby Isaac PetscInt conDof, conOff, aDof, aOff; 14100c37af3bSToby Isaac PetscInt p = closure[2*i]; 14110c37af3bSToby Isaac PetscInt o = closure[2*i+1]; 14120c37af3bSToby Isaac 14130c37af3bSToby Isaac if (p < conStart || p >= conEnd) continue; 14140c37af3bSToby Isaac if (numFields > 1) { 14150c37af3bSToby Isaac ierr = PetscSectionGetFieldDof(cSec,p,f,&conDof);CHKERRQ(ierr); 14160c37af3bSToby Isaac ierr = PetscSectionGetFieldOffset(cSec,p,f,&conOff);CHKERRQ(ierr); 14170c37af3bSToby Isaac } 14180c37af3bSToby Isaac else { 14190c37af3bSToby Isaac ierr = PetscSectionGetDof(cSec,p,&conDof);CHKERRQ(ierr); 14200c37af3bSToby Isaac ierr = PetscSectionGetOffset(cSec,p,&conOff);CHKERRQ(ierr); 14210c37af3bSToby Isaac } 14220c37af3bSToby Isaac if (!conDof) continue; 14230c37af3bSToby Isaac ierr = PetscSectionGetDof(aSec,p,&aDof);CHKERRQ(ierr); 14240c37af3bSToby Isaac ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr); 14250c37af3bSToby Isaac for (k = 0; k < aDof; k++) { 14260c37af3bSToby Isaac PetscInt a = anchors[aOff + k]; 14270c37af3bSToby Isaac PetscInt aSecDof, aSecOff; 14280c37af3bSToby Isaac 14290c37af3bSToby Isaac if (numFields > 1) { 14300c37af3bSToby Isaac ierr = PetscSectionGetFieldDof(section,a,f,&aSecDof);CHKERRQ(ierr); 14310c37af3bSToby Isaac ierr = PetscSectionGetFieldOffset(section,a,f,&aSecOff);CHKERRQ(ierr); 14320c37af3bSToby Isaac } 14330c37af3bSToby Isaac else { 14340c37af3bSToby Isaac ierr = PetscSectionGetDof(section,a,&aSecDof);CHKERRQ(ierr); 14350c37af3bSToby Isaac ierr = PetscSectionGetOffset(section,a,&aSecOff);CHKERRQ(ierr); 14360c37af3bSToby Isaac } 14370c37af3bSToby Isaac if (!aSecDof) continue; 14380c37af3bSToby Isaac 14390c37af3bSToby Isaac for (j = 0; j < closureSizeP; j++) { 14400c37af3bSToby Isaac PetscInt q = closureP[2*j]; 14410c37af3bSToby Isaac PetscInt oq = closureP[2*j+1]; 14420c37af3bSToby Isaac 14430c37af3bSToby Isaac if (q == a) { 14440c37af3bSToby Isaac PetscInt r, s, t; 14450c37af3bSToby Isaac 14460c37af3bSToby Isaac for (r = childOffsets[i]; r < childOffsets[i+1]; r++) { 14470c37af3bSToby Isaac for (s = parentOffsets[j]; s < parentOffsets[j+1]; s++) { 14480c37af3bSToby Isaac PetscScalar val; 14490c37af3bSToby Isaac PetscInt insertCol, insertRow; 14500c37af3bSToby Isaac 14510c37af3bSToby Isaac ierr = MatGetValue(Xmat,r,s,&val);CHKERRQ(ierr); 14520c37af3bSToby Isaac if (o >= 0) { 14530c37af3bSToby Isaac insertRow = conOff + fComp * (r - childOffsets[i]); 14540c37af3bSToby Isaac } 14550c37af3bSToby Isaac else { 14560c37af3bSToby Isaac insertRow = conOff + fComp * (childOffsets[i + 1] - 1 - r); 14570c37af3bSToby Isaac } 14580c37af3bSToby Isaac if (oq >= 0) { 14590c37af3bSToby Isaac insertCol = aSecOff + fComp * (s - parentOffsets[j]); 14600c37af3bSToby Isaac } 14610c37af3bSToby Isaac else { 14620c37af3bSToby Isaac insertCol = aSecOff + fComp * (parentOffsets[j + 1] - 1 - s); 14630c37af3bSToby Isaac } 14640c37af3bSToby Isaac for (t = 0; t < fComp; t++) { 14650c37af3bSToby Isaac ierr = MatSetValue(cMat,insertRow + t,insertCol + t,val,INSERT_VALUES);CHKERRQ(ierr); 14660c37af3bSToby Isaac } 14670c37af3bSToby Isaac } 14680c37af3bSToby Isaac } 14690c37af3bSToby Isaac } 14700c37af3bSToby Isaac } 14710c37af3bSToby Isaac } 14720c37af3bSToby Isaac } 14730c37af3bSToby Isaac ierr = PetscFree2(childOffsets,parentOffsets);CHKERRQ(ierr); 14740c37af3bSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,c,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 14750c37af3bSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSizeP,&closureP);CHKERRQ(ierr); 14760c37af3bSToby Isaac } 14770c37af3bSToby Isaac ierr = MatDestroy(&Amat);CHKERRQ(ierr); 14780c37af3bSToby Isaac ierr = MatDestroy(&Bmat);CHKERRQ(ierr); 14790c37af3bSToby Isaac ierr = MatDestroy(&Xmat);CHKERRQ(ierr); 14800c37af3bSToby Isaac ierr = PetscFree3(weights,pointsRef,pointsReal);CHKERRQ(ierr); 14810c37af3bSToby Isaac } 14820c37af3bSToby Isaac ierr = MatAssemblyBegin(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14830c37af3bSToby Isaac ierr = MatAssemblyEnd(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14840c37af3bSToby Isaac ierr = PetscFree6(v0,v0parent,vtmp,J,Jparent,invJparent);CHKERRQ(ierr); 14850c37af3bSToby Isaac ierr = ISRestoreIndices(aIS,&anchors);CHKERRQ(ierr); 14860c37af3bSToby Isaac 14870c37af3bSToby Isaac PetscFunctionReturn(0); 14880c37af3bSToby Isaac } 148995a0b26dSToby Isaac 149095a0b26dSToby Isaac #undef __FUNCT__ 1491f7c74593SToby Isaac #define __FUNCT__ "DMPlexComputeAnchorMatrix_Tree_FromReference" 1492f7c74593SToby Isaac static PetscErrorCode DMPlexComputeAnchorMatrix_Tree_FromReference(DM dm, PetscSection section, PetscSection conSec, Mat cMat) 149395a0b26dSToby Isaac { 149495a0b26dSToby Isaac DM refTree; 149595a0b26dSToby Isaac PetscDS ds; 1496f7c74593SToby Isaac Mat refCmat; 149795a0b26dSToby Isaac PetscInt numFields, f, pRefStart, pRefEnd, p, *rows, *cols, maxDof, maxAnDof, *perm, *iperm, pStart, pEnd, conStart, conEnd, **refPointFieldN; 149895a0b26dSToby Isaac PetscScalar ***refPointFieldMats, *pointWork; 1499f7c74593SToby Isaac PetscSection refConSec, refAnSec, anSec, refSection; 150095a0b26dSToby Isaac IS refAnIS, anIS; 150195a0b26dSToby Isaac const PetscInt *refAnchors, *anchors; 150295a0b26dSToby Isaac PetscErrorCode ierr; 150395a0b26dSToby Isaac 150495a0b26dSToby Isaac PetscFunctionBegin; 150595a0b26dSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 150695a0b26dSToby Isaac ierr = DMGetDS(dm,&ds);CHKERRQ(ierr); 150795a0b26dSToby Isaac ierr = PetscDSGetNumFields(ds,&numFields);CHKERRQ(ierr); 150895a0b26dSToby Isaac ierr = DMPlexGetReferenceTree(dm,&refTree);CHKERRQ(ierr); 150995a0b26dSToby Isaac ierr = DMSetDS(refTree,ds);CHKERRQ(ierr); 1510f7c74593SToby Isaac ierr = DMGetDefaultConstraints(refTree,&refConSec,&refCmat);CHKERRQ(ierr); 151195a0b26dSToby Isaac ierr = DMPlexGetChart(refTree,&pRefStart,&pRefEnd);CHKERRQ(ierr); 1512a17985deSToby Isaac ierr = DMPlexGetAnchors(refTree,&refAnSec,&refAnIS);CHKERRQ(ierr); 1513a17985deSToby Isaac ierr = DMPlexGetAnchors(dm,&anSec,&anIS);CHKERRQ(ierr); 151495a0b26dSToby Isaac ierr = ISGetIndices(refAnIS,&refAnchors);CHKERRQ(ierr); 151595a0b26dSToby Isaac ierr = ISGetIndices(anIS,&anchors);CHKERRQ(ierr); 151695a0b26dSToby Isaac ierr = DMGetDefaultSection(refTree,&refSection);CHKERRQ(ierr); 151795a0b26dSToby Isaac ierr = PetscSectionGetChart(refConSec,&pRefStart,&pRefEnd);CHKERRQ(ierr); 151895a0b26dSToby Isaac ierr = PetscSectionGetChart(conSec,&conStart,&conEnd);CHKERRQ(ierr); 151995a0b26dSToby Isaac ierr = PetscMalloc1(pRefEnd-pRefStart,&refPointFieldMats);CHKERRQ(ierr); 152095a0b26dSToby Isaac ierr = PetscMalloc1(pRefEnd-pRefStart,&refPointFieldN);CHKERRQ(ierr); 152195a0b26dSToby Isaac ierr = PetscSectionGetMaxDof(refConSec,&maxDof);CHKERRQ(ierr); 152295a0b26dSToby Isaac ierr = PetscSectionGetMaxDof(refAnSec,&maxAnDof);CHKERRQ(ierr); 152395a0b26dSToby Isaac ierr = PetscMalloc1(maxDof,&rows);CHKERRQ(ierr); 152495a0b26dSToby Isaac ierr = PetscMalloc1(maxDof*maxAnDof,&cols);CHKERRQ(ierr); 152595a0b26dSToby Isaac ierr = PetscMalloc1(maxDof*maxDof*maxAnDof,&pointWork);CHKERRQ(ierr); 152695a0b26dSToby Isaac 152795a0b26dSToby Isaac /* step 1: get submats for every constrained point in the reference tree */ 152895a0b26dSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 152995a0b26dSToby Isaac PetscInt parent, closureSize, *closure = NULL, pDof; 153095a0b26dSToby Isaac 153195a0b26dSToby Isaac ierr = DMPlexGetTreeParent(refTree,p,&parent,NULL);CHKERRQ(ierr); 153295a0b26dSToby Isaac ierr = PetscSectionGetDof(refConSec,p,&pDof);CHKERRQ(ierr); 153395a0b26dSToby Isaac if (!pDof || parent == p) continue; 153495a0b26dSToby Isaac 153595a0b26dSToby Isaac ierr = PetscMalloc1(numFields,&refPointFieldMats[p-pRefStart]);CHKERRQ(ierr); 153695a0b26dSToby Isaac ierr = PetscCalloc1(numFields,&refPointFieldN[p-pRefStart]);CHKERRQ(ierr); 153795a0b26dSToby Isaac ierr = DMPlexGetTransitiveClosure(refTree,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 153895a0b26dSToby Isaac for (f = 0; f < numFields; f++) { 153995a0b26dSToby Isaac PetscInt cDof, cOff, numCols, r, i, fComp; 15400dd1b1feSToby Isaac PetscObject disc; 15410dd1b1feSToby Isaac PetscClassId id; 15420dd1b1feSToby Isaac PetscFE fe = NULL; 15430dd1b1feSToby Isaac PetscFV fv = NULL; 154495a0b26dSToby Isaac 15450dd1b1feSToby Isaac ierr = PetscDSGetDiscretization(ds,f,&disc);CHKERRQ(ierr); 15460dd1b1feSToby Isaac ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 15470dd1b1feSToby Isaac if (id == PETSCFE_CLASSID) { 15480dd1b1feSToby Isaac fe = (PetscFE) disc; 154995a0b26dSToby Isaac ierr = PetscFEGetNumComponents(fe,&fComp);CHKERRQ(ierr); 15500dd1b1feSToby Isaac } 15510dd1b1feSToby Isaac else if (id == PETSCFV_CLASSID) { 15520dd1b1feSToby Isaac fv = (PetscFV) disc; 15530dd1b1feSToby Isaac ierr = PetscFVGetNumComponents(fv,&fComp);CHKERRQ(ierr); 15540dd1b1feSToby Isaac } 15550dd1b1feSToby Isaac else SETERRQ1(PetscObjectComm(disc),PETSC_ERR_ARG_UNKNOWN_TYPE, "PetscDS discretization id %d not recognized.", id); 155695a0b26dSToby Isaac 155795a0b26dSToby Isaac if (numFields > 1) { 155895a0b26dSToby Isaac ierr = PetscSectionGetFieldDof(refConSec,p,f,&cDof);CHKERRQ(ierr); 155995a0b26dSToby Isaac ierr = PetscSectionGetFieldOffset(refConSec,p,f,&cOff);CHKERRQ(ierr); 156095a0b26dSToby Isaac } 156195a0b26dSToby Isaac else { 156295a0b26dSToby Isaac ierr = PetscSectionGetDof(refConSec,p,&cDof);CHKERRQ(ierr); 156395a0b26dSToby Isaac ierr = PetscSectionGetOffset(refConSec,p,&cOff);CHKERRQ(ierr); 156495a0b26dSToby Isaac } 156595a0b26dSToby Isaac 156695a0b26dSToby Isaac if (!cDof) continue; 156795a0b26dSToby Isaac for (r = 0; r < cDof; r++) { 156895a0b26dSToby Isaac rows[r] = cOff + r; 156995a0b26dSToby Isaac } 157095a0b26dSToby Isaac numCols = 0; 157195a0b26dSToby Isaac for (i = 0; i < closureSize; i++) { 157295a0b26dSToby Isaac PetscInt q = closure[2*i]; 157395a0b26dSToby Isaac PetscInt o = closure[2*i+1]; 157495a0b26dSToby Isaac PetscInt aDof, aOff, j; 157595a0b26dSToby Isaac 157695a0b26dSToby Isaac if (numFields > 1) { 157795a0b26dSToby Isaac ierr = PetscSectionGetFieldDof(refSection,q,f,&aDof);CHKERRQ(ierr); 157895a0b26dSToby Isaac ierr = PetscSectionGetFieldOffset(refSection,q,f,&aOff);CHKERRQ(ierr); 157995a0b26dSToby Isaac } 158095a0b26dSToby Isaac else { 158195a0b26dSToby Isaac ierr = PetscSectionGetDof(refSection,q,&aDof);CHKERRQ(ierr); 158295a0b26dSToby Isaac ierr = PetscSectionGetOffset(refSection,q,&aOff);CHKERRQ(ierr); 158395a0b26dSToby Isaac } 158495a0b26dSToby Isaac 158595a0b26dSToby Isaac for (j = 0; j < aDof; j++) { 158695a0b26dSToby Isaac PetscInt node = (o >= 0) ? (j / fComp) : ((aDof - 1 - j) / fComp); 158795a0b26dSToby Isaac PetscInt comp = (j % fComp); 158895a0b26dSToby Isaac 158995a0b26dSToby Isaac cols[numCols++] = aOff + node * fComp + comp; 159095a0b26dSToby Isaac } 159195a0b26dSToby Isaac } 159295a0b26dSToby Isaac refPointFieldN[p-pRefStart][f] = numCols; 159395a0b26dSToby Isaac ierr = PetscMalloc1(cDof*numCols,&refPointFieldMats[p-pRefStart][f]);CHKERRQ(ierr); 159495a0b26dSToby Isaac ierr = MatGetValues(refCmat,cDof,rows,numCols,cols,refPointFieldMats[p-pRefStart][f]);CHKERRQ(ierr); 159595a0b26dSToby Isaac } 159695a0b26dSToby Isaac ierr = DMPlexRestoreTransitiveClosure(refTree,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 159795a0b26dSToby Isaac } 159895a0b26dSToby Isaac 159995a0b26dSToby Isaac /* step 2: compute the preorder */ 160095a0b26dSToby Isaac ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 160195a0b26dSToby Isaac ierr = PetscMalloc2(pEnd-pStart,&perm,pEnd-pStart,&iperm);CHKERRQ(ierr); 160295a0b26dSToby Isaac for (p = pStart; p < pEnd; p++) { 160395a0b26dSToby Isaac perm[p - pStart] = p; 160495a0b26dSToby Isaac iperm[p - pStart] = p-pStart; 160595a0b26dSToby Isaac } 160695a0b26dSToby Isaac for (p = 0; p < pEnd - pStart;) { 160795a0b26dSToby Isaac PetscInt point = perm[p]; 160895a0b26dSToby Isaac PetscInt parent; 160995a0b26dSToby Isaac 161095a0b26dSToby Isaac ierr = DMPlexGetTreeParent(dm,point,&parent,NULL);CHKERRQ(ierr); 161195a0b26dSToby Isaac if (parent == point) { 161295a0b26dSToby Isaac p++; 161395a0b26dSToby Isaac } 161495a0b26dSToby Isaac else { 161595a0b26dSToby Isaac PetscInt size, closureSize, *closure = NULL, i; 161695a0b26dSToby Isaac 161795a0b26dSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 161895a0b26dSToby Isaac for (i = 0; i < closureSize; i++) { 161995a0b26dSToby Isaac PetscInt q = closure[2*i]; 162095a0b26dSToby Isaac if (iperm[q-pStart] > iperm[point-pStart]) { 162195a0b26dSToby Isaac /* swap */ 162295a0b26dSToby Isaac perm[p] = q; 162395a0b26dSToby Isaac perm[iperm[q-pStart]] = point; 162495a0b26dSToby Isaac iperm[point-pStart] = iperm[q-pStart]; 162595a0b26dSToby Isaac iperm[q-pStart] = p; 162695a0b26dSToby Isaac break; 162795a0b26dSToby Isaac } 162895a0b26dSToby Isaac } 162995a0b26dSToby Isaac size = closureSize; 163095a0b26dSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 163195a0b26dSToby Isaac if (i == size) { 163295a0b26dSToby Isaac p++; 163395a0b26dSToby Isaac } 163495a0b26dSToby Isaac } 163595a0b26dSToby Isaac } 163695a0b26dSToby Isaac 163795a0b26dSToby Isaac /* step 3: fill the constraint matrix */ 163895a0b26dSToby Isaac /* we are going to use a preorder progressive fill strategy. Mat doesn't 163995a0b26dSToby Isaac * allow progressive fill without assembly, so we are going to set up the 164095a0b26dSToby Isaac * values outside of the Mat first. 164195a0b26dSToby Isaac */ 164295a0b26dSToby Isaac { 164395a0b26dSToby Isaac PetscInt nRows, row, nnz; 164495a0b26dSToby Isaac PetscBool done; 164595a0b26dSToby Isaac const PetscInt *ia, *ja; 164695a0b26dSToby Isaac PetscScalar *vals; 164795a0b26dSToby Isaac 164895a0b26dSToby Isaac ierr = MatGetRowIJ(cMat,0,PETSC_FALSE,PETSC_FALSE,&nRows,&ia,&ja,&done);CHKERRQ(ierr); 164995a0b26dSToby Isaac if (!done) SETERRQ(PetscObjectComm((PetscObject)cMat),PETSC_ERR_PLIB,"Could not get RowIJ of constraint matrix"); 165095a0b26dSToby Isaac nnz = ia[nRows]; 165195a0b26dSToby Isaac /* malloc and then zero rows right before we fill them: this way valgrind 165295a0b26dSToby Isaac * can tell if we are doing progressive fill in the wrong order */ 165395a0b26dSToby Isaac ierr = PetscMalloc1(nnz,&vals);CHKERRQ(ierr); 165495a0b26dSToby Isaac for (p = 0; p < pEnd - pStart; p++) { 165595a0b26dSToby Isaac PetscInt parent, childid, closureSize, *closure = NULL; 165695a0b26dSToby Isaac PetscInt point = perm[p], pointDof; 165795a0b26dSToby Isaac 165895a0b26dSToby Isaac ierr = DMPlexGetTreeParent(dm,point,&parent,&childid);CHKERRQ(ierr); 165995a0b26dSToby Isaac if ((point < conStart) || (point >= conEnd) || (parent == point)) continue; 166095a0b26dSToby Isaac ierr = PetscSectionGetDof(conSec,point,&pointDof);CHKERRQ(ierr); 166195a0b26dSToby Isaac if (!pointDof) continue; 166295a0b26dSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 166395a0b26dSToby Isaac for (f = 0; f < numFields; f++) { 166495a0b26dSToby Isaac PetscInt cDof, cOff, numCols, numFillCols, i, r, fComp, matOffset, offset; 166595a0b26dSToby Isaac PetscScalar *pointMat; 16660dd1b1feSToby Isaac PetscObject disc; 16670dd1b1feSToby Isaac PetscClassId id; 16680dd1b1feSToby Isaac PetscFE fe = NULL; 16690dd1b1feSToby Isaac PetscFV fv = NULL; 167095a0b26dSToby Isaac 16710dd1b1feSToby Isaac ierr = PetscDSGetDiscretization(ds,f,&disc);CHKERRQ(ierr); 16720dd1b1feSToby Isaac ierr = PetscObjectGetClassId(disc,&id);CHKERRQ(ierr); 16730dd1b1feSToby Isaac if (id == PETSCFE_CLASSID) { 16740dd1b1feSToby Isaac fe = (PetscFE) disc; 167595a0b26dSToby Isaac ierr = PetscFEGetNumComponents(fe,&fComp);CHKERRQ(ierr); 16760dd1b1feSToby Isaac } 16770dd1b1feSToby Isaac else if (id == PETSCFV_CLASSID) { 16780dd1b1feSToby Isaac fv = (PetscFV) disc; 16790dd1b1feSToby Isaac ierr = PetscFVGetNumComponents(fv,&fComp);CHKERRQ(ierr); 16800dd1b1feSToby Isaac } 168195a0b26dSToby Isaac 168295a0b26dSToby Isaac if (numFields > 1) { 168395a0b26dSToby Isaac ierr = PetscSectionGetFieldDof(conSec,point,f,&cDof);CHKERRQ(ierr); 168495a0b26dSToby Isaac ierr = PetscSectionGetFieldOffset(conSec,point,f,&cOff);CHKERRQ(ierr); 168595a0b26dSToby Isaac } 168695a0b26dSToby Isaac else { 168795a0b26dSToby Isaac ierr = PetscSectionGetDof(conSec,point,&cDof);CHKERRQ(ierr); 168895a0b26dSToby Isaac ierr = PetscSectionGetOffset(conSec,point,&cOff);CHKERRQ(ierr); 168995a0b26dSToby Isaac } 169095a0b26dSToby Isaac if (!cDof) continue; 169195a0b26dSToby Isaac 169295a0b26dSToby Isaac /* make sure that every row for this point is the same size */ 169395a0b26dSToby Isaac #if defined(PETSC_USE_DEBUG) 169495a0b26dSToby Isaac for (r = 0; r < cDof; r++) { 169595a0b26dSToby Isaac if (cDof > 1 && r) { 169695a0b26dSToby Isaac if ((ia[cOff+r+1]-ia[cOff+r]) != (ia[cOff+r]-ia[cOff+r-1])) { 169795a0b26dSToby 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])); 169895a0b26dSToby Isaac } 169995a0b26dSToby Isaac } 170095a0b26dSToby Isaac } 170195a0b26dSToby Isaac #endif 170295a0b26dSToby Isaac /* zero rows */ 170395a0b26dSToby Isaac for (i = ia[cOff] ; i< ia[cOff+cDof];i++) { 170495a0b26dSToby Isaac vals[i] = 0.; 170595a0b26dSToby Isaac } 170695a0b26dSToby Isaac matOffset = ia[cOff]; 170795a0b26dSToby Isaac numFillCols = ia[cOff+1] - matOffset; 170895a0b26dSToby Isaac pointMat = refPointFieldMats[childid-pRefStart][f]; 170995a0b26dSToby Isaac numCols = refPointFieldN[childid-pRefStart][f]; 171095a0b26dSToby Isaac offset = 0; 171195a0b26dSToby Isaac for (i = 0; i < closureSize; i++) { 171295a0b26dSToby Isaac PetscInt q = closure[2*i]; 171395a0b26dSToby Isaac PetscInt o = closure[2*i+1]; 171495a0b26dSToby Isaac PetscInt aDof, aOff, j, k, qConDof, qConOff; 171595a0b26dSToby Isaac 171695a0b26dSToby Isaac qConDof = qConOff = 0; 171795a0b26dSToby Isaac if (numFields > 1) { 171895a0b26dSToby Isaac ierr = PetscSectionGetFieldDof(section,q,f,&aDof);CHKERRQ(ierr); 171995a0b26dSToby Isaac ierr = PetscSectionGetFieldOffset(section,q,f,&aOff);CHKERRQ(ierr); 172095a0b26dSToby Isaac if (q >= conStart && q < conEnd) { 172195a0b26dSToby Isaac ierr = PetscSectionGetFieldDof(conSec,q,f,&qConDof);CHKERRQ(ierr); 172295a0b26dSToby Isaac ierr = PetscSectionGetFieldOffset(conSec,q,f,&qConOff);CHKERRQ(ierr); 172395a0b26dSToby Isaac } 172495a0b26dSToby Isaac } 172595a0b26dSToby Isaac else { 172695a0b26dSToby Isaac ierr = PetscSectionGetDof(section,q,&aDof);CHKERRQ(ierr); 172795a0b26dSToby Isaac ierr = PetscSectionGetOffset(section,q,&aOff);CHKERRQ(ierr); 172895a0b26dSToby Isaac if (q >= conStart && q < conEnd) { 172995a0b26dSToby Isaac ierr = PetscSectionGetDof(conSec,q,&qConDof);CHKERRQ(ierr); 173095a0b26dSToby Isaac ierr = PetscSectionGetOffset(conSec,q,&qConOff);CHKERRQ(ierr); 173195a0b26dSToby Isaac } 173295a0b26dSToby Isaac } 173395a0b26dSToby Isaac if (!aDof) continue; 173495a0b26dSToby Isaac if (qConDof) { 173595a0b26dSToby Isaac /* this point has anchors: its rows of the matrix should already 173695a0b26dSToby Isaac * be filled, thanks to preordering */ 173795a0b26dSToby Isaac /* first multiply into pointWork, then set in matrix */ 173895a0b26dSToby Isaac PetscInt aMatOffset = ia[qConOff]; 173995a0b26dSToby Isaac PetscInt aNumFillCols = ia[qConOff + 1] - aMatOffset; 174095a0b26dSToby Isaac for (r = 0; r < cDof; r++) { 174195a0b26dSToby Isaac for (j = 0; j < aNumFillCols; j++) { 174295a0b26dSToby Isaac PetscScalar inVal = 0; 174395a0b26dSToby Isaac for (k = 0; k < aDof; k++) { 174495a0b26dSToby Isaac PetscInt node = (o >= 0) ? (k / fComp) : ((aDof - 1 - k) / fComp); 174595a0b26dSToby Isaac PetscInt comp = (k % fComp); 174695a0b26dSToby Isaac PetscInt col = node * fComp + comp; 174795a0b26dSToby Isaac 174895a0b26dSToby Isaac inVal += pointMat[r * numCols + offset + col] * vals[aMatOffset + aNumFillCols * k + j]; 174995a0b26dSToby Isaac } 175095a0b26dSToby Isaac pointWork[r * aNumFillCols + j] = inVal; 175195a0b26dSToby Isaac } 175295a0b26dSToby Isaac } 175395a0b26dSToby Isaac /* assume that the columns are sorted, spend less time searching */ 175495a0b26dSToby Isaac for (j = 0, k = 0; j < aNumFillCols; j++) { 175595a0b26dSToby Isaac PetscInt col = ja[aMatOffset + j]; 175695a0b26dSToby Isaac for (;k < numFillCols; k++) { 175795a0b26dSToby Isaac if (ja[matOffset + k] == col) { 175895a0b26dSToby Isaac break; 175995a0b26dSToby Isaac } 176095a0b26dSToby Isaac } 176195a0b26dSToby Isaac if (k == numFillCols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"No nonzero space for (%d, %d)", cOff, col); 176295a0b26dSToby Isaac for (r = 0; r < cDof; r++) { 176395a0b26dSToby Isaac vals[matOffset + numFillCols * r + k] = pointWork[r * aNumFillCols + j]; 176495a0b26dSToby Isaac } 176595a0b26dSToby Isaac } 176695a0b26dSToby Isaac } 176795a0b26dSToby Isaac else { 176895a0b26dSToby Isaac /* find where to put this portion of pointMat into the matrix */ 176995a0b26dSToby Isaac for (k = 0; k < numFillCols; k++) { 177095a0b26dSToby Isaac if (ja[matOffset + k] == aOff) { 177195a0b26dSToby Isaac break; 177295a0b26dSToby Isaac } 177395a0b26dSToby Isaac } 177495a0b26dSToby Isaac if (k == numFillCols) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"No nonzero space for (%d, %d)", cOff, aOff); 177595a0b26dSToby Isaac for (j = 0; j < aDof; j++) { 177695a0b26dSToby Isaac PetscInt node = (o >= 0) ? (j / fComp) : ((aDof - 1 - j) / fComp); 177795a0b26dSToby Isaac PetscInt comp = (j % fComp); 177895a0b26dSToby Isaac PetscInt col = node * fComp + comp; 177995a0b26dSToby Isaac for (r = 0; r < cDof; r++) { 178095a0b26dSToby Isaac vals[matOffset + numFillCols * r + k + col] += pointMat[r * numCols + offset + col]; 178195a0b26dSToby Isaac } 178295a0b26dSToby Isaac } 178395a0b26dSToby Isaac } 178495a0b26dSToby Isaac offset += aDof; 178595a0b26dSToby Isaac } 178695a0b26dSToby Isaac } 178795a0b26dSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,parent,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 178895a0b26dSToby Isaac } 178995a0b26dSToby Isaac for (row = 0; row < nRows; row++) { 179095a0b26dSToby Isaac ierr = MatSetValues(cMat,1,&row,ia[row+1]-ia[row],&ja[ia[row]],&vals[ia[row]],INSERT_VALUES);CHKERRQ(ierr); 179195a0b26dSToby Isaac } 179295a0b26dSToby Isaac ierr = MatRestoreRowIJ(cMat,0,PETSC_FALSE,PETSC_FALSE,&nRows,&ia,&ja,&done);CHKERRQ(ierr); 179395a0b26dSToby Isaac if (!done) SETERRQ(PetscObjectComm((PetscObject)cMat),PETSC_ERR_PLIB,"Could not restore RowIJ of constraint matrix"); 179495a0b26dSToby Isaac ierr = MatAssemblyBegin(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 179595a0b26dSToby Isaac ierr = MatAssemblyEnd(cMat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 179695a0b26dSToby Isaac ierr = PetscFree(vals);CHKERRQ(ierr); 179795a0b26dSToby Isaac } 179895a0b26dSToby Isaac 179995a0b26dSToby Isaac /* clean up */ 180095a0b26dSToby Isaac ierr = ISRestoreIndices(refAnIS,&refAnchors);CHKERRQ(ierr); 180195a0b26dSToby Isaac ierr = ISRestoreIndices(anIS,&anchors);CHKERRQ(ierr); 180295a0b26dSToby Isaac ierr = PetscFree2(perm,iperm);CHKERRQ(ierr); 180395a0b26dSToby Isaac ierr = PetscFree(rows);CHKERRQ(ierr); 180495a0b26dSToby Isaac ierr = PetscFree(cols);CHKERRQ(ierr); 180595a0b26dSToby Isaac ierr = PetscFree(pointWork);CHKERRQ(ierr); 180695a0b26dSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 180795a0b26dSToby Isaac PetscInt parent, pDof; 180895a0b26dSToby Isaac 180995a0b26dSToby Isaac ierr = DMPlexGetTreeParent(refTree,p,&parent,NULL);CHKERRQ(ierr); 181095a0b26dSToby Isaac ierr = PetscSectionGetDof(refConSec,p,&pDof);CHKERRQ(ierr); 181195a0b26dSToby Isaac if (!pDof || parent == p) continue; 181295a0b26dSToby Isaac 181395a0b26dSToby Isaac for (f = 0; f < numFields; f++) { 181495a0b26dSToby Isaac PetscInt cDof; 181595a0b26dSToby Isaac 181695a0b26dSToby Isaac if (numFields > 1) { 181795a0b26dSToby Isaac ierr = PetscSectionGetFieldDof(refConSec,p,f,&cDof);CHKERRQ(ierr); 181895a0b26dSToby Isaac } 181995a0b26dSToby Isaac else { 182095a0b26dSToby Isaac ierr = PetscSectionGetDof(refConSec,p,&cDof);CHKERRQ(ierr); 182195a0b26dSToby Isaac } 182295a0b26dSToby Isaac 182395a0b26dSToby Isaac if (!cDof) continue; 182495a0b26dSToby Isaac ierr = PetscFree(refPointFieldMats[p - pRefStart][f]);CHKERRQ(ierr); 182595a0b26dSToby Isaac } 182695a0b26dSToby Isaac ierr = PetscFree(refPointFieldMats[p - pRefStart]);CHKERRQ(ierr); 182795a0b26dSToby Isaac ierr = PetscFree(refPointFieldN[p - pRefStart]);CHKERRQ(ierr); 182895a0b26dSToby Isaac } 182995a0b26dSToby Isaac ierr = PetscFree(refPointFieldMats);CHKERRQ(ierr); 183095a0b26dSToby Isaac ierr = PetscFree(refPointFieldN);CHKERRQ(ierr); 183195a0b26dSToby Isaac PetscFunctionReturn(0); 183295a0b26dSToby Isaac } 183395a0b26dSToby Isaac 18346f5f1567SToby Isaac #undef __FUNCT__ 18356f5f1567SToby Isaac #define __FUNCT__ "DMPlexTreeRefineCell" 18366f5f1567SToby Isaac /* refine a single cell on rank 0: this is not intended to provide good local refinement, only to create an example of 18376f5f1567SToby Isaac * a non-conforming mesh. Local refinement comes later */ 18386f5f1567SToby Isaac PetscErrorCode DMPlexTreeRefineCell (DM dm, PetscInt cell, DM *ncdm) 18396f5f1567SToby Isaac { 18406f5f1567SToby Isaac DM K; 1841420f55faSMatthew G. Knepley PetscMPIInt rank; 18426f5f1567SToby Isaac PetscInt dim, *pNewStart, *pNewEnd, *pNewCount, *pOldStart, *pOldEnd, offset, d, pStart, pEnd; 18436f5f1567SToby Isaac PetscInt numNewCones, *newConeSizes, *newCones, *newOrientations; 18446f5f1567SToby Isaac PetscInt *Kembedding; 18456f5f1567SToby Isaac PetscInt *cellClosure=NULL, nc; 18466f5f1567SToby Isaac PetscScalar *newVertexCoords; 18476f5f1567SToby Isaac PetscInt numPointsWithParents, *parents, *childIDs, *perm, *iperm, *preOrient, pOffset; 18486f5f1567SToby Isaac PetscSection parentSection; 18496f5f1567SToby Isaac PetscErrorCode ierr; 18506f5f1567SToby Isaac 18516f5f1567SToby Isaac PetscFunctionBegin; 18526f5f1567SToby Isaac ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 185328f4b327SMatthew G. Knepley ierr = DMGetDimension(dm,&dim);CHKERRQ(ierr); 18546f5f1567SToby Isaac ierr = DMPlexCreate(PetscObjectComm((PetscObject)dm), ncdm);CHKERRQ(ierr); 185528f4b327SMatthew G. Knepley ierr = DMSetDimension(*ncdm,dim);CHKERRQ(ierr); 18566f5f1567SToby Isaac 18576f5f1567SToby Isaac ierr = DMPlexGetChart(dm, &pStart, &pEnd);CHKERRQ(ierr); 18586f5f1567SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&parentSection);CHKERRQ(ierr); 18596f5f1567SToby Isaac ierr = DMPlexGetReferenceTree(dm,&K);CHKERRQ(ierr); 18606f5f1567SToby Isaac if (!rank) { 18616f5f1567SToby Isaac /* compute the new charts */ 18626f5f1567SToby Isaac ierr = PetscMalloc5(dim+1,&pNewCount,dim+1,&pNewStart,dim+1,&pNewEnd,dim+1,&pOldStart,dim+1,&pOldEnd);CHKERRQ(ierr); 18636f5f1567SToby Isaac offset = 0; 18646f5f1567SToby Isaac for (d = 0; d <= dim; d++) { 18656f5f1567SToby Isaac PetscInt pOldCount, kStart, kEnd, k; 18666f5f1567SToby Isaac 18676f5f1567SToby Isaac pNewStart[d] = offset; 18686f5f1567SToby Isaac ierr = DMPlexGetHeightStratum(dm,d,&pOldStart[d],&pOldEnd[d]);CHKERRQ(ierr); 18696f5f1567SToby Isaac ierr = DMPlexGetHeightStratum(K,d,&kStart,&kEnd);CHKERRQ(ierr); 18706f5f1567SToby Isaac pOldCount = pOldEnd[d] - pOldStart[d]; 18716f5f1567SToby Isaac /* adding the new points */ 18726f5f1567SToby Isaac pNewCount[d] = pOldCount + kEnd - kStart; 18736f5f1567SToby Isaac if (!d) { 18746f5f1567SToby Isaac /* removing the cell */ 18756f5f1567SToby Isaac pNewCount[d]--; 18766f5f1567SToby Isaac } 18776f5f1567SToby Isaac for (k = kStart; k < kEnd; k++) { 18786f5f1567SToby Isaac PetscInt parent; 18796f5f1567SToby Isaac ierr = DMPlexGetTreeParent(K,k,&parent,NULL);CHKERRQ(ierr); 18806f5f1567SToby Isaac if (parent == k) { 18816f5f1567SToby Isaac /* avoid double counting points that won't actually be new */ 18826f5f1567SToby Isaac pNewCount[d]--; 18836f5f1567SToby Isaac } 18846f5f1567SToby Isaac } 18856f5f1567SToby Isaac pNewEnd[d] = pNewStart[d] + pNewCount[d]; 18866f5f1567SToby Isaac offset = pNewEnd[d]; 18876f5f1567SToby Isaac 18886f5f1567SToby Isaac } 18896f5f1567SToby 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]); 18906f5f1567SToby Isaac /* get the current closure of the cell that we are removing */ 18916f5f1567SToby Isaac ierr = DMPlexGetTransitiveClosure(dm,cell,PETSC_TRUE,&nc,&cellClosure);CHKERRQ(ierr); 18926f5f1567SToby Isaac 18936f5f1567SToby Isaac ierr = PetscMalloc1(pNewEnd[dim],&newConeSizes);CHKERRQ(ierr); 18946f5f1567SToby Isaac { 18956f5f1567SToby Isaac PetscInt kStart, kEnd, k, closureSizeK, *closureK = NULL, j; 18966f5f1567SToby Isaac 18976f5f1567SToby Isaac ierr = DMPlexGetChart(K,&kStart,&kEnd);CHKERRQ(ierr); 18986f5f1567SToby Isaac ierr = PetscMalloc4(kEnd-kStart,&Kembedding,kEnd-kStart,&perm,kEnd-kStart,&iperm,kEnd-kStart,&preOrient);CHKERRQ(ierr); 18996f5f1567SToby Isaac 19006f5f1567SToby Isaac for (k = kStart; k < kEnd; k++) { 19016f5f1567SToby Isaac perm[k - kStart] = k; 19026f5f1567SToby Isaac iperm [k - kStart] = k - kStart; 19036f5f1567SToby Isaac preOrient[k - kStart] = 0; 19046f5f1567SToby Isaac } 19056f5f1567SToby Isaac 19066f5f1567SToby Isaac ierr = DMPlexGetTransitiveClosure(K,0,PETSC_TRUE,&closureSizeK,&closureK);CHKERRQ(ierr); 19076f5f1567SToby Isaac for (j = 1; j < closureSizeK; j++) { 19086f5f1567SToby Isaac PetscInt parentOrientA = closureK[2*j+1]; 19096f5f1567SToby Isaac PetscInt parentOrientB = cellClosure[2*j+1]; 19106f5f1567SToby Isaac PetscInt p, q; 19116f5f1567SToby Isaac 19126f5f1567SToby Isaac p = closureK[2*j]; 19136f5f1567SToby Isaac q = cellClosure[2*j]; 19146f5f1567SToby Isaac for (d = 0; d <= dim; d++) { 19156f5f1567SToby Isaac if (q >= pOldStart[d] && q < pOldEnd[d]) { 19166f5f1567SToby Isaac Kembedding[p] = (q - pOldStart[d]) + pNewStart[d]; 19176f5f1567SToby Isaac } 19186f5f1567SToby Isaac } 19196f5f1567SToby Isaac if (parentOrientA != parentOrientB) { 19206f5f1567SToby Isaac PetscInt numChildren, i; 19216f5f1567SToby Isaac const PetscInt *children; 19226f5f1567SToby Isaac 19236f5f1567SToby Isaac ierr = DMPlexGetTreeChildren(K,p,&numChildren,&children);CHKERRQ(ierr); 19246f5f1567SToby Isaac for (i = 0; i < numChildren; i++) { 19256f5f1567SToby Isaac PetscInt kPerm, oPerm; 19266f5f1567SToby Isaac 19276f5f1567SToby Isaac k = children[i]; 19286f5f1567SToby Isaac ierr = DMPlexReferenceTreeGetChildSymmetry(K,p,parentOrientA,0,k,parentOrientB,&oPerm,&kPerm);CHKERRQ(ierr); 19296f5f1567SToby Isaac /* perm = what refTree position I'm in */ 19306f5f1567SToby Isaac perm[kPerm-kStart] = k; 19316f5f1567SToby Isaac /* iperm = who is at this position */ 19326f5f1567SToby Isaac iperm[k-kStart] = kPerm-kStart; 19336f5f1567SToby Isaac preOrient[kPerm-kStart] = oPerm; 19346f5f1567SToby Isaac } 19356f5f1567SToby Isaac } 19366f5f1567SToby Isaac } 19376f5f1567SToby Isaac ierr = DMPlexRestoreTransitiveClosure(K,0,PETSC_TRUE,&closureSizeK,&closureK);CHKERRQ(ierr); 19386f5f1567SToby Isaac } 19396f5f1567SToby Isaac ierr = PetscSectionSetChart(parentSection,0,pNewEnd[dim]);CHKERRQ(ierr); 19406f5f1567SToby Isaac offset = 0; 19416f5f1567SToby Isaac numNewCones = 0; 19426f5f1567SToby Isaac for (d = 0; d <= dim; d++) { 19436f5f1567SToby Isaac PetscInt kStart, kEnd, k; 19446f5f1567SToby Isaac PetscInt p; 19456f5f1567SToby Isaac PetscInt size; 19466f5f1567SToby Isaac 19476f5f1567SToby Isaac for (p = pOldStart[d]; p < pOldEnd[d]; p++) { 19486f5f1567SToby Isaac /* skip cell 0 */ 19496f5f1567SToby Isaac if (p == cell) continue; 19506f5f1567SToby Isaac /* old cones to new cones */ 19516f5f1567SToby Isaac ierr = DMPlexGetConeSize(dm,p,&size);CHKERRQ(ierr); 19526f5f1567SToby Isaac newConeSizes[offset++] = size; 19536f5f1567SToby Isaac numNewCones += size; 19546f5f1567SToby Isaac } 19556f5f1567SToby Isaac 19566f5f1567SToby Isaac ierr = DMPlexGetHeightStratum(K,d,&kStart,&kEnd);CHKERRQ(ierr); 19576f5f1567SToby Isaac for (k = kStart; k < kEnd; k++) { 19586f5f1567SToby Isaac PetscInt kParent; 19596f5f1567SToby Isaac 19606f5f1567SToby Isaac ierr = DMPlexGetTreeParent(K,k,&kParent,NULL);CHKERRQ(ierr); 19616f5f1567SToby Isaac if (kParent != k) { 19626f5f1567SToby Isaac Kembedding[k] = offset; 19636f5f1567SToby Isaac ierr = DMPlexGetConeSize(K,k,&size);CHKERRQ(ierr); 19646f5f1567SToby Isaac newConeSizes[offset++] = size; 19656f5f1567SToby Isaac numNewCones += size; 19666f5f1567SToby Isaac if (kParent != 0) { 19676f5f1567SToby Isaac ierr = PetscSectionSetDof(parentSection,Kembedding[k],1);CHKERRQ(ierr); 19686f5f1567SToby Isaac } 19696f5f1567SToby Isaac } 19706f5f1567SToby Isaac } 19716f5f1567SToby Isaac } 19726f5f1567SToby Isaac 19736f5f1567SToby Isaac ierr = PetscSectionSetUp(parentSection);CHKERRQ(ierr); 19746f5f1567SToby Isaac ierr = PetscSectionGetStorageSize(parentSection,&numPointsWithParents);CHKERRQ(ierr); 19756f5f1567SToby Isaac ierr = PetscMalloc2(numNewCones,&newCones,numNewCones,&newOrientations);CHKERRQ(ierr); 19766f5f1567SToby Isaac ierr = PetscMalloc2(numPointsWithParents,&parents,numPointsWithParents,&childIDs);CHKERRQ(ierr); 19776f5f1567SToby Isaac 19786f5f1567SToby Isaac /* fill new cones */ 19796f5f1567SToby Isaac offset = 0; 19806f5f1567SToby Isaac for (d = 0; d <= dim; d++) { 19816f5f1567SToby Isaac PetscInt kStart, kEnd, k, l; 19826f5f1567SToby Isaac PetscInt p; 19836f5f1567SToby Isaac PetscInt size; 19846f5f1567SToby Isaac const PetscInt *cone, *orientation; 19856f5f1567SToby Isaac 19866f5f1567SToby Isaac for (p = pOldStart[d]; p < pOldEnd[d]; p++) { 19876f5f1567SToby Isaac /* skip cell 0 */ 19886f5f1567SToby Isaac if (p == cell) continue; 19896f5f1567SToby Isaac /* old cones to new cones */ 19906f5f1567SToby Isaac ierr = DMPlexGetConeSize(dm,p,&size);CHKERRQ(ierr); 19916f5f1567SToby Isaac ierr = DMPlexGetCone(dm,p,&cone);CHKERRQ(ierr); 19926f5f1567SToby Isaac ierr = DMPlexGetConeOrientation(dm,p,&orientation);CHKERRQ(ierr); 19936f5f1567SToby Isaac for (l = 0; l < size; l++) { 19946f5f1567SToby Isaac newCones[offset] = (cone[l] - pOldStart[d + 1]) + pNewStart[d + 1]; 19956f5f1567SToby Isaac newOrientations[offset++] = orientation[l]; 19966f5f1567SToby Isaac } 19976f5f1567SToby Isaac } 19986f5f1567SToby Isaac 19996f5f1567SToby Isaac ierr = DMPlexGetHeightStratum(K,d,&kStart,&kEnd);CHKERRQ(ierr); 20006f5f1567SToby Isaac for (k = kStart; k < kEnd; k++) { 20016f5f1567SToby Isaac PetscInt kPerm = perm[k], kParent; 20026f5f1567SToby Isaac PetscInt preO = preOrient[k]; 20036f5f1567SToby Isaac 20046f5f1567SToby Isaac ierr = DMPlexGetTreeParent(K,k,&kParent,NULL);CHKERRQ(ierr); 20056f5f1567SToby Isaac if (kParent != k) { 20066f5f1567SToby Isaac /* embed new cones */ 20076f5f1567SToby Isaac ierr = DMPlexGetConeSize(K,k,&size);CHKERRQ(ierr); 20086f5f1567SToby Isaac ierr = DMPlexGetCone(K,kPerm,&cone);CHKERRQ(ierr); 20096f5f1567SToby Isaac ierr = DMPlexGetConeOrientation(K,kPerm,&orientation);CHKERRQ(ierr); 20106f5f1567SToby Isaac for (l = 0; l < size; l++) { 20116f5f1567SToby Isaac PetscInt q, m = (preO >= 0) ? ((preO + l) % size) : ((size -(preO + 1) - l) % size); 20126f5f1567SToby Isaac PetscInt newO, lSize, oTrue; 20136f5f1567SToby Isaac 20146f5f1567SToby Isaac q = iperm[cone[m]]; 20156f5f1567SToby Isaac newCones[offset] = Kembedding[q]; 20166f5f1567SToby Isaac ierr = DMPlexGetConeSize(K,q,&lSize);CHKERRQ(ierr); 20176f5f1567SToby Isaac oTrue = orientation[m]; 20186f5f1567SToby Isaac oTrue = ((!lSize) || (preOrient[k] >= 0)) ? oTrue : -(oTrue + 2); 20196f5f1567SToby Isaac newO = DihedralCompose(lSize,oTrue,preOrient[q]); 20206f5f1567SToby Isaac newOrientations[offset++] = newO; 20216f5f1567SToby Isaac } 20226f5f1567SToby Isaac if (kParent != 0) { 20236f5f1567SToby Isaac PetscInt newPoint = Kembedding[kParent]; 20246f5f1567SToby Isaac ierr = PetscSectionGetOffset(parentSection,Kembedding[k],&pOffset);CHKERRQ(ierr); 20256f5f1567SToby Isaac parents[pOffset] = newPoint; 20266f5f1567SToby Isaac childIDs[pOffset] = k; 20276f5f1567SToby Isaac } 20286f5f1567SToby Isaac } 20296f5f1567SToby Isaac } 20306f5f1567SToby Isaac } 20316f5f1567SToby Isaac 20326f5f1567SToby Isaac ierr = PetscMalloc1(dim*(pNewEnd[dim]-pNewStart[dim]),&newVertexCoords);CHKERRQ(ierr); 20336f5f1567SToby Isaac 20346f5f1567SToby Isaac /* fill coordinates */ 20356f5f1567SToby Isaac offset = 0; 20366f5f1567SToby Isaac { 2037d90620a3SMatthew G. Knepley PetscInt kStart, kEnd, l; 20386f5f1567SToby Isaac PetscSection vSection; 20396f5f1567SToby Isaac PetscInt v; 20406f5f1567SToby Isaac Vec coords; 20416f5f1567SToby Isaac PetscScalar *coordvals; 20426f5f1567SToby Isaac PetscInt dof, off; 2043c111c6b7SMatthew G. Knepley PetscReal v0[3], J[9], detJ; 20446f5f1567SToby Isaac 20456f5f1567SToby Isaac #if defined(PETSC_USE_DEBUG) 2046d90620a3SMatthew G. Knepley { 2047d90620a3SMatthew G. Knepley PetscInt k; 20486f5f1567SToby Isaac ierr = DMPlexGetHeightStratum(K,0,&kStart,&kEnd);CHKERRQ(ierr); 20496f5f1567SToby Isaac for (k = kStart; k < kEnd; k++) { 205073a7f2aaSMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(K, k, NULL, v0, J, NULL, &detJ);CHKERRQ(ierr); 20516f5f1567SToby Isaac if (detJ <= 0.) SETERRQ1 (PETSC_COMM_SELF,PETSC_ERR_PLIB,"reference tree cell %d has bad determinant",k); 20526f5f1567SToby Isaac } 2053d90620a3SMatthew G. Knepley } 20546f5f1567SToby Isaac #endif 205573a7f2aaSMatthew G. Knepley ierr = DMPlexComputeCellGeometryFEM(dm, cell, NULL, v0, J, NULL, &detJ);CHKERRQ(ierr); 20566f5f1567SToby Isaac ierr = DMGetCoordinateSection(dm,&vSection);CHKERRQ(ierr); 20576f5f1567SToby Isaac ierr = DMGetCoordinatesLocal(dm,&coords);CHKERRQ(ierr); 20586f5f1567SToby Isaac ierr = VecGetArray(coords,&coordvals);CHKERRQ(ierr); 20596f5f1567SToby Isaac for (v = pOldStart[dim]; v < pOldEnd[dim]; v++) { 20606f5f1567SToby Isaac 20616f5f1567SToby Isaac ierr = PetscSectionGetDof(vSection,v,&dof);CHKERRQ(ierr); 20626f5f1567SToby Isaac ierr = PetscSectionGetOffset(vSection,v,&off);CHKERRQ(ierr); 20636f5f1567SToby Isaac for (l = 0; l < dof; l++) { 20646f5f1567SToby Isaac newVertexCoords[offset++] = coordvals[off + l]; 20656f5f1567SToby Isaac } 20666f5f1567SToby Isaac } 20676f5f1567SToby Isaac ierr = VecRestoreArray(coords,&coordvals);CHKERRQ(ierr); 20686f5f1567SToby Isaac 20696f5f1567SToby Isaac ierr = DMGetCoordinateSection(K,&vSection);CHKERRQ(ierr); 20706f5f1567SToby Isaac ierr = DMGetCoordinatesLocal(K,&coords);CHKERRQ(ierr); 20716f5f1567SToby Isaac ierr = VecGetArray(coords,&coordvals);CHKERRQ(ierr); 20726f5f1567SToby Isaac ierr = DMPlexGetDepthStratum(K,0,&kStart,&kEnd);CHKERRQ(ierr); 20736f5f1567SToby Isaac for (v = kStart; v < kEnd; v++) { 20749bc368c7SMatthew G. Knepley PetscReal coord[3], newCoord[3]; 20756f5f1567SToby Isaac PetscInt vPerm = perm[v]; 20766f5f1567SToby Isaac PetscInt kParent; 20776f5f1567SToby Isaac 20786f5f1567SToby Isaac ierr = DMPlexGetTreeParent(K,v,&kParent,NULL);CHKERRQ(ierr); 20796f5f1567SToby Isaac if (kParent != v) { 20806f5f1567SToby Isaac /* this is a new vertex */ 20816f5f1567SToby Isaac ierr = PetscSectionGetOffset(vSection,vPerm,&off);CHKERRQ(ierr); 20829bc368c7SMatthew G. Knepley for (l = 0; l < dim; ++l) coord[l] = PetscRealPart(coordvals[off+l]); 20839bc368c7SMatthew G. Knepley CoordinatesRefToReal(dim, dim, v0, J, coord, newCoord);CHKERRQ(ierr); 20849bc368c7SMatthew G. Knepley for (l = 0; l < dim; ++l) newVertexCoords[offset+l] = newCoord[l]; 20856f5f1567SToby Isaac offset += dim; 20866f5f1567SToby Isaac } 20876f5f1567SToby Isaac } 20886f5f1567SToby Isaac ierr = VecRestoreArray(coords,&coordvals);CHKERRQ(ierr); 20896f5f1567SToby Isaac } 20906f5f1567SToby Isaac 20916f5f1567SToby Isaac /* need to reverse the order of pNewCount: vertices first, cells last */ 20926f5f1567SToby Isaac for (d = 0; d < (dim + 1) / 2; d++) { 20936f5f1567SToby Isaac PetscInt tmp; 20946f5f1567SToby Isaac 20956f5f1567SToby Isaac tmp = pNewCount[d]; 20966f5f1567SToby Isaac pNewCount[d] = pNewCount[dim - d]; 20976f5f1567SToby Isaac pNewCount[dim - d] = tmp; 20986f5f1567SToby Isaac } 20996f5f1567SToby Isaac 21006f5f1567SToby Isaac ierr = DMPlexCreateFromDAG(*ncdm,dim,pNewCount,newConeSizes,newCones,newOrientations,newVertexCoords);CHKERRQ(ierr); 21016f5f1567SToby Isaac ierr = DMPlexSetReferenceTree(*ncdm,K);CHKERRQ(ierr); 21026f5f1567SToby Isaac ierr = DMPlexSetTree(*ncdm,parentSection,parents,childIDs);CHKERRQ(ierr); 21036f5f1567SToby Isaac 21046f5f1567SToby Isaac /* clean up */ 21056f5f1567SToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,cell,PETSC_TRUE,&nc,&cellClosure);CHKERRQ(ierr); 21066f5f1567SToby Isaac ierr = PetscFree5(pNewCount,pNewStart,pNewEnd,pOldStart,pOldEnd);CHKERRQ(ierr); 21076f5f1567SToby Isaac ierr = PetscFree(newConeSizes);CHKERRQ(ierr); 21086f5f1567SToby Isaac ierr = PetscFree2(newCones,newOrientations);CHKERRQ(ierr); 21096f5f1567SToby Isaac ierr = PetscFree(newVertexCoords);CHKERRQ(ierr); 21106f5f1567SToby Isaac ierr = PetscFree2(parents,childIDs);CHKERRQ(ierr); 21116f5f1567SToby Isaac ierr = PetscFree4(Kembedding,perm,iperm,preOrient);CHKERRQ(ierr); 21126f5f1567SToby Isaac } 21136f5f1567SToby Isaac else { 21146f5f1567SToby Isaac PetscInt p, counts[4]; 21156f5f1567SToby Isaac PetscInt *coneSizes, *cones, *orientations; 21166f5f1567SToby Isaac Vec coordVec; 21176f5f1567SToby Isaac PetscScalar *coords; 21186f5f1567SToby Isaac 21196f5f1567SToby Isaac for (d = 0; d <= dim; d++) { 21206f5f1567SToby Isaac PetscInt dStart, dEnd; 21216f5f1567SToby Isaac 21226f5f1567SToby Isaac ierr = DMPlexGetDepthStratum(dm,d,&dStart,&dEnd);CHKERRQ(ierr); 21236f5f1567SToby Isaac counts[d] = dEnd - dStart; 21246f5f1567SToby Isaac } 21256f5f1567SToby Isaac ierr = PetscMalloc1(pEnd-pStart,&coneSizes);CHKERRQ(ierr); 21266f5f1567SToby Isaac for (p = pStart; p < pEnd; p++) { 21276f5f1567SToby Isaac ierr = DMPlexGetConeSize(dm,p,&coneSizes[p-pStart]);CHKERRQ(ierr); 21286f5f1567SToby Isaac } 21296f5f1567SToby Isaac ierr = DMPlexGetCones(dm, &cones);CHKERRQ(ierr); 21306f5f1567SToby Isaac ierr = DMPlexGetConeOrientations(dm, &orientations);CHKERRQ(ierr); 21316f5f1567SToby Isaac ierr = DMGetCoordinatesLocal(dm,&coordVec);CHKERRQ(ierr); 21326f5f1567SToby Isaac ierr = VecGetArray(coordVec,&coords);CHKERRQ(ierr); 21336f5f1567SToby Isaac 21346f5f1567SToby Isaac ierr = PetscSectionSetChart(parentSection,pStart,pEnd);CHKERRQ(ierr); 21356f5f1567SToby Isaac ierr = PetscSectionSetUp(parentSection);CHKERRQ(ierr); 21366f5f1567SToby Isaac ierr = DMPlexCreateFromDAG(*ncdm,dim,counts,coneSizes,cones,orientations,NULL);CHKERRQ(ierr); 21376f5f1567SToby Isaac ierr = DMPlexSetReferenceTree(*ncdm,K);CHKERRQ(ierr); 21386f5f1567SToby Isaac ierr = DMPlexSetTree(*ncdm,parentSection,NULL,NULL);CHKERRQ(ierr); 21396f5f1567SToby Isaac ierr = VecRestoreArray(coordVec,&coords);CHKERRQ(ierr); 21406f5f1567SToby Isaac } 21416f5f1567SToby Isaac ierr = PetscSectionDestroy(&parentSection);CHKERRQ(ierr); 21426f5f1567SToby Isaac 21436f5f1567SToby Isaac PetscFunctionReturn(0); 21446f5f1567SToby Isaac } 2145