xref: /petsc/src/dm/tests/ex38.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1 
2 static char help[] = "Tests DMGlobalToLocal() for 3d DA with stencil width of 2.\n\n";
3 
4 #include <petscdm.h>
5 #include <petscdmda.h>
6 
7 int main(int argc,char **argv)
8 {
9   PetscInt         N             = 3,M=2,P=4,dof=1,rstart,rend,i;
10   PetscInt         stencil_width = 2;
11   PetscMPIInt      rank;
12   DMBoundaryType   bx           = DM_BOUNDARY_NONE,by = DM_BOUNDARY_NONE,bz = DM_BOUNDARY_NONE;
13   DMDAStencilType  stencil_type = DMDA_STENCIL_STAR;
14   DM               da;
15   Vec              global,local;
16 
17   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
18   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
19   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
20   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL));
21   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-P",&P,NULL));
22   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-dof",&dof,NULL));
23   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-stencil_width",&stencil_width,NULL));
24   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-stencil_type",(PetscInt*)&stencil_type,NULL));
25 
26   CHKERRQ(DMDACreate3d(PETSC_COMM_WORLD,bx,by,bz,stencil_type,M,N,P,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,dof,stencil_width,0,0,0,&da));
27   CHKERRQ(DMSetFromOptions(da));
28   CHKERRQ(DMSetUp(da));
29   CHKERRQ(DMCreateGlobalVector(da,&global));
30   CHKERRQ(VecGetOwnershipRange(global,&rstart,&rend));
31   for (i=rstart; i<rend; i++) CHKERRQ(VecSetValue(global,i,(PetscReal)(i + 100*rank),INSERT_VALUES));
32   CHKERRQ(VecAssemblyBegin(global));
33   CHKERRQ(VecAssemblyEnd(global));
34   CHKERRQ(DMCreateLocalVector(da,&local));
35   CHKERRQ(VecSet(local,-1));
36   CHKERRQ(DMGlobalToLocalBegin(da,global,INSERT_VALUES,local));
37   CHKERRQ(DMGlobalToLocalEnd(da,global,INSERT_VALUES,local));
38   if (rank == 0) CHKERRQ(VecView(local,0));
39   CHKERRQ(DMDestroy(&da));
40   CHKERRQ(VecDestroy(&local));
41   CHKERRQ(VecDestroy(&global));
42   CHKERRQ(PetscFinalize());
43   return 0;
44 }
45