xref: /petsc/src/dm/tests/ex19.c (revision 6a98f8dc3f2c9149905a87dc2e9d0fedaf64e09a)
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 {
13   PetscErrorCode ierr;
14   PetscInt       i,j,k,M,N,dof;
15 
16   ierr = DMDAGetInfo(da,0,&M,&N,0,0,0,0,&dof,0,0,0,0,0);CHKERRQ(ierr);
17   {
18     struct {PetscScalar inside[dof];} **mystruct;
19     ierr = 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           ierr = PetscPrintf(PETSC_COMM_WORLD,"%D %D %g\n",i,j,(double)mystruct[i][j].inside[0]);CHKERRQ(ierr);
24 
25           mystruct[i][j].inside[1] = 2.1;
26         }
27       }
28     }
29     ierr = DMDAVecRestoreArrayRead(da,global,(void*) &mystruct);CHKERRQ(ierr);
30   }
31   PetscFunctionReturn(0);
32 }
33 
34 int main(int argc,char **argv)
35 {
36   PetscInt       dof = 2,M = 3,N = 3,m = PETSC_DECIDE,n = PETSC_DECIDE;
37   PetscErrorCode ierr;
38   DM             da;
39   Vec            global,local;
40 
41   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
42   ierr = PetscOptionsGetInt(NULL,0,"-dof",&dof,0);CHKERRQ(ierr);
43   /* Create distributed array and get vectors */
44   ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,M,N,m,n,dof,1,NULL,NULL,&da);CHKERRQ(ierr);
45   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
46   ierr = DMSetUp(da);CHKERRQ(ierr);
47   ierr = DMCreateGlobalVector(da,&global);CHKERRQ(ierr);
48   ierr = DMCreateLocalVector(da,&local);CHKERRQ(ierr);
49 
50   ierr = doit(da,global);CHKERRQ(ierr);
51 
52 
53   ierr = VecView(global,0);CHKERRQ(ierr);
54 
55   /* Free memory */
56   ierr = VecDestroy(&local);CHKERRQ(ierr);
57   ierr = VecDestroy(&global);CHKERRQ(ierr);
58   ierr = DMDestroy(&da);CHKERRQ(ierr);
59   ierr = PetscFinalize();
60   return ierr;
61 }
62 
63