1*e92ee4c8SZach Atkins static char help[] = "Test degenerate near null space";
2*e92ee4c8SZach Atkins
3*e92ee4c8SZach Atkins #include <petscdmplex.h>
4*e92ee4c8SZach Atkins #include <petscsnes.h>
5*e92ee4c8SZach Atkins #include <petscds.h>
6*e92ee4c8SZach Atkins #include <petscbag.h>
7*e92ee4c8SZach Atkins #include <petscconvest.h>
8*e92ee4c8SZach Atkins
CreateMesh(MPI_Comm comm,DM * dm)9*e92ee4c8SZach Atkins static PetscErrorCode CreateMesh(MPI_Comm comm, DM *dm)
10*e92ee4c8SZach Atkins {
11*e92ee4c8SZach Atkins PetscFunctionBeginUser;
12*e92ee4c8SZach Atkins PetscCall(DMCreate(comm, dm));
13*e92ee4c8SZach Atkins PetscCall(DMSetType(*dm, DMPLEX));
14*e92ee4c8SZach Atkins PetscCall(DMSetFromOptions(*dm));
15*e92ee4c8SZach Atkins PetscCall(DMViewFromOptions(*dm, NULL, "-dm_view"));
16*e92ee4c8SZach Atkins PetscFunctionReturn(PETSC_SUCCESS);
17*e92ee4c8SZach Atkins }
18*e92ee4c8SZach Atkins
SetupBoundaries(DM dm)19*e92ee4c8SZach Atkins static PetscErrorCode SetupBoundaries(DM dm)
20*e92ee4c8SZach Atkins {
21*e92ee4c8SZach Atkins DMLabel label;
22*e92ee4c8SZach Atkins PetscInt id;
23*e92ee4c8SZach Atkins PetscInt dim;
24*e92ee4c8SZach Atkins
25*e92ee4c8SZach Atkins PetscFunctionBeginUser;
26*e92ee4c8SZach Atkins PetscCall(DMGetCoordinateDim(dm, &dim));
27*e92ee4c8SZach Atkins PetscCall(DMGetLabel(dm, "Face Sets", &label));
28*e92ee4c8SZach Atkins PetscCheck(label, PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_INCOMP, "Must have face sets label");
29*e92ee4c8SZach Atkins
30*e92ee4c8SZach Atkins if (dim == 2) {
31*e92ee4c8SZach Atkins PetscInt cmp, cmps_y[] = {0, 1};
32*e92ee4c8SZach Atkins
33*e92ee4c8SZach Atkins cmp = 0;
34*e92ee4c8SZach Atkins id = 4;
35*e92ee4c8SZach Atkins PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "left", label, 1, &id, 0, 1, &cmp, NULL, NULL, NULL, NULL));
36*e92ee4c8SZach Atkins cmp = 0;
37*e92ee4c8SZach Atkins id = 2;
38*e92ee4c8SZach Atkins PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "right", label, 1, &id, 0, 1, &cmp, NULL, NULL, NULL, NULL));
39*e92ee4c8SZach Atkins id = 1;
40*e92ee4c8SZach Atkins PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom", label, 1, &id, 0, 2, cmps_y, NULL, NULL, NULL, NULL));
41*e92ee4c8SZach Atkins } else if (dim == 3) {
42*e92ee4c8SZach Atkins PetscInt cmps_xy[] = {0, 1};
43*e92ee4c8SZach Atkins PetscInt cmps_z[] = {0, 1, 2};
44*e92ee4c8SZach Atkins
45*e92ee4c8SZach Atkins id = 6;
46*e92ee4c8SZach Atkins PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "left", label, 1, &id, 0, 2, cmps_xy, NULL, NULL, NULL, NULL));
47*e92ee4c8SZach Atkins id = 5;
48*e92ee4c8SZach Atkins PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "right", label, 1, &id, 0, 2, cmps_xy, NULL, NULL, NULL, NULL));
49*e92ee4c8SZach Atkins id = 2;
50*e92ee4c8SZach Atkins PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "bottom", label, 1, &id, 0, 3, cmps_z, NULL, NULL, NULL, NULL));
51*e92ee4c8SZach Atkins id = 3;
52*e92ee4c8SZach Atkins PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "front", label, 1, &id, 0, 2, cmps_xy, NULL, NULL, NULL, NULL));
53*e92ee4c8SZach Atkins id = 4;
54*e92ee4c8SZach Atkins PetscCall(DMAddBoundary(dm, DM_BC_ESSENTIAL, "back", label, 1, &id, 0, 2, cmps_xy, NULL, NULL, NULL, NULL));
55*e92ee4c8SZach Atkins }
56*e92ee4c8SZach Atkins PetscFunctionReturn(PETSC_SUCCESS);
57*e92ee4c8SZach Atkins }
58*e92ee4c8SZach Atkins
SetupFE(DM dm,const char name[])59*e92ee4c8SZach Atkins PetscErrorCode SetupFE(DM dm, const char name[])
60*e92ee4c8SZach Atkins {
61*e92ee4c8SZach Atkins PetscFE fe;
62*e92ee4c8SZach Atkins char prefix[PETSC_MAX_PATH_LEN];
63*e92ee4c8SZach Atkins DMPolytopeType ct;
64*e92ee4c8SZach Atkins PetscInt dim, cStart;
65*e92ee4c8SZach Atkins
66*e92ee4c8SZach Atkins PetscFunctionBegin;
67*e92ee4c8SZach Atkins /* Create finite element */
68*e92ee4c8SZach Atkins PetscCall(DMGetDimension(dm, &dim));
69*e92ee4c8SZach Atkins PetscCall(DMPlexGetHeightStratum(dm, 0, &cStart, NULL));
70*e92ee4c8SZach Atkins PetscCall(DMPlexGetCellType(dm, cStart, &ct));
71*e92ee4c8SZach Atkins PetscCall(PetscSNPrintf(prefix, PETSC_MAX_PATH_LEN, "%s_", name));
72*e92ee4c8SZach Atkins PetscCall(PetscFECreateByCell(PETSC_COMM_SELF, dim, dim, ct, name ? prefix : NULL, -1, &fe));
73*e92ee4c8SZach Atkins PetscCall(PetscObjectSetName((PetscObject)fe, name));
74*e92ee4c8SZach Atkins /* Set discretization and boundary conditions for each mesh */
75*e92ee4c8SZach Atkins PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe));
76*e92ee4c8SZach Atkins PetscCall(DMCreateDS(dm));
77*e92ee4c8SZach Atkins PetscCall(SetupBoundaries(dm));
78*e92ee4c8SZach Atkins PetscCall(PetscFEDestroy(&fe));
79*e92ee4c8SZach Atkins PetscFunctionReturn(PETSC_SUCCESS);
80*e92ee4c8SZach Atkins }
81*e92ee4c8SZach Atkins
MakeNullSpaceRigidBody(DM dm)82*e92ee4c8SZach Atkins PetscErrorCode MakeNullSpaceRigidBody(DM dm)
83*e92ee4c8SZach Atkins {
84*e92ee4c8SZach Atkins Mat A;
85*e92ee4c8SZach Atkins MatNullSpace null_space;
86*e92ee4c8SZach Atkins
87*e92ee4c8SZach Atkins PetscFunctionBegin;
88*e92ee4c8SZach Atkins /* Create null space and set onto matrix */
89*e92ee4c8SZach Atkins PetscCall(DMCreateMatrix(dm, &A));
90*e92ee4c8SZach Atkins PetscCall(DMPlexCreateRigidBody(dm, 0, &null_space));
91*e92ee4c8SZach Atkins PetscCall(MatSetNearNullSpace(A, null_space));
92*e92ee4c8SZach Atkins PetscCall(MatNullSpaceView(null_space, PETSC_VIEWER_STDOUT_WORLD));
93*e92ee4c8SZach Atkins PetscCall(MatNullSpaceDestroy(&null_space));
94*e92ee4c8SZach Atkins PetscCall(MatDestroy(&A));
95*e92ee4c8SZach Atkins PetscFunctionReturn(PETSC_SUCCESS);
96*e92ee4c8SZach Atkins }
97*e92ee4c8SZach Atkins
main(int argc,char ** argv)98*e92ee4c8SZach Atkins int main(int argc, char **argv)
99*e92ee4c8SZach Atkins {
100*e92ee4c8SZach Atkins DM dm;
101*e92ee4c8SZach Atkins
102*e92ee4c8SZach Atkins PetscFunctionBeginUser;
103*e92ee4c8SZach Atkins PetscCall(PetscInitialize(&argc, &argv, NULL, help));
104*e92ee4c8SZach Atkins PetscCall(CreateMesh(PETSC_COMM_WORLD, &dm));
105*e92ee4c8SZach Atkins PetscCall(SetupFE(dm, "displacement"));
106*e92ee4c8SZach Atkins PetscCall(MakeNullSpaceRigidBody(dm));
107*e92ee4c8SZach Atkins PetscCall(DMDestroy(&dm));
108*e92ee4c8SZach Atkins PetscCall(PetscFinalize());
109*e92ee4c8SZach Atkins return 0;
110*e92ee4c8SZach Atkins }
111*e92ee4c8SZach Atkins
112*e92ee4c8SZach Atkins /*TEST
113*e92ee4c8SZach Atkins
114*e92ee4c8SZach Atkins test:
115*e92ee4c8SZach Atkins suffix: 3d_q1
116*e92ee4c8SZach Atkins args: -dm_plex_box_faces 1,1,2 -dm_plex_simplex 0 -dm_plex_dim 3 -displacement_petscspace_degree 1
117*e92ee4c8SZach Atkins
118*e92ee4c8SZach Atkins test:
119*e92ee4c8SZach Atkins suffix: 3d_q2
120*e92ee4c8SZach Atkins args: -dm_plex_box_faces 1,1,2 -dm_plex_simplex 0 -dm_plex_dim 3 -displacement_petscspace_degree 2
121*e92ee4c8SZach Atkins
122*e92ee4c8SZach Atkins test:
123*e92ee4c8SZach Atkins suffix: 2d_q1
124*e92ee4c8SZach Atkins args: -dm_plex_box_faces 1,2 -dm_plex_simplex 0 -dm_plex_dim 2 -displacement_petscspace_degree 1
125*e92ee4c8SZach Atkins
126*e92ee4c8SZach Atkins test:
127*e92ee4c8SZach Atkins suffix: 2d_q2
128*e92ee4c8SZach Atkins args: -dm_plex_box_faces 1,2 -dm_plex_simplex 0 -dm_plex_dim 2 -displacement_petscspace_degree 2
129*e92ee4c8SZach Atkins
130*e92ee4c8SZach Atkins TEST*/
131