xref: /petsc/src/dm/tests/ex48.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[] = "Test VTK structured (.vts)  and rectilinear (.vtr) viewer support with multi-dof DMDAs.\n\
2c4762a1bSJed Brown                       Supply the -namefields flag to test with field names.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscdm.h>
5c4762a1bSJed Brown #include <petscdmda.h>
6c4762a1bSJed Brown 
7c4762a1bSJed Brown /* Helper function to name DMDA fields */
8c4762a1bSJed Brown PetscErrorCode NameFields(DM da,PetscInt dof)
9c4762a1bSJed Brown {
10c4762a1bSJed Brown   PetscInt       c;
11c4762a1bSJed Brown 
12c4762a1bSJed Brown   PetscFunctionBeginUser;
13c4762a1bSJed Brown   for (c=0; c<dof; ++c) {
14c4762a1bSJed Brown     char fieldname[256];
155f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscSNPrintf(fieldname,sizeof(fieldname),"field_%D",c));
165f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDASetFieldName(da,c,fieldname));
17c4762a1bSJed Brown   }
18c4762a1bSJed Brown   PetscFunctionReturn(0);
19c4762a1bSJed Brown }
20c4762a1bSJed Brown 
21c4762a1bSJed Brown /*
22c4762a1bSJed Brown   Write 3D DMDA vector with coordinates in VTK format
23c4762a1bSJed Brown */
24c4762a1bSJed Brown PetscErrorCode test_3d(const char filename[],PetscInt dof,PetscBool namefields)
25c4762a1bSJed Brown {
26c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
27c4762a1bSJed Brown   const PetscInt    M=10,N=15,P=30,sw=1;
28c4762a1bSJed Brown   const PetscScalar Lx=1.0,Ly=1.0,Lz=1.0;
29c4762a1bSJed Brown   DM                da;
30c4762a1bSJed Brown   Vec               v;
31c4762a1bSJed Brown   PetscViewer       view;
32c4762a1bSJed Brown   DMDALocalInfo     info;
33c4762a1bSJed Brown   PetscScalar       ****va;
34c4762a1bSJed Brown   PetscInt          i,j,k,c;
35c4762a1bSJed Brown 
365f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate3d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR, M,N,P,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,dof,sw,NULL,NULL,NULL,&da));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
395f80ce2aSJacob Faibussowitsch   if (namefields) CHKERRQ(NameFields(da,dof));
40c4762a1bSJed Brown 
415f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(da,0.0,Lx,0.0,Ly,0.0,Lz));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(da,&info));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&v));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayDOF(da,v,&va));
45c4762a1bSJed Brown   for (k=info.zs; k<info.zs+info.zm; k++) {
46c4762a1bSJed Brown     for (j=info.ys; j<info.ys+info.ym; j++) {
47c4762a1bSJed Brown       for (i=info.xs; i<info.xs+info.xm; i++) {
48c4762a1bSJed Brown         const PetscScalar x = (Lx*i)/M;
49c4762a1bSJed Brown         const PetscScalar y = (Ly*j)/N;
50c4762a1bSJed Brown         const PetscScalar z = (Lz*k)/P;
51c4762a1bSJed Brown         for (c=0; c<dof; ++ c) {
52c4762a1bSJed Brown         va[k][j][i][c] = PetscPowScalarInt(x-0.5*Lx,2)+PetscPowScalarInt(y-0.5*Ly,2)+PetscPowScalarInt(z-0.5*Lz,2) + 10.0*c;
53c4762a1bSJed Brown         }
54c4762a1bSJed Brown       }
55c4762a1bSJed Brown     }
56c4762a1bSJed Brown   }
575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayDOF(da,v,&va));
585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view));
595f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(v,view));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&view));
615f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&v));
625f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
63c4762a1bSJed Brown   return 0;
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
66c4762a1bSJed Brown /*
67c4762a1bSJed Brown   Write 2D DMDA vector with coordinates in VTK format
68c4762a1bSJed Brown */
69c4762a1bSJed Brown PetscErrorCode test_2d(const char filename[],PetscInt dof,PetscBool namefields)
70c4762a1bSJed Brown {
71c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
72c4762a1bSJed Brown   const PetscInt    M=10,N=20,sw=1;
73c4762a1bSJed Brown   const PetscScalar Lx=1.0,Ly=1.0,Lz=1.0;
74c4762a1bSJed Brown   DM                da;
75c4762a1bSJed Brown   Vec               v;
76c4762a1bSJed Brown   PetscViewer       view;
77c4762a1bSJed Brown   DMDALocalInfo     info;
78c4762a1bSJed Brown   PetscScalar       ***va;
79c4762a1bSJed Brown   PetscInt          i,j,c;
80c4762a1bSJed Brown 
815f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR, M,N,PETSC_DECIDE,PETSC_DECIDE,dof,sw,NULL,NULL,&da));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
845f80ce2aSJacob Faibussowitsch   if (namefields) CHKERRQ(NameFields(da,dof));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(da,0.0,Lx,0.0,Ly,0.0,Lz));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(da,&info));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&v));
885f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayDOF(da,v,&va));
89c4762a1bSJed Brown   for (j=info.ys; j<info.ys+info.ym; j++) {
90c4762a1bSJed Brown     for (i=info.xs; i<info.xs+info.xm; i++) {
91c4762a1bSJed Brown       const PetscScalar x = (Lx*i)/M;
92c4762a1bSJed Brown       const PetscScalar y = (Ly*j)/N;
93c4762a1bSJed Brown       for (c=0; c<dof; ++c) {
94c4762a1bSJed Brown         va[j][i][c] = PetscPowScalarInt(x-0.5*Lx,2)+PetscPowScalarInt(y-0.5*Ly,2) + 10.0*c;
95c4762a1bSJed Brown       }
96c4762a1bSJed Brown     }
97c4762a1bSJed Brown   }
985f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayDOF(da,v,&va));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(v,view));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&view));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&v));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
104c4762a1bSJed Brown   return 0;
105c4762a1bSJed Brown }
106c4762a1bSJed Brown 
107c4762a1bSJed Brown /*
108c4762a1bSJed Brown   Write a scalar and a vector field from two compatible 3d DMDAs
109c4762a1bSJed Brown */
110c4762a1bSJed Brown PetscErrorCode test_3d_compat(const char filename[],PetscInt dof,PetscBool namefields)
111c4762a1bSJed Brown {
112c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
113c4762a1bSJed Brown   const PetscInt    M=10,N=15,P=30,sw=1;
114c4762a1bSJed Brown   const PetscScalar Lx=1.0,Ly=1.0,Lz=1.0;
115c4762a1bSJed Brown   DM                da,daVector;
116c4762a1bSJed Brown   Vec               v,vVector;
117c4762a1bSJed Brown   PetscViewer       view;
118c4762a1bSJed Brown   DMDALocalInfo     info;
119c4762a1bSJed Brown   PetscScalar       ***va,****vVectora;
120c4762a1bSJed Brown   PetscInt          i,j,k,c;
121c4762a1bSJed Brown 
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate3d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR, M,N,P,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,/* dof:*/1,sw,NULL,NULL,NULL,&da));
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
1255f80ce2aSJacob Faibussowitsch   if (namefields) CHKERRQ(NameFields(da,1));
126c4762a1bSJed Brown 
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(da,0.0,Lx,0.0,Ly,0.0,Lz));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(da,&info));
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreateCompatibleDMDA(da,dof,&daVector));
1305f80ce2aSJacob Faibussowitsch   if (namefields) CHKERRQ(NameFields(daVector,dof));
1315f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&v));
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(daVector,&vVector));
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,v,&va));
1345f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayDOF(daVector,vVector,&vVectora));
135c4762a1bSJed Brown   for (k=info.zs; k<info.zs+info.zm; k++) {
136c4762a1bSJed Brown     for (j=info.ys; j<info.ys+info.ym; j++) {
137c4762a1bSJed Brown       for (i=info.xs; i<info.xs+info.xm; i++) {
138c4762a1bSJed Brown         const PetscScalar x = (Lx*i)/M;
139c4762a1bSJed Brown         const PetscScalar y = (Ly*j)/N;
140c4762a1bSJed Brown         const PetscScalar z = (Lz*k)/P;
141c4762a1bSJed Brown         va[k][j][i] = PetscPowScalarInt(x-0.5*Lx,2)+PetscPowScalarInt(y-0.5*Ly,2)+PetscPowScalarInt(z-0.5*Lz,2);
142c4762a1bSJed Brown         for (c=0; c<dof; ++c) {
143c4762a1bSJed Brown           vVectora[k][j][i][c] = PetscPowScalarInt(x-0.5*Lx,2)+PetscPowScalarInt(y-0.5*Ly,2)+PetscPowScalarInt(z-0.5*Lz,2) + 10.0*c;
144c4762a1bSJed Brown         }
145c4762a1bSJed Brown       }
146c4762a1bSJed Brown     }
147c4762a1bSJed Brown   }
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,v,&va));
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayDOF(da,v,&vVectora));
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view));
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(v,view));
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(vVector,view));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&view));
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&v));
1555f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&vVector));
1565f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
1575f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&daVector));
158c4762a1bSJed Brown   return 0;
159c4762a1bSJed Brown }
160c4762a1bSJed Brown 
161c4762a1bSJed Brown /*
162c4762a1bSJed Brown   Write a scalar and a vector field from two compatible 2d DMDAs
163c4762a1bSJed Brown */
164c4762a1bSJed Brown PetscErrorCode test_2d_compat(const char filename[],PetscInt dof,PetscBool namefields)
165c4762a1bSJed Brown {
166c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
167c4762a1bSJed Brown   const PetscInt    M=10,N=20,sw=1;
168c4762a1bSJed Brown   const PetscScalar Lx=1.0,Ly=1.0,Lz=1.0;
169c4762a1bSJed Brown   DM                da, daVector;
170c4762a1bSJed Brown   Vec               v,vVector;
171c4762a1bSJed Brown   PetscViewer       view;
172c4762a1bSJed Brown   DMDALocalInfo     info;
173c4762a1bSJed Brown   PetscScalar       **va,***vVectora;
174c4762a1bSJed Brown   PetscInt          i,j,c;
175c4762a1bSJed Brown 
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR, M,N,PETSC_DECIDE,PETSC_DECIDE,/* dof:*/ 1,sw,NULL,NULL,&da));
1775f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
1795f80ce2aSJacob Faibussowitsch   if (namefields) CHKERRQ(NameFields(da,1));
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(da,0.0,Lx,0.0,Ly,0.0,Lz));
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreateCompatibleDMDA(da,dof,&daVector));
1825f80ce2aSJacob Faibussowitsch   if (namefields) CHKERRQ(NameFields(daVector,dof));
1835f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(da,&info));
1845f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&v));
1855f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(daVector,&vVector));
1865f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,v,&va));
1875f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayDOF(daVector,vVector,&vVectora));
188c4762a1bSJed Brown   for (j=info.ys; j<info.ys+info.ym; j++) {
189c4762a1bSJed Brown     for (i=info.xs; i<info.xs+info.xm; i++) {
190c4762a1bSJed Brown       const PetscScalar x = (Lx*i)/M;
191c4762a1bSJed Brown       const PetscScalar y = (Ly*j)/N;
192c4762a1bSJed Brown       va[j][i] = PetscPowScalarInt(x-0.5*Lx,2)+PetscPowScalarInt(y-0.5*Ly,2);
193c4762a1bSJed Brown       for (c=0; c<dof; ++c) {
194c4762a1bSJed Brown         vVectora[j][i][c] = PetscPowScalarInt(x-0.5*Lx,2)+PetscPowScalarInt(y-0.5*Ly,2) + 10.0*c;
195c4762a1bSJed Brown       }
196c4762a1bSJed Brown     }
197c4762a1bSJed Brown   }
1985f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,v,&va));
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayDOF(daVector,vVector,&vVectora));
2005f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view));
2015f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(v,view));
2025f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(vVector,view));
2035f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&view));
2045f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&v));
2055f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&vVector));
2065f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
2075f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&daVector));
208c4762a1bSJed Brown   return 0;
209c4762a1bSJed Brown }
210c4762a1bSJed Brown 
211c4762a1bSJed Brown int main(int argc, char *argv[])
212c4762a1bSJed Brown {
213c4762a1bSJed Brown   PetscInt       dof;
214c4762a1bSJed Brown   PetscBool      namefields;
215c4762a1bSJed Brown 
216*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,0,help));
217c4762a1bSJed Brown   dof = 2;
2185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-dof",&dof,NULL));
219c4762a1bSJed Brown   namefields = PETSC_FALSE;
2205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-namefields",&namefields,NULL));
2215f80ce2aSJacob Faibussowitsch   CHKERRQ(test_3d("3d.vtr",dof,namefields));
2225f80ce2aSJacob Faibussowitsch   CHKERRQ(test_2d("2d.vtr",dof,namefields));
2235f80ce2aSJacob Faibussowitsch   CHKERRQ(test_3d_compat("3d_compat.vtr",dof,namefields));
2245f80ce2aSJacob Faibussowitsch   CHKERRQ(test_2d_compat("2d_compat.vtr",dof,namefields));
2255f80ce2aSJacob Faibussowitsch   CHKERRQ(test_3d("3d.vts",dof,namefields));
2265f80ce2aSJacob Faibussowitsch   CHKERRQ(test_2d("2d.vts",dof,namefields));
2275f80ce2aSJacob Faibussowitsch   CHKERRQ(test_3d_compat("3d_compat.vts",dof,namefields));
2285f80ce2aSJacob Faibussowitsch   CHKERRQ(test_2d_compat("2d_compat.vts",dof,namefields));
229*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
230*b122ec5aSJacob Faibussowitsch   return 0;
231c4762a1bSJed Brown }
232c4762a1bSJed Brown 
233c4762a1bSJed Brown /*TEST
234c4762a1bSJed Brown 
235c4762a1bSJed Brown    build:
236c4762a1bSJed Brown       requires: !complex
237c4762a1bSJed Brown 
238c4762a1bSJed Brown    test:
239c4762a1bSJed Brown       suffix: 1
240c4762a1bSJed Brown       nsize: 2
241c4762a1bSJed Brown       args: -dof 2
242c4762a1bSJed Brown 
243c4762a1bSJed Brown    test:
244c4762a1bSJed Brown       suffix: 2
245c4762a1bSJed Brown       nsize: 2
246c4762a1bSJed Brown       args: -dof 2
247c4762a1bSJed Brown 
248c4762a1bSJed Brown    test:
249c4762a1bSJed Brown       suffix: 3
250c4762a1bSJed Brown       nsize: 2
251c4762a1bSJed Brown       args: -dof 3
252c4762a1bSJed Brown 
253c4762a1bSJed Brown    test:
254c4762a1bSJed Brown       suffix: 4
255c4762a1bSJed Brown       nsize: 1
256c4762a1bSJed Brown       args: -dof 2 -namefields
257c4762a1bSJed Brown 
258c4762a1bSJed Brown    test:
259c4762a1bSJed Brown       suffix: 5
260c4762a1bSJed Brown       nsize: 2
261c4762a1bSJed Brown       args: -dof 4 -namefields
262c4762a1bSJed Brown 
263c4762a1bSJed Brown TEST*/
264