xref: /petsc/src/dm/tests/ex4.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests various 2-dimensional DMDA routines.\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscdm.h>
5*c4762a1bSJed Brown #include <petscdmda.h>
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown int main(int argc,char **argv)
8*c4762a1bSJed Brown {
9*c4762a1bSJed Brown   PetscMPIInt      rank;
10*c4762a1bSJed Brown   PetscErrorCode   ierr;
11*c4762a1bSJed Brown   PetscInt         M = 10,N = 8,m = PETSC_DECIDE;
12*c4762a1bSJed Brown   PetscInt         s =2,w=2,n = PETSC_DECIDE,nloc,l,i,j,kk;
13*c4762a1bSJed Brown   PetscInt         Xs,Xm,Ys,Ym,iloc,*iglobal;
14*c4762a1bSJed Brown   const PetscInt   *ltog;
15*c4762a1bSJed Brown   PetscInt         *lx       = NULL,*ly = NULL;
16*c4762a1bSJed Brown   PetscBool        testorder = PETSC_FALSE,flg;
17*c4762a1bSJed Brown   DMBoundaryType   bx        = DM_BOUNDARY_NONE,by= DM_BOUNDARY_NONE;
18*c4762a1bSJed Brown   DM               da;
19*c4762a1bSJed Brown   PetscViewer      viewer;
20*c4762a1bSJed Brown   Vec              local,global;
21*c4762a1bSJed Brown   PetscScalar      value;
22*c4762a1bSJed Brown   DMDAStencilType  st = DMDA_STENCIL_BOX;
23*c4762a1bSJed Brown   AO               ao;
24*c4762a1bSJed Brown 
25*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
26*c4762a1bSJed Brown   ierr = PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",300,0,400,400,&viewer);CHKERRQ(ierr);
27*c4762a1bSJed Brown 
28*c4762a1bSJed Brown   /* Readoptions */
29*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-NX",&M,NULL);CHKERRQ(ierr);
30*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-NY",&N,NULL);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
32*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
33*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-s",&s,NULL);CHKERRQ(ierr);
34*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-w",&w,NULL);CHKERRQ(ierr);
35*c4762a1bSJed Brown 
36*c4762a1bSJed Brown   flg  = PETSC_FALSE;
37*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-xperiodic",&flg,NULL);CHKERRQ(ierr); if (flg) bx = DM_BOUNDARY_PERIODIC;
38*c4762a1bSJed Brown   flg  = PETSC_FALSE;
39*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-yperiodic",&flg,NULL);CHKERRQ(ierr); if (flg) by = DM_BOUNDARY_PERIODIC;
40*c4762a1bSJed Brown   flg  = PETSC_FALSE;
41*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-xghosted",&flg,NULL);CHKERRQ(ierr); if (flg) bx = DM_BOUNDARY_GHOSTED;
42*c4762a1bSJed Brown   flg  = PETSC_FALSE;
43*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-yghosted",&flg,NULL);CHKERRQ(ierr); if (flg) by = DM_BOUNDARY_GHOSTED;
44*c4762a1bSJed Brown   flg  = PETSC_FALSE;
45*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-star",&flg,NULL);CHKERRQ(ierr); if (flg) st = DMDA_STENCIL_STAR;
46*c4762a1bSJed Brown   flg  = PETSC_FALSE;
47*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-box",&flg,NULL);CHKERRQ(ierr); if (flg) st = DMDA_STENCIL_BOX;
48*c4762a1bSJed Brown   flg  = PETSC_FALSE;
49*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-testorder",&testorder,NULL);CHKERRQ(ierr);
50*c4762a1bSJed Brown   /*
51*c4762a1bSJed Brown       Test putting two nodes in x and y on each processor, exact last processor
52*c4762a1bSJed Brown       in x and y gets the rest.
53*c4762a1bSJed Brown   */
54*c4762a1bSJed Brown   flg  = PETSC_FALSE;
55*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-distribute",&flg,NULL);CHKERRQ(ierr);
56*c4762a1bSJed Brown   if (flg) {
57*c4762a1bSJed Brown     if (m == PETSC_DECIDE) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Must set -m option with -distribute option");
58*c4762a1bSJed Brown     ierr = PetscMalloc1(m,&lx);CHKERRQ(ierr);
59*c4762a1bSJed Brown     for (i=0; i<m-1; i++) { lx[i] = 4;}
60*c4762a1bSJed Brown     lx[m-1] = M - 4*(m-1);
61*c4762a1bSJed Brown     if (n == PETSC_DECIDE) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Must set -n option with -distribute option");
62*c4762a1bSJed Brown     ierr = PetscMalloc1(n,&ly);CHKERRQ(ierr);
63*c4762a1bSJed Brown     for (i=0; i<n-1; i++) { ly[i] = 2;}
64*c4762a1bSJed Brown     ly[n-1] = N - 2*(n-1);
65*c4762a1bSJed Brown   }
66*c4762a1bSJed Brown 
67*c4762a1bSJed Brown 
68*c4762a1bSJed Brown   /* Create distributed array and get vectors */
69*c4762a1bSJed Brown   ierr = DMDACreate2d(PETSC_COMM_WORLD,bx,by,st,M,N,m,n,w,s,lx,ly,&da);CHKERRQ(ierr);
70*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
71*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
72*c4762a1bSJed Brown   ierr = PetscFree(lx);CHKERRQ(ierr);
73*c4762a1bSJed Brown   ierr = PetscFree(ly);CHKERRQ(ierr);
74*c4762a1bSJed Brown 
75*c4762a1bSJed Brown   ierr = DMView(da,viewer);CHKERRQ(ierr);
76*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&global);CHKERRQ(ierr);
77*c4762a1bSJed Brown   ierr = DMCreateLocalVector(da,&local);CHKERRQ(ierr);
78*c4762a1bSJed Brown 
79*c4762a1bSJed Brown   /* Set global vector; send ghost points to local vectors */
80*c4762a1bSJed Brown   value = 1;
81*c4762a1bSJed Brown   ierr = VecSet(global,value);CHKERRQ(ierr);
82*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);CHKERRQ(ierr);
83*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);CHKERRQ(ierr);
84*c4762a1bSJed Brown 
85*c4762a1bSJed Brown   /* Scale local vectors according to processor rank; pass to global vector */
86*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
87*c4762a1bSJed Brown   value = rank;
88*c4762a1bSJed Brown   ierr = VecScale(local,value);CHKERRQ(ierr);
89*c4762a1bSJed Brown   ierr = DMLocalToGlobalBegin(da,local,INSERT_VALUES,global);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = DMLocalToGlobalEnd(da,local,INSERT_VALUES,global);CHKERRQ(ierr);
91*c4762a1bSJed Brown 
92*c4762a1bSJed Brown   if (!testorder) { /* turn off printing when testing ordering mappings */
93*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"\nGlobal Vectors:\n");CHKERRQ(ierr);
94*c4762a1bSJed Brown     ierr = VecView(global,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
95*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"\n\n");CHKERRQ(ierr);
96*c4762a1bSJed Brown   }
97*c4762a1bSJed Brown 
98*c4762a1bSJed Brown   /* Send ghost points to local vectors */
99*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,global,INSERT_VALUES,local);CHKERRQ(ierr);
100*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,global,INSERT_VALUES,local);CHKERRQ(ierr);
101*c4762a1bSJed Brown 
102*c4762a1bSJed Brown   flg  = PETSC_FALSE;
103*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-local_print",&flg,NULL);CHKERRQ(ierr);
104*c4762a1bSJed Brown   if (flg) {
105*c4762a1bSJed Brown     PetscViewer sviewer;
106*c4762a1bSJed Brown 
107*c4762a1bSJed Brown     ierr = PetscViewerASCIIPushSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
108*c4762a1bSJed Brown     ierr = PetscSynchronizedPrintf(PETSC_COMM_WORLD,"\nLocal Vector: processor %d\n",rank);CHKERRQ(ierr);
109*c4762a1bSJed Brown     ierr = PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
110*c4762a1bSJed Brown     ierr = VecView(local,sviewer);CHKERRQ(ierr);
111*c4762a1bSJed Brown     ierr = PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr);
112*c4762a1bSJed Brown     ierr = PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
113*c4762a1bSJed Brown     ierr = PetscViewerASCIIPopSynchronized(PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
114*c4762a1bSJed Brown   }
115*c4762a1bSJed Brown 
116*c4762a1bSJed Brown   /* Tests mappings betweeen application/PETSc orderings */
117*c4762a1bSJed Brown   if (testorder) {
118*c4762a1bSJed Brown     ISLocalToGlobalMapping ltogm;
119*c4762a1bSJed Brown 
120*c4762a1bSJed Brown     ierr = DMGetLocalToGlobalMapping(da,&ltogm);CHKERRQ(ierr);
121*c4762a1bSJed Brown     ierr = ISLocalToGlobalMappingGetSize(ltogm,&nloc);CHKERRQ(ierr);
122*c4762a1bSJed Brown     ierr = ISLocalToGlobalMappingGetIndices(ltogm,&ltog);CHKERRQ(ierr);
123*c4762a1bSJed Brown     ierr = DMDAGetGhostCorners(da,&Xs,&Ys,NULL,&Xm,&Ym,NULL);CHKERRQ(ierr);
124*c4762a1bSJed Brown     ierr = DMDAGetAO(da,&ao);CHKERRQ(ierr);
125*c4762a1bSJed Brown     ierr = PetscMalloc1(nloc,&iglobal);CHKERRQ(ierr);
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown     /* Set iglobal to be global indices for each processor's local and ghost nodes,
128*c4762a1bSJed Brown        using the DMDA ordering of grid points */
129*c4762a1bSJed Brown     kk = 0;
130*c4762a1bSJed Brown     for (j=Ys; j<Ys+Ym; j++) {
131*c4762a1bSJed Brown       for (i=Xs; i<Xs+Xm; i++) {
132*c4762a1bSJed Brown         iloc = w*((j-Ys)*Xm + i-Xs);
133*c4762a1bSJed Brown         for (l=0; l<w; l++) {
134*c4762a1bSJed Brown           iglobal[kk++] = ltog[iloc+l];
135*c4762a1bSJed Brown         }
136*c4762a1bSJed Brown       }
137*c4762a1bSJed Brown     }
138*c4762a1bSJed Brown 
139*c4762a1bSJed Brown     /* Map this to the application ordering (which for DMDAs is just the natural ordering
140*c4762a1bSJed Brown        that would be used for 1 processor, numbering most rapidly by x, then y) */
141*c4762a1bSJed Brown     ierr = AOPetscToApplication(ao,nloc,iglobal);CHKERRQ(ierr);
142*c4762a1bSJed Brown 
143*c4762a1bSJed Brown     /* Then map the application ordering back to the PETSc DMDA ordering */
144*c4762a1bSJed Brown     ierr = AOApplicationToPetsc(ao,nloc,iglobal);CHKERRQ(ierr);
145*c4762a1bSJed Brown 
146*c4762a1bSJed Brown     /* Verify the mappings */
147*c4762a1bSJed Brown     kk=0;
148*c4762a1bSJed Brown     for (j=Ys; j<Ys+Ym; j++) {
149*c4762a1bSJed Brown       for (i=Xs; i<Xs+Xm; i++) {
150*c4762a1bSJed Brown         iloc = w*((j-Ys)*Xm + i-Xs);
151*c4762a1bSJed Brown         for (l=0; l<w; l++) {
152*c4762a1bSJed Brown           if (iglobal[kk] != ltog[iloc+l]) {
153*c4762a1bSJed Brown             ierr = PetscFPrintf(PETSC_COMM_SELF,stdout,"[%d] Problem with mapping: j=%D, i=%D, l=%D, petsc1=%D, petsc2=%D\n",rank,j,i,l,ltog[iloc+l],iglobal[kk]);CHKERRQ(ierr);
154*c4762a1bSJed Brown           }
155*c4762a1bSJed Brown           kk++;
156*c4762a1bSJed Brown         }
157*c4762a1bSJed Brown       }
158*c4762a1bSJed Brown     }
159*c4762a1bSJed Brown     ierr = PetscFree(iglobal);CHKERRQ(ierr);
160*c4762a1bSJed Brown     ierr = ISLocalToGlobalMappingRestoreIndices(ltogm,&ltog);CHKERRQ(ierr);
161*c4762a1bSJed Brown   }
162*c4762a1bSJed Brown 
163*c4762a1bSJed Brown   /* Free memory */
164*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
165*c4762a1bSJed Brown   ierr = VecDestroy(&local);CHKERRQ(ierr);
166*c4762a1bSJed Brown   ierr = VecDestroy(&global);CHKERRQ(ierr);
167*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
168*c4762a1bSJed Brown 
169*c4762a1bSJed Brown   ierr = PetscFinalize();
170*c4762a1bSJed Brown   return ierr;
171*c4762a1bSJed Brown }
172*c4762a1bSJed Brown 
173*c4762a1bSJed Brown 
174*c4762a1bSJed Brown /*TEST
175*c4762a1bSJed Brown 
176*c4762a1bSJed Brown    test:
177*c4762a1bSJed Brown       nsize: 4
178*c4762a1bSJed Brown       args: -nox
179*c4762a1bSJed Brown       filter: grep -v -i Object
180*c4762a1bSJed Brown       requires: x
181*c4762a1bSJed Brown 
182*c4762a1bSJed Brown    test:
183*c4762a1bSJed Brown       suffix: 2
184*c4762a1bSJed Brown       args: -testorder -nox
185*c4762a1bSJed Brown       requires: x
186*c4762a1bSJed Brown 
187*c4762a1bSJed Brown TEST*/
188