xref: /petsc/src/dm/tests/ex48.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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];
1563a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(fieldname,sizeof(fieldname),"field_%" PetscInt_FMT,c));
169566063dSJacob Faibussowitsch     PetscCall(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 
369566063dSJacob Faibussowitsch   PetscCall(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));
379566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
389566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
399566063dSJacob Faibussowitsch   if (namefields) PetscCall(NameFields(da,dof));
40c4762a1bSJed Brown 
419566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da,0.0,Lx,0.0,Ly,0.0,Lz));
429566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
439566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da,&v));
449566063dSJacob Faibussowitsch   PetscCall(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   }
579566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOF(da,v,&va));
589566063dSJacob Faibussowitsch   PetscCall(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view));
599566063dSJacob Faibussowitsch   PetscCall(VecView(v,view));
609566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&view));
619566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v));
629566063dSJacob Faibussowitsch   PetscCall(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 
819566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR, M,N,PETSC_DECIDE,PETSC_DECIDE,dof,sw,NULL,NULL,&da));
829566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
839566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
849566063dSJacob Faibussowitsch   if (namefields) PetscCall(NameFields(da,dof));
859566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da,0.0,Lx,0.0,Ly,0.0,Lz));
869566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
879566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da,&v));
889566063dSJacob Faibussowitsch   PetscCall(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   }
989566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOF(da,v,&va));
999566063dSJacob Faibussowitsch   PetscCall(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view));
1009566063dSJacob Faibussowitsch   PetscCall(VecView(v,view));
1019566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&view));
1029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v));
1039566063dSJacob Faibussowitsch   PetscCall(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 
1229566063dSJacob Faibussowitsch   PetscCall(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));
1239566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
1249566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
1259566063dSJacob Faibussowitsch   if (namefields) PetscCall(NameFields(da,1));
126c4762a1bSJed Brown 
1279566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da,0.0,Lx,0.0,Ly,0.0,Lz));
1289566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
1299566063dSJacob Faibussowitsch   PetscCall(DMDACreateCompatibleDMDA(da,dof,&daVector));
1309566063dSJacob Faibussowitsch   if (namefields) PetscCall(NameFields(daVector,dof));
1319566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da,&v));
1329566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(daVector,&vVector));
1339566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,v,&va));
1349566063dSJacob Faibussowitsch   PetscCall(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   }
1489566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,v,&va));
1499566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOF(da,v,&vVectora));
1509566063dSJacob Faibussowitsch   PetscCall(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view));
1519566063dSJacob Faibussowitsch   PetscCall(VecView(v,view));
1529566063dSJacob Faibussowitsch   PetscCall(VecView(vVector,view));
1539566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&view));
1549566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v));
1559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vVector));
1569566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1579566063dSJacob Faibussowitsch   PetscCall(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 
1769566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR, M,N,PETSC_DECIDE,PETSC_DECIDE,/* dof:*/ 1,sw,NULL,NULL,&da));
1779566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
1789566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
1799566063dSJacob Faibussowitsch   if (namefields) PetscCall(NameFields(da,1));
1809566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da,0.0,Lx,0.0,Ly,0.0,Lz));
1819566063dSJacob Faibussowitsch   PetscCall(DMDACreateCompatibleDMDA(da,dof,&daVector));
1829566063dSJacob Faibussowitsch   if (namefields) PetscCall(NameFields(daVector,dof));
1839566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da,&info));
1849566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da,&v));
1859566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(daVector,&vVector));
1869566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da,v,&va));
1879566063dSJacob Faibussowitsch   PetscCall(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   }
1989566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da,v,&va));
1999566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOF(daVector,vVector,&vVectora));
2009566063dSJacob Faibussowitsch   PetscCall(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view));
2019566063dSJacob Faibussowitsch   PetscCall(VecView(v,view));
2029566063dSJacob Faibussowitsch   PetscCall(VecView(vVector,view));
2039566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&view));
2049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v));
2059566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vVector));
2069566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
2079566063dSJacob Faibussowitsch   PetscCall(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*327415f7SBarry Smith   PetscFunctionBeginUser;
2179566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,0,help));
218c4762a1bSJed Brown   dof = 2;
2199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL,"-dof",&dof,NULL));
220c4762a1bSJed Brown   namefields = PETSC_FALSE;
2219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-namefields",&namefields,NULL));
2229566063dSJacob Faibussowitsch   PetscCall(test_3d("3d.vtr",dof,namefields));
2239566063dSJacob Faibussowitsch   PetscCall(test_2d("2d.vtr",dof,namefields));
2249566063dSJacob Faibussowitsch   PetscCall(test_3d_compat("3d_compat.vtr",dof,namefields));
2259566063dSJacob Faibussowitsch   PetscCall(test_2d_compat("2d_compat.vtr",dof,namefields));
2269566063dSJacob Faibussowitsch   PetscCall(test_3d("3d.vts",dof,namefields));
2279566063dSJacob Faibussowitsch   PetscCall(test_2d("2d.vts",dof,namefields));
2289566063dSJacob Faibussowitsch   PetscCall(test_3d_compat("3d_compat.vts",dof,namefields));
2299566063dSJacob Faibussowitsch   PetscCall(test_2d_compat("2d_compat.vts",dof,namefields));
2309566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
231b122ec5aSJacob Faibussowitsch   return 0;
232c4762a1bSJed Brown }
233c4762a1bSJed Brown 
234c4762a1bSJed Brown /*TEST
235c4762a1bSJed Brown 
236c4762a1bSJed Brown    build:
237c4762a1bSJed Brown       requires: !complex
238c4762a1bSJed Brown 
239c4762a1bSJed Brown    test:
240c4762a1bSJed Brown       suffix: 1
241c4762a1bSJed Brown       nsize: 2
242c4762a1bSJed Brown       args: -dof 2
243c4762a1bSJed Brown 
244c4762a1bSJed Brown    test:
245c4762a1bSJed Brown       suffix: 2
246c4762a1bSJed Brown       nsize: 2
247c4762a1bSJed Brown       args: -dof 2
248c4762a1bSJed Brown 
249c4762a1bSJed Brown    test:
250c4762a1bSJed Brown       suffix: 3
251c4762a1bSJed Brown       nsize: 2
252c4762a1bSJed Brown       args: -dof 3
253c4762a1bSJed Brown 
254c4762a1bSJed Brown    test:
255c4762a1bSJed Brown       suffix: 4
256c4762a1bSJed Brown       nsize: 1
257c4762a1bSJed Brown       args: -dof 2 -namefields
258c4762a1bSJed Brown 
259c4762a1bSJed Brown    test:
260c4762a1bSJed Brown       suffix: 5
261c4762a1bSJed Brown       nsize: 2
262c4762a1bSJed Brown       args: -dof 4 -namefields
263c4762a1bSJed Brown 
264c4762a1bSJed Brown TEST*/
265