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