xref: /petsc/src/dm/tests/ex48.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Test VTK structured (.vts)  and rectilinear (.vtr) viewer support with multi-dof DMDAs.\n\
2*c4762a1bSJed Brown                       Supply the -namefields flag to test with field names.\n\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown #include <petscdm.h>
5*c4762a1bSJed Brown #include <petscdmda.h>
6*c4762a1bSJed Brown 
7*c4762a1bSJed Brown /* Helper function to name DMDA fields */
8*c4762a1bSJed Brown PetscErrorCode NameFields(DM da,PetscInt dof)
9*c4762a1bSJed Brown {
10*c4762a1bSJed Brown   PetscErrorCode ierr;
11*c4762a1bSJed Brown   PetscInt       c;
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown   PetscFunctionBeginUser;
14*c4762a1bSJed Brown   for (c=0; c<dof; ++c) {
15*c4762a1bSJed Brown     char fieldname[256];
16*c4762a1bSJed Brown     ierr = PetscSNPrintf(fieldname,sizeof(fieldname),"field_%D",c);CHKERRQ(ierr);
17*c4762a1bSJed Brown     ierr = DMDASetFieldName(da,c,fieldname);CHKERRQ(ierr);
18*c4762a1bSJed Brown   }
19*c4762a1bSJed Brown   PetscFunctionReturn(0);
20*c4762a1bSJed Brown }
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown /*
23*c4762a1bSJed Brown   Write 3D DMDA vector with coordinates in VTK format
24*c4762a1bSJed Brown */
25*c4762a1bSJed Brown PetscErrorCode test_3d(const char filename[],PetscInt dof,PetscBool namefields)
26*c4762a1bSJed Brown {
27*c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
28*c4762a1bSJed Brown   const PetscInt    M=10,N=15,P=30,sw=1;
29*c4762a1bSJed Brown   const PetscScalar Lx=1.0,Ly=1.0,Lz=1.0;
30*c4762a1bSJed Brown   DM                da;
31*c4762a1bSJed Brown   Vec               v;
32*c4762a1bSJed Brown   PetscViewer       view;
33*c4762a1bSJed Brown   DMDALocalInfo     info;
34*c4762a1bSJed Brown   PetscScalar       ****va;
35*c4762a1bSJed Brown   PetscInt          i,j,k,c;
36*c4762a1bSJed Brown   PetscErrorCode    ierr;
37*c4762a1bSJed Brown 
38*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);
39*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
40*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
41*c4762a1bSJed Brown   if (namefields) {ierr = NameFields(da,dof);CHKERRQ(ierr);}
42*c4762a1bSJed Brown 
43*c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(da,0.0,Lx,0.0,Ly,0.0,Lz);CHKERRQ(ierr);
44*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
45*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&v);CHKERRQ(ierr);
46*c4762a1bSJed Brown   ierr = DMDAVecGetArrayDOF(da,v,&va);CHKERRQ(ierr);
47*c4762a1bSJed Brown   for (k=info.zs; k<info.zs+info.zm; k++) {
48*c4762a1bSJed Brown     for (j=info.ys; j<info.ys+info.ym; j++) {
49*c4762a1bSJed Brown       for (i=info.xs; i<info.xs+info.xm; i++) {
50*c4762a1bSJed Brown         const PetscScalar x = (Lx*i)/M;
51*c4762a1bSJed Brown         const PetscScalar y = (Ly*j)/N;
52*c4762a1bSJed Brown         const PetscScalar z = (Lz*k)/P;
53*c4762a1bSJed Brown         for (c=0; c<dof; ++ c) {
54*c4762a1bSJed 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;
55*c4762a1bSJed Brown         }
56*c4762a1bSJed Brown       }
57*c4762a1bSJed Brown     }
58*c4762a1bSJed Brown   }
59*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayDOF(da,v,&va);CHKERRQ(ierr);
60*c4762a1bSJed Brown   ierr = PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view);CHKERRQ(ierr);
61*c4762a1bSJed Brown   ierr = VecView(v,view);CHKERRQ(ierr);
62*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&view);CHKERRQ(ierr);
63*c4762a1bSJed Brown   ierr = VecDestroy(&v);CHKERRQ(ierr);
64*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
65*c4762a1bSJed Brown   return 0;
66*c4762a1bSJed Brown }
67*c4762a1bSJed Brown 
68*c4762a1bSJed Brown 
69*c4762a1bSJed Brown /*
70*c4762a1bSJed Brown   Write 2D DMDA vector with coordinates in VTK format
71*c4762a1bSJed Brown */
72*c4762a1bSJed Brown PetscErrorCode test_2d(const char filename[],PetscInt dof,PetscBool namefields)
73*c4762a1bSJed Brown {
74*c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
75*c4762a1bSJed Brown   const PetscInt    M=10,N=20,sw=1;
76*c4762a1bSJed Brown   const PetscScalar Lx=1.0,Ly=1.0,Lz=1.0;
77*c4762a1bSJed Brown   DM                da;
78*c4762a1bSJed Brown   Vec               v;
79*c4762a1bSJed Brown   PetscViewer       view;
80*c4762a1bSJed Brown   DMDALocalInfo     info;
81*c4762a1bSJed Brown   PetscScalar       ***va;
82*c4762a1bSJed Brown   PetscInt          i,j,c;
83*c4762a1bSJed Brown   PetscErrorCode    ierr;
84*c4762a1bSJed Brown 
85*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);
86*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
87*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
88*c4762a1bSJed Brown   if (namefields) {ierr = NameFields(da,dof);CHKERRQ(ierr);}
89*c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(da,0.0,Lx,0.0,Ly,0.0,Lz);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&v);CHKERRQ(ierr);
92*c4762a1bSJed Brown   ierr = DMDAVecGetArrayDOF(da,v,&va);CHKERRQ(ierr);
93*c4762a1bSJed Brown   for (j=info.ys; j<info.ys+info.ym; j++) {
94*c4762a1bSJed Brown     for (i=info.xs; i<info.xs+info.xm; i++) {
95*c4762a1bSJed Brown       const PetscScalar x = (Lx*i)/M;
96*c4762a1bSJed Brown       const PetscScalar y = (Ly*j)/N;
97*c4762a1bSJed Brown       for (c=0; c<dof; ++c) {
98*c4762a1bSJed Brown         va[j][i][c] = PetscPowScalarInt(x-0.5*Lx,2)+PetscPowScalarInt(y-0.5*Ly,2) + 10.0*c;
99*c4762a1bSJed Brown       }
100*c4762a1bSJed Brown     }
101*c4762a1bSJed Brown   }
102*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayDOF(da,v,&va);CHKERRQ(ierr);
103*c4762a1bSJed Brown   ierr = PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view);CHKERRQ(ierr);
104*c4762a1bSJed Brown   ierr = VecView(v,view);CHKERRQ(ierr);
105*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&view);CHKERRQ(ierr);
106*c4762a1bSJed Brown   ierr = VecDestroy(&v);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
108*c4762a1bSJed Brown   return 0;
109*c4762a1bSJed Brown }
110*c4762a1bSJed Brown 
111*c4762a1bSJed Brown /*
112*c4762a1bSJed Brown   Write a scalar and a vector field from two compatible 3d DMDAs
113*c4762a1bSJed Brown */
114*c4762a1bSJed Brown PetscErrorCode test_3d_compat(const char filename[],PetscInt dof,PetscBool namefields)
115*c4762a1bSJed Brown {
116*c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
117*c4762a1bSJed Brown   const PetscInt    M=10,N=15,P=30,sw=1;
118*c4762a1bSJed Brown   const PetscScalar Lx=1.0,Ly=1.0,Lz=1.0;
119*c4762a1bSJed Brown   DM                da,daVector;
120*c4762a1bSJed Brown   Vec               v,vVector;
121*c4762a1bSJed Brown   PetscViewer       view;
122*c4762a1bSJed Brown   DMDALocalInfo     info;
123*c4762a1bSJed Brown   PetscScalar       ***va,****vVectora;
124*c4762a1bSJed Brown   PetscInt          i,j,k,c;
125*c4762a1bSJed Brown   PetscErrorCode    ierr;
126*c4762a1bSJed Brown 
127*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:*/1,sw,NULL,NULL,NULL,&da);CHKERRQ(ierr);
128*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
129*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
130*c4762a1bSJed Brown   if (namefields) {ierr = NameFields(da,1);CHKERRQ(ierr);}
131*c4762a1bSJed Brown 
132*c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(da,0.0,Lx,0.0,Ly,0.0,Lz);CHKERRQ(ierr);
133*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
134*c4762a1bSJed Brown   ierr = DMDACreateCompatibleDMDA(da,dof,&daVector);CHKERRQ(ierr);
135*c4762a1bSJed Brown   if (namefields) {ierr = NameFields(daVector,dof);CHKERRQ(ierr);}
136*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&v);CHKERRQ(ierr);
137*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(daVector,&vVector);CHKERRQ(ierr);
138*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,v,&va);CHKERRQ(ierr);
139*c4762a1bSJed Brown   ierr = DMDAVecGetArrayDOF(daVector,vVector,&vVectora);CHKERRQ(ierr);
140*c4762a1bSJed Brown   for (k=info.zs; k<info.zs+info.zm; k++) {
141*c4762a1bSJed Brown     for (j=info.ys; j<info.ys+info.ym; j++) {
142*c4762a1bSJed Brown       for (i=info.xs; i<info.xs+info.xm; i++) {
143*c4762a1bSJed Brown         const PetscScalar x = (Lx*i)/M;
144*c4762a1bSJed Brown         const PetscScalar y = (Ly*j)/N;
145*c4762a1bSJed Brown         const PetscScalar z = (Lz*k)/P;
146*c4762a1bSJed Brown         va[k][j][i] = PetscPowScalarInt(x-0.5*Lx,2)+PetscPowScalarInt(y-0.5*Ly,2)+PetscPowScalarInt(z-0.5*Lz,2);
147*c4762a1bSJed Brown         for (c=0; c<dof; ++c) {
148*c4762a1bSJed 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;
149*c4762a1bSJed Brown         }
150*c4762a1bSJed Brown       }
151*c4762a1bSJed Brown     }
152*c4762a1bSJed Brown   }
153*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,v,&va);CHKERRQ(ierr);
154*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayDOF(da,v,&vVectora);CHKERRQ(ierr);
155*c4762a1bSJed Brown   ierr = PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view);CHKERRQ(ierr);
156*c4762a1bSJed Brown   ierr = VecView(v,view);CHKERRQ(ierr);
157*c4762a1bSJed Brown   ierr = VecView(vVector,view);CHKERRQ(ierr);
158*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&view);CHKERRQ(ierr);
159*c4762a1bSJed Brown   ierr = VecDestroy(&v);CHKERRQ(ierr);
160*c4762a1bSJed Brown   ierr = VecDestroy(&vVector);CHKERRQ(ierr);
161*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
162*c4762a1bSJed Brown   ierr = DMDestroy(&daVector);CHKERRQ(ierr);
163*c4762a1bSJed Brown   return 0;
164*c4762a1bSJed Brown }
165*c4762a1bSJed Brown 
166*c4762a1bSJed Brown /*
167*c4762a1bSJed Brown   Write a scalar and a vector field from two compatible 2d DMDAs
168*c4762a1bSJed Brown */
169*c4762a1bSJed Brown PetscErrorCode test_2d_compat(const char filename[],PetscInt dof,PetscBool namefields)
170*c4762a1bSJed Brown {
171*c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
172*c4762a1bSJed Brown   const PetscInt    M=10,N=20,sw=1;
173*c4762a1bSJed Brown   const PetscScalar Lx=1.0,Ly=1.0,Lz=1.0;
174*c4762a1bSJed Brown   DM                da, daVector;
175*c4762a1bSJed Brown   Vec               v,vVector;
176*c4762a1bSJed Brown   PetscViewer       view;
177*c4762a1bSJed Brown   DMDALocalInfo     info;
178*c4762a1bSJed Brown   PetscScalar       **va,***vVectora;
179*c4762a1bSJed Brown   PetscInt          i,j,c;
180*c4762a1bSJed Brown   PetscErrorCode    ierr;
181*c4762a1bSJed Brown 
182*c4762a1bSJed Brown   ierr = DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR, M,N,PETSC_DECIDE,PETSC_DECIDE,/* dof:*/ 1,sw,NULL,NULL,&da);CHKERRQ(ierr);
183*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
184*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
185*c4762a1bSJed Brown   if (namefields) {ierr = NameFields(da,1);CHKERRQ(ierr);}
186*c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(da,0.0,Lx,0.0,Ly,0.0,Lz);CHKERRQ(ierr);
187*c4762a1bSJed Brown   ierr = DMDACreateCompatibleDMDA(da,dof,&daVector);CHKERRQ(ierr);
188*c4762a1bSJed Brown   if (namefields) {ierr = NameFields(daVector,dof);CHKERRQ(ierr);}
189*c4762a1bSJed Brown   ierr = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
190*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&v);CHKERRQ(ierr);
191*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(daVector,&vVector);CHKERRQ(ierr);
192*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,v,&va);CHKERRQ(ierr);
193*c4762a1bSJed Brown   ierr = DMDAVecGetArrayDOF(daVector,vVector,&vVectora);CHKERRQ(ierr);
194*c4762a1bSJed Brown   for (j=info.ys; j<info.ys+info.ym; j++) {
195*c4762a1bSJed Brown     for (i=info.xs; i<info.xs+info.xm; i++) {
196*c4762a1bSJed Brown       const PetscScalar x = (Lx*i)/M;
197*c4762a1bSJed Brown       const PetscScalar y = (Ly*j)/N;
198*c4762a1bSJed Brown       va[j][i] = PetscPowScalarInt(x-0.5*Lx,2)+PetscPowScalarInt(y-0.5*Ly,2);
199*c4762a1bSJed Brown       for (c=0; c<dof; ++c) {
200*c4762a1bSJed Brown         vVectora[j][i][c] = PetscPowScalarInt(x-0.5*Lx,2)+PetscPowScalarInt(y-0.5*Ly,2) + 10.0*c;
201*c4762a1bSJed Brown       }
202*c4762a1bSJed Brown     }
203*c4762a1bSJed Brown   }
204*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,v,&va);CHKERRQ(ierr);
205*c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayDOF(daVector,vVector,&vVectora);CHKERRQ(ierr);
206*c4762a1bSJed Brown   ierr = PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view);CHKERRQ(ierr);
207*c4762a1bSJed Brown   ierr = VecView(v,view);CHKERRQ(ierr);
208*c4762a1bSJed Brown   ierr = VecView(vVector,view);CHKERRQ(ierr);
209*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&view);CHKERRQ(ierr);
210*c4762a1bSJed Brown   ierr = VecDestroy(&v);CHKERRQ(ierr);
211*c4762a1bSJed Brown   ierr = VecDestroy(&vVector);CHKERRQ(ierr);
212*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
213*c4762a1bSJed Brown   ierr = DMDestroy(&daVector);CHKERRQ(ierr);
214*c4762a1bSJed Brown   return 0;
215*c4762a1bSJed Brown }
216*c4762a1bSJed Brown 
217*c4762a1bSJed Brown int main(int argc, char *argv[])
218*c4762a1bSJed Brown {
219*c4762a1bSJed Brown   PetscErrorCode ierr;
220*c4762a1bSJed Brown   PetscInt       dof;
221*c4762a1bSJed Brown   PetscBool      namefields;
222*c4762a1bSJed Brown 
223*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr;
224*c4762a1bSJed Brown   dof = 2;
225*c4762a1bSJed Brown   ierr = PetscOptionsGetInt(NULL,NULL,"-dof",&dof,NULL);CHKERRQ(ierr);
226*c4762a1bSJed Brown   namefields = PETSC_FALSE;
227*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-namefields",&namefields,NULL);CHKERRQ(ierr);
228*c4762a1bSJed Brown   ierr = test_3d("3d.vtr",dof,namefields);CHKERRQ(ierr);
229*c4762a1bSJed Brown   ierr = test_2d("2d.vtr",dof,namefields);CHKERRQ(ierr);
230*c4762a1bSJed Brown   ierr = test_3d_compat("3d_compat.vtr",dof,namefields);CHKERRQ(ierr);
231*c4762a1bSJed Brown   ierr = test_2d_compat("2d_compat.vtr",dof,namefields);CHKERRQ(ierr);
232*c4762a1bSJed Brown   ierr = test_3d("3d.vts",dof,namefields);CHKERRQ(ierr);
233*c4762a1bSJed Brown   ierr = test_2d("2d.vts",dof,namefields);CHKERRQ(ierr);
234*c4762a1bSJed Brown   ierr = test_3d_compat("3d_compat.vts",dof,namefields);CHKERRQ(ierr);
235*c4762a1bSJed Brown   ierr = test_2d_compat("2d_compat.vts",dof,namefields);CHKERRQ(ierr);
236*c4762a1bSJed Brown   ierr = PetscFinalize();
237*c4762a1bSJed Brown   return ierr;
238*c4762a1bSJed Brown }
239*c4762a1bSJed Brown 
240*c4762a1bSJed Brown /*TEST
241*c4762a1bSJed Brown 
242*c4762a1bSJed Brown    build:
243*c4762a1bSJed Brown       requires: !complex
244*c4762a1bSJed Brown 
245*c4762a1bSJed Brown    test:
246*c4762a1bSJed Brown       suffix: 1
247*c4762a1bSJed Brown       nsize: 2
248*c4762a1bSJed Brown       args: -dof 2
249*c4762a1bSJed Brown 
250*c4762a1bSJed Brown    test:
251*c4762a1bSJed Brown       suffix: 2
252*c4762a1bSJed Brown       nsize: 2
253*c4762a1bSJed Brown       args: -dof 2
254*c4762a1bSJed Brown 
255*c4762a1bSJed Brown    test:
256*c4762a1bSJed Brown       suffix: 3
257*c4762a1bSJed Brown       nsize: 2
258*c4762a1bSJed Brown       args: -dof 3
259*c4762a1bSJed Brown 
260*c4762a1bSJed Brown    test:
261*c4762a1bSJed Brown       suffix: 4
262*c4762a1bSJed Brown       nsize: 1
263*c4762a1bSJed Brown       args: -dof 2 -namefields
264*c4762a1bSJed Brown 
265*c4762a1bSJed Brown    test:
266*c4762a1bSJed Brown       suffix: 5
267*c4762a1bSJed Brown       nsize: 2
268*c4762a1bSJed Brown       args: -dof 4 -namefields
269*c4762a1bSJed Brown 
270*c4762a1bSJed Brown TEST*/
271