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