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