xref: /petsc/src/dm/tutorials/ex3.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Tests DMCreateInterpolation() for nonuniform DMDA coordinates.\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscdm.h>
5*c4762a1bSJed Brown #include <petscdmda.h>
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown PetscErrorCode SetCoordinates1d(DM da)
8*c4762a1bSJed Brown {
9*c4762a1bSJed Brown   PetscErrorCode ierr;
10*c4762a1bSJed Brown   PetscInt       i,start,m;
11*c4762a1bSJed Brown   Vec            local,global;
12*c4762a1bSJed Brown   PetscScalar    *coors,*coorslocal;
13*c4762a1bSJed Brown   DM             cda;
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown   PetscFunctionBeginUser;
16*c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
17*c4762a1bSJed Brown   ierr = DMGetCoordinateDM(da,&cda);CHKERRQ(ierr);
18*c4762a1bSJed Brown   ierr = DMGetCoordinates(da,&global);CHKERRQ(ierr);
19*c4762a1bSJed Brown   ierr = DMGetCoordinatesLocal(da,&local);CHKERRQ(ierr);
20*c4762a1bSJed Brown   ierr = DMDAVecGetArray(cda,global,&coors);CHKERRQ(ierr);
21*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(cda,local,&coorslocal);CHKERRQ(ierr);
22*c4762a1bSJed Brown   ierr = DMDAGetCorners(cda,&start,0,0,&m,0,0);CHKERRQ(ierr);
23*c4762a1bSJed Brown   for (i=start; i<start+m; i++) {
24*c4762a1bSJed Brown     if (i % 2) {
25*c4762a1bSJed Brown       coors[i] = coorslocal[i-1] + .1*(coorslocal[i+1] - coorslocal[i-1]);
26*c4762a1bSJed Brown     }
27*c4762a1bSJed Brown   }
28*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(cda,global,&coors);CHKERRQ(ierr);
29*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(cda,local,&coorslocal);CHKERRQ(ierr);
30*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(cda,global,INSERT_VALUES,local);CHKERRQ(ierr);
31*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(cda,global,INSERT_VALUES,local);CHKERRQ(ierr);
32*c4762a1bSJed Brown   PetscFunctionReturn(0);
33*c4762a1bSJed Brown }
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown PetscErrorCode SetCoordinates2d(DM da)
36*c4762a1bSJed Brown {
37*c4762a1bSJed Brown   PetscErrorCode ierr;
38*c4762a1bSJed Brown   PetscInt       i,j,mstart,m,nstart,n;
39*c4762a1bSJed Brown   Vec            local,global;
40*c4762a1bSJed Brown   DMDACoor2d     **coors,**coorslocal;
41*c4762a1bSJed Brown   DM             cda;
42*c4762a1bSJed Brown 
43*c4762a1bSJed Brown   PetscFunctionBeginUser;
44*c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
45*c4762a1bSJed Brown   ierr = DMGetCoordinateDM(da,&cda);CHKERRQ(ierr);
46*c4762a1bSJed Brown   ierr = DMGetCoordinates(da,&global);CHKERRQ(ierr);
47*c4762a1bSJed Brown   ierr = DMGetCoordinatesLocal(da,&local);CHKERRQ(ierr);
48*c4762a1bSJed Brown   ierr = DMDAVecGetArray(cda,global,&coors);CHKERRQ(ierr);
49*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(cda,local,&coorslocal);CHKERRQ(ierr);
50*c4762a1bSJed Brown   ierr = DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0);CHKERRQ(ierr);
51*c4762a1bSJed Brown   for (i=mstart; i<mstart+m; i++) {
52*c4762a1bSJed Brown     for (j=nstart; j<nstart+n; j++) {
53*c4762a1bSJed Brown       if (i % 2) {
54*c4762a1bSJed Brown         coors[j][i].x = coorslocal[j][i-1].x + .1*(coorslocal[j][i+1].x - coorslocal[j][i-1].x);
55*c4762a1bSJed Brown       }
56*c4762a1bSJed Brown       if (j % 2) {
57*c4762a1bSJed Brown         coors[j][i].y = coorslocal[j-1][i].y + .3*(coorslocal[j+1][i].y - coorslocal[j-1][i].y);
58*c4762a1bSJed Brown       }
59*c4762a1bSJed Brown     }
60*c4762a1bSJed Brown   }
61*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(cda,global,&coors);CHKERRQ(ierr);
62*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(cda,local,&coorslocal);CHKERRQ(ierr);
63*c4762a1bSJed Brown 
64*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(cda,global,INSERT_VALUES,local);CHKERRQ(ierr);
65*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(cda,global,INSERT_VALUES,local);CHKERRQ(ierr);
66*c4762a1bSJed Brown   PetscFunctionReturn(0);
67*c4762a1bSJed Brown }
68*c4762a1bSJed Brown 
69*c4762a1bSJed Brown PetscErrorCode SetCoordinates3d(DM da)
70*c4762a1bSJed Brown {
71*c4762a1bSJed Brown   PetscErrorCode ierr;
72*c4762a1bSJed Brown   PetscInt       i,j,mstart,m,nstart,n,pstart,p,k;
73*c4762a1bSJed Brown   Vec            local,global;
74*c4762a1bSJed Brown   DMDACoor3d     ***coors,***coorslocal;
75*c4762a1bSJed Brown   DM             cda;
76*c4762a1bSJed Brown 
77*c4762a1bSJed Brown   PetscFunctionBeginUser;
78*c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
79*c4762a1bSJed Brown   ierr = DMGetCoordinateDM(da,&cda);CHKERRQ(ierr);
80*c4762a1bSJed Brown   ierr = DMGetCoordinates(da,&global);CHKERRQ(ierr);
81*c4762a1bSJed Brown   ierr = DMGetCoordinatesLocal(da,&local);CHKERRQ(ierr);
82*c4762a1bSJed Brown   ierr = DMDAVecGetArray(cda,global,&coors);CHKERRQ(ierr);
83*c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(cda,local,&coorslocal);CHKERRQ(ierr);
84*c4762a1bSJed Brown   ierr = DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p);CHKERRQ(ierr);
85*c4762a1bSJed Brown   for (i=mstart; i<mstart+m; i++) {
86*c4762a1bSJed Brown     for (j=nstart; j<nstart+n; j++) {
87*c4762a1bSJed Brown       for (k=pstart; k<pstart+p; k++) {
88*c4762a1bSJed Brown         if (i % 2) {
89*c4762a1bSJed Brown           coors[k][j][i].x = coorslocal[k][j][i-1].x + .1*(coorslocal[k][j][i+1].x - coorslocal[k][j][i-1].x);
90*c4762a1bSJed Brown         }
91*c4762a1bSJed Brown         if (j % 2) {
92*c4762a1bSJed Brown           coors[k][j][i].y = coorslocal[k][j-1][i].y + .3*(coorslocal[k][j+1][i].y - coorslocal[k][j-1][i].y);
93*c4762a1bSJed Brown         }
94*c4762a1bSJed Brown         if (k % 2) {
95*c4762a1bSJed Brown           coors[k][j][i].z = coorslocal[k-1][j][i].z + .4*(coorslocal[k+1][j][i].z - coorslocal[k-1][j][i].z);
96*c4762a1bSJed Brown         }
97*c4762a1bSJed Brown       }
98*c4762a1bSJed Brown     }
99*c4762a1bSJed Brown   }
100*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(cda,global,&coors);CHKERRQ(ierr);
101*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(cda,local,&coorslocal);CHKERRQ(ierr);
102*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(cda,global,INSERT_VALUES,local);CHKERRQ(ierr);
103*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(cda,global,INSERT_VALUES,local);CHKERRQ(ierr);
104*c4762a1bSJed Brown   PetscFunctionReturn(0);
105*c4762a1bSJed Brown }
106*c4762a1bSJed Brown 
107*c4762a1bSJed Brown int main(int argc,char **argv)
108*c4762a1bSJed Brown {
109*c4762a1bSJed Brown   PetscInt         M = 5,N = 4,P = 3, m = PETSC_DECIDE,n = PETSC_DECIDE,p = PETSC_DECIDE,dim = 1;
110*c4762a1bSJed Brown   PetscErrorCode   ierr;
111*c4762a1bSJed Brown   DM               dac,daf;
112*c4762a1bSJed Brown   DMBoundaryType   bx    = DM_BOUNDARY_NONE,by=DM_BOUNDARY_NONE,bz=DM_BOUNDARY_NONE;
113*c4762a1bSJed Brown   DMDAStencilType  stype = DMDA_STENCIL_BOX;
114*c4762a1bSJed Brown   Mat              A;
115*c4762a1bSJed Brown 
116*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
117*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL);CHKERRQ(ierr);
118*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL);CHKERRQ(ierr);
119*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-P",&P,NULL);CHKERRQ(ierr);
120*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);CHKERRQ(ierr);
121*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL);CHKERRQ(ierr);
123*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL);CHKERRQ(ierr);
124*c4762a1bSJed Brown 
125*c4762a1bSJed Brown   /* Create distributed array and get vectors */
126*c4762a1bSJed Brown   if (dim == 1) {
127*c4762a1bSJed Brown     ierr = DMDACreate1d(PETSC_COMM_WORLD,bx,M,1,1,NULL,&dac);CHKERRQ(ierr);
128*c4762a1bSJed Brown   } else if (dim == 2) {
129*c4762a1bSJed Brown     ierr = DMDACreate2d(PETSC_COMM_WORLD,bx,by,stype,M,N,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&dac);CHKERRQ(ierr);
130*c4762a1bSJed Brown   } else if (dim == 3) {
131*c4762a1bSJed Brown     ierr = DMDACreate3d(PETSC_COMM_WORLD,bx,by,bz,stype,M,N,P,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,NULL,&dac);CHKERRQ(ierr);
132*c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"dim must be 1,2, or 3");
133*c4762a1bSJed Brown   ierr = DMSetFromOptions(dac);CHKERRQ(ierr);
134*c4762a1bSJed Brown   ierr = DMSetUp(dac);CHKERRQ(ierr);
135*c4762a1bSJed Brown 
136*c4762a1bSJed Brown   ierr = DMRefine(dac,PETSC_COMM_WORLD,&daf);CHKERRQ(ierr);
137*c4762a1bSJed Brown 
138*c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(dac,0.0,1.0,0.0,1.0,0.0,1.0);CHKERRQ(ierr);
139*c4762a1bSJed Brown   if (dim == 1) {
140*c4762a1bSJed Brown     ierr = SetCoordinates1d(daf);CHKERRQ(ierr);
141*c4762a1bSJed Brown   } else if (dim == 2) {
142*c4762a1bSJed Brown     ierr = SetCoordinates2d(daf);CHKERRQ(ierr);
143*c4762a1bSJed Brown   } else if (dim == 3) {
144*c4762a1bSJed Brown     ierr = SetCoordinates3d(daf);CHKERRQ(ierr);
145*c4762a1bSJed Brown   }
146*c4762a1bSJed Brown   ierr = DMCreateInterpolation(dac,daf,&A,0);CHKERRQ(ierr);
147*c4762a1bSJed Brown 
148*c4762a1bSJed Brown   /* Free memory */
149*c4762a1bSJed Brown   ierr = DMDestroy(&dac);CHKERRQ(ierr);
150*c4762a1bSJed Brown   ierr = DMDestroy(&daf);CHKERRQ(ierr);
151*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
152*c4762a1bSJed Brown   ierr = PetscFinalize();
153*c4762a1bSJed Brown   return ierr;
154*c4762a1bSJed Brown }
155*c4762a1bSJed Brown 
156*c4762a1bSJed Brown 
157*c4762a1bSJed Brown 
158*c4762a1bSJed Brown /*TEST
159*c4762a1bSJed Brown 
160*c4762a1bSJed Brown    test:
161*c4762a1bSJed Brown       nsize: 3
162*c4762a1bSJed Brown       args: -mat_view
163*c4762a1bSJed Brown 
164*c4762a1bSJed Brown    test:
165*c4762a1bSJed Brown       suffix: 2
166*c4762a1bSJed Brown       nsize: 3
167*c4762a1bSJed Brown       args: -mat_view -dim 2
168*c4762a1bSJed Brown 
169*c4762a1bSJed Brown    test:
170*c4762a1bSJed Brown       suffix: 3
171*c4762a1bSJed Brown       nsize: 3
172*c4762a1bSJed Brown       args: -mat_view -dim 3
173*c4762a1bSJed Brown 
174*c4762a1bSJed Brown TEST*/
175