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 69f3102b2SMatthew G. Knepley #include <petscksp.h> 79f3102b2SMatthew G. Knepley 8919ab0a2SMatthew G. Knepley #undef __FUNCT__ 9919ab0a2SMatthew G. Knepley #define __FUNCT__ "DMCoarsen_Plex" 10919ab0a2SMatthew G. Knepley PetscErrorCode DMCoarsen_Plex(DM dm, MPI_Comm comm, DM *dmCoarsened) 11919ab0a2SMatthew G. Knepley { 12919ab0a2SMatthew G. Knepley DM_Plex *mesh = (DM_Plex*) dm->data; 13919ab0a2SMatthew G. Knepley DM udm, coordDM; 14919ab0a2SMatthew G. Knepley DMLabel bd; 159f3102b2SMatthew G. Knepley Mat A; 169f3102b2SMatthew G. Knepley Vec coordinates, mb, mx; 17919ab0a2SMatthew G. Knepley PetscSection coordSection; 18919ab0a2SMatthew G. Knepley const PetscScalar *coords; 19919ab0a2SMatthew G. Knepley double *coarseCoords; 20919ab0a2SMatthew G. Knepley IS bdIS; 219f3102b2SMatthew G. Knepley PetscReal *x, *y, *z, *eqns, *metric; 229f3102b2SMatthew G. Knepley PetscReal coarseRatio = PetscSqr(0.5); 23919ab0a2SMatthew G. Knepley const PetscInt *faces; 24919ab0a2SMatthew G. Knepley PetscInt *cells, *bdFaces, *bdFaceIds; 259f3102b2SMatthew G. Knepley PetscInt dim, numCorners, cStart, cEnd, numCells, numCoarseCells, c, vStart, vEnd, numVertices, numCoarseVertices, v, numBdFaces, numCoarseBdFaces, f, maxConeSize, size, bdSize, coff; 26919ab0a2SMatthew G. Knepley PetscErrorCode ierr; 27919ab0a2SMatthew G. Knepley 28919ab0a2SMatthew G. Knepley PetscFunctionBegin; 29919ab0a2SMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC 30919ab0a2SMatthew G. Knepley if (!mesh->coarseMesh) { 31919ab0a2SMatthew G. Knepley ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr); 32919ab0a2SMatthew G. Knepley ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr); 33919ab0a2SMatthew G. Knepley ierr = DMGetDefaultSection(coordDM, &coordSection);CHKERRQ(ierr); 34919ab0a2SMatthew G. Knepley ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr); 35919ab0a2SMatthew G. Knepley ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr); 36919ab0a2SMatthew G. Knepley ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr); 37919ab0a2SMatthew G. Knepley ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr); 38919ab0a2SMatthew G. Knepley ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr); 39919ab0a2SMatthew G. Knepley numCells = cEnd - cStart; 40919ab0a2SMatthew G. Knepley numVertices = vEnd - vStart; 41d2d4c474SMatthew G. Knepley ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices*PetscSqr(dim), &metric, numCells*maxConeSize, &cells);CHKERRQ(ierr); 42919ab0a2SMatthew G. Knepley ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr); 43d2d4c474SMatthew G. Knepley for (v = vStart; v < vEnd; ++v) { 44919ab0a2SMatthew G. Knepley PetscInt off; 45919ab0a2SMatthew G. Knepley 46919ab0a2SMatthew G. Knepley ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr); 47d2d4c474SMatthew G. Knepley x[v-vStart] = coords[off+0]; 48d2d4c474SMatthew G. Knepley y[v-vStart] = coords[off+1]; 49d2d4c474SMatthew G. Knepley if (dim > 2) z[v-vStart] = coords[off+2]; 50919ab0a2SMatthew G. Knepley } 51919ab0a2SMatthew G. Knepley ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr); 52919ab0a2SMatthew G. Knepley for (c = 0, coff = 0; c < numCells; ++c) { 53919ab0a2SMatthew G. Knepley const PetscInt *cone; 54919ab0a2SMatthew G. Knepley PetscInt coneSize, cl; 55919ab0a2SMatthew G. Knepley 56919ab0a2SMatthew G. Knepley ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr); 57919ab0a2SMatthew G. Knepley ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr); 58d2d4c474SMatthew G. Knepley for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart; 59919ab0a2SMatthew G. Knepley } 60919ab0a2SMatthew G. Knepley switch (dim) { 61919ab0a2SMatthew G. Knepley case 2: 62919ab0a2SMatthew G. Knepley pragmatic_2d_init(&numVertices, &numCells, cells, x, y); 63919ab0a2SMatthew G. Knepley break; 64919ab0a2SMatthew G. Knepley case 3: 65919ab0a2SMatthew G. Knepley pragmatic_3d_init(&numVertices, &numCells, cells, x, y, z); 66919ab0a2SMatthew G. Knepley break; 67919ab0a2SMatthew G. Knepley default: 68919ab0a2SMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic coarsening defined for dimension %d", dim); 69919ab0a2SMatthew G. Knepley } 70919ab0a2SMatthew G. Knepley /* Create boundary mesh */ 71919ab0a2SMatthew G. Knepley ierr = DMLabelCreate("boundary", &bd);CHKERRQ(ierr); 72919ab0a2SMatthew G. Knepley ierr = DMPlexMarkBoundaryFaces(dm, bd);CHKERRQ(ierr); 73919ab0a2SMatthew G. Knepley ierr = DMLabelGetStratumIS(bd, 1, &bdIS);CHKERRQ(ierr); 74919ab0a2SMatthew G. Knepley ierr = DMLabelGetStratumSize(bd, 1, &numBdFaces);CHKERRQ(ierr); 75919ab0a2SMatthew G. Knepley ierr = ISGetIndices(bdIS, &faces);CHKERRQ(ierr); 76919ab0a2SMatthew G. Knepley for (f = 0, bdSize = 0; f < numBdFaces; ++f) { 77919ab0a2SMatthew G. Knepley PetscInt *closure = NULL; 78919ab0a2SMatthew G. Knepley PetscInt closureSize, cl; 79919ab0a2SMatthew G. Knepley 80d2d4c474SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, faces[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 81919ab0a2SMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 82919ab0a2SMatthew G. Knepley if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize; 83919ab0a2SMatthew G. Knepley } 84919ab0a2SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 85919ab0a2SMatthew G. Knepley } 86919ab0a2SMatthew G. Knepley ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr); 87919ab0a2SMatthew G. Knepley for (f = 0, bdSize = 0; f < numBdFaces; ++f) { 88919ab0a2SMatthew G. Knepley PetscInt *closure = NULL; 89919ab0a2SMatthew G. Knepley PetscInt closureSize, cl; 90919ab0a2SMatthew G. Knepley 91d2d4c474SMatthew G. Knepley ierr = DMPlexGetTransitiveClosure(dm, faces[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 92919ab0a2SMatthew G. Knepley for (cl = 0; cl < closureSize*2; cl += 2) { 93d2d4c474SMatthew G. Knepley if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart; 94919ab0a2SMatthew G. Knepley } 95919ab0a2SMatthew G. Knepley /* TODO Fix */ 96919ab0a2SMatthew G. Knepley bdFaceIds[f] = 1; 97919ab0a2SMatthew G. Knepley ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr); 98919ab0a2SMatthew G. Knepley } 99919ab0a2SMatthew G. Knepley ierr = ISDestroy(&bdIS);CHKERRQ(ierr); 100d2d4c474SMatthew G. Knepley ierr = DMLabelDestroy(&bd);CHKERRQ(ierr); 101919ab0a2SMatthew G. Knepley pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds); 102919ab0a2SMatthew G. Knepley /* Create metric */ 1039f3102b2SMatthew G. Knepley size = (dim*(dim+1))/2; 1049f3102b2SMatthew G. Knepley ierr = PetscMalloc1(PetscSqr(size), &eqns);CHKERRQ(ierr); 1059f3102b2SMatthew G. Knepley ierr = MatCreateSeqDense(PETSC_COMM_SELF, size, size, eqns, &A);CHKERRQ(ierr); 1069f3102b2SMatthew G. Knepley ierr = MatCreateVecs(A, &mx, &mb);CHKERRQ(ierr); 1079f3102b2SMatthew G. Knepley ierr = VecSet(mb, 1.0);CHKERRQ(ierr); 1089f3102b2SMatthew G. Knepley for (c = 0; c < numCells; ++c) { 1099f3102b2SMatthew G. Knepley const PetscScalar *sol; 1109f3102b2SMatthew G. Knepley PetscScalar *cellCoords = NULL; 1119f3102b2SMatthew G. Knepley PetscReal e[3], vol; 1129f3102b2SMatthew G. Knepley const PetscInt *cone; 1139f3102b2SMatthew G. Knepley PetscInt coneSize, cl, i, j, d, r; 1149f3102b2SMatthew G. Knepley 1159f3102b2SMatthew G. Knepley ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, c, NULL, &cellCoords);CHKERRQ(ierr); 1169f3102b2SMatthew G. Knepley /* Only works for simplices */ 1179f3102b2SMatthew G. Knepley for (i = 0, r = 0; i < dim+1; ++i) { 1189f3102b2SMatthew G. Knepley for (j = 0; j < i; ++j, ++r) { 1199f3102b2SMatthew G. Knepley for (d = 0; d < dim; ++d) e[d] = cellCoords[i*dim+d] - cellCoords[j*dim+d]; 1209f3102b2SMatthew G. Knepley /* FORTRAN ORDERING */ 1219f3102b2SMatthew G. Knepley if (dim == 2) { 1229f3102b2SMatthew G. Knepley eqns[0*size+r] = PetscSqr(e[0]); 1239f3102b2SMatthew G. Knepley eqns[1*size+r] = 2.0*e[0]*e[1]; 1249f3102b2SMatthew G. Knepley eqns[2*size+r] = PetscSqr(e[1]); 1259f3102b2SMatthew G. Knepley } else { 1269f3102b2SMatthew G. Knepley eqns[0*size+r] = PetscSqr(e[0]); 1279f3102b2SMatthew G. Knepley eqns[1*size+r] = 2.0*e[0]*e[1]; 1289f3102b2SMatthew G. Knepley eqns[2*size+r] = 2.0*e[0]*e[2]; 1299f3102b2SMatthew G. Knepley eqns[3*size+r] = PetscSqr(e[1]); 1309f3102b2SMatthew G. Knepley eqns[4*size+r] = 2.0*e[1]*e[2]; 1319f3102b2SMatthew G. Knepley eqns[5*size+r] = PetscSqr(e[2]); 1329f3102b2SMatthew G. Knepley } 1339f3102b2SMatthew G. Knepley } 1349f3102b2SMatthew G. Knepley } 1359f3102b2SMatthew G. Knepley ierr = MatSetUnfactored(A);CHKERRQ(ierr); 1369f3102b2SMatthew G. Knepley ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, c, NULL, &cellCoords);CHKERRQ(ierr); 1379f3102b2SMatthew G. Knepley ierr = MatLUFactor(A, NULL, NULL, NULL);CHKERRQ(ierr); 1389f3102b2SMatthew G. Knepley ierr = MatSolve(A, mb, mx);CHKERRQ(ierr); 1399f3102b2SMatthew G. Knepley ierr = VecGetArrayRead(mx, &sol);CHKERRQ(ierr); 1409f3102b2SMatthew G. Knepley ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL);CHKERRQ(ierr); 1419f3102b2SMatthew G. Knepley ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr); 1429f3102b2SMatthew G. Knepley ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr); 1439f3102b2SMatthew G. Knepley for (cl = 0; cl < coneSize; ++cl) { 1449f3102b2SMatthew G. Knepley const PetscInt v = cone[cl] - vStart; 1459f3102b2SMatthew G. Knepley 1469f3102b2SMatthew G. Knepley if (dim == 2) { 1479f3102b2SMatthew G. Knepley metric[v*4+0] += vol*coarseRatio*sol[0]; 1489f3102b2SMatthew G. Knepley metric[v*4+1] += vol*coarseRatio*sol[1]; 1499f3102b2SMatthew G. Knepley metric[v*4+2] += vol*coarseRatio*sol[1]; 1509f3102b2SMatthew G. Knepley metric[v*4+3] += vol*coarseRatio*sol[2]; 1519f3102b2SMatthew G. Knepley } else { 1529f3102b2SMatthew G. Knepley metric[v*9+0] += vol*coarseRatio*sol[0]; 1539f3102b2SMatthew G. Knepley metric[v*9+1] += vol*coarseRatio*sol[1]; 1549f3102b2SMatthew G. Knepley metric[v*9+3] += vol*coarseRatio*sol[1]; 1559f3102b2SMatthew G. Knepley metric[v*9+2] += vol*coarseRatio*sol[2]; 1569f3102b2SMatthew G. Knepley metric[v*9+6] += vol*coarseRatio*sol[2]; 1579f3102b2SMatthew G. Knepley metric[v*9+4] += vol*coarseRatio*sol[3]; 1589f3102b2SMatthew G. Knepley metric[v*9+5] += vol*coarseRatio*sol[4]; 1599f3102b2SMatthew G. Knepley metric[v*9+7] += vol*coarseRatio*sol[4]; 1609f3102b2SMatthew G. Knepley metric[v*9+8] += vol*coarseRatio*sol[5]; 1619f3102b2SMatthew G. Knepley } 1629f3102b2SMatthew G. Knepley } 1639f3102b2SMatthew G. Knepley ierr = VecRestoreArrayRead(mx, &sol);CHKERRQ(ierr); 1649f3102b2SMatthew G. Knepley } 1659f3102b2SMatthew G. Knepley for (v = 0; v < numVertices; ++v) { 1669f3102b2SMatthew G. Knepley const PetscInt *support; 1679f3102b2SMatthew G. Knepley PetscInt supportSize, s; 1689f3102b2SMatthew G. Knepley PetscReal vol, totVol = 0.0; 1699f3102b2SMatthew G. Knepley 1709f3102b2SMatthew G. Knepley ierr = DMPlexGetSupport(udm, v+vStart, &support);CHKERRQ(ierr); 1719f3102b2SMatthew G. Knepley ierr = DMPlexGetSupportSize(udm, v+vStart, &supportSize);CHKERRQ(ierr); 1729f3102b2SMatthew G. Knepley for (s = 0; s < supportSize; ++s) {ierr = DMPlexComputeCellGeometryFVM(dm, support[s], &vol, NULL, NULL);CHKERRQ(ierr); totVol += vol;} 1739f3102b2SMatthew G. Knepley for (s = 0; s < PetscSqr(dim); ++s) metric[v*PetscSqr(dim)+s] /= totVol; 1749f3102b2SMatthew G. Knepley } 1759f3102b2SMatthew G. Knepley ierr = VecDestroy(&mx);CHKERRQ(ierr); 1769f3102b2SMatthew G. Knepley ierr = VecDestroy(&mb);CHKERRQ(ierr); 1779f3102b2SMatthew G. Knepley ierr = MatDestroy(&A);CHKERRQ(ierr); 1789f3102b2SMatthew G. Knepley ierr = DMDestroy(&udm);CHKERRQ(ierr); 1799f3102b2SMatthew G. Knepley ierr = PetscFree(eqns);CHKERRQ(ierr); 180919ab0a2SMatthew G. Knepley pragmatic_set_metric(metric); 181919ab0a2SMatthew G. Knepley pragmatic_adapt(); 182919ab0a2SMatthew G. Knepley /* Read out mesh */ 183919ab0a2SMatthew G. Knepley pragmatic_get_info(&numCoarseVertices, &numCoarseCells); 184d2d4c474SMatthew G. Knepley ierr = PetscMalloc1(numCoarseVertices*dim, &coarseCoords);CHKERRQ(ierr); 185919ab0a2SMatthew G. Knepley switch (dim) { 186919ab0a2SMatthew G. Knepley case 2: 187919ab0a2SMatthew G. Knepley pragmatic_get_coords_2d(x, y); 188919ab0a2SMatthew G. Knepley numCorners = 3; 189919ab0a2SMatthew G. Knepley for (v = 0; v < numCoarseVertices; ++v) {coarseCoords[v*2+0] = x[v]; coarseCoords[v*2+1] = y[v];} 190919ab0a2SMatthew G. Knepley break; 191919ab0a2SMatthew G. Knepley case 3: 192919ab0a2SMatthew G. Knepley pragmatic_get_coords_3d(x, y, z); 193919ab0a2SMatthew G. Knepley numCorners = 4; 194919ab0a2SMatthew 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];} 195919ab0a2SMatthew G. Knepley break; 196919ab0a2SMatthew G. Knepley default: 197919ab0a2SMatthew G. Knepley SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic coarsening defined for dimension %d", dim); 198919ab0a2SMatthew G. Knepley } 199919ab0a2SMatthew G. Knepley pragmatic_get_elements(cells); 200919ab0a2SMatthew G. Knepley /* TODO Read out markers for boundary */ 201919ab0a2SMatthew G. Knepley ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject) dm), dim, numCoarseCells, numCoarseVertices, numCorners, PETSC_TRUE, cells, dim, coarseCoords, &mesh->coarseMesh);CHKERRQ(ierr); 202919ab0a2SMatthew G. Knepley pragmatic_finalize(); 203d2d4c474SMatthew G. Knepley ierr = PetscFree5(x, y, z, metric, cells);CHKERRQ(ierr); 204919ab0a2SMatthew G. Knepley ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr); 205919ab0a2SMatthew G. Knepley ierr = PetscFree(coarseCoords);CHKERRQ(ierr); 206919ab0a2SMatthew G. Knepley } 207919ab0a2SMatthew G. Knepley #endif 208919ab0a2SMatthew G. Knepley ierr = PetscObjectReference((PetscObject) mesh->coarseMesh);CHKERRQ(ierr); 209919ab0a2SMatthew G. Knepley *dmCoarsened = mesh->coarseMesh; 210919ab0a2SMatthew G. Knepley PetscFunctionReturn(0); 211919ab0a2SMatthew G. Knepley } 212*b653a561SMatthew G. Knepley 213*b653a561SMatthew G. Knepley #undef __FUNCT__ 214*b653a561SMatthew G. Knepley #define __FUNCT__ "DMCoarsenHierarchy_Plex" 215*b653a561SMatthew G. Knepley PetscErrorCode DMCoarsenHierarchy_Plex(DM dm, PetscInt nlevels, DM dmCoarsened[]) 216*b653a561SMatthew G. Knepley { 217*b653a561SMatthew G. Knepley DM rdm = dm; 218*b653a561SMatthew G. Knepley PetscInt c; 219*b653a561SMatthew G. Knepley PetscErrorCode ierr; 220*b653a561SMatthew G. Knepley 221*b653a561SMatthew G. Knepley PetscFunctionBegin; 222*b653a561SMatthew G. Knepley for (c = nlevels-1; c >= 0; --c) { 223*b653a561SMatthew G. Knepley CellRefiner cellRefiner; 224*b653a561SMatthew G. Knepley 225*b653a561SMatthew G. Knepley ierr = DMCoarsen(rdm, PetscObjectComm((PetscObject) dm), &dmCoarsened[c]);CHKERRQ(ierr); 226*b653a561SMatthew G. Knepley ierr = DMPlexCopyBoundary(rdm, dmCoarsened[c]);CHKERRQ(ierr); 227*b653a561SMatthew G. Knepley ierr = DMPlexSetCoarseDM(rdm, dmCoarsened[c]);CHKERRQ(ierr); 228*b653a561SMatthew G. Knepley rdm = dmCoarsened[c]; 229*b653a561SMatthew G. Knepley } 230*b653a561SMatthew G. Knepley PetscFunctionReturn(0); 231*b653a561SMatthew G. Knepley } 232