xref: /petsc/src/dm/impls/plex/tests/ex27.c (revision 6a98f8dc3f2c9149905a87dc2e9d0fedaf64e09a)
1 static char help[] = "Test section ordering for FEM discretizations\n\n";
2 
3 #include <petscdmplex.h>
4 #include <petscds.h>
5 
6 typedef struct {
7   PetscInt  dim;                          /* The topological mesh dimension */
8   PetscInt  faces[3];                     /* Number of faces per dimension */
9   PetscBool simplex;                      /* Use simplices or hexes */
10 } AppCtx;
11 
12 static PetscErrorCode ProcessOptions(MPI_Comm comm, AppCtx *options)
13 {
14   PetscInt       dim;
15   PetscErrorCode ierr;
16 
17   PetscFunctionBegin;
18   options->dim      = 1;
19   options->faces[0] = 1;
20   options->faces[1] = 1;
21   options->faces[2] = 1;
22   options->simplex  = PETSC_TRUE;
23 
24   ierr = PetscOptionsBegin(comm, "", "Plex ex27", "DMPLEX");CHKERRQ(ierr);
25   ierr = PetscOptionsRangeInt("-dim", "The topological mesh dimension", "ex27.c", options->dim, &options->dim, NULL,1,3);CHKERRQ(ierr);
26   ierr = PetscOptionsBool("-simplex", "Use simplices if true, otherwise tensor product cells", "ex27.c", options->simplex, &options->simplex, NULL);CHKERRQ(ierr);
27   dim  = options->dim;
28   ierr = PetscOptionsIntArray("-faces", "Number of faces per dimension", "ex27.c", options->faces, &dim, NULL);CHKERRQ(ierr);
29   if (dim) options->dim = dim;
30   ierr = PetscOptionsEnd();
31 
32   if (options->dim > 3) SETERRQ1(comm, PETSC_ERR_ARG_OUTOFRANGE, "Dimension %d must in [0, 4)", options->dim);
33   PetscFunctionReturn(0);
34 }
35 
36 static PetscErrorCode CreateMesh(MPI_Comm comm, AppCtx *options, DM *dm)
37 {
38   PetscInt       dim     = options->dim;
39   PetscInt      *faces   = options->faces;
40   PetscBool      simplex = options->simplex;
41   PetscErrorCode ierr;
42 
43   PetscFunctionBegin;
44   ierr = DMPlexCreateBoxMesh(comm, dim, simplex, faces, NULL, NULL, NULL, PETSC_TRUE, dm);CHKERRQ(ierr);
45   ierr = PetscObjectSetName((PetscObject) *dm, "Serial Mesh");CHKERRQ(ierr);
46   ierr = DMViewFromOptions(*dm, NULL, "-dm_view");CHKERRQ(ierr);
47   PetscFunctionReturn(0);
48 }
49 
50 static PetscErrorCode TestLocalDofOrder(DM dm, AppCtx *ctx)
51 {
52   PetscFE        fe[3];
53   PetscSection   s;
54   PetscInt       dim, Nf, f;
55   PetscErrorCode ierr;
56 
57   PetscFunctionBegin;
58   ierr = DMGetDimension(dm, &dim);CHKERRQ(ierr);
59   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, dim, ctx->simplex, "field0_", -1, &fe[0]);CHKERRQ(ierr);
60   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1,   ctx->simplex, "field1_", -1, &fe[1]);CHKERRQ(ierr);
61   ierr = PetscFECreateDefault(PetscObjectComm((PetscObject) dm), dim, 1,   ctx->simplex, "field2_", -1, &fe[2]);CHKERRQ(ierr);
62 
63   ierr = DMSetField(dm, 0, NULL, (PetscObject) fe[0]);CHKERRQ(ierr);
64   ierr = DMSetField(dm, 1, NULL, (PetscObject) fe[1]);CHKERRQ(ierr);
65   ierr = DMSetField(dm, 2, NULL, (PetscObject) fe[2]);CHKERRQ(ierr);
66   ierr = DMCreateDS(dm);CHKERRQ(ierr);
67   ierr = DMGetLocalSection(dm, &s);CHKERRQ(ierr);
68   ierr = PetscObjectViewFromOptions((PetscObject) s, NULL, "-dof_view");CHKERRQ(ierr);
69 
70   ierr = DMGetNumFields(dm, &Nf);CHKERRQ(ierr);
71   for (f = 0; f < Nf; ++f) {ierr = PetscFEDestroy(&fe[f]);CHKERRQ(ierr);}
72   PetscFunctionReturn(0);
73 }
74 
75 int main(int argc, char **argv)
76 {
77   MPI_Comm       comm;
78   DM             dm;
79   AppCtx         ctx;
80   PetscErrorCode ierr;
81 
82   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
83   comm = PETSC_COMM_WORLD;
84   ierr = ProcessOptions(comm, &ctx);CHKERRQ(ierr);
85   ierr = CreateMesh(comm, &ctx, &dm);CHKERRQ(ierr);
86   ierr = TestLocalDofOrder(dm, &ctx);CHKERRQ(ierr);
87   ierr = DMDestroy(&dm);CHKERRQ(ierr);
88   ierr = PetscFinalize();
89   return ierr;
90 }
91 
92 /*TEST
93 
94   test:
95     suffix: tri_pm
96     requires: triangle
97     args: -dim 2 -field0_petscspace_degree 2 -field1_petscspace_degree 1 -field2_petscspace_degree 1 -dm_view -dof_view
98 
99   test:
100     suffix: quad_pm
101     requires:
102     args: -dim 2 -simplex 0 -field0_petscspace_degree 2 -field1_petscspace_degree 1 -field2_petscspace_degree 1 -dm_view -dof_view
103 
104   test:
105     suffix: tri_fm
106     requires: triangle
107     args: -dim 2 -field0_petscspace_degree 2 -field1_petscspace_degree 1 -field2_petscspace_degree 1 -petscsection_point_major 0 -dm_view -dof_view
108 
109   test:
110     suffix: quad_fm
111     requires:
112     args: -dim 2 -simplex 0 -field0_petscspace_degree 2 -field1_petscspace_degree 1 -field2_petscspace_degree 1 -petscsection_point_major 0 -dm_view -dof_view
113 
114 TEST*/
115