xref: /petsc/src/dm/impls/plex/plexadapt.c (revision 4190e864e78a3951651949aa7edac2b157a361bb)
10bbe5d31Sbarral #include <petsc/private/dmpleximpl.h>   /*I      "petscdmplex.h"   I*/
20bbe5d31Sbarral #ifdef PETSC_HAVE_PRAGMATIC
30bbe5d31Sbarral #include <pragmatic/cpragmatic.h>
40bbe5d31Sbarral #endif
50bbe5d31Sbarral 
60bbe5d31Sbarral 
70bbe5d31Sbarral 
80bbe5d31Sbarral #undef __FUNCT__
90bbe5d31Sbarral #define __FUNCT__ "DMPlexAdaptdm->coarseMesh"
10*4190e864Sbarral PetscErrorCode DMPlexAdapt(DM dm, Vec metric, const char bdyLabelName[], DM *dmCoarsened)
110bbe5d31Sbarral {
120bbe5d31Sbarral #ifdef PETSC_HAVE_PRAGMATIC
130bbe5d31Sbarral   DM             udm, coordDM;
14*4190e864Sbarral   DMLabel        bd, bdtags;
150bbe5d31Sbarral   Vec            coordinates;
160bbe5d31Sbarral   PetscSection   coordSection;
170bbe5d31Sbarral   const PetscScalar *coords;
180bbe5d31Sbarral   double        *coarseCoords;
190bbe5d31Sbarral   IS             bdIS;
200bbe5d31Sbarral   PetscReal     *x, *y, *z;
210bbe5d31Sbarral   const PetscInt *faces;
220bbe5d31Sbarral   PetscInt      *cells, *ccells, *bdFaces, *bdFaceIds;
230bbe5d31Sbarral   PetscInt       dim, numCorners, cStart, cEnd, numCells, numCoarseCells, c, vStart, vEnd, numVertices, numCoarseVertices, v, numBdFaces, f, maxConeSize, bdSize, coff;
240bbe5d31Sbarral   const PetscScalar   *metricArray;
250bbe5d31Sbarral   PetscReal     *met;
260bbe5d31Sbarral #endif
270bbe5d31Sbarral   PetscErrorCode ierr;
280bbe5d31Sbarral   MPI_Comm       comm;
290bbe5d31Sbarral 
300bbe5d31Sbarral   PetscFunctionBegin;
310bbe5d31Sbarral   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
320bbe5d31Sbarral #ifdef PETSC_HAVE_PRAGMATIC
330bbe5d31Sbarral   if (!dm->coarseMesh) {
340bbe5d31Sbarral     ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
350bbe5d31Sbarral     ierr = DMGetCoordinateDM(dm, &coordDM);CHKERRQ(ierr);
360bbe5d31Sbarral     ierr = DMGetDefaultSection(coordDM, &coordSection);CHKERRQ(ierr);
370bbe5d31Sbarral     ierr = DMGetCoordinatesLocal(dm, &coordinates);CHKERRQ(ierr);
380bbe5d31Sbarral     ierr = DMPlexGetHeightStratum(dm, 0, &cStart, &cEnd);CHKERRQ(ierr);
390bbe5d31Sbarral     ierr = DMPlexGetDepthStratum(dm, 0, &vStart, &vEnd);CHKERRQ(ierr);
400bbe5d31Sbarral     ierr = DMPlexUninterpolate(dm, &udm);CHKERRQ(ierr);
410bbe5d31Sbarral     ierr = DMPlexGetMaxSizes(udm, &maxConeSize, NULL);CHKERRQ(ierr);
420bbe5d31Sbarral     numCells    = cEnd - cStart;
430bbe5d31Sbarral     numVertices = vEnd - vStart;
440bbe5d31Sbarral     ierr = PetscCalloc5(numVertices, &x, numVertices, &y, numVertices, &z, numCells*maxConeSize, &cells, dim*dim*numVertices, &met);CHKERRQ(ierr);
450bbe5d31Sbarral     ierr = VecGetArrayRead(coordinates, &coords);CHKERRQ(ierr);
460bbe5d31Sbarral     ierr = VecGetArrayRead(metric, &metricArray);CHKERRQ(ierr);
470bbe5d31Sbarral     for (v = vStart; v < vEnd; ++v) {
480bbe5d31Sbarral       PetscInt off;
490bbe5d31Sbarral 
500bbe5d31Sbarral       ierr = PetscSectionGetOffset(coordSection, v, &off);CHKERRQ(ierr);
510bbe5d31Sbarral       x[v-vStart] = coords[off+0];
520bbe5d31Sbarral       y[v-vStart] = coords[off+1];
530bbe5d31Sbarral       if (dim > 2) z[v-vStart] = coords[off+2];
540bbe5d31Sbarral 
550bbe5d31Sbarral       ierr = PetscMemcpy(&met[dim*dim*(v-vStart)], &metricArray[dim*off], dim*dim*sizeof(PetscScalar));CHKERRQ(ierr);
560bbe5d31Sbarral     }
570bbe5d31Sbarral     ierr = VecRestoreArrayRead(coordinates, &coords);CHKERRQ(ierr);
580bbe5d31Sbarral     ierr = VecRestoreArrayRead(metric, &metricArray);CHKERRQ(ierr);
590bbe5d31Sbarral     for (c = 0, coff = 0; c < numCells; ++c) {
600bbe5d31Sbarral       const PetscInt *cone;
610bbe5d31Sbarral       PetscInt        coneSize, cl;
620bbe5d31Sbarral 
630bbe5d31Sbarral       ierr = DMPlexGetConeSize(udm, c, &coneSize);CHKERRQ(ierr);
640bbe5d31Sbarral       ierr = DMPlexGetCone(udm, c, &cone);CHKERRQ(ierr);
650bbe5d31Sbarral       for (cl = 0; cl < coneSize; ++cl) cells[coff++] = cone[cl] - vStart;
660bbe5d31Sbarral     }
670bbe5d31Sbarral     switch (dim) {
680bbe5d31Sbarral     case 2:
690bbe5d31Sbarral       pragmatic_2d_init(&numVertices, &numCells, cells, x, y);
700bbe5d31Sbarral       break;
710bbe5d31Sbarral     case 3:
720bbe5d31Sbarral       pragmatic_3d_init(&numVertices, &numCells, cells, x, y, z);
730bbe5d31Sbarral       break;
740bbe5d31Sbarral     default:
750bbe5d31Sbarral       SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic adaptation defined for dimension %d", dim);
760bbe5d31Sbarral     }
770bbe5d31Sbarral     /* Create boundary mesh */
78*4190e864Sbarral     if ( bdyLabelName ) {
79*4190e864Sbarral       ierr = DMGetLabel(dm, bdyLabelName, &bdtags);CHKERRQ(ierr);
80*4190e864Sbarral     }
81*4190e864Sbarral     // TODO fix if bdyLabelName = "boundary"
82*4190e864Sbarral     // TODO a way to use only bdyLabelName would be to loop over its strata, if I can't get all the strata at once
830bbe5d31Sbarral     ierr = DMLabelCreate("boundary", &bd);CHKERRQ(ierr);
840bbe5d31Sbarral     ierr = DMPlexMarkBoundaryFaces(dm, bd);CHKERRQ(ierr);
850bbe5d31Sbarral     ierr = DMLabelGetStratumIS(bd, 1, &bdIS);CHKERRQ(ierr);
860bbe5d31Sbarral     ierr = DMLabelGetStratumSize(bd, 1, &numBdFaces);CHKERRQ(ierr);
870bbe5d31Sbarral     ierr = ISGetIndices(bdIS, &faces);CHKERRQ(ierr);
880bbe5d31Sbarral     for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
890bbe5d31Sbarral       PetscInt *closure = NULL;
900bbe5d31Sbarral       PetscInt  closureSize, cl;
910bbe5d31Sbarral 
920bbe5d31Sbarral       ierr = DMPlexGetTransitiveClosure(dm, faces[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
930bbe5d31Sbarral       for (cl = 0; cl < closureSize*2; cl += 2) {
940bbe5d31Sbarral         if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) ++bdSize;
950bbe5d31Sbarral       }
960bbe5d31Sbarral       ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
970bbe5d31Sbarral     }
980bbe5d31Sbarral     ierr = PetscMalloc2(bdSize, &bdFaces, numBdFaces, &bdFaceIds);CHKERRQ(ierr);
990bbe5d31Sbarral     for (f = 0, bdSize = 0; f < numBdFaces; ++f) {
1000bbe5d31Sbarral       PetscInt *closure = NULL;
1010bbe5d31Sbarral       PetscInt  closureSize, cl;
1020bbe5d31Sbarral 
1030bbe5d31Sbarral       ierr = DMPlexGetTransitiveClosure(dm, faces[f], PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1040bbe5d31Sbarral       for (cl = 0; cl < closureSize*2; cl += 2) {
1050bbe5d31Sbarral         if ((closure[cl] >= vStart) && (closure[cl] < vEnd)) bdFaces[bdSize++] = closure[cl] - vStart;
1060bbe5d31Sbarral       }
107*4190e864Sbarral       if ( bdyLabelName )
108*4190e864Sbarral         DMLabelGetValue(bdtags, faces[f], &bdFaceIds[f]);
109*4190e864Sbarral       else
1100bbe5d31Sbarral         bdFaceIds[f] = 1;
1110bbe5d31Sbarral       ierr = DMPlexRestoreTransitiveClosure(dm, f, PETSC_TRUE, &closureSize, &closure);CHKERRQ(ierr);
1120bbe5d31Sbarral     }
1130bbe5d31Sbarral     ierr = ISDestroy(&bdIS);CHKERRQ(ierr);
1140bbe5d31Sbarral     ierr = DMLabelDestroy(&bd);CHKERRQ(ierr);
1150bbe5d31Sbarral     pragmatic_set_boundary(&numBdFaces, bdFaces, bdFaceIds);
1160bbe5d31Sbarral     pragmatic_set_metric(met);
1170bbe5d31Sbarral     pragmatic_adapt();
1180bbe5d31Sbarral     /* Read out mesh */
1190bbe5d31Sbarral     pragmatic_get_info(&numCoarseVertices, &numCoarseCells);
1200bbe5d31Sbarral     ierr = PetscMalloc1(numCoarseVertices*dim, &coarseCoords);CHKERRQ(ierr);
1210bbe5d31Sbarral     switch (dim) {
1220bbe5d31Sbarral     case 2:
1230bbe5d31Sbarral       ierr = PetscMalloc2(numCoarseVertices, &x, numCoarseVertices, &y);CHKERRQ(ierr);
1240bbe5d31Sbarral       pragmatic_get_coords_2d(x, y);
1250bbe5d31Sbarral       numCorners = 3;
1260bbe5d31Sbarral       for (v = 0; v < numCoarseVertices; ++v) {coarseCoords[v*2+0] = x[v]; coarseCoords[v*2+1] = y[v];}
1270bbe5d31Sbarral       break;
1280bbe5d31Sbarral     case 3:
1290bbe5d31Sbarral       ierr = PetscMalloc3(numCoarseVertices, &x, numCoarseVertices, &y, numCoarseVertices, &z);CHKERRQ(ierr);
1300bbe5d31Sbarral       pragmatic_get_coords_3d(x, y, z);
1310bbe5d31Sbarral       numCorners = 4;
1320bbe5d31Sbarral       for (v = 0; v < numCoarseVertices; ++v) {coarseCoords[v*3+0] = x[v]; coarseCoords[v*3+1] = y[v]; coarseCoords[v*3+2] = z[v];}
1330bbe5d31Sbarral       break;
1340bbe5d31Sbarral     default:
1350bbe5d31Sbarral       SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_ARG_OUTOFRANGE, "No Pragmatic coarsening defined for dimension %d", dim);
1360bbe5d31Sbarral     }
1370bbe5d31Sbarral     ierr = PetscMalloc1(numCoarseCells*(dim+1), &ccells);CHKERRQ(ierr); // only for simplicial meshes
1380bbe5d31Sbarral     pragmatic_get_elements(ccells);
1390bbe5d31Sbarral     /* TODO Read out markers for boundary */
1400bbe5d31Sbarral     ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject) dm), dim, numCoarseCells, numCoarseVertices, numCorners, PETSC_TRUE, ccells, dim, coarseCoords, &dm->coarseMesh);CHKERRQ(ierr);
1410bbe5d31Sbarral     pragmatic_finalize();
1420bbe5d31Sbarral     ierr = PetscFree4(x, y, z, cells);CHKERRQ(ierr);
1430bbe5d31Sbarral     ierr = PetscFree2(bdFaces, bdFaceIds);CHKERRQ(ierr);
1440bbe5d31Sbarral     ierr = PetscFree(coarseCoords);CHKERRQ(ierr);
1450bbe5d31Sbarral   }
1460bbe5d31Sbarral #endif
1470bbe5d31Sbarral   ierr = PetscObjectReference((PetscObject) dm->coarseMesh);CHKERRQ(ierr);
1480bbe5d31Sbarral   *dmCoarsened = dm->coarseMesh;
1490bbe5d31Sbarral   PetscFunctionReturn(0);
1500bbe5d31Sbarral }
151