1*919ab0a2SMatthew G. Knepley #include <petsc/private/dmpleximpl.h> /*I "petscdmplex.h" I*/ 2*919ab0a2SMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC 3*919ab0a2SMatthew G. Knepley #include <pragmatic/cpragmatic.h> 4*919ab0a2SMatthew G. Knepley #endif 5*919ab0a2SMatthew G. Knepley 6*919ab0a2SMatthew G. Knepley #undef __FUNCT__ 7*919ab0a2SMatthew G. Knepley #define __FUNCT__ "DMCoarsen_Plex" 8*919ab0a2SMatthew G. Knepley PetscErrorCode DMCoarsen_Plex(DM dm, MPI_Comm comm, DM *dmCoarsened) 9*919ab0a2SMatthew G. Knepley { 10*919ab0a2SMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 11*919ab0a2SMatthew G. Knepley DM udm, coordDM; 12*919ab0a2SMatthew G. Knepley DMLabel bd; 13*919ab0a2SMatthew G. Knepley Vec coordinates; 14*919ab0a2SMatthew G. Knepley PetscSection coordSection; 15*919ab0a2SMatthew G. Knepley const PetscScalar *coords; 16*919ab0a2SMatthew G. Knepley double *coarseCoords; 17*919ab0a2SMatthew G. Knepley IS bdIS; 18*919ab0a2SMatthew G. Knepley PetscReal *x, *y, *z, *metric; 19*919ab0a2SMatthew G. Knepley const PetscInt *faces; 20*919ab0a2SMatthew G. Knepley PetscInt *cells, *bdFaces, *bdFaceIds; 21*919ab0a2SMatthew G. Knepley PetscInt dim, numCorners, cStart, cEnd, numCells, numCoarseCells, c, vStart, vEnd, numVertices, numCoarseVertices, v, numBdFaces, numCoarseBdFaces, f, maxConeSize, bdSize, coff; 22*919ab0a2SMatthew G. Knepley PetscErrorCode ierr; 23*919ab0a2SMatthew G. Knepley 24*919ab0a2SMatthew G. Knepley PetscFunctionBegin; 25*919ab0a2SMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC 26*919ab0a2SMatthew G. Knepley if (!mesh->coarseMesh) { 27*919ab0a2SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 28*919ab0a2SMatthew G. Knepley ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr); 29*919ab0a2SMatthew G. Knepley ierr = DMGetDefaultSection(coordDM, &coordSection);CHKERRQ(ierr); 30*919ab0a2SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 31*919ab0a2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 32*919ab0a2SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 33*919ab0a2SMatthew G. Knepley ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr); 34*919ab0a2SMatthew G. Knepley ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr); 35*919ab0a2SMatthew G. Knepley numCells = cEnd - cStart; 36*919ab0a2SMatthew G. Knepley numVertices = vEnd - vStart; 37*919ab0a2SMatthew G. Knepley ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices, &metric, numCells*maxConeSize, &cells);CHKERRQ(ierr); 38*919ab0a2SMatthew G. Knepley ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 39*919ab0a2SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 40*919ab0a2SMatthew G. Knepley PetscInt off; 41*919ab0a2SMatthew G. Knepley 42*919ab0a2SMatthew G. Knepley ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 43*919ab0a2SMatthew G. Knepley x[v] = coords[off+0]; 44*919ab0a2SMatthew G. Knepley y[v] = coords[off+1]; 45*919ab0a2SMatthew G. Knepley if (dim > 2) z[v] = coords[off+2]; 46*919ab0a2SMatthew G. Knepley } 47*919ab0a2SMatthew G. Knepley ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 48*919ab0a2SMatthew G. Knepley for (c = 0, coff = 0; c < numCells; ++c) { 49*919ab0a2SMatthew G. Knepley const PetscInt *cone; 50*919ab0a2SMatthew G. Knepley PetscInt coneSize, cl; 51*919ab0a2SMatthew G. Knepley 52*919ab0a2SMatthew G. Knepley ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr); 53*919ab0a2SMatthew G. Knepley ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr); 54*919ab0a2SMatthew G. Knepley for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl]; 55*919ab0a2SMatthew G. Knepley } 56*919ab0a2SMatthew G. Knepley switch (dim) { 57*919ab0a2SMatthew G. Knepley case 2: 58*919ab0a2SMatthew G. Knepley pragmatic_2d_init(&numVertices, &numCells, cells, x, y); 59*919ab0a2SMatthew G. Knepley break; 60*919ab0a2SMatthew G. Knepley case 3: 61*919ab0a2SMatthew G. Knepley pragmatic_3d_init(&numVertices, &numCells, cells, x, y, z); 62*919ab0a2SMatthew G. Knepley break; 63*919ab0a2SMatthew G. Knepley default: 64*919ab0a2SMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic coarsening defined for dimension %d", dim); 65*919ab0a2SMatthew G. Knepley } 66*919ab0a2SMatthew G. Knepley ierr = DMDestroy(&udm);CHKERRQ(ierr); 67*919ab0a2SMatthew G. Knepley /* Create boundary mesh */ 68*919ab0a2SMatthew G. Knepley ierr = DMLabelCreate("boundary", &bd);CHKERRQ(ierr); 69*919ab0a2SMatthew G. Knepley ierr = DMPlexMarkBoundaryFaces(dm, bd);CHKERRQ(ierr); 70*919ab0a2SMatthew G. Knepley ierr = DMLabelGetStratumIS(bd, 1, &bdIS);CHKERRQ(ierr); 71*919ab0a2SMatthew G. Knepley ierr = DMLabelGetStratumSize(bd, 1, &numBdFaces);CHKERRQ(ierr); 72*919ab0a2SMatthew G. Knepley ierr = ISGetIndices(bdIS, &faces);CHKERRQ(ierr); 73*919ab0a2SMatthew G. Knepley for (f = 0, bdSize = 0; f < numBdFaces; ++f) { 74*919ab0a2SMatthew G. Knepley PetscInt *closure = NULL; 75*919ab0a2SMatthew G. Knepley PetscInt closureSize, cl; 76*919ab0a2SMatthew G. Knepley 77*919ab0a2SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 78*919ab0a2SMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 79*919ab0a2SMatthew G. Knepley if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize; 80*919ab0a2SMatthew G. Knepley } 81*919ab0a2SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 82*919ab0a2SMatthew G. Knepley } 83*919ab0a2SMatthew G. Knepley ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr); 84*919ab0a2SMatthew G. Knepley for (f = 0, bdSize = 0; f < numBdFaces; ++f) { 85*919ab0a2SMatthew G. Knepley PetscInt *closure = NULL; 86*919ab0a2SMatthew G. Knepley PetscInt closureSize, cl; 87*919ab0a2SMatthew G. Knepley 88*919ab0a2SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 89*919ab0a2SMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 90*919ab0a2SMatthew G. Knepley if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl]; 91*919ab0a2SMatthew G. Knepley } 92*919ab0a2SMatthew G. Knepley /* TODO Fix */ 93*919ab0a2SMatthew G. Knepley bdFaceIds[f] = 1; 94*919ab0a2SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 95*919ab0a2SMatthew G. Knepley } 96*919ab0a2SMatthew G. Knepley ierr = ISDestroy(&bdIS);CHKERRQ(ierr); 97*919ab0a2SMatthew G. Knepley pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds); 98*919ab0a2SMatthew G. Knepley /* Create metric */ 99*919ab0a2SMatthew G. Knepley if (dim == 2) {for (v = 0; v < numVertices; ++v) metric[v*4+0] = metric[v*4+3] = 1.0;} 100*919ab0a2SMatthew G. Knepley else {for (v = 0; v < numVertices; ++v) metric[v*9+0] = metric[v*9+4] = metric[v*9+8] = 1.0;} 101*919ab0a2SMatthew G. Knepley pragmatic_set_metric(metric); 102*919ab0a2SMatthew G. Knepley pragmatic_adapt(); 103*919ab0a2SMatthew G. Knepley /* Read out mesh */ 104*919ab0a2SMatthew G. Knepley pragmatic_get_info(&numCoarseVertices, &numCoarseCells); 105*919ab0a2SMatthew G. Knepley switch (dim) { 106*919ab0a2SMatthew G. Knepley case 2: 107*919ab0a2SMatthew G. Knepley pragmatic_get_coords_2d(x, y); 108*919ab0a2SMatthew G. Knepley numCorners = 3; 109*919ab0a2SMatthew G. Knepley for (v = 0; v < numCoarseVertices; ++v) {coarseCoords[v*2+0] = x[v]; coarseCoords[v*2+1] = y[v];} 110*919ab0a2SMatthew G. Knepley break; 111*919ab0a2SMatthew G. Knepley case 3: 112*919ab0a2SMatthew G. Knepley pragmatic_get_coords_3d(x, y, z); 113*919ab0a2SMatthew G. Knepley numCorners = 4; 114*919ab0a2SMatthew 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];} 115*919ab0a2SMatthew G. Knepley break; 116*919ab0a2SMatthew G. Knepley default: 117*919ab0a2SMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic coarsening defined for dimension %d", dim); 118*919ab0a2SMatthew G. Knepley } 119*919ab0a2SMatthew G. Knepley pragmatic_get_elements(cells); 120*919ab0a2SMatthew G. Knepley /* TODO Read out markers for boundary */ 121*919ab0a2SMatthew G. Knepley ierr = PetscMalloc1(numCoarseVertices*dim, &coarseCoords);CHKERRQ(ierr); 122*919ab0a2SMatthew G. Knepley ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject) dm), dim, numCoarseCells, numCoarseVertices, numCorners, PETSC_TRUE, cells, dim, coarseCoords, &mesh->coarseMesh);CHKERRQ(ierr); 123*919ab0a2SMatthew G. Knepley pragmatic_finalize(); 124*919ab0a2SMatthew G. Knepley ierr = PetscFree4(x, y, z, metric);CHKERRQ(ierr); 125*919ab0a2SMatthew G. Knepley ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr); 126*919ab0a2SMatthew G. Knepley ierr = PetscFree(coarseCoords);CHKERRQ(ierr); 127*919ab0a2SMatthew G. Knepley } 128*919ab0a2SMatthew G. Knepley #endif 129*919ab0a2SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) mesh->coarseMesh);CHKERRQ(ierr); 130*919ab0a2SMatthew G. Knepley *dmCoarsened = mesh->coarseMesh; 131*919ab0a2SMatthew G. Knepley PetscFunctionReturn(0); 132*919ab0a2SMatthew G. Knepley } 133