xref: /petsc/src/dm/dt/fe/tests/ex4.c (revision 1bd0d7ef42bb7aa3da4c59e50d8c0998a28fb2f9)
1*1bd0d7efSStefano Zampini static const char help[] = "Tests DMCreateMassMatrix and DMCreateMassMatrixLumped";
2*1bd0d7efSStefano Zampini 
3*1bd0d7efSStefano Zampini #include <petscdmplex.h>
4*1bd0d7efSStefano Zampini #include <petscfe.h>
5*1bd0d7efSStefano Zampini #include <petscds.h>
6*1bd0d7efSStefano Zampini 
7*1bd0d7efSStefano Zampini static void volume(PetscInt dim, PetscInt Nf, PetscInt NfAux, const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[], const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[], PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar obj[])
8*1bd0d7efSStefano Zampini {
9*1bd0d7efSStefano Zampini   obj[0] = 1.0;
10*1bd0d7efSStefano Zampini }
11*1bd0d7efSStefano Zampini 
12*1bd0d7efSStefano Zampini static PetscErrorCode SetupDiscretization(DM dm)
13*1bd0d7efSStefano Zampini {
14*1bd0d7efSStefano Zampini   PetscFE   fe0, fe1, fe2;
15*1bd0d7efSStefano Zampini   DM        plex;
16*1bd0d7efSStefano Zampini   PetscInt  dim;
17*1bd0d7efSStefano Zampini   PetscBool simplex;
18*1bd0d7efSStefano Zampini   PetscDS   ds;
19*1bd0d7efSStefano Zampini 
20*1bd0d7efSStefano Zampini   PetscFunctionBeginUser;
21*1bd0d7efSStefano Zampini   PetscCall(DMGetDimension(dm, &dim));
22*1bd0d7efSStefano Zampini   PetscCall(DMConvert(dm, DMPLEX, &plex));
23*1bd0d7efSStefano Zampini   PetscCall(DMPlexIsSimplex(plex, &simplex));
24*1bd0d7efSStefano Zampini   PetscCall(MPIU_Allreduce(MPI_IN_PLACE, &simplex, 1, MPIU_BOOL, MPI_LOR, PetscObjectComm((PetscObject)dm)));
25*1bd0d7efSStefano Zampini   PetscCall(DMDestroy(&plex));
26*1bd0d7efSStefano Zampini 
27*1bd0d7efSStefano Zampini   /* Create finite element */
28*1bd0d7efSStefano Zampini   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 1, simplex, "f0_", -1, &fe0));
29*1bd0d7efSStefano Zampini   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, 2, simplex, "f1_", -1, &fe1));
30*1bd0d7efSStefano Zampini   PetscCall(PetscFECreateDefault(PetscObjectComm((PetscObject)dm), dim, dim, simplex, "f2_", -1, &fe2));
31*1bd0d7efSStefano Zampini   PetscCall(DMSetField(dm, 0, NULL, (PetscObject)fe0));
32*1bd0d7efSStefano Zampini   PetscCall(DMSetField(dm, 1, NULL, (PetscObject)fe1));
33*1bd0d7efSStefano Zampini   PetscCall(DMSetField(dm, 2, NULL, (PetscObject)fe2));
34*1bd0d7efSStefano Zampini   PetscCall(DMCreateDS(dm));
35*1bd0d7efSStefano Zampini   PetscCall(DMGetDS(dm, &ds));
36*1bd0d7efSStefano Zampini   PetscCall(PetscDSSetObjective(ds, 0, volume));
37*1bd0d7efSStefano Zampini   PetscCall(PetscDSSetObjective(ds, 1, volume));
38*1bd0d7efSStefano Zampini   PetscCall(PetscDSSetObjective(ds, 2, volume));
39*1bd0d7efSStefano Zampini   PetscCall(PetscFEDestroy(&fe0));
40*1bd0d7efSStefano Zampini   PetscCall(PetscFEDestroy(&fe1));
41*1bd0d7efSStefano Zampini   PetscCall(PetscFEDestroy(&fe2));
42*1bd0d7efSStefano Zampini   PetscFunctionReturn(PETSC_SUCCESS);
43*1bd0d7efSStefano Zampini }
44*1bd0d7efSStefano Zampini 
45*1bd0d7efSStefano Zampini #define CheckVals(a, b, rtol, atol, msg) \
46*1bd0d7efSStefano Zampini   do { \
47*1bd0d7efSStefano Zampini     if (!PetscIsCloseAtTolScalar(a, b, rtol, atol)) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s: %g (%g - %g)\n", msg, (double)PetscAbsScalar(a - b), (double)PetscAbsScalar(a), (double)PetscAbsScalar(b))); \
48*1bd0d7efSStefano Zampini   } while (0)
49*1bd0d7efSStefano Zampini 
50*1bd0d7efSStefano Zampini int main(int argc, char **argv)
51*1bd0d7efSStefano Zampini {
52*1bd0d7efSStefano Zampini   DM          dm;
53*1bd0d7efSStefano Zampini   Mat         M;
54*1bd0d7efSStefano Zampini   Vec         llM, lM, ones, work;
55*1bd0d7efSStefano Zampini   PetscReal   rtol = PETSC_SMALL, atol = 0.0;
56*1bd0d7efSStefano Zampini   PetscScalar vals[6];
57*1bd0d7efSStefano Zampini   PetscInt    dim;
58*1bd0d7efSStefano Zampini   PetscBool   amr = PETSC_FALSE;
59*1bd0d7efSStefano Zampini 
60*1bd0d7efSStefano Zampini   PetscFunctionBeginUser;
61*1bd0d7efSStefano Zampini   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
62*1bd0d7efSStefano Zampini   PetscCall(PetscOptionsGetBool(NULL, NULL, "-amr", &amr, NULL));
63*1bd0d7efSStefano Zampini   PetscCall(DMCreate(PETSC_COMM_WORLD, &dm));
64*1bd0d7efSStefano Zampini   PetscCall(DMSetType(dm, DMPLEX));
65*1bd0d7efSStefano Zampini   PetscCall(DMSetFromOptions(dm));
66*1bd0d7efSStefano Zampini   PetscCall(DMGetDimension(dm, &dim));
67*1bd0d7efSStefano Zampini   if (amr) {
68*1bd0d7efSStefano Zampini     DM dmConv;
69*1bd0d7efSStefano Zampini 
70*1bd0d7efSStefano Zampini     PetscCall(DMConvert(dm, dim == 2 ? DMP4EST : DMP8EST, &dmConv));
71*1bd0d7efSStefano Zampini     PetscCall(DMDestroy(&dm));
72*1bd0d7efSStefano Zampini     dm = dmConv;
73*1bd0d7efSStefano Zampini     PetscCall(DMSetFromOptions(dm));
74*1bd0d7efSStefano Zampini     PetscCall(DMSetUp(dm));
75*1bd0d7efSStefano Zampini   }
76*1bd0d7efSStefano Zampini   PetscCall(DMViewFromOptions(dm, NULL, "-mesh_view"));
77*1bd0d7efSStefano Zampini   PetscCall(SetupDiscretization(dm));
78*1bd0d7efSStefano Zampini   PetscCall(DMGetGlobalVector(dm, &ones));
79*1bd0d7efSStefano Zampini   PetscCall(DMGetGlobalVector(dm, &work));
80*1bd0d7efSStefano Zampini   PetscCall(VecSet(ones, 1.0));
81*1bd0d7efSStefano Zampini   PetscCall(DMCreateMassMatrixLumped(dm, &llM, &lM));
82*1bd0d7efSStefano Zampini   PetscCall(VecViewFromOptions(llM, NULL, "-local_lumped_mass_view"));
83*1bd0d7efSStefano Zampini   PetscCall(VecViewFromOptions(lM, NULL, "-lumped_mass_view"));
84*1bd0d7efSStefano Zampini   PetscCall(DMCreateMassMatrix(dm, NULL, &M));
85*1bd0d7efSStefano Zampini   PetscCall(MatViewFromOptions(M, NULL, "-mass_view"));
86*1bd0d7efSStefano Zampini   PetscCall(MatMult(M, ones, work));
87*1bd0d7efSStefano Zampini   PetscCall(VecViewFromOptions(work, NULL, "-mass_rowsum_view"));
88*1bd0d7efSStefano Zampini   PetscCall(VecSum(work, &vals[3]));
89*1bd0d7efSStefano Zampini   PetscCall(VecSet(work, 0.0));
90*1bd0d7efSStefano Zampini   PetscCall(DMLocalToGlobal(dm, llM, ADD_VALUES, work));
91*1bd0d7efSStefano Zampini   PetscCall(VecSum(work, &vals[4]));
92*1bd0d7efSStefano Zampini   PetscCall(VecSum(lM, &vals[5]));
93*1bd0d7efSStefano Zampini   PetscCall(DMPlexComputeIntegralFEM(dm, ones, vals, NULL));
94*1bd0d7efSStefano Zampini   CheckVals(vals[0], vals[1], rtol, atol, "Error volume");
95*1bd0d7efSStefano Zampini   CheckVals((3 + dim) * vals[0], vals[3], rtol, atol, "Error mass");
96*1bd0d7efSStefano Zampini   CheckVals((3 + dim) * vals[0], vals[4], rtol, atol, "Error local lumped mass");
97*1bd0d7efSStefano Zampini   CheckVals((3 + dim) * vals[0], vals[5], rtol, atol, "Error lumped mass");
98*1bd0d7efSStefano Zampini   PetscCall(MatDestroy(&M));
99*1bd0d7efSStefano Zampini   PetscCall(VecDestroy(&lM));
100*1bd0d7efSStefano Zampini   PetscCall(VecDestroy(&llM));
101*1bd0d7efSStefano Zampini   PetscCall(DMRestoreGlobalVector(dm, &work));
102*1bd0d7efSStefano Zampini   PetscCall(DMRestoreGlobalVector(dm, &ones));
103*1bd0d7efSStefano Zampini   PetscCall(DMDestroy(&dm));
104*1bd0d7efSStefano Zampini   PetscCall(PetscFinalize());
105*1bd0d7efSStefano Zampini   return 0;
106*1bd0d7efSStefano Zampini }
107*1bd0d7efSStefano Zampini 
108*1bd0d7efSStefano Zampini /*TEST
109*1bd0d7efSStefano Zampini 
110*1bd0d7efSStefano Zampini   test:
111*1bd0d7efSStefano Zampini     suffix: 0
112*1bd0d7efSStefano Zampini     output_file: output/ex4_ok.out
113*1bd0d7efSStefano Zampini     nsize: {{1 2}}
114*1bd0d7efSStefano Zampini     args: -dm_plex_box_faces 3,3,3 -dm_plex_dim {{2 3}} -dm_plex_simplex 0 -f0_petscspace_degree {{0 1}} -f1_petscspace_degree {{0 1}} -f2_petscspace_degree {{0 1}} -petscpartitioner_type simple
115*1bd0d7efSStefano Zampini 
116*1bd0d7efSStefano Zampini   test:
117*1bd0d7efSStefano Zampini     TODO: broken
118*1bd0d7efSStefano Zampini     suffix: 0_amr
119*1bd0d7efSStefano Zampini     output_file: output/ex4_ok.out
120*1bd0d7efSStefano Zampini     requires: p4est
121*1bd0d7efSStefano Zampini     nsize: {{1 2}}
122*1bd0d7efSStefano Zampini     args: -dm_plex_dim {{2 3}} -dm_plex_simplex 0 -f0_petscspace_degree {{0 1}} -f1_petscspace_degree {{0 1}} -f2_petscspace_degree {{0 1}} -amr -dm_p4est_refine_pattern hash -dm_forest_initial_refinement 0 -dm_forest_maximum_refinement 2
123*1bd0d7efSStefano Zampini 
124*1bd0d7efSStefano Zampini   test:
125*1bd0d7efSStefano Zampini     suffix: 1
126*1bd0d7efSStefano Zampini     output_file: output/ex4_ok.out
127*1bd0d7efSStefano Zampini     requires: triangle
128*1bd0d7efSStefano Zampini     nsize: {{1 2}}
129*1bd0d7efSStefano Zampini     args: -dm_plex_box_faces 3,3 -dm_plex_dim 2 -dm_plex_simplex 1 -f0_petscspace_degree {{0 1}} -f1_petscspace_degree {{0 1}} -f2_petscspace_degree {{0 1}} -petscpartitioner_type simple
130*1bd0d7efSStefano Zampini 
131*1bd0d7efSStefano Zampini   test:
132*1bd0d7efSStefano Zampini     suffix: 2
133*1bd0d7efSStefano Zampini     output_file: output/ex4_ok.out
134*1bd0d7efSStefano Zampini     requires: ctetgen
135*1bd0d7efSStefano Zampini     nsize: {{1 2}}
136*1bd0d7efSStefano Zampini     args: -dm_plex_box_faces 3,3,3 -dm_plex_dim 3 -dm_plex_simplex 1 -f0_petscspace_degree {{0 1}} -f1_petscspace_degree {{0 1}} -f2_petscspace_degree {{0 1}} -petscpartitioner_type simple
137*1bd0d7efSStefano Zampini 
138*1bd0d7efSStefano Zampini TEST*/
139