1d6a7ad0dSToby Isaac #include <petsc-private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2d6a7ad0dSToby Isaac #include <../src/sys/utils/hash.h> 3d6a7ad0dSToby Isaac #include <petsc-private/isimpl.h> 4d6a7ad0dSToby Isaac #include <petscsf.h> 5d6a7ad0dSToby Isaac 6d6a7ad0dSToby Isaac /** hierarchy routines */ 7d6a7ad0dSToby Isaac 8d6a7ad0dSToby Isaac #undef __FUNCT__ 9d6a7ad0dSToby Isaac #define __FUNCT__ "DMPlexSetReferenceTree" 10d6a7ad0dSToby Isaac /*@ 11d6a7ad0dSToby Isaac DMPlexSetReferenceTree - set the reference tree for hierarchically non-conforming meshes. 12d6a7ad0dSToby Isaac 13d6a7ad0dSToby Isaac Not collective 14d6a7ad0dSToby Isaac 15d6a7ad0dSToby Isaac Input Parameters: 16d6a7ad0dSToby Isaac + dm - The DMPlex object 17d6a7ad0dSToby Isaac - ref - The reference tree DMPlex object 18d6a7ad0dSToby Isaac 190b7167a0SToby Isaac Level: intermediate 20d6a7ad0dSToby Isaac 21da43764aSToby Isaac .seealso: DMPlexGetReferenceTree(), DMPlexCreateDefaultReferenceTree() 22d6a7ad0dSToby Isaac @*/ 23d6a7ad0dSToby Isaac PetscErrorCode DMPlexSetReferenceTree(DM dm, DM ref) 24d6a7ad0dSToby Isaac { 25d6a7ad0dSToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 26d6a7ad0dSToby Isaac PetscErrorCode ierr; 27d6a7ad0dSToby Isaac 28d6a7ad0dSToby Isaac PetscFunctionBegin; 29d6a7ad0dSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 30d6a7ad0dSToby Isaac PetscValidHeaderSpecific(ref, DM_CLASSID, 2); 31d6a7ad0dSToby Isaac ierr = PetscObjectReference((PetscObject)ref);CHKERRQ(ierr); 32d6a7ad0dSToby Isaac ierr = DMDestroy(&mesh->referenceTree);CHKERRQ(ierr); 33d6a7ad0dSToby Isaac mesh->referenceTree = ref; 34d6a7ad0dSToby Isaac PetscFunctionReturn(0); 35d6a7ad0dSToby Isaac } 36d6a7ad0dSToby Isaac 37d6a7ad0dSToby Isaac #undef __FUNCT__ 38d6a7ad0dSToby Isaac #define __FUNCT__ "DMPlexGetReferenceTree" 39d6a7ad0dSToby Isaac /*@ 40d6a7ad0dSToby Isaac DMPlexGetReferenceTree - get the reference tree for hierarchically non-conforming meshes. 41d6a7ad0dSToby Isaac 42d6a7ad0dSToby Isaac Not collective 43d6a7ad0dSToby Isaac 44d6a7ad0dSToby Isaac Input Parameters: 45d6a7ad0dSToby Isaac . dm - The DMPlex object 46d6a7ad0dSToby Isaac 47d6a7ad0dSToby Isaac Output Parameters 48d6a7ad0dSToby Isaac . ref - The reference tree DMPlex object 49d6a7ad0dSToby Isaac 500b7167a0SToby Isaac Level: intermediate 51d6a7ad0dSToby Isaac 52da43764aSToby Isaac .seealso: DMPlexSetReferenceTree(), DMPlexCreateDefaultReferenceTree() 53d6a7ad0dSToby Isaac @*/ 54d6a7ad0dSToby Isaac PetscErrorCode DMPlexGetReferenceTree(DM dm, DM *ref) 55d6a7ad0dSToby Isaac { 56d6a7ad0dSToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 57d6a7ad0dSToby Isaac 58d6a7ad0dSToby Isaac PetscFunctionBegin; 59d6a7ad0dSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 60d6a7ad0dSToby Isaac PetscValidPointer(ref,2); 61d6a7ad0dSToby Isaac *ref = mesh->referenceTree; 62d6a7ad0dSToby Isaac PetscFunctionReturn(0); 63d6a7ad0dSToby Isaac } 64d6a7ad0dSToby Isaac 65da43764aSToby Isaac #undef __FUNCT__ 66da43764aSToby Isaac #define __FUNCT__ "DMPlexCreateDefaultReferenceTree" 67da43764aSToby Isaac /*@ 68da43764aSToby Isaac DMPlexCreateDefaultReferenceTree - create a reference tree for isotropic hierarchical mesh refinement. 69da43764aSToby Isaac 70da43764aSToby Isaac Collective on comm 71da43764aSToby Isaac 72da43764aSToby Isaac Input Parameters: 73da43764aSToby Isaac + comm - the MPI communicator 74da43764aSToby Isaac . dim - the spatial dimension 75da43764aSToby Isaac - simplex - Flag for simplex, otherwise use a tensor-product cell 76da43764aSToby Isaac 77da43764aSToby Isaac Output Parameters: 78da43764aSToby Isaac . ref - the reference tree DMPlex object 79da43764aSToby Isaac 800b7167a0SToby Isaac Level: intermediate 81da43764aSToby Isaac 82da43764aSToby Isaac .keywords: reference cell 83da43764aSToby Isaac .seealso: DMPlexSetReferenceTree(), DMPlexGetReferenceTree() 84da43764aSToby Isaac @*/ 85da43764aSToby Isaac PetscErrorCode DMPlexCreateDefaultReferenceTree(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *ref) 86da43764aSToby Isaac { 87da43764aSToby Isaac DM K, Kref; 8810f7e118SToby Isaac PetscInt p, pStart, pEnd, pRefStart, pRefEnd, d, offset, parentSize, *parents, *childIDs; 89da43764aSToby Isaac PetscInt *permvals, *unionCones, *coneSizes, *unionOrientations, numUnionPoints, *numDimPoints, numCones, numVerts; 90da43764aSToby Isaac DMLabel identity, identityRef; 9110f7e118SToby Isaac PetscSection unionSection, unionConeSection, parentSection; 92da43764aSToby Isaac PetscScalar *unionCoords; 93da43764aSToby Isaac IS perm; 94da43764aSToby Isaac PetscErrorCode ierr; 95da43764aSToby Isaac 96da43764aSToby Isaac PetscFunctionBegin; 97da43764aSToby Isaac /* create a reference element */ 98da43764aSToby Isaac ierr = DMPlexCreateReferenceCell(comm, dim, simplex, &K);CHKERRQ(ierr); 99da43764aSToby Isaac ierr = DMPlexCreateLabel(K, "identity");CHKERRQ(ierr); 100da43764aSToby Isaac ierr = DMPlexGetLabel(K, "identity", &identity);CHKERRQ(ierr); 101da43764aSToby Isaac ierr = DMPlexGetChart(K, &pStart, &pEnd);CHKERRQ(ierr); 102da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 103da43764aSToby Isaac ierr = DMLabelSetValue(identity, p, p);CHKERRQ(ierr); 104da43764aSToby Isaac } 105da43764aSToby Isaac /* refine it */ 106da43764aSToby Isaac ierr = DMRefine(K,comm,&Kref);CHKERRQ(ierr); 107da43764aSToby Isaac 108da43764aSToby Isaac /* the reference tree is the union of these two, without duplicating 109da43764aSToby Isaac * points that appear in both */ 110da43764aSToby Isaac ierr = DMPlexGetLabel(Kref, "identity", &identityRef);CHKERRQ(ierr); 111da43764aSToby Isaac ierr = DMPlexGetChart(Kref, &pRefStart, &pRefEnd);CHKERRQ(ierr); 112da43764aSToby Isaac ierr = PetscSectionCreate(comm, &unionSection);CHKERRQ(ierr); 113da43764aSToby Isaac ierr = PetscSectionSetChart(unionSection, 0, (pEnd - pStart) + (pRefEnd - pRefStart));CHKERRQ(ierr); 114da43764aSToby Isaac /* count points that will go in the union */ 115da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 116da43764aSToby Isaac ierr = PetscSectionSetDof(unionSection, p - pStart, 1);CHKERRQ(ierr); 117da43764aSToby Isaac } 118da43764aSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 119da43764aSToby Isaac PetscInt q, qSize; 120da43764aSToby Isaac ierr = DMLabelGetValue(identityRef, p, &q);CHKERRQ(ierr); 121da43764aSToby Isaac ierr = DMLabelGetStratumSize(identityRef, q, &qSize);CHKERRQ(ierr); 122da43764aSToby Isaac if (qSize > 1) { 123da43764aSToby Isaac ierr = PetscSectionSetDof(unionSection, p - pRefStart + (pEnd - pStart), 1);CHKERRQ(ierr); 124da43764aSToby Isaac } 125da43764aSToby Isaac } 126da43764aSToby Isaac ierr = PetscMalloc1((pEnd - pStart) + (pRefEnd - pRefStart),&permvals);CHKERRQ(ierr); 127da43764aSToby Isaac offset = 0; 128da43764aSToby Isaac /* stratify points in the union by topological dimension */ 129da43764aSToby Isaac for (d = 0; d <= dim; d++) { 130da43764aSToby Isaac PetscInt cStart, cEnd, c; 131da43764aSToby Isaac 132da43764aSToby Isaac ierr = DMPlexGetHeightStratum(K, d, &cStart, &cEnd);CHKERRQ(ierr); 133da43764aSToby Isaac for (c = cStart; c < cEnd; c++) { 134da43764aSToby Isaac permvals[offset++] = c; 135da43764aSToby Isaac } 136da43764aSToby Isaac 137da43764aSToby Isaac ierr = DMPlexGetHeightStratum(Kref, d, &cStart, &cEnd);CHKERRQ(ierr); 138da43764aSToby Isaac for (c = cStart; c < cEnd; c++) { 139da43764aSToby Isaac permvals[offset++] = c + (pEnd - pStart); 140da43764aSToby Isaac } 141da43764aSToby Isaac } 142da43764aSToby Isaac ierr = ISCreateGeneral(comm, (pEnd - pStart) + (pRefEnd - pRefStart), permvals, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr); 143da43764aSToby Isaac ierr = PetscSectionSetPermutation(unionSection,perm);CHKERRQ(ierr); 144da43764aSToby Isaac ierr = PetscSectionSetUp(unionSection);CHKERRQ(ierr); 145da43764aSToby Isaac ierr = PetscSectionGetStorageSize(unionSection,&numUnionPoints);CHKERRQ(ierr); 146da43764aSToby Isaac ierr = PetscMalloc2(numUnionPoints,&coneSizes,dim+1,&numDimPoints);CHKERRQ(ierr); 147da43764aSToby Isaac /* count dimension points */ 148da43764aSToby Isaac for (d = 0; d <= dim; d++) { 149da43764aSToby Isaac PetscInt cStart, cOff, cOff2; 150da43764aSToby Isaac ierr = DMPlexGetHeightStratum(K,d,&cStart,NULL);CHKERRQ(ierr); 151da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,cStart-pStart,&cOff);CHKERRQ(ierr); 152da43764aSToby Isaac if (d < dim) { 153da43764aSToby Isaac ierr = DMPlexGetHeightStratum(K,d+1,&cStart,NULL);CHKERRQ(ierr); 154da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,cStart-pStart,&cOff2);CHKERRQ(ierr); 155da43764aSToby Isaac } 156da43764aSToby Isaac else { 157da43764aSToby Isaac cOff2 = numUnionPoints; 158da43764aSToby Isaac } 159da43764aSToby Isaac numDimPoints[dim - d] = cOff2 - cOff; 160da43764aSToby Isaac } 161da43764aSToby Isaac ierr = PetscSectionCreate(comm, &unionConeSection);CHKERRQ(ierr); 162da43764aSToby Isaac ierr = PetscSectionSetChart(unionConeSection, 0, numUnionPoints);CHKERRQ(ierr); 163da43764aSToby Isaac /* count the cones in the union */ 164da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 165da43764aSToby Isaac PetscInt dof, uOff; 166da43764aSToby Isaac 167da43764aSToby Isaac ierr = DMPlexGetConeSize(K, p, &dof);CHKERRQ(ierr); 168da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pStart,&uOff);CHKERRQ(ierr); 169da43764aSToby Isaac ierr = PetscSectionSetDof(unionConeSection, uOff, dof);CHKERRQ(ierr); 170da43764aSToby Isaac coneSizes[uOff] = dof; 171da43764aSToby Isaac } 172da43764aSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 173da43764aSToby Isaac PetscInt dof, uDof, uOff; 174da43764aSToby Isaac 175da43764aSToby Isaac ierr = DMPlexGetConeSize(Kref, p, &dof);CHKERRQ(ierr); 176da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 177da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 178da43764aSToby Isaac if (uDof) { 179da43764aSToby Isaac ierr = PetscSectionSetDof(unionConeSection, uOff, dof);CHKERRQ(ierr); 180da43764aSToby Isaac coneSizes[uOff] = dof; 181da43764aSToby Isaac } 182da43764aSToby Isaac } 183da43764aSToby Isaac ierr = PetscSectionSetUp(unionConeSection);CHKERRQ(ierr); 184da43764aSToby Isaac ierr = PetscSectionGetStorageSize(unionConeSection,&numCones);CHKERRQ(ierr); 185da43764aSToby Isaac ierr = PetscMalloc2(numCones,&unionCones,numCones,&unionOrientations);CHKERRQ(ierr); 186da43764aSToby Isaac /* write the cones in the union */ 187da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 188da43764aSToby Isaac PetscInt dof, uOff, c, cOff; 189da43764aSToby Isaac const PetscInt *cone, *orientation; 190da43764aSToby Isaac 191da43764aSToby Isaac ierr = DMPlexGetConeSize(K, p, &dof);CHKERRQ(ierr); 192da43764aSToby Isaac ierr = DMPlexGetCone(K, p, &cone);CHKERRQ(ierr); 193da43764aSToby Isaac ierr = DMPlexGetConeOrientation(K, p, &orientation);CHKERRQ(ierr); 194da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pStart,&uOff);CHKERRQ(ierr); 195da43764aSToby Isaac ierr = PetscSectionGetOffset(unionConeSection,uOff,&cOff);CHKERRQ(ierr); 196da43764aSToby Isaac for (c = 0; c < dof; c++) { 197da43764aSToby Isaac PetscInt e, eOff; 198da43764aSToby Isaac e = cone[c]; 199da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, e - pStart, &eOff);CHKERRQ(ierr); 200da43764aSToby Isaac unionCones[cOff + c] = eOff; 201da43764aSToby Isaac unionOrientations[cOff + c] = orientation[c]; 202da43764aSToby Isaac } 203da43764aSToby Isaac } 204da43764aSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 205da43764aSToby Isaac PetscInt dof, uDof, uOff, c, cOff; 206da43764aSToby Isaac const PetscInt *cone, *orientation; 207da43764aSToby Isaac 208da43764aSToby Isaac ierr = DMPlexGetConeSize(Kref, p, &dof);CHKERRQ(ierr); 209da43764aSToby Isaac ierr = DMPlexGetCone(Kref, p, &cone);CHKERRQ(ierr); 210da43764aSToby Isaac ierr = DMPlexGetConeOrientation(Kref, p, &orientation);CHKERRQ(ierr); 211da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 212da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 213da43764aSToby Isaac if (uDof) { 214da43764aSToby Isaac ierr = PetscSectionGetOffset(unionConeSection,uOff,&cOff);CHKERRQ(ierr); 215da43764aSToby Isaac for (c = 0; c < dof; c++) { 216da43764aSToby Isaac PetscInt e, eOff, eDof; 217da43764aSToby Isaac 218da43764aSToby Isaac e = cone[c]; 219da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection, e - pRefStart + (pEnd - pStart),&eDof);CHKERRQ(ierr); 220da43764aSToby Isaac if (eDof) { 221da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, e - pRefStart + (pEnd - pStart), &eOff);CHKERRQ(ierr); 222da43764aSToby Isaac } 223da43764aSToby Isaac else { 224da43764aSToby Isaac ierr = DMLabelGetValue(identityRef, e, &e);CHKERRQ(ierr); 225da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, e - pStart, &eOff);CHKERRQ(ierr); 226da43764aSToby Isaac } 227da43764aSToby Isaac unionCones[cOff + c] = eOff; 228da43764aSToby Isaac unionOrientations[cOff + c] = orientation[c]; 229da43764aSToby Isaac } 230da43764aSToby Isaac } 231da43764aSToby Isaac } 232da43764aSToby Isaac /* get the coordinates */ 233da43764aSToby Isaac { 234da43764aSToby Isaac PetscInt vStart, vEnd, vRefStart, vRefEnd, v, vDof, vOff; 235da43764aSToby Isaac PetscSection KcoordsSec, KrefCoordsSec; 236da43764aSToby Isaac Vec KcoordsVec, KrefCoordsVec; 237da43764aSToby Isaac PetscScalar *Kcoords; 238da43764aSToby Isaac 239da43764aSToby Isaac DMGetCoordinateSection(K, &KcoordsSec);CHKERRQ(ierr); 240da43764aSToby Isaac DMGetCoordinatesLocal(K, &KcoordsVec);CHKERRQ(ierr); 241da43764aSToby Isaac DMGetCoordinateSection(Kref, &KrefCoordsSec);CHKERRQ(ierr); 242da43764aSToby Isaac DMGetCoordinatesLocal(Kref, &KrefCoordsVec);CHKERRQ(ierr); 243da43764aSToby Isaac 244da43764aSToby Isaac numVerts = numDimPoints[0]; 245da43764aSToby Isaac ierr = PetscMalloc1(numVerts * dim,&unionCoords);CHKERRQ(ierr); 246da43764aSToby Isaac ierr = DMPlexGetDepthStratum(K,0,&vStart,&vEnd);CHKERRQ(ierr); 247da43764aSToby Isaac 248da43764aSToby Isaac offset = 0; 249da43764aSToby Isaac for (v = vStart; v < vEnd; v++) { 250da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,v - pStart,&vOff);CHKERRQ(ierr); 251da43764aSToby Isaac ierr = VecGetValuesSection(KcoordsVec, KcoordsSec, v, &Kcoords);CHKERRQ(ierr); 252da43764aSToby Isaac for (d = 0; d < dim; d++) { 253da43764aSToby Isaac unionCoords[offset * dim + d] = Kcoords[d]; 254da43764aSToby Isaac } 255da43764aSToby Isaac offset++; 256da43764aSToby Isaac } 257da43764aSToby Isaac ierr = DMPlexGetDepthStratum(Kref,0,&vRefStart,&vRefEnd);CHKERRQ(ierr); 258da43764aSToby Isaac for (v = vRefStart; v < vRefEnd; v++) { 259da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection,v - pRefStart + (pEnd - pStart),&vDof);CHKERRQ(ierr); 260da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,v - pRefStart + (pEnd - pStart),&vOff);CHKERRQ(ierr); 261da43764aSToby Isaac ierr = VecGetValuesSection(KrefCoordsVec, KrefCoordsSec, v, &Kcoords);CHKERRQ(ierr); 262da43764aSToby Isaac if (vDof) { 263da43764aSToby Isaac for (d = 0; d < dim; d++) { 264da43764aSToby Isaac unionCoords[offset * dim + d] = Kcoords[d]; 265da43764aSToby Isaac } 266da43764aSToby Isaac offset++; 267da43764aSToby Isaac } 268da43764aSToby Isaac } 269da43764aSToby Isaac } 270da43764aSToby Isaac ierr = DMCreate(comm,ref);CHKERRQ(ierr); 271da43764aSToby Isaac ierr = DMSetType(*ref,DMPLEX);CHKERRQ(ierr); 272da43764aSToby Isaac ierr = DMPlexSetDimension(*ref,dim);CHKERRQ(ierr); 273da43764aSToby Isaac ierr = DMPlexCreateFromDAG(*ref,dim,numDimPoints,coneSizes,unionCones,unionOrientations,unionCoords);CHKERRQ(ierr); 27410f7e118SToby Isaac /* set the tree */ 27510f7e118SToby Isaac ierr = PetscSectionCreate(comm,&parentSection);CHKERRQ(ierr); 27610f7e118SToby Isaac ierr = PetscSectionSetChart(parentSection,0,numUnionPoints);CHKERRQ(ierr); 27710f7e118SToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 27810f7e118SToby Isaac PetscInt uDof, uOff; 27910f7e118SToby Isaac 28010f7e118SToby Isaac ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 28110f7e118SToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 28210f7e118SToby Isaac if (uDof) { 28310f7e118SToby Isaac PetscSectionSetDof(parentSection,uOff,1);CHKERRQ(ierr); 28410f7e118SToby Isaac } 28510f7e118SToby Isaac } 28610f7e118SToby Isaac ierr = PetscSectionSetUp(parentSection);CHKERRQ(ierr); 28710f7e118SToby Isaac ierr = PetscSectionGetStorageSize(parentSection,&parentSize);CHKERRQ(ierr); 28810f7e118SToby Isaac ierr = PetscMalloc2(parentSize,&parents,parentSize,&childIDs);CHKERRQ(ierr); 28910f7e118SToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 29010f7e118SToby Isaac PetscInt uDof, uOff; 29110f7e118SToby Isaac 29210f7e118SToby Isaac ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 29310f7e118SToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 29410f7e118SToby Isaac if (uDof) { 29510f7e118SToby Isaac PetscInt pOff, parent, parentU; 29610f7e118SToby Isaac PetscSectionGetOffset(parentSection,uOff,&pOff);CHKERRQ(ierr); 29710f7e118SToby Isaac DMLabelGetValue(identityRef,p,&parent);CHKERRQ(ierr); 29810f7e118SToby Isaac ierr = PetscSectionGetOffset(unionSection, parent - pStart,&parentU);CHKERRQ(ierr); 29910f7e118SToby Isaac parents[pOff] = parentU; 30010f7e118SToby Isaac childIDs[pOff] = uOff; 30110f7e118SToby Isaac } 30210f7e118SToby Isaac } 30310f7e118SToby Isaac ierr = DMPlexSetTree(*ref,parentSection,parents,childIDs);CHKERRQ(ierr); 30410f7e118SToby Isaac ierr = PetscSectionDestroy(&parentSection);CHKERRQ(ierr); 30510f7e118SToby Isaac ierr = PetscFree2(parents,childIDs);CHKERRQ(ierr); 30610f7e118SToby Isaac 307da43764aSToby Isaac /* clean up */ 308da43764aSToby Isaac ierr = PetscSectionDestroy(&unionSection);CHKERRQ(ierr); 309da43764aSToby Isaac ierr = PetscSectionDestroy(&unionConeSection);CHKERRQ(ierr); 310da43764aSToby Isaac ierr = ISDestroy(&perm);CHKERRQ(ierr); 311da43764aSToby Isaac ierr = PetscFree(unionCoords);CHKERRQ(ierr); 312da43764aSToby Isaac ierr = PetscFree2(unionCones,unionOrientations);CHKERRQ(ierr); 313da43764aSToby Isaac ierr = PetscFree2(coneSizes,numDimPoints);CHKERRQ(ierr); 314da43764aSToby Isaac ierr = DMDestroy(&K);CHKERRQ(ierr); 315da43764aSToby Isaac ierr = DMDestroy(&Kref);CHKERRQ(ierr); 316da43764aSToby Isaac PetscFunctionReturn(0); 317da43764aSToby Isaac } 318da43764aSToby Isaac 319d961a43aSToby Isaac #undef __FUNCT__ 320878b19aaSToby Isaac #define __FUNCT__ "DMPlexTreeSymmetrize" 321878b19aaSToby Isaac static PetscErrorCode DMPlexTreeSymmetrize(DM dm) 322878b19aaSToby Isaac { 323878b19aaSToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 324878b19aaSToby Isaac PetscSection childSec, pSec; 325878b19aaSToby Isaac PetscInt p, pSize, cSize, parMax = PETSC_MIN_INT, parMin = PETSC_MAX_INT; 326878b19aaSToby Isaac PetscInt *offsets, *children, pStart, pEnd; 327878b19aaSToby Isaac PetscErrorCode ierr; 328878b19aaSToby Isaac 329878b19aaSToby Isaac PetscFunctionBegin; 330878b19aaSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 331878b19aaSToby Isaac ierr = PetscSectionDestroy(&mesh->childSection);CHKERRQ(ierr); 332878b19aaSToby Isaac ierr = PetscFree(mesh->children);CHKERRQ(ierr); 333878b19aaSToby Isaac pSec = mesh->parentSection; 334878b19aaSToby Isaac if (!pSec) PetscFunctionReturn(0); 335878b19aaSToby Isaac ierr = PetscSectionGetStorageSize(pSec,&pSize);CHKERRQ(ierr); 336878b19aaSToby Isaac for (p = 0; p < pSize; p++) { 337878b19aaSToby Isaac PetscInt par = mesh->parents[p]; 338878b19aaSToby Isaac 339878b19aaSToby Isaac parMax = PetscMax(parMax,par+1); 340878b19aaSToby Isaac parMin = PetscMin(parMin,par); 341878b19aaSToby Isaac } 342878b19aaSToby Isaac if (parMin > parMax) { 343878b19aaSToby Isaac parMin = -1; 344878b19aaSToby Isaac parMax = -1; 345878b19aaSToby Isaac } 346878b19aaSToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)pSec),&childSec);CHKERRQ(ierr); 347878b19aaSToby Isaac ierr = PetscSectionSetChart(childSec,parMin,parMax);CHKERRQ(ierr); 348878b19aaSToby Isaac for (p = 0; p < pSize; p++) { 349878b19aaSToby Isaac PetscInt par = mesh->parents[p]; 350878b19aaSToby Isaac 351878b19aaSToby Isaac ierr = PetscSectionAddDof(childSec,par,1);CHKERRQ(ierr); 352878b19aaSToby Isaac } 353878b19aaSToby Isaac ierr = PetscSectionSetUp(childSec);CHKERRQ(ierr); 354878b19aaSToby Isaac ierr = PetscSectionGetStorageSize(childSec,&cSize);CHKERRQ(ierr); 355878b19aaSToby Isaac ierr = PetscMalloc1(cSize,&children);CHKERRQ(ierr); 356878b19aaSToby Isaac ierr = PetscCalloc1(parMax-parMin,&offsets);CHKERRQ(ierr); 357878b19aaSToby Isaac ierr = PetscSectionGetChart(pSec,&pStart,&pEnd);CHKERRQ(ierr); 358878b19aaSToby Isaac for (p = pStart; p < pEnd; p++) { 359878b19aaSToby Isaac PetscInt dof, off, i; 360878b19aaSToby Isaac 361878b19aaSToby Isaac ierr = PetscSectionGetDof(pSec,p,&dof);CHKERRQ(ierr); 362878b19aaSToby Isaac ierr = PetscSectionGetOffset(pSec,p,&off);CHKERRQ(ierr); 363878b19aaSToby Isaac for (i = 0; i < dof; i++) { 364878b19aaSToby Isaac PetscInt par = mesh->parents[off + i], cOff; 365878b19aaSToby Isaac 366878b19aaSToby Isaac ierr = PetscSectionGetOffset(childSec,par,&cOff);CHKERRQ(ierr); 367878b19aaSToby Isaac children[cOff + offsets[par-parMin]++] = p; 368878b19aaSToby Isaac } 369878b19aaSToby Isaac } 370878b19aaSToby Isaac mesh->childSection = childSec; 371878b19aaSToby Isaac mesh->children = children; 372878b19aaSToby Isaac ierr = PetscFree(offsets);CHKERRQ(ierr); 373878b19aaSToby Isaac PetscFunctionReturn(0); 374878b19aaSToby Isaac } 375878b19aaSToby Isaac 376878b19aaSToby Isaac #undef __FUNCT__ 377*6dd5a8c8SToby Isaac #define __FUNCT__ "AnchorsFlatten" 378*6dd5a8c8SToby Isaac static PetscErrorCode AnchorsFlatten (PetscSection section, IS is, PetscSection *sectionNew, IS *isNew) 379*6dd5a8c8SToby Isaac { 380*6dd5a8c8SToby Isaac PetscInt pStart, pEnd, size, sizeNew, i, p, *valsNew = NULL; 381*6dd5a8c8SToby Isaac const PetscInt *vals; 382*6dd5a8c8SToby Isaac PetscSection secNew; 383*6dd5a8c8SToby Isaac PetscBool anyNew, globalAnyNew; 384*6dd5a8c8SToby Isaac PetscBool compress; 385*6dd5a8c8SToby Isaac PetscErrorCode ierr; 386*6dd5a8c8SToby Isaac 387*6dd5a8c8SToby Isaac PetscFunctionBegin; 388*6dd5a8c8SToby Isaac ierr = PetscSectionGetChart(section,&pStart,&pEnd);CHKERRQ(ierr); 389*6dd5a8c8SToby Isaac ierr = ISGetLocalSize(is,&size);CHKERRQ(ierr); 390*6dd5a8c8SToby Isaac ierr = ISGetIndices(is,&vals);CHKERRQ(ierr); 391*6dd5a8c8SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)section),&secNew);CHKERRQ(ierr); 392*6dd5a8c8SToby Isaac ierr = PetscSectionSetChart(secNew,pStart,pEnd);CHKERRQ(ierr); 393*6dd5a8c8SToby Isaac for (i = 0; i < size; i++) { 394*6dd5a8c8SToby Isaac PetscInt dof; 395*6dd5a8c8SToby Isaac 396*6dd5a8c8SToby Isaac p = vals[i]; 397*6dd5a8c8SToby Isaac if (p < pStart || p >= pEnd) continue; 398*6dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 399*6dd5a8c8SToby Isaac if (dof) break; 400*6dd5a8c8SToby Isaac } 401*6dd5a8c8SToby Isaac if (i == size) { 402*6dd5a8c8SToby Isaac ierr = PetscSectionSetUp(secNew);CHKERRQ(ierr); 403*6dd5a8c8SToby Isaac anyNew = PETSC_FALSE; 404*6dd5a8c8SToby Isaac compress = PETSC_FALSE; 405*6dd5a8c8SToby Isaac sizeNew = 0; 406*6dd5a8c8SToby Isaac } 407*6dd5a8c8SToby Isaac else { 408*6dd5a8c8SToby Isaac anyNew = PETSC_TRUE; 409*6dd5a8c8SToby Isaac for (p = pStart; p < pEnd; p++) { 410*6dd5a8c8SToby Isaac PetscInt dof, off; 411*6dd5a8c8SToby Isaac 412*6dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 413*6dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 414*6dd5a8c8SToby Isaac for (i = 0; i < dof; i++) { 415*6dd5a8c8SToby Isaac PetscInt q = vals[off + i], qDof = 0; 416*6dd5a8c8SToby Isaac 417*6dd5a8c8SToby Isaac if (q >= pStart && q < pEnd) { 418*6dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, q, &qDof);CHKERRQ(ierr); 419*6dd5a8c8SToby Isaac } 420*6dd5a8c8SToby Isaac if (qDof) { 421*6dd5a8c8SToby Isaac ierr = PetscSectionAddDof(secNew, p, qDof);CHKERRQ(ierr); 422*6dd5a8c8SToby Isaac } 423*6dd5a8c8SToby Isaac else { 424*6dd5a8c8SToby Isaac ierr = PetscSectionAddDof(secNew, p, 1);CHKERRQ(ierr); 425*6dd5a8c8SToby Isaac } 426*6dd5a8c8SToby Isaac } 427*6dd5a8c8SToby Isaac } 428*6dd5a8c8SToby Isaac ierr = PetscSectionSetUp(secNew);CHKERRQ(ierr); 429*6dd5a8c8SToby Isaac ierr = PetscSectionGetStorageSize(secNew,&sizeNew);CHKERRQ(ierr); 430*6dd5a8c8SToby Isaac ierr = PetscMalloc1(sizeNew,&valsNew);CHKERRQ(ierr); 431*6dd5a8c8SToby Isaac compress = PETSC_FALSE; 432*6dd5a8c8SToby Isaac for (p = pStart; p < pEnd; p++) { 433*6dd5a8c8SToby Isaac PetscInt dof, off, count, offNew, dofNew; 434*6dd5a8c8SToby Isaac 435*6dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, p, &dof);CHKERRQ(ierr); 436*6dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(section, p, &off);CHKERRQ(ierr); 437*6dd5a8c8SToby Isaac ierr = PetscSectionGetDof(secNew, p, &dofNew);CHKERRQ(ierr); 438*6dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(secNew, p, &offNew);CHKERRQ(ierr); 439*6dd5a8c8SToby Isaac count = 0; 440*6dd5a8c8SToby Isaac for (i = 0; i < dof; i++) { 441*6dd5a8c8SToby Isaac PetscInt q = vals[off + i], qDof = 0, qOff = 0, j; 442*6dd5a8c8SToby Isaac 443*6dd5a8c8SToby Isaac if (q >= pStart && q < pEnd) { 444*6dd5a8c8SToby Isaac ierr = PetscSectionGetDof(section, q, &qDof);CHKERRQ(ierr); 445*6dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(section, q, &qOff);CHKERRQ(ierr); 446*6dd5a8c8SToby Isaac } 447*6dd5a8c8SToby Isaac if (qDof) { 448*6dd5a8c8SToby Isaac PetscInt oldCount = count; 449*6dd5a8c8SToby Isaac 450*6dd5a8c8SToby Isaac for (j = 0; j < qDof; j++) { 451*6dd5a8c8SToby Isaac PetscInt k, r = vals[qOff + j]; 452*6dd5a8c8SToby Isaac 453*6dd5a8c8SToby Isaac for (k = 0; k < oldCount; k++) { 454*6dd5a8c8SToby Isaac if (valsNew[offNew + k] == r) { 455*6dd5a8c8SToby Isaac break; 456*6dd5a8c8SToby Isaac } 457*6dd5a8c8SToby Isaac } 458*6dd5a8c8SToby Isaac if (k == oldCount) { 459*6dd5a8c8SToby Isaac valsNew[offNew + count++] = r; 460*6dd5a8c8SToby Isaac } 461*6dd5a8c8SToby Isaac } 462*6dd5a8c8SToby Isaac } 463*6dd5a8c8SToby Isaac else { 464*6dd5a8c8SToby Isaac PetscInt k, oldCount = count; 465*6dd5a8c8SToby Isaac 466*6dd5a8c8SToby Isaac for (k = 0; k < oldCount; k++) { 467*6dd5a8c8SToby Isaac if (valsNew[offNew + k] == q) { 468*6dd5a8c8SToby Isaac break; 469*6dd5a8c8SToby Isaac } 470*6dd5a8c8SToby Isaac } 471*6dd5a8c8SToby Isaac if (k == oldCount) { 472*6dd5a8c8SToby Isaac valsNew[offNew + count++] = q; 473*6dd5a8c8SToby Isaac } 474*6dd5a8c8SToby Isaac } 475*6dd5a8c8SToby Isaac } 476*6dd5a8c8SToby Isaac if (count < dofNew) { 477*6dd5a8c8SToby Isaac ierr = PetscSectionSetDof(secNew, p, count);CHKERRQ(ierr); 478*6dd5a8c8SToby Isaac compress = PETSC_TRUE; 479*6dd5a8c8SToby Isaac } 480*6dd5a8c8SToby Isaac } 481*6dd5a8c8SToby Isaac } 482*6dd5a8c8SToby Isaac ierr = ISRestoreIndices(is,&vals);CHKERRQ(ierr); 483*6dd5a8c8SToby Isaac ierr = MPI_Allreduce(&anyNew,&globalAnyNew,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)secNew));CHKERRQ(ierr); 484*6dd5a8c8SToby Isaac if (!globalAnyNew) { 485*6dd5a8c8SToby Isaac ierr = PetscSectionDestroy(&secNew);CHKERRQ(ierr); 486*6dd5a8c8SToby Isaac *sectionNew = NULL; 487*6dd5a8c8SToby Isaac *isNew = NULL; 488*6dd5a8c8SToby Isaac } 489*6dd5a8c8SToby Isaac else { 490*6dd5a8c8SToby Isaac PetscBool globalCompress; 491*6dd5a8c8SToby Isaac 492*6dd5a8c8SToby Isaac ierr = MPI_Allreduce(&compress,&globalCompress,1,MPIU_BOOL,MPI_LOR,PetscObjectComm((PetscObject)secNew));CHKERRQ(ierr); 493*6dd5a8c8SToby Isaac if (compress) { 494*6dd5a8c8SToby Isaac PetscSection secComp; 495*6dd5a8c8SToby Isaac PetscInt *valsComp = NULL; 496*6dd5a8c8SToby Isaac 497*6dd5a8c8SToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)section),&secComp);CHKERRQ(ierr); 498*6dd5a8c8SToby Isaac ierr = PetscSectionSetChart(secComp,pStart,pEnd);CHKERRQ(ierr); 499*6dd5a8c8SToby Isaac for (p = pStart; p < pEnd; p++) { 500*6dd5a8c8SToby Isaac PetscInt dof; 501*6dd5a8c8SToby Isaac 502*6dd5a8c8SToby Isaac ierr = PetscSectionGetDof(secNew, p, &dof);CHKERRQ(ierr); 503*6dd5a8c8SToby Isaac ierr = PetscSectionSetDof(secComp, p, dof);CHKERRQ(ierr); 504*6dd5a8c8SToby Isaac } 505*6dd5a8c8SToby Isaac ierr = PetscSectionSetUp(secComp);CHKERRQ(ierr); 506*6dd5a8c8SToby Isaac ierr = PetscSectionGetStorageSize(secComp,&sizeNew);CHKERRQ(ierr); 507*6dd5a8c8SToby Isaac ierr = PetscMalloc1(sizeNew,&valsComp);CHKERRQ(ierr); 508*6dd5a8c8SToby Isaac for (p = pStart; p < pEnd; p++) { 509*6dd5a8c8SToby Isaac PetscInt dof, off, offNew, j; 510*6dd5a8c8SToby Isaac 511*6dd5a8c8SToby Isaac ierr = PetscSectionGetDof(secNew, p, &dof);CHKERRQ(ierr); 512*6dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(secNew, p, &off);CHKERRQ(ierr); 513*6dd5a8c8SToby Isaac ierr = PetscSectionGetOffset(secComp, p, &offNew);CHKERRQ(ierr); 514*6dd5a8c8SToby Isaac for (j = 0; j < dof; j++) { 515*6dd5a8c8SToby Isaac valsComp[offNew + j] = valsNew[off + j]; 516*6dd5a8c8SToby Isaac } 517*6dd5a8c8SToby Isaac } 518*6dd5a8c8SToby Isaac ierr = PetscSectionDestroy(&secNew);CHKERRQ(ierr); 519*6dd5a8c8SToby Isaac secNew = secComp; 520*6dd5a8c8SToby Isaac ierr = PetscFree(valsNew);CHKERRQ(ierr); 521*6dd5a8c8SToby Isaac valsNew = valsComp; 522*6dd5a8c8SToby Isaac } 523*6dd5a8c8SToby Isaac ierr = ISCreateGeneral(PetscObjectComm((PetscObject)is),sizeNew,valsNew,PETSC_OWN_POINTER,isNew);CHKERRQ(ierr); 524*6dd5a8c8SToby Isaac } 525*6dd5a8c8SToby Isaac PetscFunctionReturn(0); 526*6dd5a8c8SToby Isaac } 527*6dd5a8c8SToby Isaac 528*6dd5a8c8SToby Isaac #undef __FUNCT__ 52966af876cSToby Isaac #define __FUNCT__ "DMPlexComputeConstraints_Tree" 53066af876cSToby Isaac static PetscErrorCode DMPlexComputeConstraints_Tree(DM dm) 53166af876cSToby Isaac { 53266af876cSToby Isaac PetscInt p, pStart, pEnd, *anchors, size; 53366af876cSToby Isaac PetscInt aMin = PETSC_MAX_INT, aMax = PETSC_MIN_INT; 53466af876cSToby Isaac PetscSection aSec; 53566af876cSToby Isaac IS aIS; 53666af876cSToby Isaac PetscErrorCode ierr; 53766af876cSToby Isaac 53866af876cSToby Isaac PetscFunctionBegin; 53966af876cSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 54066af876cSToby Isaac ierr = DMPlexGetChart(dm,&pStart,&pEnd);CHKERRQ(ierr); 54166af876cSToby Isaac for (p = pStart; p < pEnd; p++) { 54266af876cSToby Isaac PetscInt parent; 54366af876cSToby Isaac 54466af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 54566af876cSToby Isaac if (parent != p) { 54666af876cSToby Isaac aMin = PetscMin(aMin,p); 54766af876cSToby Isaac aMax = PetscMax(aMax,p+1); 54866af876cSToby Isaac } 54966af876cSToby Isaac } 55066af876cSToby Isaac if (aMin > aMax) { 55166af876cSToby Isaac aMin = -1; 55266af876cSToby Isaac aMax = -1; 55366af876cSToby Isaac } 55466af876cSToby Isaac ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&aSec);CHKERRQ(ierr); 55566af876cSToby Isaac ierr = PetscSectionSetChart(aSec,aMin,aMax);CHKERRQ(ierr); 55666af876cSToby Isaac for (p = aMin; p < aMax; p++) { 55766af876cSToby Isaac PetscInt parent, ancestor = p; 55866af876cSToby Isaac 55966af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 56066af876cSToby Isaac while (parent != ancestor) { 56166af876cSToby Isaac ancestor = parent; 56266af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,ancestor,&parent,NULL);CHKERRQ(ierr); 56366af876cSToby Isaac } 56466af876cSToby Isaac if (ancestor != p) { 56566af876cSToby Isaac PetscInt closureSize, *closure = NULL; 56666af876cSToby Isaac 56766af876cSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 56866af876cSToby Isaac ierr = PetscSectionSetDof(aSec,p,closureSize);CHKERRQ(ierr); 56966af876cSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 57066af876cSToby Isaac } 57166af876cSToby Isaac } 57266af876cSToby Isaac ierr = PetscSectionSetUp(aSec);CHKERRQ(ierr); 57366af876cSToby Isaac ierr = PetscSectionGetStorageSize(aSec,&size);CHKERRQ(ierr); 57466af876cSToby Isaac ierr = PetscMalloc1(size,&anchors);CHKERRQ(ierr); 57566af876cSToby Isaac for (p = aMin; p < aMax; p++) { 57666af876cSToby Isaac PetscInt parent, ancestor = p; 57766af876cSToby Isaac 57866af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,p,&parent,NULL);CHKERRQ(ierr); 57966af876cSToby Isaac while (parent != ancestor) { 58066af876cSToby Isaac ancestor = parent; 58166af876cSToby Isaac ierr = DMPlexGetTreeParent(dm,ancestor,&parent,NULL);CHKERRQ(ierr); 58266af876cSToby Isaac } 58366af876cSToby Isaac if (ancestor != p) { 58466af876cSToby Isaac PetscInt j, closureSize, *closure = NULL, aOff; 58566af876cSToby Isaac 58666af876cSToby Isaac ierr = PetscSectionGetOffset(aSec,p,&aOff);CHKERRQ(ierr); 58766af876cSToby Isaac 58866af876cSToby Isaac ierr = DMPlexGetTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 58966af876cSToby Isaac for (j = 0; j < closureSize; j++) { 59066af876cSToby Isaac anchors[aOff + j] = closure[2*j]; 59166af876cSToby Isaac } 59266af876cSToby Isaac ierr = DMPlexRestoreTransitiveClosure(dm,ancestor,PETSC_TRUE,&closureSize,&closure);CHKERRQ(ierr); 59366af876cSToby Isaac } 59466af876cSToby Isaac } 59566af876cSToby Isaac ierr = ISCreateGeneral(PetscObjectComm((PetscObject)dm),size,anchors,PETSC_OWN_POINTER,&aIS);CHKERRQ(ierr); 596*6dd5a8c8SToby Isaac { 597*6dd5a8c8SToby Isaac PetscSection aSecNew = aSec; 598*6dd5a8c8SToby Isaac IS aISNew = aIS; 599*6dd5a8c8SToby Isaac 600*6dd5a8c8SToby Isaac ierr = PetscObjectReference((PetscObject)aSec);CHKERRQ(ierr); 601*6dd5a8c8SToby Isaac ierr = PetscObjectReference((PetscObject)aIS);CHKERRQ(ierr); 602*6dd5a8c8SToby Isaac while (aSecNew) { 603*6dd5a8c8SToby Isaac ierr = PetscSectionDestroy(&aSec);CHKERRQ(ierr); 604*6dd5a8c8SToby Isaac ierr = ISDestroy(&aIS);CHKERRQ(ierr); 605*6dd5a8c8SToby Isaac aSec = aSecNew; 606*6dd5a8c8SToby Isaac aIS = aISNew; 607*6dd5a8c8SToby Isaac aSecNew = NULL; 608*6dd5a8c8SToby Isaac aISNew = NULL; 609*6dd5a8c8SToby Isaac ierr = AnchorsFlatten(aSec,aIS,&aSecNew,&aISNew);CHKERRQ(ierr); 610*6dd5a8c8SToby Isaac } 611*6dd5a8c8SToby Isaac } 61266af876cSToby Isaac ierr = DMPlexSetConstraints(dm,aSec,aIS);CHKERRQ(ierr); 61366af876cSToby Isaac ierr = PetscSectionDestroy(&aSec);CHKERRQ(ierr); 61466af876cSToby Isaac ierr = ISDestroy(&aIS);CHKERRQ(ierr); 61566af876cSToby Isaac PetscFunctionReturn(0); 61666af876cSToby Isaac } 61766af876cSToby Isaac 61866af876cSToby Isaac #undef __FUNCT__ 6190b7167a0SToby Isaac #define __FUNCT__ "DMPlexSetTree" 6200b7167a0SToby Isaac /*@ 6210b7167a0SToby Isaac DMPlexSetTree - set the tree that describes the hierarchy of non-conforming mesh points. This routine also creates 6220b7167a0SToby Isaac the point-to-point constraints determined by the tree: a point is constained to the points in the closure of its 6230b7167a0SToby Isaac tree root. 6240b7167a0SToby Isaac 6250b7167a0SToby Isaac Collective on dm 6260b7167a0SToby Isaac 6270b7167a0SToby Isaac Input Parameters: 6280b7167a0SToby Isaac + dm - the DMPlex object 6290b7167a0SToby Isaac . parentSection - a section describing the tree: a point has a parent if it has 1 dof in the section; the section 6300b7167a0SToby Isaac offset indexes the parent and childID list; the reference count of parentSection is incremented 6310b7167a0SToby Isaac . parents - a list of the point parents; copied, can be destroyed 6320b7167a0SToby Isaac - childIDs - identifies the relationship of the child point to the parent point; if there is a reference tree, then 6330b7167a0SToby Isaac the child corresponds to the point in the reference tree with index childIDs; copied, can be destroyed 6340b7167a0SToby Isaac 6350b7167a0SToby Isaac Level: intermediate 6360b7167a0SToby Isaac 637b2f41788SToby Isaac .seealso: DMPlexGetTree(), DMPlexSetReferenceTree(), DMPlexSetConstraints(), DMPlexGetTreeParent(), DMPlexGetTreeChildren() 6380b7167a0SToby Isaac @*/ 639b2f41788SToby Isaac PetscErrorCode DMPlexSetTree(DM dm, PetscSection parentSection, PetscInt parents[], PetscInt childIDs[]) 6400b7167a0SToby Isaac { 6410b7167a0SToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 6420b7167a0SToby Isaac PetscInt size; 6430b7167a0SToby Isaac PetscErrorCode ierr; 6440b7167a0SToby Isaac 6450b7167a0SToby Isaac PetscFunctionBegin; 6460b7167a0SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 6470b7167a0SToby Isaac PetscValidHeaderSpecific(parentSection, PETSC_SECTION_CLASSID, 2); 6480b7167a0SToby Isaac ierr = PetscObjectReference((PetscObject)parentSection);CHKERRQ(ierr); 6490b7167a0SToby Isaac ierr = PetscSectionDestroy(&mesh->parentSection);CHKERRQ(ierr); 6500b7167a0SToby Isaac mesh->parentSection = parentSection; 6510b7167a0SToby Isaac ierr = PetscSectionGetStorageSize(parentSection,&size);CHKERRQ(ierr); 6520b7167a0SToby Isaac if (parents != mesh->parents) { 6530b7167a0SToby Isaac ierr = PetscFree(mesh->parents);CHKERRQ(ierr); 6540b7167a0SToby Isaac ierr = PetscMalloc1(size,&mesh->parents);CHKERRQ(ierr); 6550b7167a0SToby Isaac ierr = PetscMemcpy(mesh->parents, parents, size * sizeof(*parents));CHKERRQ(ierr); 6560b7167a0SToby Isaac } 6570b7167a0SToby Isaac if (childIDs != mesh->childIDs) { 6580b7167a0SToby Isaac ierr = PetscFree(mesh->childIDs);CHKERRQ(ierr); 6590b7167a0SToby Isaac ierr = PetscMalloc1(size,&mesh->childIDs);CHKERRQ(ierr); 6600b7167a0SToby Isaac ierr = PetscMemcpy(mesh->childIDs, childIDs, size * sizeof(*childIDs));CHKERRQ(ierr); 6610b7167a0SToby Isaac } 662878b19aaSToby Isaac ierr = DMPlexTreeSymmetrize(dm);CHKERRQ(ierr); 66366af876cSToby Isaac ierr = DMPlexComputeConstraints_Tree(dm);CHKERRQ(ierr); 6640b7167a0SToby Isaac PetscFunctionReturn(0); 6650b7167a0SToby Isaac } 6660b7167a0SToby Isaac 6670b7167a0SToby Isaac #undef __FUNCT__ 668b2f41788SToby Isaac #define __FUNCT__ "DMPlexGetTree" 669b2f41788SToby Isaac /*@ 670b2f41788SToby Isaac DMPlexGetTree - get the tree that describes the hierarchy of non-conforming mesh points. 671b2f41788SToby Isaac Collective on dm 672b2f41788SToby Isaac 673b2f41788SToby Isaac Input Parameters: 674b2f41788SToby Isaac . dm - the DMPlex object 675b2f41788SToby Isaac 676b2f41788SToby Isaac Output Parameters: 677b2f41788SToby Isaac + parentSection - a section describing the tree: a point has a parent if it has 1 dof in the section; the section 678b2f41788SToby Isaac offset indexes the parent and childID list 679b2f41788SToby Isaac . parents - a list of the point parents 680b2f41788SToby Isaac . childIDs - identifies the relationship of the child point to the parent point; if there is a reference tree, then 681b2f41788SToby Isaac the child corresponds to the point in the reference tree with index childID 682b2f41788SToby Isaac . childSection - the inverse of the parent section 683b2f41788SToby Isaac - children - a list of the point children 684b2f41788SToby Isaac 685b2f41788SToby Isaac Level: intermediate 686b2f41788SToby Isaac 687b2f41788SToby Isaac .seealso: DMPlexSetTree(), DMPlexSetReferenceTree(), DMPlexSetConstraints(), DMPlexGetTreeParent(), DMPlexGetTreeChildren() 688b2f41788SToby Isaac @*/ 689b2f41788SToby Isaac PetscErrorCode DMPlexGetTree(DM dm, PetscSection *parentSection, PetscInt *parents[], PetscInt *childIDs[], PetscSection *childSection, PetscInt *children[]) 690b2f41788SToby Isaac { 691b2f41788SToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 692b2f41788SToby Isaac 693b2f41788SToby Isaac PetscFunctionBegin; 694b2f41788SToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 695b2f41788SToby Isaac if (parentSection) *parentSection = mesh->parentSection; 696b2f41788SToby Isaac if (parents) *parents = mesh->parents; 697b2f41788SToby Isaac if (childIDs) *childIDs = mesh->childIDs; 698b2f41788SToby Isaac if (childSection) *childSection = mesh->childSection; 699b2f41788SToby Isaac if (children) *children = mesh->children; 700b2f41788SToby Isaac PetscFunctionReturn(0); 701b2f41788SToby Isaac } 702b2f41788SToby Isaac 703b2f41788SToby Isaac #undef __FUNCT__ 704d961a43aSToby Isaac #define __FUNCT__ "DMPlexGetTreeParent" 705d961a43aSToby Isaac /*@ 706d961a43aSToby Isaac DMPlexGetTreeParent - get the parent of a point in the tree describing the point hierarchy (not the Sieve DAG) 707d961a43aSToby Isaac 708d961a43aSToby Isaac Input Parameters: 709d961a43aSToby Isaac + dm - the DMPlex object 710d961a43aSToby Isaac - point - the query point 711d961a43aSToby Isaac 712d961a43aSToby Isaac Output Parameters: 713d961a43aSToby Isaac + parent - if not NULL, set to the parent of the point, or the point itself if the point does not have a parent 714d961a43aSToby Isaac - childID - if not NULL, set to the child ID of the point with respect to its parent, or 0 if the point 715d961a43aSToby Isaac does not have a parent 716d961a43aSToby Isaac 717d961a43aSToby Isaac Level: intermediate 718d961a43aSToby Isaac 719d961a43aSToby Isaac .seealso: DMPlexSetTree(), DMPlexGetTree(), DMPlexGetTreeChildren() 720d961a43aSToby Isaac @*/ 721d961a43aSToby Isaac PetscErrorCode DMPlexGetTreeParent(DM dm, PetscInt point, PetscInt *parent, PetscInt *childID) 722d961a43aSToby Isaac { 723d961a43aSToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 724d961a43aSToby Isaac PetscSection pSec; 725d961a43aSToby Isaac PetscErrorCode ierr; 726d961a43aSToby Isaac 727d961a43aSToby Isaac PetscFunctionBegin; 728d961a43aSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 729d961a43aSToby Isaac pSec = mesh->parentSection; 730d961a43aSToby Isaac if (pSec && point >= pSec->pStart && point < pSec->pEnd) { 731d961a43aSToby Isaac PetscInt dof; 732d961a43aSToby Isaac 733d961a43aSToby Isaac ierr = PetscSectionGetDof (pSec, point, &dof);CHKERRQ(ierr); 734d961a43aSToby Isaac if (dof) { 735d961a43aSToby Isaac PetscInt off; 736d961a43aSToby Isaac 737d961a43aSToby Isaac ierr = PetscSectionGetOffset (pSec, point, &off);CHKERRQ(ierr); 738d961a43aSToby Isaac if (parent) *parent = mesh->parents[off]; 739d961a43aSToby Isaac if (childID) *childID = mesh->childIDs[off]; 740d961a43aSToby Isaac PetscFunctionReturn(0); 741d961a43aSToby Isaac } 742d961a43aSToby Isaac } 743d961a43aSToby Isaac if (parent) { 744d961a43aSToby Isaac *parent = point; 745d961a43aSToby Isaac } 746d961a43aSToby Isaac if (childID) { 747d961a43aSToby Isaac *childID = 0; 748d961a43aSToby Isaac } 749d961a43aSToby Isaac PetscFunctionReturn(0); 750d961a43aSToby Isaac } 751d961a43aSToby Isaac 752d961a43aSToby Isaac #undef __FUNCT__ 753d961a43aSToby Isaac #define __FUNCT__ "DMPlexGetTreeChildren" 754d961a43aSToby Isaac /*@C 755d961a43aSToby Isaac DMPlexGetTreeChildren - get the children of a point in the tree describing the point hierarchy (not the Sieve DAG) 756d961a43aSToby Isaac 757d961a43aSToby Isaac Input Parameters: 758d961a43aSToby Isaac + dm - the DMPlex object 759d961a43aSToby Isaac - point - the query point 760d961a43aSToby Isaac 761d961a43aSToby Isaac Output Parameters: 762d961a43aSToby Isaac + numChildren - if not NULL, set to the number of children 763d961a43aSToby Isaac - children - if not NULL, set to a list children, or set to NULL if the point has no children 764d961a43aSToby Isaac 765d961a43aSToby Isaac Level: intermediate 766d961a43aSToby Isaac 767d961a43aSToby Isaac Fortran Notes: 768d961a43aSToby Isaac Since it returns an array, this routine is only available in Fortran 90, and you must 769d961a43aSToby Isaac include petsc.h90 in your code. 770d961a43aSToby Isaac 771d961a43aSToby Isaac .seealso: DMPlexSetTree(), DMPlexGetTree(), DMPlexGetTreeParent() 772d961a43aSToby Isaac @*/ 773d961a43aSToby Isaac PetscErrorCode DMPlexGetTreeChildren(DM dm, PetscInt point, PetscInt *numChildren, const PetscInt *children[]) 774d961a43aSToby Isaac { 775d961a43aSToby Isaac DM_Plex *mesh = (DM_Plex *)dm->data; 776d961a43aSToby Isaac PetscSection childSec; 777d961a43aSToby Isaac PetscInt dof = 0; 778d961a43aSToby Isaac PetscErrorCode ierr; 779d961a43aSToby Isaac 780d961a43aSToby Isaac PetscFunctionBegin; 781d961a43aSToby Isaac PetscValidHeaderSpecific(dm, DM_CLASSID, 1); 782d961a43aSToby Isaac childSec = mesh->childSection; 783d961a43aSToby Isaac if (childSec && point >= childSec->pStart && point < childSec->pEnd) { 784d961a43aSToby Isaac ierr = PetscSectionGetDof (childSec, point, &dof);CHKERRQ(ierr); 785d961a43aSToby Isaac } 786d961a43aSToby Isaac if (numChildren) *numChildren = dof; 787d961a43aSToby Isaac if (children) { 788d961a43aSToby Isaac if (dof) { 789d961a43aSToby Isaac PetscInt off; 790d961a43aSToby Isaac 791d961a43aSToby Isaac ierr = PetscSectionGetOffset (childSec, point, &off);CHKERRQ(ierr); 792d961a43aSToby Isaac *children = &mesh->children[off]; 793d961a43aSToby Isaac } 794d961a43aSToby Isaac else { 795d961a43aSToby Isaac *children = NULL; 796d961a43aSToby Isaac } 797d961a43aSToby Isaac } 798d961a43aSToby Isaac PetscFunctionReturn(0); 799d961a43aSToby Isaac } 800