1 static const char help[] = "Tests for mesh extrusion"; 2 3 #include <petscdmplex.h> 4 5 typedef struct { 6 char bdLabel[PETSC_MAX_PATH_LEN]; /* The boundary label name */ 7 PetscInt Nbd; /* The number of boundary markers to extrude, 0 for all */ 8 PetscInt bd[64]; /* The boundary markers to be extruded */ 9 } AppCtx; 10 11 PETSC_EXTERN PetscErrorCode pyramidNormal(PetscInt, PetscReal, const PetscReal[], PetscInt, PetscScalar[], void *); 12 13 /* The pyramid apex is at (0.5, 0.5, -1) */ 14 PetscErrorCode pyramidNormal(PetscInt dim, PetscReal time, const PetscReal x[], PetscInt r, PetscScalar u[], void *ctx) 15 { 16 PetscReal apex[3] = {0.5, 0.5, -1.0}; 17 PetscInt d; 18 19 for (d = 0; d < dim; ++d) u[d] = x[d] - apex[d]; 20 for (d = dim; d < 3; ++d) u[d] = 0.0 - apex[d]; 21 return PETSC_SUCCESS; 22 } 23 24 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options) 25 { 26 PetscInt n = 64; 27 PetscBool flg; 28 29 PetscFunctionBeginUser; 30 PetscCall(PetscStrncpy(options->bdLabel, "marker", sizeof(options->bdLabel))); 31 PetscOptionsBegin(comm, "", "Parallel Mesh Adaptation Options", "DMPLEX"); 32 PetscCall(PetscOptionsString("-label", "The boundary label name", "ex44.c", options->bdLabel, options->bdLabel, sizeof(options->bdLabel), NULL)); 33 PetscCall(PetscOptionsIntArray("-bd", "The boundaries to be extruded", "ex44.c", options->bd, &n, &flg)); 34 options->Nbd = flg ? n : 0; 35 PetscOptionsEnd(); 36 PetscFunctionReturn(PETSC_SUCCESS); 37 } 38 39 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm) 40 { 41 PetscFunctionBegin; 42 PetscCall(DMCreate(comm, dm)); 43 PetscCall(DMSetType(*dm, DMPLEX)); 44 PetscCall(DMSetFromOptions(*dm)); 45 PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view")); 46 PetscFunctionReturn(PETSC_SUCCESS); 47 } 48 49 static PetscErrorCode CreateAdaptLabel(DM dm, AppCtx *ctx, DMLabel *adaptLabel) 50 { 51 DMLabel label; 52 PetscInt b; 53 54 PetscFunctionBegin; 55 if (!ctx->Nbd) { 56 *adaptLabel = NULL; 57 PetscFunctionReturn(PETSC_SUCCESS); 58 } 59 PetscCall(DMGetLabel(dm, ctx->bdLabel, &label)); 60 PetscCall(DMLabelCreate(PETSC_COMM_SELF, "Adaptation Label", adaptLabel)); 61 for (b = 0; b < ctx->Nbd; ++b) { 62 IS bdIS; 63 const PetscInt *points; 64 PetscInt n, i; 65 66 PetscCall(DMLabelGetStratumIS(label, ctx->bd[b], &bdIS)); 67 if (!bdIS) continue; 68 PetscCall(ISGetLocalSize(bdIS, &n)); 69 PetscCall(ISGetIndices(bdIS, &points)); 70 for (i = 0; i < n; ++i) PetscCall(DMLabelSetValue(*adaptLabel, points[i], DM_ADAPT_REFINE)); 71 PetscCall(ISRestoreIndices(bdIS, &points)); 72 PetscCall(ISDestroy(&bdIS)); 73 } 74 PetscFunctionReturn(PETSC_SUCCESS); 75 } 76 77 int main(int argc, char **argv) 78 { 79 DM dm, dma; 80 DMLabel adaptLabel; 81 AppCtx ctx; 82 83 PetscFunctionBeginUser; 84 PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 85 PetscCall(ProcessOptions(PETSC_COMM_WORLD, &ctx)); 86 PetscCall(CreateMesh(PETSC_COMM_WORLD, &ctx, &dm)); 87 PetscCall(CreateAdaptLabel(dm, &ctx, &adaptLabel)); 88 if (adaptLabel) { 89 PetscCall(DMAdaptLabel(dm, adaptLabel, &dma)); 90 } else { 91 PetscCall(DMExtrude(dm, 3, &dma)); 92 } 93 PetscCall(PetscObjectSetName((PetscObject)dma, "Adapted Mesh")); 94 PetscCall(DMLabelDestroy(&adaptLabel)); 95 PetscCall(DMDestroy(&dm)); 96 PetscCall(DMViewFromOptions(dma, NULL, "-adapt_dm_view")); 97 PetscCall(DMDestroy(&dma)); 98 PetscCall(PetscFinalize()); 99 return 0; 100 } 101 102 /*TEST 103 104 test: 105 suffix: seg_periodic_0 106 args: -dm_plex_dim 1 -dm_plex_box_faces 3 -dm_plex_transform_extrude_periodic -dm_plex_transform_extrude_use_tensor 0 \ 107 -dm_view -adapt_dm_view -dm_plex_check_all 108 109 test: 110 suffix: tri_tensor_0 111 requires: triangle 112 args: -dm_plex_transform_extrude_use_tensor {{0 1}separate output} \ 113 -dm_view -adapt_dm_view -dm_plex_check_all 114 115 test: 116 suffix: quad_tensor_0 117 args: -dm_plex_simplex 0 -dm_plex_transform_extrude_use_tensor {{0 1}separate output} \ 118 -dm_view -adapt_dm_view -dm_plex_check_all 119 120 test: 121 suffix: quad_normal_0 122 args: -dm_plex_simplex 0 -dm_plex_transform_extrude_normal 0,1,1 \ 123 -dm_view -adapt_dm_view -dm_plex_check_all 124 125 test: 126 suffix: quad_normal_1 127 args: -dm_plex_simplex 0 -dm_plex_transform_extrude_normal_function pyramidNormal \ 128 -dm_view -adapt_dm_view -dm_plex_check_all 129 130 test: 131 suffix: quad_symmetric_0 132 args: -dm_plex_simplex 0 -dm_plex_transform_extrude_symmetric \ 133 -dm_view -adapt_dm_view -dm_plex_check_all 134 135 test: 136 suffix: quad_label 137 args: -dm_plex_simplex 0 -dm_plex_transform_label_replica_inc {{0 100}separate output} \ 138 -dm_view -adapt_dm_view -dm_plex_check_all 139 140 test: 141 suffix: quad_periodic_0 142 args: -dm_plex_simplex 0 -dm_plex_transform_extrude_periodic -dm_plex_transform_extrude_use_tensor 0 \ 143 -dm_view -adapt_dm_view -dm_plex_check_all 144 145 testset: 146 args: -dm_adaptor cellrefiner -dm_plex_transform_type extrude \ 147 -dm_view -adapt_dm_view 148 149 test: 150 suffix: tri_adapt_0 151 requires: triangle 152 args: -dm_plex_box_faces 2,2 -dm_plex_separate_marker -bd 1,3 \ 153 -dm_plex_transform_extrude_thickness 0.5 154 155 test: 156 suffix: tri_adapt_1 157 requires: triangle 158 nsize: 2 159 args: -dm_plex_box_faces 2,2 -bd 1 \ 160 -dm_plex_transform_extrude_thickness 0.5 -petscpartitioner_type simple 161 162 test: 163 suffix: quad_adapt_0 164 args: -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -dm_plex_separate_marker -bd 1,3 \ 165 -dm_plex_transform_extrude_thickness 0.5 166 167 test: 168 suffix: quad_adapt_1 169 nsize: 2 170 args: -dm_plex_simplex 0 -dm_plex_box_faces 2,2 -bd 1 \ 171 -dm_plex_transform_extrude_thickness 0.5 -petscpartitioner_type simple 172 173 test: 174 suffix: tet_adapt_0 175 requires: ctetgen 176 args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_separate_marker -bd 1,3 \ 177 -dm_plex_transform_extrude_thickness 0.5 178 179 test: 180 suffix: tet_adapt_1 181 requires: ctetgen 182 nsize: 2 183 args: -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -bd 1 \ 184 -dm_plex_transform_extrude_thickness 0.5 -petscpartitioner_type simple 185 186 test: 187 suffix: hex_adapt_0 188 args: -dm_plex_simplex 0 -dm_plex_dim 3 -dm_plex_box_faces 2,2,2 -dm_plex_separate_marker -bd 1,3 \ 189 -dm_plex_transform_extrude_thickness 0.5 190 191 TEST*/ 192