xref: /petsc/src/dm/tests/ex47.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown static char help[] = "Test VTK structured grid (.vts) viewer support\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscdm.h>
4c4762a1bSJed Brown #include <petscdmda.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown /*
7c4762a1bSJed Brown   Write 3D DMDA vector with coordinates in VTK .vts format
8c4762a1bSJed Brown 
9c4762a1bSJed Brown */
10c4762a1bSJed Brown PetscErrorCode test_3d(const char filename[])
11c4762a1bSJed Brown {
12c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
13c4762a1bSJed Brown   const PetscInt    M=10,N=15,P=30,dof=1,sw=1;
14c4762a1bSJed Brown   const PetscScalar Lx=1.0,Ly=1.0,Lz=1.0;
15c4762a1bSJed Brown   DM                da;
16c4762a1bSJed Brown   Vec               v;
17c4762a1bSJed Brown   PetscViewer       view;
18c4762a1bSJed Brown   DMDALocalInfo     info;
19c4762a1bSJed Brown   PetscScalar       ***va;
20c4762a1bSJed Brown   PetscInt          i,j,k;
21c4762a1bSJed Brown 
22*5f80ce2aSJacob 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));
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
25c4762a1bSJed Brown 
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(da,0.0,Lx,0.0,Ly,0.0,Lz));
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(da,&info));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&v));
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,v,&va));
30c4762a1bSJed Brown   for (k=info.zs; k<info.zs+info.zm; k++) {
31c4762a1bSJed Brown     for (j=info.ys; j<info.ys+info.ym; j++) {
32c4762a1bSJed Brown       for (i=info.xs; i<info.xs+info.xm; i++) {
33c4762a1bSJed Brown         PetscScalar x = (Lx*i)/M;
34c4762a1bSJed Brown         PetscScalar y = (Ly*j)/N;
35c4762a1bSJed Brown         PetscScalar z = (Lz*k)/P;
36c4762a1bSJed Brown         va[k][j][i] = PetscPowScalarInt(x-0.5*Lx,2)+PetscPowScalarInt(y-0.5*Ly,2)+PetscPowScalarInt(z-0.5*Lz,2);
37c4762a1bSJed Brown       }
38c4762a1bSJed Brown     }
39c4762a1bSJed Brown   }
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,v,&va));
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(v,view));
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&view));
44*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&v));
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
46c4762a1bSJed Brown   return 0;
47c4762a1bSJed Brown }
48c4762a1bSJed Brown 
49c4762a1bSJed Brown /*
50c4762a1bSJed Brown   Write 2D DMDA vector with coordinates in VTK .vts format
51c4762a1bSJed Brown 
52c4762a1bSJed Brown */
53c4762a1bSJed Brown PetscErrorCode test_2d(const char filename[])
54c4762a1bSJed Brown {
55c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
56c4762a1bSJed Brown   const PetscInt    M=10,N=20,dof=1,sw=1;
57c4762a1bSJed Brown   const PetscScalar Lx=1.0,Ly=1.0,Lz=1.0;
58c4762a1bSJed Brown   DM                da;
59c4762a1bSJed Brown   Vec               v;
60c4762a1bSJed Brown   PetscViewer       view;
61c4762a1bSJed Brown   DMDALocalInfo     info;
62c4762a1bSJed Brown   PetscScalar       **va;
63c4762a1bSJed Brown   PetscInt          i,j;
64c4762a1bSJed Brown 
65*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR, M,N,PETSC_DECIDE,PETSC_DECIDE,dof,sw,NULL,NULL,&da));
66*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
67*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
68*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(da,0.0,Lx,0.0,Ly,0.0,Lz));
69*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(da,&info));
70*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&v));
71*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,v,&va));
72c4762a1bSJed Brown   for (j=info.ys; j<info.ys+info.ym; j++) {
73c4762a1bSJed Brown     for (i=info.xs; i<info.xs+info.xm; i++) {
74c4762a1bSJed Brown       PetscScalar x = (Lx*i)/M;
75c4762a1bSJed Brown       PetscScalar y = (Ly*j)/N;
76c4762a1bSJed Brown       va[j][i] = PetscPowScalarInt(x-0.5*Lx,2)+PetscPowScalarInt(y-0.5*Ly,2);
77c4762a1bSJed Brown     }
78c4762a1bSJed Brown   }
79*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,v,&va));
80*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view));
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(v,view));
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&view));
83*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&v));
84*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
85c4762a1bSJed Brown   return 0;
86c4762a1bSJed Brown }
87c4762a1bSJed Brown 
88c4762a1bSJed Brown /*
89c4762a1bSJed Brown   Write 2D DMDA vector without coordinates in VTK .vts format
90c4762a1bSJed Brown 
91c4762a1bSJed Brown */
92c4762a1bSJed Brown PetscErrorCode test_2d_nocoord(const char filename[])
93c4762a1bSJed Brown {
94c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
95c4762a1bSJed Brown   const PetscInt    M=10,N=20,dof=1,sw=1;
96c4762a1bSJed Brown   const PetscScalar Lx=1.0,Ly=1.0;
97c4762a1bSJed Brown   DM                da;
98c4762a1bSJed Brown   Vec               v;
99c4762a1bSJed Brown   PetscViewer       view;
100c4762a1bSJed Brown   DMDALocalInfo     info;
101c4762a1bSJed Brown   PetscScalar       **va;
102c4762a1bSJed Brown   PetscInt          i,j;
103c4762a1bSJed Brown 
104*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,PETSC_DECIDE,PETSC_DECIDE,dof,sw,NULL,NULL,&da));
105*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
106*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
107*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(da,&info));
108*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&v));
109*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,v,&va));
110c4762a1bSJed Brown   for (j=info.ys; j<info.ys+info.ym; j++) {
111c4762a1bSJed Brown     for (i=info.xs; i<info.xs+info.xm; i++) {
112c4762a1bSJed Brown       PetscScalar x = (Lx*i)/M;
113c4762a1bSJed Brown       PetscScalar y = (Ly*j)/N;
114c4762a1bSJed Brown       va[j][i] = PetscPowScalarInt(x-0.5*Lx,2)+PetscPowScalarInt(y-0.5*Ly,2);
115c4762a1bSJed Brown     }
116c4762a1bSJed Brown   }
117*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,v,&va));
118*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view));
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(v,view));
120*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&view));
121*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&v));
122*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
123c4762a1bSJed Brown   return 0;
124c4762a1bSJed Brown }
125c4762a1bSJed Brown 
126c4762a1bSJed Brown /*
127c4762a1bSJed Brown   Write 3D DMDA vector without coordinates in VTK .vts format
128c4762a1bSJed Brown 
129c4762a1bSJed Brown */
130c4762a1bSJed Brown PetscErrorCode test_3d_nocoord(const char filename[])
131c4762a1bSJed Brown {
132c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
133c4762a1bSJed Brown   const PetscInt    M=10,N=20,P=30,dof=1,sw=1;
134c4762a1bSJed Brown   const PetscScalar Lx=1.0,Ly=1.0,Lz=1.0;
135c4762a1bSJed Brown   DM                da;
136c4762a1bSJed Brown   Vec               v;
137c4762a1bSJed Brown   PetscViewer       view;
138c4762a1bSJed Brown   DMDALocalInfo     info;
139c4762a1bSJed Brown   PetscScalar       ***va;
140c4762a1bSJed Brown   PetscInt          i,j,k;
141c4762a1bSJed Brown 
142*5f80ce2aSJacob 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));
143*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
144*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
145c4762a1bSJed Brown 
146*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(da,&info));
147*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&v));
148*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,v,&va));
149c4762a1bSJed Brown   for (k=info.zs; k<info.zs+info.zm; k++) {
150c4762a1bSJed Brown     for (j=info.ys; j<info.ys+info.ym; j++) {
151c4762a1bSJed Brown       for (i=info.xs; i<info.xs+info.xm; i++) {
152c4762a1bSJed Brown         PetscScalar x = (Lx*i)/M;
153c4762a1bSJed Brown         PetscScalar y = (Ly*j)/N;
154c4762a1bSJed Brown         PetscScalar z = (Lz*k)/P;
155c4762a1bSJed Brown         va[k][j][i] = PetscPowScalarInt(x-0.5*Lx,2)+PetscPowScalarInt(y-0.5*Ly,2)+PetscPowScalarInt(z-0.5*Lz,2);
156c4762a1bSJed Brown       }
157c4762a1bSJed Brown     }
158c4762a1bSJed Brown   }
159*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,v,&va));
160*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view));
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(v,view));
162*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&view));
163*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&v));
164*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
165c4762a1bSJed Brown   return 0;
166c4762a1bSJed Brown }
167c4762a1bSJed Brown 
168c4762a1bSJed Brown int main(int argc, char *argv[])
169c4762a1bSJed Brown {
170c4762a1bSJed Brown   PetscErrorCode ierr;
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
173*5f80ce2aSJacob Faibussowitsch   CHKERRQ(test_3d("3d.vts"));
174*5f80ce2aSJacob Faibussowitsch   CHKERRQ(test_2d("2d.vts"));
175*5f80ce2aSJacob Faibussowitsch   CHKERRQ(test_2d_nocoord("2d_nocoord.vts"));
176*5f80ce2aSJacob Faibussowitsch   CHKERRQ(test_3d_nocoord("3d_nocoord.vts"));
177c4762a1bSJed Brown   ierr = PetscFinalize();
178c4762a1bSJed Brown   return ierr;
179c4762a1bSJed Brown }
180c4762a1bSJed Brown 
181c4762a1bSJed Brown /*TEST
182c4762a1bSJed Brown 
183c4762a1bSJed Brown    build:
184c4762a1bSJed Brown       requires: !complex
185c4762a1bSJed Brown 
186c4762a1bSJed Brown    test:
187c4762a1bSJed Brown       nsize: 2
188c4762a1bSJed Brown 
189c4762a1bSJed Brown TEST*/
190