125538a83SMatthew G. Knepley static char help[] = "Tests for ephemeral meshes.\n"; 225538a83SMatthew G. Knepley 325538a83SMatthew G. Knepley #include <petscdmplex.h> 425538a83SMatthew G. Knepley #include <petscdmplextransform.h> 5*0fcd6f6bSMatthew G. Knepley #include <petscdmlabelephemeral.h> 625538a83SMatthew G. Knepley 725538a83SMatthew G. Knepley /* 825538a83SMatthew G. Knepley Use 925538a83SMatthew G. Knepley 1025538a83SMatthew G. Knepley -ref_dm_view -eph_dm_view 1125538a83SMatthew G. Knepley 1225538a83SMatthew G. Knepley to view the concrete and ephemeral meshes from the first transformation, and 1325538a83SMatthew G. Knepley 1425538a83SMatthew G. Knepley -ref_dm_sec_view -eph_dm_sec_view 1525538a83SMatthew G. Knepley 1625538a83SMatthew G. Knepley for the second. 1725538a83SMatthew G. Knepley */ 1825538a83SMatthew G. Knepley 1925538a83SMatthew G. Knepley // Should remove when I have an API for everything 2025538a83SMatthew G. Knepley #include <petsc/private/dmplextransformimpl.h> 2125538a83SMatthew G. Knepley 2225538a83SMatthew G. Knepley typedef struct { 2325538a83SMatthew G. Knepley DMLabel active; /* Label for transform */ 2425538a83SMatthew G. Knepley PetscBool second; /* Flag to execute a second transformation */ 2525538a83SMatthew G. Knepley PetscBool concrete; /* Flag to use the concrete mesh for the second transformation */ 2625538a83SMatthew G. Knepley } AppCtx; 2725538a83SMatthew G. Knepley 2825538a83SMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 2925538a83SMatthew G. Knepley { 3025538a83SMatthew G. Knepley PetscInt cells[1024], Nc = 1024; 3125538a83SMatthew G. Knepley PetscBool flg; 3225538a83SMatthew G. Knepley 3325538a83SMatthew G. Knepley PetscFunctionBeginUser; 3425538a83SMatthew G. Knepley options->active = NULL; 3525538a83SMatthew G. Knepley options->second = PETSC_FALSE; 3625538a83SMatthew G. Knepley options->concrete = PETSC_FALSE; 3725538a83SMatthew G. Knepley 3825538a83SMatthew G. Knepley PetscOptionsBegin(comm, "", "Ephemeral Meshing Options", "DMPLEX"); 3925538a83SMatthew G. Knepley PetscCall(PetscOptionsIntArray("-cells", "Cells to mark for transformation", "ex57.c", cells, &Nc, &flg)); 4025538a83SMatthew G. Knepley if (flg) { 4125538a83SMatthew G. Knepley PetscCall(DMLabelCreate(comm, "active", &options->active)); 4225538a83SMatthew G. Knepley for (PetscInt c = 0; c < Nc; ++c) PetscCall(DMLabelSetValue(options->active, cells[c], DM_ADAPT_REFINE)); 4325538a83SMatthew G. Knepley } 4425538a83SMatthew G. Knepley PetscCall(PetscOptionsBool("-second", "Use a second transformation", "ex57.c", options->second, &options->second, &flg)); 4525538a83SMatthew G. Knepley PetscCall(PetscOptionsBool("-concrete", "Use concrete mesh for the second transformation", "ex57.c", options->concrete, &options->concrete, &flg)); 4625538a83SMatthew G. Knepley PetscOptionsEnd(); 4725538a83SMatthew G. Knepley PetscFunctionReturn(0); 4825538a83SMatthew G. Knepley } 4925538a83SMatthew G. Knepley 5025538a83SMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm) 5125538a83SMatthew G. Knepley { 5225538a83SMatthew G. Knepley PetscFunctionBeginUser; 5325538a83SMatthew G. Knepley PetscCall(DMCreate(comm, dm)); 5425538a83SMatthew G. Knepley PetscCall(DMSetType(*dm, DMPLEX)); 5525538a83SMatthew G. Knepley PetscCall(DMSetFromOptions(*dm)); 5625538a83SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)*dm, "Mesh")); 5725538a83SMatthew G. Knepley PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 5825538a83SMatthew G. Knepley PetscFunctionReturn(0); 5925538a83SMatthew G. Knepley } 6025538a83SMatthew G. Knepley 6125538a83SMatthew G. Knepley static PetscErrorCode CreateTransform(DM dm, DMLabel active, const char prefix[], DMPlexTransform *tr) 6225538a83SMatthew G. Knepley { 6325538a83SMatthew G. Knepley PetscFunctionBeginUser; 6425538a83SMatthew G. Knepley PetscCall(DMPlexTransformCreate(PetscObjectComm((PetscObject)dm), tr)); 6525538a83SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)*tr, "Transform")); 6625538a83SMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*tr, prefix)); 6725538a83SMatthew G. Knepley PetscCall(DMPlexTransformSetDM(*tr, dm)); 6825538a83SMatthew G. Knepley PetscCall(DMPlexTransformSetActive(*tr, active)); 6925538a83SMatthew G. Knepley PetscCall(DMPlexTransformSetFromOptions(*tr)); 7025538a83SMatthew G. Knepley PetscCall(DMPlexTransformSetUp(*tr)); 7125538a83SMatthew G. Knepley 7225538a83SMatthew G. Knepley PetscCall(DMSetApplicationContext(dm, *tr)); 7325538a83SMatthew G. Knepley PetscCall(PetscObjectViewFromOptions((PetscObject)*tr, NULL, "-dm_plex_transform_view")); 7425538a83SMatthew G. Knepley PetscFunctionReturn(0); 7525538a83SMatthew G. Knepley } 7625538a83SMatthew G. Knepley 7725538a83SMatthew G. Knepley static PetscErrorCode CreateEphemeralMesh(DMPlexTransform tr, DM *tdm) 7825538a83SMatthew G. Knepley { 7925538a83SMatthew G. Knepley PetscFunctionBegin; 80*0fcd6f6bSMatthew G. Knepley PetscCall(DMPlexCreateEphemeral(tr, tdm)); 8125538a83SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)*tdm, "Ephemeral Mesh")); 8225538a83SMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*tdm, "eph_")); 8325538a83SMatthew G. Knepley PetscFunctionReturn(0); 8425538a83SMatthew G. Knepley } 8525538a83SMatthew G. Knepley 8625538a83SMatthew G. Knepley static PetscErrorCode CreateConcreteMesh(DMPlexTransform tr, DM *rdm) 8725538a83SMatthew G. Knepley { 8825538a83SMatthew G. Knepley DM cdm, codm, rcodm; 8925538a83SMatthew G. Knepley 9025538a83SMatthew G. Knepley PetscFunctionBegin; 9125538a83SMatthew G. Knepley PetscCall(DMPlexTransformGetDM(tr, &cdm)); 9225538a83SMatthew G. Knepley PetscCall(DMPlexTransformApply(tr, cdm, rdm)); 9325538a83SMatthew G. Knepley PetscCall(DMSetCoarsenLevel(*rdm, cdm->leveldown)); 9425538a83SMatthew G. Knepley PetscCall(DMSetRefineLevel(*rdm, cdm->levelup + 1)); 9525538a83SMatthew G. Knepley PetscCall(DMCopyDisc(cdm, *rdm)); 9625538a83SMatthew G. Knepley PetscCall(DMGetCoordinateDM(cdm, &codm)); 9725538a83SMatthew G. Knepley PetscCall(DMGetCoordinateDM(*rdm, &rcodm)); 9825538a83SMatthew G. Knepley PetscCall(DMCopyDisc(codm, rcodm)); 9925538a83SMatthew G. Knepley PetscCall(DMPlexTransformCreateDiscLabels(tr, *rdm)); 10025538a83SMatthew G. Knepley PetscCall(DMSetCoarseDM(*rdm, cdm)); 10125538a83SMatthew G. Knepley PetscCall(DMPlexSetRegularRefinement(*rdm, PETSC_TRUE)); 10225538a83SMatthew G. Knepley if (rdm) { 10325538a83SMatthew G. Knepley ((DM_Plex *)(*rdm)->data)->printFEM = ((DM_Plex *)cdm->data)->printFEM; 10425538a83SMatthew G. Knepley ((DM_Plex *)(*rdm)->data)->printL2 = ((DM_Plex *)cdm->data)->printL2; 10525538a83SMatthew G. Knepley } 10625538a83SMatthew G. Knepley PetscCall(PetscObjectSetName((PetscObject)*rdm, "Concrete Mesh")); 10725538a83SMatthew G. Knepley PetscCall(PetscObjectSetOptionsPrefix((PetscObject)*rdm, "ref_")); 10825538a83SMatthew G. Knepley PetscFunctionReturn(0); 10925538a83SMatthew G. Knepley } 11025538a83SMatthew G. Knepley 11125538a83SMatthew G. Knepley static PetscErrorCode CompareMeshes(DM dmA, DM dmB, DM dm) 11225538a83SMatthew G. Knepley { 11325538a83SMatthew G. Knepley PetscInt dim, dimB, pStart, pEnd, pStartB, pEndB; 11425538a83SMatthew G. Knepley MPI_Comm comm; 11525538a83SMatthew G. Knepley 11625538a83SMatthew G. Knepley PetscFunctionBegin; 11725538a83SMatthew G. Knepley PetscCall(PetscObjectGetComm((PetscObject)dmA, &comm)); 11825538a83SMatthew G. Knepley PetscCall(DMGetDimension(dmA, &dim)); 11925538a83SMatthew G. Knepley PetscCall(DMGetDimension(dmB, &dimB)); 12025538a83SMatthew G. Knepley PetscCheck(dim == dimB, comm, PETSC_ERR_ARG_INCOMP, "Dimension from dmA %" PetscInt_FMT " != %" PetscInt_FMT " from dmB", dim, dimB); 12125538a83SMatthew G. Knepley PetscCall(DMPlexGetChart(dmA, &pStart, &pEnd)); 12225538a83SMatthew G. Knepley PetscCall(DMPlexGetChart(dmB, &pStartB, &pEndB)); 12325538a83SMatthew G. Knepley PetscCheck(pStart == pStartB && pEnd == pEndB, comm, PETSC_ERR_ARG_INCOMP, "Chart from dmA (%" PetscInt_FMT ", %" PetscInt_FMT ") does not match (%" PetscInt_FMT ", %" PetscInt_FMT ") for dmB", pStart, pEnd, pStartB, pEndB); 12425538a83SMatthew G. Knepley for (PetscInt p = pStart; p < pEnd; ++p) { 12525538a83SMatthew G. Knepley const PetscInt *cone, *ornt, *coneB, *orntB; 12625538a83SMatthew G. Knepley PetscInt coneSize, coneSizeB; 12725538a83SMatthew G. Knepley 12825538a83SMatthew G. Knepley PetscCall(DMPlexGetConeSize(dmA, p, &coneSize)); 12925538a83SMatthew G. Knepley PetscCall(DMPlexGetConeSize(dmB, p, &coneSizeB)); 13025538a83SMatthew G. Knepley PetscCheck(coneSize == coneSizeB, comm, PETSC_ERR_ARG_INCOMP, "Cone size for %" PetscInt_FMT " from dmA %" PetscInt_FMT " does not match %" PetscInt_FMT " for dmB", p, coneSize, coneSizeB); 13125538a83SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dmA, p, &cone, &ornt)); 13225538a83SMatthew G. Knepley PetscCall(DMPlexGetOrientedCone(dmB, p, &coneB, &orntB)); 13325538a83SMatthew G. Knepley for (PetscInt c = 0; c < coneSize; ++c) { 13425538a83SMatthew G. Knepley PetscCheck(cone[c] == coneB[c] && ornt[c] == orntB[c], comm, PETSC_ERR_ARG_INCOMP, "Cone point %" PetscInt_FMT " for point %" PetscInt_FMT " from dmA (%" PetscInt_FMT ", %" PetscInt_FMT ") does not match (%" PetscInt_FMT ", %" PetscInt_FMT ") for dmB", c, p, cone[c], ornt[c], coneB[c], orntB[c]); 13525538a83SMatthew G. Knepley } 13625538a83SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dmA, p, &cone, &ornt)); 13725538a83SMatthew G. Knepley PetscCall(DMPlexRestoreOrientedCone(dmB, p, &coneB, &orntB)); 13825538a83SMatthew G. Knepley } 13925538a83SMatthew G. Knepley PetscFunctionReturn(0); 14025538a83SMatthew G. Knepley } 14125538a83SMatthew G. Knepley 14225538a83SMatthew G. Knepley int main(int argc, char *argv[]) 14325538a83SMatthew G. Knepley { 14425538a83SMatthew G. Knepley DM dm, tdm, rdm; 14525538a83SMatthew G. Knepley DMPlexTransform tr; 14625538a83SMatthew G. Knepley AppCtx user; 14725538a83SMatthew G. Knepley MPI_Comm comm; 14825538a83SMatthew G. Knepley 14925538a83SMatthew G. Knepley PetscFunctionBeginUser; 15025538a83SMatthew G. Knepley PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 15125538a83SMatthew G. Knepley comm = PETSC_COMM_WORLD; 15225538a83SMatthew G. Knepley PetscCall(ProcessOptions(comm, &user)); 15325538a83SMatthew G. Knepley PetscCall(CreateMesh(comm, &dm)); 15425538a83SMatthew G. Knepley PetscCall(CreateTransform(dm, user.active, "first_", &tr)); 15525538a83SMatthew G. Knepley PetscCall(CreateEphemeralMesh(tr, &tdm)); 15625538a83SMatthew G. Knepley PetscCall(CreateConcreteMesh(tr, &rdm)); 15725538a83SMatthew G. Knepley if (user.second) { 15825538a83SMatthew G. Knepley DMPlexTransform tr2; 15925538a83SMatthew G. Knepley DM tdm2, rdm2; 16025538a83SMatthew G. Knepley 16125538a83SMatthew G. Knepley PetscCall(DMViewFromOptions(rdm, NULL, "-dm_sec_view")); 16225538a83SMatthew G. Knepley PetscCall(DMViewFromOptions(tdm, NULL, "-dm_sec_view")); 16325538a83SMatthew G. Knepley if (user.concrete) PetscCall(CreateTransform(rdm, user.active, "second_", &tr2)); 16425538a83SMatthew G. Knepley else PetscCall(CreateTransform(tdm, user.active, "second_", &tr2)); 16525538a83SMatthew G. Knepley PetscCall(CreateEphemeralMesh(tr2, &tdm2)); 16625538a83SMatthew G. Knepley PetscCall(CreateConcreteMesh(tr2, &rdm2)); 16725538a83SMatthew G. Knepley PetscCall(DMDestroy(&tdm)); 16825538a83SMatthew G. Knepley PetscCall(DMDestroy(&rdm)); 16925538a83SMatthew G. Knepley PetscCall(DMPlexTransformDestroy(&tr2)); 17025538a83SMatthew G. Knepley tdm = tdm2; 17125538a83SMatthew G. Knepley rdm = rdm2; 17225538a83SMatthew G. Knepley } 17325538a83SMatthew G. Knepley PetscCall(DMViewFromOptions(tdm, NULL, "-dm_view")); 17425538a83SMatthew G. Knepley PetscCall(DMViewFromOptions(rdm, NULL, "-dm_view")); 17525538a83SMatthew G. Knepley PetscCall(CompareMeshes(rdm, tdm, dm)); 17625538a83SMatthew G. Knepley PetscCall(DMPlexTransformDestroy(&tr)); 17725538a83SMatthew G. Knepley PetscCall(DMDestroy(&tdm)); 17825538a83SMatthew G. Knepley PetscCall(DMDestroy(&rdm)); 17925538a83SMatthew G. Knepley PetscCall(DMDestroy(&dm)); 18025538a83SMatthew G. Knepley PetscCall(DMLabelDestroy(&user.active)); 18125538a83SMatthew G. Knepley PetscCall(PetscFinalize()); 18225538a83SMatthew G. Knepley return 0; 18325538a83SMatthew G. Knepley } 18425538a83SMatthew G. Knepley 18525538a83SMatthew G. Knepley /*TEST 18625538a83SMatthew G. Knepley 18725538a83SMatthew G. Knepley # Tests for regular refinement of whole meshes 18825538a83SMatthew G. Knepley test: 18925538a83SMatthew G. Knepley suffix: tri 19025538a83SMatthew G. Knepley requires: triangle 19125538a83SMatthew G. Knepley args: -first_dm_plex_transform_view ::ascii_info_detail 19225538a83SMatthew G. Knepley 19325538a83SMatthew G. Knepley test: 19425538a83SMatthew G. Knepley suffix: quad 19525538a83SMatthew G. Knepley args: -dm_plex_simplex 0 19625538a83SMatthew G. Knepley 197*0fcd6f6bSMatthew G. Knepley # Here I am checking that the 'marker' label is correct for the ephemeral mesh 19825538a83SMatthew G. Knepley test: 19925538a83SMatthew G. Knepley suffix: tet 20025538a83SMatthew G. Knepley requires: ctetgen 201*0fcd6f6bSMatthew G. Knepley args: -dm_plex_dim 3 -ref_dm_view -eph_dm_view -eph_dm_plex_view_labels marker 20225538a83SMatthew G. Knepley 20325538a83SMatthew G. Knepley test: 20425538a83SMatthew G. Knepley suffix: hex 20525538a83SMatthew G. Knepley args: -dm_plex_dim 3 -dm_plex_simplex 0 20625538a83SMatthew G. Knepley 20725538a83SMatthew G. Knepley # Tests for filter patches 20825538a83SMatthew G. Knepley testset: 20925538a83SMatthew G. Knepley args: -first_dm_plex_transform_type transform_filter -ref_dm_view 21025538a83SMatthew G. Knepley 21125538a83SMatthew G. Knepley test: 21225538a83SMatthew G. Knepley suffix: tri_patch 21325538a83SMatthew G. Knepley requires: triangle 21425538a83SMatthew G. Knepley args: -cells 0,1,2,4 21525538a83SMatthew G. Knepley 21625538a83SMatthew G. Knepley test: 21725538a83SMatthew G. Knepley suffix: quad_patch 21825538a83SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_box_faces 3,3 -cells 0,1,3,4 21925538a83SMatthew G. Knepley 22025538a83SMatthew G. Knepley # Tests for refined filter patches 22125538a83SMatthew G. Knepley testset: 222*0fcd6f6bSMatthew G. Knepley args: -first_dm_plex_transform_type transform_filter -ref_dm_view -eph_dm_view -second 22325538a83SMatthew G. Knepley 22425538a83SMatthew G. Knepley test: 22525538a83SMatthew G. Knepley suffix: tri_patch_ref 22625538a83SMatthew G. Knepley requires: triangle 22725538a83SMatthew G. Knepley args: -cells 0,1,2,4 22825538a83SMatthew G. Knepley 22925538a83SMatthew G. Knepley test: 23025538a83SMatthew G. Knepley suffix: tri_patch_ref_concrete 23125538a83SMatthew G. Knepley requires: triangle 23225538a83SMatthew G. Knepley args: -cells 0,1,2,4 -concrete -first_dm_plex_transform_view ::ascii_info_detail 23325538a83SMatthew G. Knepley 23425538a83SMatthew G. Knepley # Tests for boundary layer refinement 23525538a83SMatthew G. Knepley test: 23625538a83SMatthew G. Knepley suffix: quad_bl 23725538a83SMatthew G. Knepley args: -dm_plex_simplex 0 -dm_plex_dim 1 -dm_plex_box_faces 5 -dm_extrude 2 -cells 0,2,4,6,8 \ 23825538a83SMatthew G. Knepley -first_dm_plex_transform_type refine_boundary_layer -first_dm_plex_transform_bl_splits 4 \ 23925538a83SMatthew G. Knepley -ref_dm_view 24025538a83SMatthew G. Knepley 24125538a83SMatthew G. Knepley TEST*/ 242