xref: /petsc/src/dm/impls/plex/tests/ex41.c (revision e600fa544e2bb197ca2af9b6e65ea465976dec56)
1568a717dSMatthew G. Knepley static const char help[] = "Tests for adaptive refinement";
2568a717dSMatthew G. Knepley 
3568a717dSMatthew G. Knepley #include <petscdmplex.h>
4568a717dSMatthew G. Knepley 
5568a717dSMatthew G. Knepley typedef struct {
6568a717dSMatthew G. Knepley   PetscBool metric;  /* Flag to use metric adaptation, instead of tagging */
7568a717dSMatthew G. Knepley   PetscInt *refcell; /* A cell to be refined on each process */
8568a717dSMatthew G. Knepley } AppCtx;
9568a717dSMatthew G. Knepley 
10568a717dSMatthew G. Knepley static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
11568a717dSMatthew G. Knepley {
12568a717dSMatthew G. Knepley   PetscMPIInt    size;
13568a717dSMatthew G. Knepley   PetscInt       n;
14568a717dSMatthew G. Knepley   PetscErrorCode ierr;
15568a717dSMatthew G. Knepley 
16568a717dSMatthew G. Knepley   PetscFunctionBeginUser;
17568a717dSMatthew G. Knepley   options->metric  = PETSC_FALSE;
1855b25c41SPierre Jolivet   ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr);
19568a717dSMatthew G. Knepley   ierr = PetscCalloc1(size, &options->refcell);CHKERRQ(ierr);
20568a717dSMatthew G. Knepley   n    = size;
21568a717dSMatthew G. Knepley 
22568a717dSMatthew G. Knepley   ierr = PetscOptionsBegin(comm, "", "Parallel Mesh Adaptation Options", "DMPLEX");CHKERRQ(ierr);
23568a717dSMatthew G. Knepley   ierr = PetscOptionsBool("-metric", "Flag for metric refinement", "ex41.c", options->metric, &options->metric, NULL);CHKERRQ(ierr);
24568a717dSMatthew G. Knepley   ierr = PetscOptionsIntArray("-refcell", "The cell to be refined", "ex41.c", options->refcell, &n, NULL);CHKERRQ(ierr);
252c71b3e2SJacob Faibussowitsch   PetscCheckFalse(n && n != size,comm, PETSC_ERR_ARG_SIZ, "Only gave %D cells to refine, must give one for all %D processes", n, size);
261e1ea65dSPierre Jolivet   ierr = PetscOptionsEnd();CHKERRQ(ierr);
27568a717dSMatthew G. Knepley   PetscFunctionReturn(0);
28568a717dSMatthew G. Knepley }
29568a717dSMatthew G. Knepley 
30568a717dSMatthew G. Knepley static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *ctx, DM *dm)
31568a717dSMatthew G. Knepley {
32568a717dSMatthew G. Knepley   PetscErrorCode ierr;
33568a717dSMatthew G. Knepley 
34568a717dSMatthew G. Knepley   PetscFunctionBegin;
3530602db0SMatthew G. Knepley   ierr = DMCreate(comm, dm);CHKERRQ(ierr);
3630602db0SMatthew G. Knepley   ierr = DMSetType(*dm, DMPLEX);CHKERRQ(ierr);
37568a717dSMatthew G. Knepley   ierr = DMSetFromOptions(*dm);CHKERRQ(ierr);
38568a717dSMatthew G. Knepley   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
39568a717dSMatthew G. Knepley   PetscFunctionReturn(0);
40568a717dSMatthew G. Knepley }
41568a717dSMatthew G. Knepley 
42568a717dSMatthew G. Knepley static PetscErrorCode CreateAdaptLabel(DM dm, AppCtx *ctx, DMLabel *adaptLabel)
43568a717dSMatthew G. Knepley {
44568a717dSMatthew G. Knepley   PetscMPIInt    rank;
45568a717dSMatthew G. Knepley   PetscErrorCode ierr;
46568a717dSMatthew G. Knepley 
47568a717dSMatthew G. Knepley   PetscFunctionBegin;
4855b25c41SPierre Jolivet   ierr = MPI_Comm_rank(PetscObjectComm((PetscObject) dm), &rank);CHKERRMPI(ierr);
49568a717dSMatthew G. Knepley   ierr = DMLabelCreate(PETSC_COMM_SELF, "Adaptation Label", adaptLabel);CHKERRQ(ierr);
50568a717dSMatthew G. Knepley   if (ctx->refcell[rank] >= 0) {ierr = DMLabelSetValue(*adaptLabel, ctx->refcell[rank], DM_ADAPT_REFINE);CHKERRQ(ierr);}
51568a717dSMatthew G. Knepley   PetscFunctionReturn(0);
52568a717dSMatthew G. Knepley }
53568a717dSMatthew G. Knepley 
54568a717dSMatthew G. Knepley int main(int argc, char **argv)
55568a717dSMatthew G. Knepley {
56568a717dSMatthew G. Knepley   DM             dm, dma;
57568a717dSMatthew G. Knepley   DMLabel        adaptLabel;
58568a717dSMatthew G. Knepley   AppCtx         ctx;
59568a717dSMatthew G. Knepley   PetscErrorCode ierr;
60568a717dSMatthew G. Knepley 
61568a717dSMatthew G. Knepley   ierr = PetscInitialize(&argc, &argv, NULL, help); if (ierr) return ierr;
62568a717dSMatthew G. Knepley   ierr = ProcessOptions(PETSC_COMM_WORLD, &ctx);CHKERRQ(ierr);
63568a717dSMatthew G. Knepley   ierr = CreateMesh(PETSC_COMM_WORLD, &ctx, &dm);CHKERRQ(ierr);
64568a717dSMatthew G. Knepley   ierr = CreateAdaptLabel(dm, &ctx, &adaptLabel);CHKERRQ(ierr);
65568a717dSMatthew G. Knepley   ierr = DMAdaptLabel(dm, adaptLabel, &dma);CHKERRQ(ierr);
66568a717dSMatthew G. Knepley   ierr = PetscObjectSetName((PetscObject) dma, "Adapted Mesh");CHKERRQ(ierr);
67568a717dSMatthew G. Knepley   ierr = DMLabelDestroy(&adaptLabel);CHKERRQ(ierr);
68568a717dSMatthew G. Knepley   ierr = DMDestroy(&dm);CHKERRQ(ierr);
69568a717dSMatthew G. Knepley   ierr = DMViewFromOptions(dma, NULL, "-adapt_dm_view");CHKERRQ(ierr);
70568a717dSMatthew G. Knepley   ierr = DMDestroy(&dma);CHKERRQ(ierr);
71568a717dSMatthew G. Knepley   ierr = PetscFree(ctx.refcell);CHKERRQ(ierr);
72568a717dSMatthew G. Knepley   ierr = PetscFinalize();
73568a717dSMatthew G. Knepley   return ierr;
74568a717dSMatthew G. Knepley }
75568a717dSMatthew G. Knepley 
76568a717dSMatthew G. Knepley /*TEST
77568a717dSMatthew G. Knepley 
78012bc364SMatthew G. Knepley   testset:
79c0517cd5SMatthew G. Knepley     args: -dm_adaptor cellrefiner -dm_plex_transform_type refine_sbr
80012bc364SMatthew G. Knepley 
81568a717dSMatthew G. Knepley     test:
82568a717dSMatthew G. Knepley       suffix: 0
83568a717dSMatthew G. Knepley       requires: triangle
84012bc364SMatthew G. Knepley       args: -dm_view -adapt_dm_view
85568a717dSMatthew G. Knepley 
86568a717dSMatthew G. Knepley     test:
87568a717dSMatthew G. Knepley       suffix: 1
88568a717dSMatthew G. Knepley       requires: triangle
89012bc364SMatthew G. Knepley       args: -dm_coord_space 0 -refcell 2 -dm_view ::ascii_info_detail -adapt_dm_view ::ascii_info_detail
90568a717dSMatthew G. Knepley 
91568a717dSMatthew G. Knepley     test:
92568a717dSMatthew G. Knepley       suffix: 2
93568a717dSMatthew G. Knepley       requires: triangle
94568a717dSMatthew G. Knepley       nsize: 2
95*e600fa54SMatthew G. Knepley       args: -refcell 2,-1 -petscpartitioner_type simple -dm_view -adapt_dm_view
96568a717dSMatthew G. Knepley 
97568a717dSMatthew G. Knepley TEST*/
98