xref: /petsc/src/dm/tests/ex19.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
1c4762a1bSJed Brown static char help[] = "Tests DMDA with variable multiple degrees of freedom per node.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown    This code only compiles with gcc, since it is not ANSI C
5c4762a1bSJed Brown */
6c4762a1bSJed Brown 
7c4762a1bSJed Brown #include <petscdm.h>
8c4762a1bSJed Brown #include <petscdmda.h>
9c4762a1bSJed Brown 
doit(DM da,Vec global)10d71ae5a4SJacob Faibussowitsch PetscErrorCode doit(DM da, Vec global)
11d71ae5a4SJacob Faibussowitsch {
12c4762a1bSJed Brown   PetscInt i, j, k, M, N, dof;
13c4762a1bSJed Brown 
149566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &M, &N, 0, 0, 0, 0, &dof, 0, 0, 0, 0, 0));
15c4762a1bSJed Brown   {
169371c9d4SSatish Balay     struct {
179371c9d4SSatish Balay       PetscScalar inside[dof];
189371c9d4SSatish Balay     } **mystruct;
199566063dSJacob Faibussowitsch     PetscCall(DMDAVecGetArrayRead(da, global, (void *)&mystruct));
20c4762a1bSJed Brown     for (i = 0; i < N; i++) {
21c4762a1bSJed Brown       for (j = 0; j < M; j++) {
22c4762a1bSJed Brown         for (k = 0; k < dof; k++) {
2363a3b9bcSJacob Faibussowitsch           PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%" PetscInt_FMT " %" PetscInt_FMT " %g\n", i, j, (double)mystruct[i][j].inside[0]));
24c4762a1bSJed Brown 
25c4762a1bSJed Brown           mystruct[i][j].inside[1] = 2.1;
26c4762a1bSJed Brown         }
27c4762a1bSJed Brown       }
28c4762a1bSJed Brown     }
299566063dSJacob Faibussowitsch     PetscCall(DMDAVecRestoreArrayRead(da, global, (void *)&mystruct));
30c4762a1bSJed Brown   }
313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
32c4762a1bSJed Brown }
33c4762a1bSJed Brown 
main(int argc,char ** argv)34d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
35d71ae5a4SJacob Faibussowitsch {
36c4762a1bSJed Brown   PetscInt dof = 2, M = 3, N = 3, m = PETSC_DECIDE, n = PETSC_DECIDE;
37c4762a1bSJed Brown   DM       da;
38c4762a1bSJed Brown   Vec      global, local;
39c4762a1bSJed Brown 
40327415f7SBarry Smith   PetscFunctionBeginUser;
41*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
429566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, 0, "-dof", &dof, 0));
43c4762a1bSJed Brown   /* Create distributed array and get vectors */
449566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, M, N, m, n, dof, 1, NULL, NULL, &da));
459566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
469566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
479566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &global));
489566063dSJacob Faibussowitsch   PetscCall(DMCreateLocalVector(da, &local));
49c4762a1bSJed Brown 
509566063dSJacob Faibussowitsch   PetscCall(doit(da, global));
51c4762a1bSJed Brown 
529566063dSJacob Faibussowitsch   PetscCall(VecView(global, 0));
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   /* Free memory */
559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&local));
569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&global));
579566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
589566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
59b122ec5aSJacob Faibussowitsch   return 0;
60c4762a1bSJed Brown }
61