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