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