1919ab0a2SMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2919ab0a2SMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC 3919ab0a2SMatthew G. Knepley #include <pragmatic/cpragmatic.h> 4919ab0a2SMatthew G. Knepley #endif 5919ab0a2SMatthew G. Knepley 6919ab0a2SMatthew G. Knepley #undef __FUNCT__ 7919ab0a2SMatthew G. Knepley #define __FUNCT__ "DMCoarsen_Plex" 8919ab0a2SMatthew G. Knepley PetscErrorCode DMCoarsen_Plex(DM dm, MPI_Comm comm, DM *dmCoarsened) 9919ab0a2SMatthew G. Knepley { 10919ab0a2SMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 11919ab0a2SMatthew G. Knepley DM udm, coordDM; 12919ab0a2SMatthew G. Knepley DMLabel bd; 13919ab0a2SMatthew G. Knepley Vec coordinates; 14919ab0a2SMatthew G. Knepley PetscSection coordSection; 15919ab0a2SMatthew G. Knepley const PetscScalar *coords; 16919ab0a2SMatthew G. Knepley double *coarseCoords; 17919ab0a2SMatthew G. Knepley IS bdIS; 18919ab0a2SMatthew G. Knepley PetscReal *x, *y, *z, *metric; 19919ab0a2SMatthew G. Knepley const PetscInt *faces; 20919ab0a2SMatthew G. Knepley PetscInt *cells, *bdFaces, *bdFaceIds; 21919ab0a2SMatthew G. Knepley PetscInt dim, numCorners, cStart, cEnd, numCells, numCoarseCells, c, vStart, vEnd, numVertices, numCoarseVertices, v, numBdFaces, numCoarseBdFaces, f, maxConeSize, bdSize, coff; 22919ab0a2SMatthew G. Knepley PetscErrorCode ierr; 23919ab0a2SMatthew G. Knepley 24919ab0a2SMatthew G. Knepley PetscFunctionBegin; 25919ab0a2SMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC 26919ab0a2SMatthew G. Knepley if (!mesh->coarseMesh) { 27919ab0a2SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 28919ab0a2SMatthew G. Knepley ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr); 29919ab0a2SMatthew G. Knepley ierr = DMGetDefaultSection(coordDM, &coordSection);CHKERRQ(ierr); 30919ab0a2SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 31919ab0a2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 32919ab0a2SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 33919ab0a2SMatthew G. Knepley ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr); 34919ab0a2SMatthew G. Knepley ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr); 35919ab0a2SMatthew G. Knepley numCells = cEnd - cStart; 36919ab0a2SMatthew G. Knepley numVertices = vEnd - vStart; 37*d2d4c474SMatthew G. Knepley ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices*PetscSqr(dim), &metric, numCells*maxConeSize, &cells);CHKERRQ(ierr); 38919ab0a2SMatthew G. Knepley ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 39*d2d4c474SMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 40919ab0a2SMatthew G. Knepley PetscInt off; 41919ab0a2SMatthew G. Knepley 42919ab0a2SMatthew G. Knepley ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 43*d2d4c474SMatthew G. Knepley x[v-vStart] = coords[off+0]; 44*d2d4c474SMatthew G. Knepley y[v-vStart] = coords[off+1]; 45*d2d4c474SMatthew G. Knepley if (dim > 2) z[v-vStart] = coords[off+2]; 46919ab0a2SMatthew G. Knepley } 47919ab0a2SMatthew G. Knepley ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 48919ab0a2SMatthew G. Knepley for (c = 0, coff = 0; c < numCells; ++c) { 49919ab0a2SMatthew G. Knepley const PetscInt *cone; 50919ab0a2SMatthew G. Knepley PetscInt coneSize, cl; 51919ab0a2SMatthew G. Knepley 52919ab0a2SMatthew G. Knepley ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr); 53919ab0a2SMatthew G. Knepley ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr); 54*d2d4c474SMatthew G. Knepley for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart; 55919ab0a2SMatthew G. Knepley } 56919ab0a2SMatthew G. Knepley switch (dim) { 57919ab0a2SMatthew G. Knepley case 2: 58919ab0a2SMatthew G. Knepley pragmatic_2d_init(&numVertices, &numCells, cells, x, y); 59919ab0a2SMatthew G. Knepley break; 60919ab0a2SMatthew G. Knepley case 3: 61919ab0a2SMatthew G. Knepley pragmatic_3d_init(&numVertices, &numCells, cells, x, y, z); 62919ab0a2SMatthew G. Knepley break; 63919ab0a2SMatthew G. Knepley default: 64919ab0a2SMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic coarsening defined for dimension %d", dim); 65919ab0a2SMatthew G. Knepley } 66919ab0a2SMatthew G. Knepley ierr = DMDestroy(&udm);CHKERRQ(ierr); 67919ab0a2SMatthew G. Knepley /* Create boundary mesh */ 68919ab0a2SMatthew G. Knepley ierr = DMLabelCreate("boundary", &bd);CHKERRQ(ierr); 69919ab0a2SMatthew G. Knepley ierr = DMPlexMarkBoundaryFaces(dm, bd);CHKERRQ(ierr); 70919ab0a2SMatthew G. Knepley ierr = DMLabelGetStratumIS(bd, 1, &bdIS);CHKERRQ(ierr); 71919ab0a2SMatthew G. Knepley ierr = DMLabelGetStratumSize(bd, 1, &numBdFaces);CHKERRQ(ierr); 72919ab0a2SMatthew G. Knepley ierr = ISGetIndices(bdIS, &faces);CHKERRQ(ierr); 73919ab0a2SMatthew G. Knepley for (f = 0, bdSize = 0; f < numBdFaces; ++f) { 74919ab0a2SMatthew G. Knepley PetscInt *closure = NULL; 75919ab0a2SMatthew G. Knepley PetscInt closureSize, cl; 76919ab0a2SMatthew G. Knepley 77*d2d4c474SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, faces[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 78919ab0a2SMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 79919ab0a2SMatthew G. Knepley if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize; 80919ab0a2SMatthew G. Knepley } 81919ab0a2SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 82919ab0a2SMatthew G. Knepley } 83919ab0a2SMatthew G. Knepley ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr); 84919ab0a2SMatthew G. Knepley for (f = 0, bdSize = 0; f < numBdFaces; ++f) { 85919ab0a2SMatthew G. Knepley PetscInt *closure = NULL; 86919ab0a2SMatthew G. Knepley PetscInt closureSize, cl; 87919ab0a2SMatthew G. Knepley 88*d2d4c474SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, faces[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 89919ab0a2SMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 90*d2d4c474SMatthew G. Knepley if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart; 91919ab0a2SMatthew G. Knepley } 92919ab0a2SMatthew G. Knepley /* TODO Fix */ 93919ab0a2SMatthew G. Knepley bdFaceIds[f] = 1; 94919ab0a2SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 95919ab0a2SMatthew G. Knepley } 96919ab0a2SMatthew G. Knepley ierr = ISDestroy(&bdIS);CHKERRQ(ierr); 97*d2d4c474SMatthew G. Knepley ierr = DMLabelDestroy(&bd);CHKERRQ(ierr); 98919ab0a2SMatthew G. Knepley pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds); 99919ab0a2SMatthew G. Knepley /* Create metric */ 100*d2d4c474SMatthew G. Knepley if (dim == 2) {for (v = 0; v < numVertices; ++v) metric[v*4+0] = metric[v*4+3] = 0.5;} 101*d2d4c474SMatthew G. Knepley else {for (v = 0; v < numVertices; ++v) metric[v*9+0] = metric[v*9+4] = metric[v*9+8] = 0.5;} 102919ab0a2SMatthew G. Knepley pragmatic_set_metric(metric); 103919ab0a2SMatthew G. Knepley pragmatic_adapt(); 104919ab0a2SMatthew G. Knepley /* Read out mesh */ 105919ab0a2SMatthew G. Knepley pragmatic_get_info(&numCoarseVertices, &numCoarseCells); 106*d2d4c474SMatthew G. Knepley ierr = PetscMalloc1(numCoarseVertices*dim, &coarseCoords);CHKERRQ(ierr); 107919ab0a2SMatthew G. Knepley switch (dim) { 108919ab0a2SMatthew G. Knepley case 2: 109919ab0a2SMatthew G. Knepley pragmatic_get_coords_2d(x, y); 110919ab0a2SMatthew G. Knepley numCorners = 3; 111919ab0a2SMatthew G. Knepley for (v = 0; v < numCoarseVertices; ++v) {coarseCoords[v*2+0] = x[v]; coarseCoords[v*2+1] = y[v];} 112919ab0a2SMatthew G. Knepley break; 113919ab0a2SMatthew G. Knepley case 3: 114919ab0a2SMatthew G. Knepley pragmatic_get_coords_3d(x, y, z); 115919ab0a2SMatthew G. Knepley numCorners = 4; 116919ab0a2SMatthew G. Knepley for (v = 0; v < numCoarseVertices; ++v) {coarseCoords[v*3+0] = x[v]; coarseCoords[v*3+1] = y[v]; coarseCoords[v*3+2] = z[v];} 117919ab0a2SMatthew G. Knepley break; 118919ab0a2SMatthew G. Knepley default: 119919ab0a2SMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic coarsening defined for dimension %d", dim); 120919ab0a2SMatthew G. Knepley } 121919ab0a2SMatthew G. Knepley pragmatic_get_elements(cells); 122919ab0a2SMatthew G. Knepley /* TODO Read out markers for boundary */ 123919ab0a2SMatthew G. Knepley ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject) dm), dim, numCoarseCells, numCoarseVertices, numCorners, PETSC_TRUE, cells, dim, coarseCoords, &mesh->coarseMesh);CHKERRQ(ierr); 124919ab0a2SMatthew G. Knepley pragmatic_finalize(); 125*d2d4c474SMatthew G. Knepley ierr = PetscFree5(x, y, z, metric, cells);CHKERRQ(ierr); 126919ab0a2SMatthew G. Knepley ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr); 127919ab0a2SMatthew G. Knepley ierr = PetscFree(coarseCoords);CHKERRQ(ierr); 128919ab0a2SMatthew G. Knepley } 129919ab0a2SMatthew G. Knepley #endif 130919ab0a2SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) mesh->coarseMesh);CHKERRQ(ierr); 131919ab0a2SMatthew G. Knepley *dmCoarsened = mesh->coarseMesh; 132919ab0a2SMatthew G. Knepley PetscFunctionReturn(0); 133919ab0a2SMatthew G. Knepley } 134