xref: /petsc/src/dm/tests/ex51.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1 
2 static char help[] = "Tests DMDAGlobalToNaturalAllCreate() using contour plotting for 2d DMDAs.\n\n";
3 
4 #include <petscdm.h>
5 #include <petscdmda.h>
6 #include <petscdraw.h>
7 
8 int main(int argc,char **argv)
9 {
10   PetscInt         i,j,M = 10,N = 8,m = PETSC_DECIDE,n = PETSC_DECIDE;
11   PetscMPIInt      rank;
12   PetscBool        flg = PETSC_FALSE;
13   DM               da;
14   PetscViewer      viewer;
15   Vec              localall,global;
16   PetscScalar      value,*vlocal;
17   DMBoundaryType   bx    = DM_BOUNDARY_NONE,by = DM_BOUNDARY_NONE;
18   DMDAStencilType  stype = DMDA_STENCIL_BOX;
19   VecScatter       tolocalall,fromlocalall;
20   PetscInt         start,end;
21 
22   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
23   CHKERRQ(PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",300,0,300,300,&viewer));
24 
25   /* Read options */
26   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
27   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL));
28   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
29   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
30   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-star_stencil",&flg,NULL));
31   if (flg) stype = DMDA_STENCIL_STAR;
32 
33   /* Create distributed array and get vectors */
34   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,bx,by,stype,M,N,m,n,1,1,NULL,NULL,&da));
35   CHKERRQ(DMSetFromOptions(da));
36   CHKERRQ(DMSetUp(da));
37 
38   CHKERRQ(DMCreateGlobalVector(da,&global));
39   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,M*N,&localall));
40 
41   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
42   CHKERRQ(VecGetOwnershipRange(global,&start,&end));
43   for (i=start; i<end; i++) {
44     value = 5.0*rank;
45     CHKERRQ(VecSetValues(global,1,&i,&value,INSERT_VALUES));
46   }
47   CHKERRQ(VecView(global,viewer));
48 
49   /*
50      Create Scatter from global DMDA parallel vector to local vector that
51    contains all entries
52   */
53   CHKERRQ(DMDAGlobalToNaturalAllCreate(da,&tolocalall));
54   CHKERRQ(DMDANaturalAllToGlobalCreate(da,&fromlocalall));
55 
56   CHKERRQ(VecScatterBegin(tolocalall,global,localall,INSERT_VALUES,SCATTER_FORWARD));
57   CHKERRQ(VecScatterEnd(tolocalall,global,localall,INSERT_VALUES,SCATTER_FORWARD));
58 
59   CHKERRQ(VecGetArray(localall,&vlocal));
60   for (j=0; j<N; j++) {
61     for (i=0; i<M; i++) {
62       *vlocal++ += i + j*M;
63     }
64   }
65   CHKERRQ(VecRestoreArray(localall,&vlocal));
66 
67   /* scatter back to global vector */
68   CHKERRQ(VecScatterBegin(fromlocalall,localall,global,INSERT_VALUES,SCATTER_FORWARD));
69   CHKERRQ(VecScatterEnd(fromlocalall,localall,global,INSERT_VALUES,SCATTER_FORWARD));
70 
71   CHKERRQ(VecView(global,viewer));
72 
73   /* Free memory */
74   CHKERRQ(VecScatterDestroy(&tolocalall));
75   CHKERRQ(VecScatterDestroy(&fromlocalall));
76   CHKERRQ(PetscViewerDestroy(&viewer));
77   CHKERRQ(VecDestroy(&localall));
78   CHKERRQ(VecDestroy(&global));
79   CHKERRQ(DMDestroy(&da));
80   CHKERRQ(PetscFinalize());
81   return 0;
82 }
83 
84 /*TEST
85 
86    build:
87      requires: !complex
88 
89    test:
90       nsize: 3
91 
92 TEST*/
93