xref: /petsc/src/dm/impls/stag/tests/ex51.c (revision 609caa7c8c030312b00807b4f015fd827bb80932)
1 static char help[] = "DMStag slot test (to excerpt for manual)\n\n";
2 
3 #include <petscdmstag.h>
4 
main(int argc,char ** argv)5 int main(int argc, char **argv)
6 {
7   DM              dm;
8   Vec             x;
9   PetscInt        s_x, s_y, s_z, n_x, n_y, n_z, n_e_x, n_e_y, n_e_z, slot_vertex_2;
10   PetscScalar ****x_array;
11 
12   const DMStagStencilLocation location_vertex = DMSTAG_BACK_DOWN_LEFT;
13   const PetscInt              dof0 = 2, dof1 = 2, dof2 = 2, dof3 = 2, N_x = 3, N_y = 3, N_z = 3;
14 
15   PetscFunctionBeginUser;
16   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
17 
18   PetscCall(DMStagCreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, N_x, N_y, N_z, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, dof0, dof1, dof2, dof3, DMSTAG_STENCIL_BOX, 1, NULL, NULL, NULL, &dm));
19   PetscCall(DMSetFromOptions(dm));
20   PetscCall(DMSetUp(dm));
21 
22   PetscCall(DMCreateLocalVector(dm, &x));
23   PetscCall(VecZeroEntries(x));
24 
25   /* Set the second component of all vertex dof to 2.0 */
26   PetscCall(DMStagGetCorners(dm, &s_x, &s_y, &s_z, &n_x, &n_y, &n_z, &n_e_x, &n_e_y, &n_e_z));
27   PetscCall(DMStagGetLocationSlot(dm, location_vertex, 1, &slot_vertex_2));
28   PetscCall(DMStagVecGetArray(dm, x, &x_array));
29   for (PetscInt k = s_z; k < s_z + n_z + n_e_z; ++k) {
30     for (PetscInt j = s_y; j < s_y + n_y + n_e_y; ++j) {
31       for (PetscInt i = s_x; i < s_x + n_x + n_e_x; ++i) x_array[k][j][i][slot_vertex_2] = 2.0;
32     }
33   }
34   PetscCall(DMStagVecRestoreArray(dm, x, &x_array));
35   PetscCall(VecDestroy(&x));
36   PetscCall(DMDestroy(&dm));
37   PetscCall(PetscFinalize());
38   return 0;
39 }
40 
41 /*TEST
42 
43    test:
44      output_file: output/empty.out
45 
46 TEST*/
47