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 19*da43764aSToby Isaac Level: beginner 20d6a7ad0dSToby Isaac 21*da43764aSToby 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 50*da43764aSToby Isaac Level: beginner 51d6a7ad0dSToby Isaac 52*da43764aSToby 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 65*da43764aSToby Isaac #undef __FUNCT__ 66*da43764aSToby Isaac #define __FUNCT__ "DMPlexCreateDefaultReferenceTree" 67*da43764aSToby Isaac /*@ 68*da43764aSToby Isaac DMPlexCreateDefaultReferenceTree - create a reference tree for isotropic hierarchical mesh refinement. 69*da43764aSToby Isaac 70*da43764aSToby Isaac Collective on comm 71*da43764aSToby Isaac 72*da43764aSToby Isaac Input Parameters: 73*da43764aSToby Isaac + comm - the MPI communicator 74*da43764aSToby Isaac . dim - the spatial dimension 75*da43764aSToby Isaac - simplex - Flag for simplex, otherwise use a tensor-product cell 76*da43764aSToby Isaac 77*da43764aSToby Isaac Output Parameters: 78*da43764aSToby Isaac . ref - the reference tree DMPlex object 79*da43764aSToby Isaac 80*da43764aSToby Isaac Level: beginner 81*da43764aSToby Isaac 82*da43764aSToby Isaac .keywords: reference cell 83*da43764aSToby Isaac .seealso: DMPlexSetReferenceTree(), DMPlexGetReferenceTree() 84*da43764aSToby Isaac @*/ 85*da43764aSToby Isaac PetscErrorCode DMPlexCreateDefaultReferenceTree(MPI_Comm comm, PetscInt dim, PetscBool simplex, DM *ref) 86*da43764aSToby Isaac { 87*da43764aSToby Isaac DM K, Kref; 88*da43764aSToby Isaac PetscInt p, pStart, pEnd, pRefStart, pRefEnd, d, offset; 89*da43764aSToby Isaac PetscInt *permvals, *unionCones, *coneSizes, *unionOrientations, numUnionPoints, *numDimPoints, numCones, numVerts; 90*da43764aSToby Isaac DMLabel identity, identityRef; 91*da43764aSToby Isaac PetscSection unionSection, unionConeSection; 92*da43764aSToby Isaac PetscScalar *unionCoords; 93*da43764aSToby Isaac IS perm; 94*da43764aSToby Isaac PetscErrorCode ierr; 95*da43764aSToby Isaac 96*da43764aSToby Isaac PetscFunctionBegin; 97*da43764aSToby Isaac /* create a reference element */ 98*da43764aSToby Isaac ierr = DMPlexCreateReferenceCell(comm, dim, simplex, &K);CHKERRQ(ierr); 99*da43764aSToby Isaac ierr = DMPlexCreateLabel(K, "identity");CHKERRQ(ierr); 100*da43764aSToby Isaac ierr = DMPlexGetLabel(K, "identity", &identity);CHKERRQ(ierr); 101*da43764aSToby Isaac ierr = DMPlexGetChart(K, &pStart, &pEnd);CHKERRQ(ierr); 102*da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 103*da43764aSToby Isaac ierr = DMLabelSetValue(identity, p, p);CHKERRQ(ierr); 104*da43764aSToby Isaac } 105*da43764aSToby Isaac /* refine it */ 106*da43764aSToby Isaac ierr = DMRefine(K,comm,&Kref);CHKERRQ(ierr); 107*da43764aSToby Isaac 108*da43764aSToby Isaac /* the reference tree is the union of these two, without duplicating 109*da43764aSToby Isaac * points that appear in both */ 110*da43764aSToby Isaac ierr = DMPlexGetLabel(Kref, "identity", &identityRef);CHKERRQ(ierr); 111*da43764aSToby Isaac ierr = DMPlexGetChart(Kref, &pRefStart, &pRefEnd);CHKERRQ(ierr); 112*da43764aSToby Isaac ierr = PetscSectionCreate(comm, &unionSection);CHKERRQ(ierr); 113*da43764aSToby Isaac ierr = PetscSectionSetChart(unionSection, 0, (pEnd - pStart) + (pRefEnd - pRefStart));CHKERRQ(ierr); 114*da43764aSToby Isaac /* count points that will go in the union */ 115*da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 116*da43764aSToby Isaac ierr = PetscSectionSetDof(unionSection, p - pStart, 1);CHKERRQ(ierr); 117*da43764aSToby Isaac } 118*da43764aSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 119*da43764aSToby Isaac PetscInt q, qSize; 120*da43764aSToby Isaac ierr = DMLabelGetValue(identityRef, p, &q);CHKERRQ(ierr); 121*da43764aSToby Isaac ierr = DMLabelGetStratumSize(identityRef, q, &qSize);CHKERRQ(ierr); 122*da43764aSToby Isaac if (qSize > 1) { 123*da43764aSToby Isaac ierr = PetscSectionSetDof(unionSection, p - pRefStart + (pEnd - pStart), 1);CHKERRQ(ierr); 124*da43764aSToby Isaac } 125*da43764aSToby Isaac } 126*da43764aSToby Isaac ierr = PetscMalloc1((pEnd - pStart) + (pRefEnd - pRefStart),&permvals);CHKERRQ(ierr); 127*da43764aSToby Isaac offset = 0; 128*da43764aSToby Isaac /* stratify points in the union by topological dimension */ 129*da43764aSToby Isaac for (d = 0; d <= dim; d++) { 130*da43764aSToby Isaac PetscInt cStart, cEnd, c; 131*da43764aSToby Isaac 132*da43764aSToby Isaac ierr = DMPlexGetHeightStratum(K, d, &cStart, &cEnd);CHKERRQ(ierr); 133*da43764aSToby Isaac for (c = cStart; c < cEnd; c++) { 134*da43764aSToby Isaac permvals[offset++] = c; 135*da43764aSToby Isaac } 136*da43764aSToby Isaac 137*da43764aSToby Isaac ierr = DMPlexGetHeightStratum(Kref, d, &cStart, &cEnd);CHKERRQ(ierr); 138*da43764aSToby Isaac for (c = cStart; c < cEnd; c++) { 139*da43764aSToby Isaac permvals[offset++] = c + (pEnd - pStart); 140*da43764aSToby Isaac } 141*da43764aSToby Isaac } 142*da43764aSToby Isaac ierr = ISCreateGeneral(comm, (pEnd - pStart) + (pRefEnd - pRefStart), permvals, PETSC_OWN_POINTER, &perm);CHKERRQ(ierr); 143*da43764aSToby Isaac ierr = PetscSectionSetPermutation(unionSection,perm);CHKERRQ(ierr); 144*da43764aSToby Isaac ierr = PetscSectionSetUp(unionSection);CHKERRQ(ierr); 145*da43764aSToby Isaac ierr = PetscSectionGetStorageSize(unionSection,&numUnionPoints);CHKERRQ(ierr); 146*da43764aSToby Isaac ierr = PetscMalloc2(numUnionPoints,&coneSizes,dim+1,&numDimPoints);CHKERRQ(ierr); 147*da43764aSToby Isaac /* count dimension points */ 148*da43764aSToby Isaac for (d = 0; d <= dim; d++) { 149*da43764aSToby Isaac PetscInt cStart, cOff, cOff2; 150*da43764aSToby Isaac ierr = DMPlexGetHeightStratum(K,d,&cStart,NULL);CHKERRQ(ierr); 151*da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,cStart-pStart,&cOff);CHKERRQ(ierr); 152*da43764aSToby Isaac if (d < dim) { 153*da43764aSToby Isaac ierr = DMPlexGetHeightStratum(K,d+1,&cStart,NULL);CHKERRQ(ierr); 154*da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,cStart-pStart,&cOff2);CHKERRQ(ierr); 155*da43764aSToby Isaac } 156*da43764aSToby Isaac else { 157*da43764aSToby Isaac cOff2 = numUnionPoints; 158*da43764aSToby Isaac } 159*da43764aSToby Isaac numDimPoints[dim - d] = cOff2 - cOff; 160*da43764aSToby Isaac } 161*da43764aSToby Isaac ierr = PetscSectionCreate(comm, &unionConeSection);CHKERRQ(ierr); 162*da43764aSToby Isaac ierr = PetscSectionSetChart(unionConeSection, 0, numUnionPoints);CHKERRQ(ierr); 163*da43764aSToby Isaac /* count the cones in the union */ 164*da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 165*da43764aSToby Isaac PetscInt dof, uOff; 166*da43764aSToby Isaac 167*da43764aSToby Isaac ierr = DMPlexGetConeSize(K, p, &dof);CHKERRQ(ierr); 168*da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pStart,&uOff);CHKERRQ(ierr); 169*da43764aSToby Isaac ierr = PetscSectionSetDof(unionConeSection, uOff, dof);CHKERRQ(ierr); 170*da43764aSToby Isaac coneSizes[uOff] = dof; 171*da43764aSToby Isaac } 172*da43764aSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 173*da43764aSToby Isaac PetscInt dof, uDof, uOff; 174*da43764aSToby Isaac 175*da43764aSToby Isaac ierr = DMPlexGetConeSize(Kref, p, &dof);CHKERRQ(ierr); 176*da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 177*da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 178*da43764aSToby Isaac if (uDof) { 179*da43764aSToby Isaac ierr = PetscSectionSetDof(unionConeSection, uOff, dof);CHKERRQ(ierr); 180*da43764aSToby Isaac coneSizes[uOff] = dof; 181*da43764aSToby Isaac } 182*da43764aSToby Isaac } 183*da43764aSToby Isaac ierr = PetscSectionSetUp(unionConeSection);CHKERRQ(ierr); 184*da43764aSToby Isaac ierr = PetscSectionGetStorageSize(unionConeSection,&numCones);CHKERRQ(ierr); 185*da43764aSToby Isaac ierr = PetscMalloc2(numCones,&unionCones,numCones,&unionOrientations);CHKERRQ(ierr); 186*da43764aSToby Isaac /* write the cones in the union */ 187*da43764aSToby Isaac for (p = pStart; p < pEnd; p++) { 188*da43764aSToby Isaac PetscInt dof, uOff, c, cOff; 189*da43764aSToby Isaac const PetscInt *cone, *orientation; 190*da43764aSToby Isaac 191*da43764aSToby Isaac ierr = DMPlexGetConeSize(K, p, &dof);CHKERRQ(ierr); 192*da43764aSToby Isaac ierr = DMPlexGetCone(K, p, &cone);CHKERRQ(ierr); 193*da43764aSToby Isaac ierr = DMPlexGetConeOrientation(K, p, &orientation);CHKERRQ(ierr); 194*da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pStart,&uOff);CHKERRQ(ierr); 195*da43764aSToby Isaac ierr = PetscSectionGetOffset(unionConeSection,uOff,&cOff);CHKERRQ(ierr); 196*da43764aSToby Isaac for (c = 0; c < dof; c++) { 197*da43764aSToby Isaac PetscInt e, eOff; 198*da43764aSToby Isaac e = cone[c]; 199*da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, e - pStart, &eOff);CHKERRQ(ierr); 200*da43764aSToby Isaac unionCones[cOff + c] = eOff; 201*da43764aSToby Isaac unionOrientations[cOff + c] = orientation[c]; 202*da43764aSToby Isaac } 203*da43764aSToby Isaac } 204*da43764aSToby Isaac for (p = pRefStart; p < pRefEnd; p++) { 205*da43764aSToby Isaac PetscInt dof, uDof, uOff, c, cOff; 206*da43764aSToby Isaac const PetscInt *cone, *orientation; 207*da43764aSToby Isaac 208*da43764aSToby Isaac ierr = DMPlexGetConeSize(Kref, p, &dof);CHKERRQ(ierr); 209*da43764aSToby Isaac ierr = DMPlexGetCone(Kref, p, &cone);CHKERRQ(ierr); 210*da43764aSToby Isaac ierr = DMPlexGetConeOrientation(Kref, p, &orientation);CHKERRQ(ierr); 211*da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection, p - pRefStart + (pEnd - pStart),&uDof);CHKERRQ(ierr); 212*da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, p - pRefStart + (pEnd - pStart),&uOff);CHKERRQ(ierr); 213*da43764aSToby Isaac if (uDof) { 214*da43764aSToby Isaac ierr = PetscSectionGetOffset(unionConeSection,uOff,&cOff);CHKERRQ(ierr); 215*da43764aSToby Isaac for (c = 0; c < dof; c++) { 216*da43764aSToby Isaac PetscInt e, eOff, eDof; 217*da43764aSToby Isaac 218*da43764aSToby Isaac e = cone[c]; 219*da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection, e - pRefStart + (pEnd - pStart),&eDof);CHKERRQ(ierr); 220*da43764aSToby Isaac if (eDof) { 221*da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, e - pRefStart + (pEnd - pStart), &eOff);CHKERRQ(ierr); 222*da43764aSToby Isaac } 223*da43764aSToby Isaac else { 224*da43764aSToby Isaac ierr = DMLabelGetValue(identityRef, e, &e);CHKERRQ(ierr); 225*da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection, e - pStart, &eOff);CHKERRQ(ierr); 226*da43764aSToby Isaac } 227*da43764aSToby Isaac unionCones[cOff + c] = eOff; 228*da43764aSToby Isaac unionOrientations[cOff + c] = orientation[c]; 229*da43764aSToby Isaac } 230*da43764aSToby Isaac } 231*da43764aSToby Isaac } 232*da43764aSToby Isaac /* get the coordinates */ 233*da43764aSToby Isaac { 234*da43764aSToby Isaac PetscInt vStart, vEnd, vRefStart, vRefEnd, v, vDof, vOff; 235*da43764aSToby Isaac PetscSection KcoordsSec, KrefCoordsSec; 236*da43764aSToby Isaac Vec KcoordsVec, KrefCoordsVec; 237*da43764aSToby Isaac PetscScalar *Kcoords; 238*da43764aSToby Isaac 239*da43764aSToby Isaac DMGetCoordinateSection(K, &KcoordsSec);CHKERRQ(ierr); 240*da43764aSToby Isaac DMGetCoordinatesLocal(K, &KcoordsVec);CHKERRQ(ierr); 241*da43764aSToby Isaac DMGetCoordinateSection(Kref, &KrefCoordsSec);CHKERRQ(ierr); 242*da43764aSToby Isaac DMGetCoordinatesLocal(Kref, &KrefCoordsVec);CHKERRQ(ierr); 243*da43764aSToby Isaac 244*da43764aSToby Isaac numVerts = numDimPoints[0]; 245*da43764aSToby Isaac ierr = PetscMalloc1(numVerts * dim,&unionCoords);CHKERRQ(ierr); 246*da43764aSToby Isaac ierr = DMPlexGetDepthStratum(K,0,&vStart,&vEnd);CHKERRQ(ierr); 247*da43764aSToby Isaac 248*da43764aSToby Isaac offset = 0; 249*da43764aSToby Isaac for (v = vStart; v < vEnd; v++) { 250*da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,v - pStart,&vOff);CHKERRQ(ierr); 251*da43764aSToby Isaac ierr = VecGetValuesSection(KcoordsVec, KcoordsSec, v, &Kcoords);CHKERRQ(ierr); 252*da43764aSToby Isaac for (d = 0; d < dim; d++) { 253*da43764aSToby Isaac unionCoords[offset * dim + d] = Kcoords[d]; 254*da43764aSToby Isaac } 255*da43764aSToby Isaac offset++; 256*da43764aSToby Isaac } 257*da43764aSToby Isaac ierr = DMPlexGetDepthStratum(Kref,0,&vRefStart,&vRefEnd);CHKERRQ(ierr); 258*da43764aSToby Isaac for (v = vRefStart; v < vRefEnd; v++) { 259*da43764aSToby Isaac ierr = PetscSectionGetDof(unionSection,v - pRefStart + (pEnd - pStart),&vDof);CHKERRQ(ierr); 260*da43764aSToby Isaac ierr = PetscSectionGetOffset(unionSection,v - pRefStart + (pEnd - pStart),&vOff);CHKERRQ(ierr); 261*da43764aSToby Isaac ierr = VecGetValuesSection(KrefCoordsVec, KrefCoordsSec, v, &Kcoords);CHKERRQ(ierr); 262*da43764aSToby Isaac if (vDof) { 263*da43764aSToby Isaac for (d = 0; d < dim; d++) { 264*da43764aSToby Isaac unionCoords[offset * dim + d] = Kcoords[d]; 265*da43764aSToby Isaac } 266*da43764aSToby Isaac offset++; 267*da43764aSToby Isaac } 268*da43764aSToby Isaac } 269*da43764aSToby Isaac } 270*da43764aSToby Isaac ierr = DMCreate(comm,ref);CHKERRQ(ierr); 271*da43764aSToby Isaac ierr = DMSetType(*ref,DMPLEX);CHKERRQ(ierr); 272*da43764aSToby Isaac ierr = DMPlexSetDimension(*ref,dim);CHKERRQ(ierr); 273*da43764aSToby Isaac ierr = DMPlexCreateFromDAG(*ref,dim,numDimPoints,coneSizes,unionCones,unionOrientations,unionCoords);CHKERRQ(ierr); 274*da43764aSToby Isaac /* clean up */ 275*da43764aSToby Isaac ierr = PetscSectionDestroy(&unionSection);CHKERRQ(ierr); 276*da43764aSToby Isaac ierr = PetscSectionDestroy(&unionConeSection);CHKERRQ(ierr); 277*da43764aSToby Isaac ierr = ISDestroy(&perm);CHKERRQ(ierr); 278*da43764aSToby Isaac ierr = PetscFree(unionCoords);CHKERRQ(ierr); 279*da43764aSToby Isaac ierr = PetscFree2(unionCones,unionOrientations);CHKERRQ(ierr); 280*da43764aSToby Isaac ierr = PetscFree2(coneSizes,numDimPoints);CHKERRQ(ierr); 281*da43764aSToby Isaac ierr = DMDestroy(&K);CHKERRQ(ierr); 282*da43764aSToby Isaac ierr = DMDestroy(&Kref);CHKERRQ(ierr); 283*da43764aSToby Isaac PetscFunctionReturn(0); 284*da43764aSToby Isaac } 285*da43764aSToby Isaac 286