xref: /petsc/src/dm/impls/plex/plexcoarsen.c (revision a8fb8f2974ee66bb003afa65900fab794a3643ec)
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 {
1217f047d8SMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC
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;
2577902bb5SPatrick Farrell   PetscInt       dim, numCorners, cStart, cEnd, numCells, numCoarseCells, c, vStart, vEnd, numVertices, numCoarseVertices, v, numBdFaces, f, maxConeSize, size, bdSize, coff;
2617f047d8SMatthew G. Knepley #endif
27919ab0a2SMatthew G. Knepley   PetscErrorCode ierr;
28919ab0a2SMatthew G. Knepley 
29919ab0a2SMatthew G. Knepley   PetscFunctionBegin;
30919ab0a2SMatthew G. Knepley #ifdef PETSC_HAVE_PRAGMATIC
31*a8fb8f29SToby Isaac   if (!dm->coarseMesh) {
32919ab0a2SMatthew G. Knepley     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
33919ab0a2SMatthew G. Knepley     ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr);
34919ab0a2SMatthew G. Knepley     ierr = DMGetDefaultSection(coordDM, &coordSection);CHKERRQ(ierr);
35919ab0a2SMatthew G. Knepley     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
36919ab0a2SMatthew G. Knepley     ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
37919ab0a2SMatthew G. Knepley     ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
38919ab0a2SMatthew G. Knepley     ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
39919ab0a2SMatthew G. Knepley     ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr);
40919ab0a2SMatthew G. Knepley     numCells    = cEnd - cStart;
41919ab0a2SMatthew G. Knepley     numVertices = vEnd - vStart;
42d2d4c474SMatthew G. Knepley     ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numVertices*PetscSqr(dim), &metric, numCells*maxConeSize, &cells);CHKERRQ(ierr);
43919ab0a2SMatthew G. Knepley     ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
44d2d4c474SMatthew G. Knepley     for (v = vStart; v < vEnd; ++v) {
45919ab0a2SMatthew G. Knepley       PetscInt off;
46919ab0a2SMatthew G. Knepley 
47919ab0a2SMatthew G. Knepley       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
48d2d4c474SMatthew G. Knepley       x[v-vStart] = coords[off+0];
49d2d4c474SMatthew G. Knepley       y[v-vStart] = coords[off+1];
50d2d4c474SMatthew G. Knepley       if (dim > 2) z[v-vStart] = coords[off+2];
51919ab0a2SMatthew G. Knepley     }
52919ab0a2SMatthew G. Knepley     ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
53919ab0a2SMatthew G. Knepley     for (c = 0, coff = 0; c < numCells; ++c) {
54919ab0a2SMatthew G. Knepley       const PetscInt *cone;
55919ab0a2SMatthew G. Knepley       PetscInt        coneSize, cl;
56919ab0a2SMatthew G. Knepley 
57919ab0a2SMatthew G. Knepley       ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
58919ab0a2SMatthew G. Knepley       ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
59d2d4c474SMatthew G. Knepley       for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart;
60919ab0a2SMatthew G. Knepley     }
61919ab0a2SMatthew G. Knepley     switch (dim) {
62919ab0a2SMatthew G. Knepley     case 2:
63919ab0a2SMatthew G. Knepley       pragmatic_2d_init(&numVertices, &numCells, cells, x, y);
64919ab0a2SMatthew G. Knepley       break;
65919ab0a2SMatthew G. Knepley     case 3:
66919ab0a2SMatthew G. Knepley       pragmatic_3d_init(&numVertices, &numCells, cells, x, y, z);
67919ab0a2SMatthew G. Knepley       break;
68919ab0a2SMatthew G. Knepley     default:
69919ab0a2SMatthew G. Knepley       SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic coarsening defined for dimension %d", dim);
70919ab0a2SMatthew G. Knepley     }
71919ab0a2SMatthew G. Knepley     /* Create boundary mesh */
72919ab0a2SMatthew G. Knepley     ierr = DMLabelCreate("boundary", &bd);CHKERRQ(ierr);
73919ab0a2SMatthew G. Knepley     ierr = DMPlexMarkBoundaryFaces(dm, bd);CHKERRQ(ierr);
74919ab0a2SMatthew G. Knepley     ierr = DMLabelGetStratumIS(bd, 1, &bdIS);CHKERRQ(ierr);
75919ab0a2SMatthew G. Knepley     ierr = DMLabelGetStratumSize(bd, 1, &numBdFaces);CHKERRQ(ierr);
76919ab0a2SMatthew G. Knepley     ierr = ISGetIndices(bdIS, &faces);CHKERRQ(ierr);
77919ab0a2SMatthew G. Knepley     for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
78919ab0a2SMatthew G. Knepley       PetscInt *closure = NULL;
79919ab0a2SMatthew G. Knepley       PetscInt  closureSize, cl;
80919ab0a2SMatthew G. Knepley 
81d2d4c474SMatthew G. Knepley       ierr = DMPlexGetTransitiveClosure(dm, faces[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
82919ab0a2SMatthew G. Knepley       for (cl = 0; cl < closureSize*2; cl += 2) {
83919ab0a2SMatthew G. Knepley         if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
84919ab0a2SMatthew G. Knepley       }
85919ab0a2SMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
86919ab0a2SMatthew G. Knepley     }
87919ab0a2SMatthew G. Knepley     ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr);
88919ab0a2SMatthew G. Knepley     for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
89919ab0a2SMatthew G. Knepley       PetscInt *closure = NULL;
90919ab0a2SMatthew G. Knepley       PetscInt  closureSize, cl;
91919ab0a2SMatthew G. Knepley 
92d2d4c474SMatthew G. Knepley       ierr = DMPlexGetTransitiveClosure(dm, faces[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
93919ab0a2SMatthew G. Knepley       for (cl = 0; cl < closureSize*2; cl += 2) {
94d2d4c474SMatthew G. Knepley         if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart;
95919ab0a2SMatthew G. Knepley       }
96919ab0a2SMatthew G. Knepley       /* TODO Fix */
97919ab0a2SMatthew G. Knepley       bdFaceIds[f] = 1;
98919ab0a2SMatthew G. Knepley       ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
99919ab0a2SMatthew G. Knepley     }
100919ab0a2SMatthew G. Knepley     ierr = ISDestroy(&bdIS);CHKERRQ(ierr);
101d2d4c474SMatthew G. Knepley     ierr = DMLabelDestroy(&bd);CHKERRQ(ierr);
102919ab0a2SMatthew G. Knepley     pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds);
103919ab0a2SMatthew G. Knepley     /* Create metric */
1049f3102b2SMatthew G. Knepley     size = (dim*(dim+1))/2;
1059f3102b2SMatthew G. Knepley     ierr = PetscMalloc1(PetscSqr(size), &eqns);CHKERRQ(ierr);
1069f3102b2SMatthew G. Knepley     ierr = MatCreateSeqDense(PETSC_COMM_SELF, size, size, eqns, &A);CHKERRQ(ierr);
1079f3102b2SMatthew G. Knepley     ierr = MatCreateVecs(A, &mx, &mb);CHKERRQ(ierr);
1089f3102b2SMatthew G. Knepley     ierr = VecSet(mb, 1.0);CHKERRQ(ierr);
1099f3102b2SMatthew G. Knepley     for (c = 0; c < numCells; ++c) {
1109f3102b2SMatthew G. Knepley       const PetscScalar *sol;
1119f3102b2SMatthew G. Knepley       PetscScalar       *cellCoords = NULL;
1129f3102b2SMatthew G. Knepley       PetscReal          e[3], vol;
1139f3102b2SMatthew G. Knepley       const PetscInt    *cone;
1149f3102b2SMatthew G. Knepley       PetscInt           coneSize, cl, i, j, d, r;
1159f3102b2SMatthew G. Knepley 
1169f3102b2SMatthew G. Knepley       ierr = DMPlexVecGetClosure(dm, coordSection, coordinates, c, NULL, &cellCoords);CHKERRQ(ierr);
1179f3102b2SMatthew G. Knepley       /* Only works for simplices */
1189f3102b2SMatthew G. Knepley       for (i = 0, r = 0; i < dim+1; ++i) {
1199f3102b2SMatthew G. Knepley         for (j = 0; j < i; ++j, ++r) {
1209f3102b2SMatthew G. Knepley           for (d = 0; d < dim; ++d) e[d] = cellCoords[i*dim+d] - cellCoords[j*dim+d];
1219f3102b2SMatthew G. Knepley           /* FORTRAN ORDERING */
1229f3102b2SMatthew G. Knepley           if (dim == 2) {
1239f3102b2SMatthew G. Knepley             eqns[0*size+r] = PetscSqr(e[0]);
1249f3102b2SMatthew G. Knepley             eqns[1*size+r] = 2.0*e[0]*e[1];
1259f3102b2SMatthew G. Knepley             eqns[2*size+r] = PetscSqr(e[1]);
1269f3102b2SMatthew G. Knepley           } else {
1279f3102b2SMatthew G. Knepley             eqns[0*size+r] = PetscSqr(e[0]);
1289f3102b2SMatthew G. Knepley             eqns[1*size+r] = 2.0*e[0]*e[1];
1299f3102b2SMatthew G. Knepley             eqns[2*size+r] = 2.0*e[0]*e[2];
1309f3102b2SMatthew G. Knepley             eqns[3*size+r] = PetscSqr(e[1]);
1319f3102b2SMatthew G. Knepley             eqns[4*size+r] = 2.0*e[1]*e[2];
1329f3102b2SMatthew G. Knepley             eqns[5*size+r] = PetscSqr(e[2]);
1339f3102b2SMatthew G. Knepley           }
1349f3102b2SMatthew G. Knepley         }
1359f3102b2SMatthew G. Knepley       }
1369f3102b2SMatthew G. Knepley       ierr = MatSetUnfactored(A);CHKERRQ(ierr);
1379f3102b2SMatthew G. Knepley       ierr = DMPlexVecRestoreClosure(dm, coordSection, coordinates, c, NULL, &cellCoords);CHKERRQ(ierr);
1389f3102b2SMatthew G. Knepley       ierr = MatLUFactor(A, NULL, NULL, NULL);CHKERRQ(ierr);
1399f3102b2SMatthew G. Knepley       ierr = MatSolve(A, mb, mx);CHKERRQ(ierr);
1409f3102b2SMatthew G. Knepley       ierr = VecGetArrayRead(mx, &sol);CHKERRQ(ierr);
1419f3102b2SMatthew G. Knepley       ierr = DMPlexComputeCellGeometryFVM(dm, c, &vol, NULL, NULL);CHKERRQ(ierr);
1429f3102b2SMatthew G. Knepley       ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
1439f3102b2SMatthew G. Knepley       ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
1449f3102b2SMatthew G. Knepley       for (cl = 0; cl < coneSize; ++cl) {
1459f3102b2SMatthew G. Knepley         const PetscInt v = cone[cl] - vStart;
1469f3102b2SMatthew G. Knepley 
1479f3102b2SMatthew G. Knepley         if (dim == 2) {
1489f3102b2SMatthew G. Knepley           metric[v*4+0] += vol*coarseRatio*sol[0];
1499f3102b2SMatthew G. Knepley           metric[v*4+1] += vol*coarseRatio*sol[1];
1509f3102b2SMatthew G. Knepley           metric[v*4+2] += vol*coarseRatio*sol[1];
1519f3102b2SMatthew G. Knepley           metric[v*4+3] += vol*coarseRatio*sol[2];
1529f3102b2SMatthew G. Knepley         } else {
1539f3102b2SMatthew G. Knepley           metric[v*9+0] += vol*coarseRatio*sol[0];
1549f3102b2SMatthew G. Knepley           metric[v*9+1] += vol*coarseRatio*sol[1];
1559f3102b2SMatthew G. Knepley           metric[v*9+3] += vol*coarseRatio*sol[1];
1569f3102b2SMatthew G. Knepley           metric[v*9+2] += vol*coarseRatio*sol[2];
1579f3102b2SMatthew G. Knepley           metric[v*9+6] += vol*coarseRatio*sol[2];
1589f3102b2SMatthew G. Knepley           metric[v*9+4] += vol*coarseRatio*sol[3];
1599f3102b2SMatthew G. Knepley           metric[v*9+5] += vol*coarseRatio*sol[4];
1609f3102b2SMatthew G. Knepley           metric[v*9+7] += vol*coarseRatio*sol[4];
1619f3102b2SMatthew G. Knepley           metric[v*9+8] += vol*coarseRatio*sol[5];
1629f3102b2SMatthew G. Knepley         }
1639f3102b2SMatthew G. Knepley       }
1649f3102b2SMatthew G. Knepley       ierr = VecRestoreArrayRead(mx, &sol);CHKERRQ(ierr);
1659f3102b2SMatthew G. Knepley     }
1669f3102b2SMatthew G. Knepley     for (v = 0; v < numVertices; ++v) {
1679f3102b2SMatthew G. Knepley       const PetscInt *support;
1689f3102b2SMatthew G. Knepley       PetscInt        supportSize, s;
1699f3102b2SMatthew G. Knepley       PetscReal       vol, totVol = 0.0;
1709f3102b2SMatthew G. Knepley 
1719f3102b2SMatthew G. Knepley       ierr = DMPlexGetSupport(udm, v+vStart, &support);CHKERRQ(ierr);
1729f3102b2SMatthew G. Knepley       ierr = DMPlexGetSupportSize(udm, v+vStart, &supportSize);CHKERRQ(ierr);
1739f3102b2SMatthew G. Knepley       for (s = 0; s < supportSize; ++s) {ierr = DMPlexComputeCellGeometryFVM(dm, support[s], &vol, NULL, NULL);CHKERRQ(ierr); totVol += vol;}
1749f3102b2SMatthew G. Knepley       for (s = 0; s < PetscSqr(dim); ++s) metric[v*PetscSqr(dim)+s] /= totVol;
1759f3102b2SMatthew G. Knepley     }
1769f3102b2SMatthew G. Knepley     ierr = VecDestroy(&mx);CHKERRQ(ierr);
1779f3102b2SMatthew G. Knepley     ierr = VecDestroy(&mb);CHKERRQ(ierr);
1789f3102b2SMatthew G. Knepley     ierr = MatDestroy(&A);CHKERRQ(ierr);
1799f3102b2SMatthew G. Knepley     ierr = DMDestroy(&udm);CHKERRQ(ierr);
1809f3102b2SMatthew G. Knepley     ierr = PetscFree(eqns);CHKERRQ(ierr);
181919ab0a2SMatthew G. Knepley     pragmatic_set_metric(metric);
182919ab0a2SMatthew G. Knepley     pragmatic_adapt();
183919ab0a2SMatthew G. Knepley     /* Read out mesh */
184919ab0a2SMatthew G. Knepley     pragmatic_get_info(&numCoarseVertices, &numCoarseCells);
185d2d4c474SMatthew G. Knepley     ierr = PetscMalloc1(numCoarseVertices*dim, &coarseCoords);CHKERRQ(ierr);
186919ab0a2SMatthew G. Knepley     switch (dim) {
187919ab0a2SMatthew G. Knepley     case 2:
188919ab0a2SMatthew G. Knepley       pragmatic_get_coords_2d(x, y);
189919ab0a2SMatthew G. Knepley       numCorners = 3;
190919ab0a2SMatthew G. Knepley       for (v = 0; v < numCoarseVertices; ++v) {coarseCoords[v*2+0] = x[v]; coarseCoords[v*2+1] = y[v];}
191919ab0a2SMatthew G. Knepley       break;
192919ab0a2SMatthew G. Knepley     case 3:
193919ab0a2SMatthew G. Knepley       pragmatic_get_coords_3d(x, y, z);
194919ab0a2SMatthew G. Knepley       numCorners = 4;
195919ab0a2SMatthew 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];}
196919ab0a2SMatthew G. Knepley       break;
197919ab0a2SMatthew G. Knepley     default:
198919ab0a2SMatthew G. Knepley       SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic coarsening defined for dimension %d", dim);
199919ab0a2SMatthew G. Knepley     }
200919ab0a2SMatthew G. Knepley     pragmatic_get_elements(cells);
201919ab0a2SMatthew G. Knepley     /* TODO Read out markers for boundary */
202*a8fb8f29SToby Isaac     ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject) dm), dim, numCoarseCells, numCoarseVertices, numCorners, PETSC_TRUE, cells, dim, coarseCoords, &dm->coarseMesh);CHKERRQ(ierr);
203919ab0a2SMatthew G. Knepley     pragmatic_finalize();
204d2d4c474SMatthew G. Knepley     ierr = PetscFree5(x, y, z, metric, cells);CHKERRQ(ierr);
205919ab0a2SMatthew G. Knepley     ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr);
206919ab0a2SMatthew G. Knepley     ierr = PetscFree(coarseCoords);CHKERRQ(ierr);
207919ab0a2SMatthew G. Knepley   }
208919ab0a2SMatthew G. Knepley #endif
209*a8fb8f29SToby Isaac   ierr = PetscObjectReference((PetscObject) dm->coarseMesh);CHKERRQ(ierr);
210*a8fb8f29SToby Isaac   *dmCoarsened = dm->coarseMesh;
211919ab0a2SMatthew G. Knepley   PetscFunctionReturn(0);
212919ab0a2SMatthew G. Knepley }
213b653a561SMatthew G. Knepley 
214b653a561SMatthew G. Knepley #undef __FUNCT__
215b653a561SMatthew G. Knepley #define __FUNCT__ "DMCoarsenHierarchy_Plex"
216b653a561SMatthew G. Knepley PetscErrorCode DMCoarsenHierarchy_Plex(DM dm, PetscInt nlevels, DM dmCoarsened[])
217b653a561SMatthew G. Knepley {
218b653a561SMatthew G. Knepley   DM             rdm = dm;
219b653a561SMatthew G. Knepley   PetscInt       c;
220b653a561SMatthew G. Knepley   PetscErrorCode ierr;
221b653a561SMatthew G. Knepley 
222b653a561SMatthew G. Knepley   PetscFunctionBegin;
223b653a561SMatthew G. Knepley   for (c = nlevels-1; c >= 0; --c) {
224b653a561SMatthew G. Knepley     ierr = DMCoarsen(rdm, PetscObjectComm((PetscObject) dm), &dmCoarsened[c]);CHKERRQ(ierr);
225b653a561SMatthew G. Knepley     ierr = DMPlexCopyBoundary(rdm, dmCoarsened[c]);CHKERRQ(ierr);
226*a8fb8f29SToby Isaac     ierr = DMSetCoarseDM(rdm, dmCoarsened[c]);CHKERRQ(ierr);
227b653a561SMatthew G. Knepley     rdm  = dmCoarsened[c];
228b653a561SMatthew G. Knepley   }
229b653a561SMatthew G. Knepley   PetscFunctionReturn(0);
230b653a561SMatthew G. Knepley }
231