xref: /petsc/src/dm/tutorials/ex3.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests DMCreateInterpolation() for nonuniform DMDA coordinates.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscdm.h>
5c4762a1bSJed Brown #include <petscdmda.h>
6c4762a1bSJed Brown 
7c4762a1bSJed Brown PetscErrorCode SetCoordinates1d(DM da)
8c4762a1bSJed Brown {
9c4762a1bSJed Brown   PetscInt       i,start,m;
10c4762a1bSJed Brown   Vec            local,global;
11c4762a1bSJed Brown   PetscScalar    *coors,*coorslocal;
12c4762a1bSJed Brown   DM             cda;
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   PetscFunctionBeginUser;
155f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0));
165f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(da,&cda));
175f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinates(da,&global));
185f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(da,&local));
195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(cda,global,&coors));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(cda,local,&coorslocal));
215f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(cda,&start,0,0,&m,0,0));
22c4762a1bSJed Brown   for (i=start; i<start+m; i++) {
23c4762a1bSJed Brown     if (i % 2) {
24c4762a1bSJed Brown       coors[i] = coorslocal[i-1] + .1*(coorslocal[i+1] - coorslocal[i-1]);
25c4762a1bSJed Brown     }
26c4762a1bSJed Brown   }
275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(cda,global,&coors));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(cda,local,&coorslocal));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(cda,global,INSERT_VALUES,local));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(cda,global,INSERT_VALUES,local));
31c4762a1bSJed Brown   PetscFunctionReturn(0);
32c4762a1bSJed Brown }
33c4762a1bSJed Brown 
34c4762a1bSJed Brown PetscErrorCode SetCoordinates2d(DM da)
35c4762a1bSJed Brown {
36c4762a1bSJed Brown   PetscInt       i,j,mstart,m,nstart,n;
37c4762a1bSJed Brown   Vec            local,global;
38c4762a1bSJed Brown   DMDACoor2d     **coors,**coorslocal;
39c4762a1bSJed Brown   DM             cda;
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   PetscFunctionBeginUser;
425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(da,&cda));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinates(da,&global));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(da,&local));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(cda,global,&coors));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(cda,local,&coorslocal));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(cda,&mstart,&nstart,0,&m,&n,0));
49c4762a1bSJed Brown   for (i=mstart; i<mstart+m; i++) {
50c4762a1bSJed Brown     for (j=nstart; j<nstart+n; j++) {
51c4762a1bSJed Brown       if (i % 2) {
52c4762a1bSJed Brown         coors[j][i].x = coorslocal[j][i-1].x + .1*(coorslocal[j][i+1].x - coorslocal[j][i-1].x);
53c4762a1bSJed Brown       }
54c4762a1bSJed Brown       if (j % 2) {
55c4762a1bSJed Brown         coors[j][i].y = coorslocal[j-1][i].y + .3*(coorslocal[j+1][i].y - coorslocal[j-1][i].y);
56c4762a1bSJed Brown       }
57c4762a1bSJed Brown     }
58c4762a1bSJed Brown   }
595f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(cda,global,&coors));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(cda,local,&coorslocal));
61c4762a1bSJed Brown 
625f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(cda,global,INSERT_VALUES,local));
635f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(cda,global,INSERT_VALUES,local));
64c4762a1bSJed Brown   PetscFunctionReturn(0);
65c4762a1bSJed Brown }
66c4762a1bSJed Brown 
67c4762a1bSJed Brown PetscErrorCode SetCoordinates3d(DM da)
68c4762a1bSJed Brown {
69c4762a1bSJed Brown   PetscInt       i,j,mstart,m,nstart,n,pstart,p,k;
70c4762a1bSJed Brown   Vec            local,global;
71c4762a1bSJed Brown   DMDACoor3d     ***coors,***coorslocal;
72c4762a1bSJed Brown   DM             cda;
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   PetscFunctionBeginUser;
755f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(da,0.0,1.0,0.0,1.0,0.0,1.0));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinateDM(da,&cda));
775f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinates(da,&global));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetCoordinatesLocal(da,&local));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(cda,global,&coors));
805f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(cda,local,&coorslocal));
815f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(cda,&mstart,&nstart,&pstart,&m,&n,&p));
82c4762a1bSJed Brown   for (i=mstart; i<mstart+m; i++) {
83c4762a1bSJed Brown     for (j=nstart; j<nstart+n; j++) {
84c4762a1bSJed Brown       for (k=pstart; k<pstart+p; k++) {
85c4762a1bSJed Brown         if (i % 2) {
86c4762a1bSJed 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);
87c4762a1bSJed Brown         }
88c4762a1bSJed Brown         if (j % 2) {
89c4762a1bSJed 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);
90c4762a1bSJed Brown         }
91c4762a1bSJed Brown         if (k % 2) {
92c4762a1bSJed 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);
93c4762a1bSJed Brown         }
94c4762a1bSJed Brown       }
95c4762a1bSJed Brown     }
96c4762a1bSJed Brown   }
975f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(cda,global,&coors));
985f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(cda,local,&coorslocal));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(cda,global,INSERT_VALUES,local));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(cda,global,INSERT_VALUES,local));
101c4762a1bSJed Brown   PetscFunctionReturn(0);
102c4762a1bSJed Brown }
103c4762a1bSJed Brown 
104c4762a1bSJed Brown int main(int argc,char **argv)
105c4762a1bSJed Brown {
106c4762a1bSJed Brown   PetscInt         M = 5,N = 4,P = 3, m = PETSC_DECIDE,n = PETSC_DECIDE,p = PETSC_DECIDE,dim = 1;
107c4762a1bSJed Brown   DM               dac,daf;
108c4762a1bSJed Brown   DMBoundaryType   bx    = DM_BOUNDARY_NONE,by=DM_BOUNDARY_NONE,bz=DM_BOUNDARY_NONE;
109c4762a1bSJed Brown   DMDAStencilType  stype = DMDA_STENCIL_BOX;
110c4762a1bSJed Brown   Mat              A;
111c4762a1bSJed Brown 
112*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-M",&M,NULL));
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-N",&N,NULL));
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-P",&P,NULL));
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-p",&p,NULL));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-dim",&dim,NULL));
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   /* Create distributed array and get vectors */
122c4762a1bSJed Brown   if (dim == 1) {
1235f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDACreate1d(PETSC_COMM_WORLD,bx,M,1,1,NULL,&dac));
124c4762a1bSJed Brown   } else if (dim == 2) {
1255f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,bx,by,stype,M,N,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&dac));
126c4762a1bSJed Brown   } else if (dim == 3) {
1275f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDACreate3d(PETSC_COMM_WORLD,bx,by,bz,stype,M,N,P,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,NULL,&dac));
128c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"dim must be 1,2, or 3");
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(dac));
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(dac));
131c4762a1bSJed Brown 
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRefine(dac,PETSC_COMM_WORLD,&daf));
133c4762a1bSJed Brown 
1345f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(dac,0.0,1.0,0.0,1.0,0.0,1.0));
135c4762a1bSJed Brown   if (dim == 1) {
1365f80ce2aSJacob Faibussowitsch     CHKERRQ(SetCoordinates1d(daf));
137c4762a1bSJed Brown   } else if (dim == 2) {
1385f80ce2aSJacob Faibussowitsch     CHKERRQ(SetCoordinates2d(daf));
139c4762a1bSJed Brown   } else if (dim == 3) {
1405f80ce2aSJacob Faibussowitsch     CHKERRQ(SetCoordinates3d(daf));
141c4762a1bSJed Brown   }
1425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateInterpolation(dac,daf,&A,0));
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* Free memory */
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&dac));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&daf));
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
148*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
149*b122ec5aSJacob Faibussowitsch   return 0;
150c4762a1bSJed Brown }
151c4762a1bSJed Brown 
152c4762a1bSJed Brown /*TEST
153c4762a1bSJed Brown 
154c4762a1bSJed Brown    test:
155c4762a1bSJed Brown       nsize: 3
156c4762a1bSJed Brown       args: -mat_view
157c4762a1bSJed Brown 
158c4762a1bSJed Brown    test:
159c4762a1bSJed Brown       suffix: 2
160c4762a1bSJed Brown       nsize: 3
161c4762a1bSJed Brown       args: -mat_view -dim 2
162c4762a1bSJed Brown 
163c4762a1bSJed Brown    test:
164c4762a1bSJed Brown       suffix: 3
165c4762a1bSJed Brown       nsize: 3
166c4762a1bSJed Brown       args: -mat_view -dim 3
167c4762a1bSJed Brown 
168c4762a1bSJed Brown TEST*/
169