xref: /petsc/src/dm/tests/ex19.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 
2 static char help[] = "Tests DMDA with variable multiple degrees of freedom per node.\n\n";
3 
4 /*
5    This code only compiles with gcc, since it is not ANSI C
6 */
7 
8 #include <petscdm.h>
9 #include <petscdmda.h>
10 
11 PetscErrorCode doit(DM da, Vec global) {
12   PetscInt i, j, k, M, N, dof;
13 
14   PetscCall(DMDAGetInfo(da, 0, &M, &N, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
15   {
16     struct {
17       PetscScalar inside[dof];
18     } * *mystruct;
19     PetscCall(DMDAVecGetArrayRead(da, global, (void *)&mystruct));
20     for (i = 0; i < N; i++) {
21       for (j = 0; j < M; j++) {
22         for (k = 0; k < dof; k++) {
23           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT " %" PetscInt_FMT " %g\n", i, j, (double)mystruct[i][j].inside[0]));
24 
25           mystruct[i][j].inside[1] = 2.1;
26         }
27       }
28     }
29     PetscCall(DMDAVecRestoreArrayRead(da, global, (void *)&mystruct));
30   }
31   PetscFunctionReturn(0);
32 }
33 
34 int main(int argc, char **argv) {
35   PetscInt dof = 2, M = 3, N = 3, m = PETSC_DECIDE, n = PETSC_DECIDE;
36   DM       da;
37   Vec      global, local;
38 
39   PetscFunctionBeginUser;
40   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
41   PetscCall(PetscOptionsGetInt(NULL, 0, "-dof", &dof, 0));
42   /* Create distributed array and get vectors */
43   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, m, n, dof, 1, NULL, NULL, &da));
44   PetscCall(DMSetFromOptions(da));
45   PetscCall(DMSetUp(da));
46   PetscCall(DMCreateGlobalVector(da, &global));
47   PetscCall(DMCreateLocalVector(da, &local));
48 
49   PetscCall(doit(da, global));
50 
51   PetscCall(VecView(global, 0));
52 
53   /* Free memory */
54   PetscCall(VecDestroy(&local));
55   PetscCall(VecDestroy(&global));
56   PetscCall(DMDestroy(&da));
57   PetscCall(PetscFinalize());
58   return 0;
59 }
60