xref: /petsc/src/dm/tutorials/ex12.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests DMGetGlobalVector() and DMRestoreGlobalVector().\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown /*
5*c4762a1bSJed Brown Use the options
6*c4762a1bSJed Brown      -da_grid_x <nx> - number of grid points in x direction, if M < 0
7*c4762a1bSJed Brown      -da_grid_y <ny> - number of grid points in y direction, if N < 0
8*c4762a1bSJed Brown      -da_processors_x <MX> number of processors in x directio
9*c4762a1bSJed Brown      -da_processors_y <MY> number of processors in x direction
10*c4762a1bSJed Brown */
11*c4762a1bSJed Brown 
12*c4762a1bSJed Brown #include <petscdm.h>
13*c4762a1bSJed Brown #include <petscdmda.h>
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown int main(int argc,char **argv)
16*c4762a1bSJed Brown {
17*c4762a1bSJed Brown   PetscInt         M = 10,N = 8;
18*c4762a1bSJed Brown   PetscErrorCode   ierr;
19*c4762a1bSJed Brown   PetscBool        flg = PETSC_FALSE;
20*c4762a1bSJed Brown   DM               da;
21*c4762a1bSJed Brown   Vec              global1,global2,global3;
22*c4762a1bSJed Brown   DMBoundaryType   bx    = DM_BOUNDARY_NONE,by = DM_BOUNDARY_NONE;
23*c4762a1bSJed Brown   DMDAStencilType  stype = DMDA_STENCIL_BOX;
24*c4762a1bSJed Brown 
25*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
26*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-star_stencil",&flg,NULL);CHKERRQ(ierr);
27*c4762a1bSJed Brown   if (flg) stype = DMDA_STENCIL_STAR;
28*c4762a1bSJed Brown 
29*c4762a1bSJed Brown   /* Create distributed array and get vectors */
30*c4762a1bSJed Brown   ierr = DMDACreate2d(PETSC_COMM_WORLD,bx,by,stype,M,N,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
33*c4762a1bSJed Brown   ierr = DMGetGlobalVector(da,&global1);CHKERRQ(ierr);
34*c4762a1bSJed Brown   ierr = DMGetGlobalVector(da,&global2);CHKERRQ(ierr);
35*c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(da,&global1);CHKERRQ(ierr);
36*c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(da,&global2);CHKERRQ(ierr);
37*c4762a1bSJed Brown   ierr = DMGetGlobalVector(da,&global1);CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = DMGetGlobalVector(da,&global3);CHKERRQ(ierr);
39*c4762a1bSJed Brown   ierr = DMGetGlobalVector(da,&global2);CHKERRQ(ierr);
40*c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(da,&global1);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(da,&global3);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(da,&global2);CHKERRQ(ierr);
43*c4762a1bSJed Brown   ierr = DMGetGlobalVector(da,&global1);CHKERRQ(ierr);
44*c4762a1bSJed Brown   ierr = DMGetGlobalVector(da,&global3);CHKERRQ(ierr);
45*c4762a1bSJed Brown   ierr = DMGetGlobalVector(da,&global2);CHKERRQ(ierr);
46*c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(da,&global1);CHKERRQ(ierr);
47*c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(da,&global3);CHKERRQ(ierr);
48*c4762a1bSJed Brown   ierr = DMRestoreGlobalVector(da,&global2);CHKERRQ(ierr);
49*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
50*c4762a1bSJed Brown   ierr = PetscFinalize();
51*c4762a1bSJed Brown   return ierr;
52*c4762a1bSJed Brown }
53*c4762a1bSJed Brown 
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown 
56*c4762a1bSJed Brown /*TEST
57*c4762a1bSJed Brown 
58*c4762a1bSJed Brown    test:
59*c4762a1bSJed Brown 
60*c4762a1bSJed Brown TEST*/
61