xref: /petsc/src/dm/tests/ex42.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown /* -*- Mode: C++; c-basic-offset:2 ; indent-tabs-mode:nil ; -*- */
2c4762a1bSJed Brown 
3c4762a1bSJed Brown static char help[] = "Test VTK Rectilinear grid (.vtr) viewer support\n\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown #include <petscdm.h>
6c4762a1bSJed Brown #include <petscdmda.h>
7c4762a1bSJed Brown 
8c4762a1bSJed Brown /*
9c4762a1bSJed Brown   Write 3D DMDA vector with coordinates in VTK VTR format
10c4762a1bSJed Brown 
11c4762a1bSJed Brown */
12c4762a1bSJed Brown PetscErrorCode test_3d(const char filename[])
13c4762a1bSJed Brown {
14c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
15c4762a1bSJed Brown   const PetscInt    M=10, N=15, P=30, dof=1, sw=1;
16c4762a1bSJed Brown   const PetscScalar Lx=1.0, Ly=1.0, Lz=1.0;
17c4762a1bSJed Brown   DM                da;
18c4762a1bSJed Brown   Vec               v;
19c4762a1bSJed Brown   PetscViewer       view;
20c4762a1bSJed Brown   DMDALocalInfo     info;
21c4762a1bSJed Brown   PetscScalar       ***va;
22c4762a1bSJed Brown   PetscInt          i,j,k;
23c4762a1bSJed Brown 
245f80ce2aSJacob 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));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
27c4762a1bSJed Brown 
285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(da,0.0,Lx,0.0,Ly,0.0,Lz));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(da,&info));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&v));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,v,&va));
32c4762a1bSJed Brown   for (k=info.zs; k<info.zs+info.zm; k++) {
33c4762a1bSJed Brown     for (j=info.ys; j<info.ys+info.ym; j++) {
34c4762a1bSJed Brown       for (i=info.xs; i<info.xs+info.xm; i++) {
35c4762a1bSJed Brown         PetscScalar x = (Lx*i)/M;
36c4762a1bSJed Brown         PetscScalar y = (Ly*j)/N;
37c4762a1bSJed Brown         PetscScalar z = (Lz*k)/P;
38c4762a1bSJed Brown         va[k][j][i] = PetscPowScalarInt(x-0.5*Lx,2)+PetscPowScalarInt(y-0.5*Ly,2)+PetscPowScalarInt(z-0.5*Lz,2);
39c4762a1bSJed Brown       }
40c4762a1bSJed Brown     }
41c4762a1bSJed Brown   }
425f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,v,&va));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(v,view));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&view));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&v));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
48c4762a1bSJed Brown   return 0;
49c4762a1bSJed Brown }
50c4762a1bSJed Brown 
51c4762a1bSJed Brown /*
52c4762a1bSJed Brown   Write 2D DMDA vector with coordinates in VTK VTR format
53c4762a1bSJed Brown 
54c4762a1bSJed Brown */
55c4762a1bSJed Brown PetscErrorCode test_2d(const char filename[])
56c4762a1bSJed Brown {
57c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
58c4762a1bSJed Brown   const PetscInt    M=10, N=20, dof=1, sw=1;
59c4762a1bSJed Brown   const PetscScalar Lx=1.0, Ly=1.0, Lz=1.0;
60c4762a1bSJed Brown   DM                da;
61c4762a1bSJed Brown   Vec               v;
62c4762a1bSJed Brown   PetscViewer       view;
63c4762a1bSJed Brown   DMDALocalInfo     info;
64c4762a1bSJed Brown   PetscScalar       **va;
65c4762a1bSJed Brown   PetscInt          i,j;
66c4762a1bSJed Brown 
675f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR, M,N,PETSC_DECIDE,PETSC_DECIDE,dof,sw,NULL,NULL,&da));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
705f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(da,0.0,Lx,0.0,Ly,0.0,Lz));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(da,&info));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&v));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,v,&va));
74c4762a1bSJed Brown   for (j=info.ys; j<info.ys+info.ym; j++) {
75c4762a1bSJed Brown     for (i=info.xs; i<info.xs+info.xm; i++) {
76c4762a1bSJed Brown       PetscScalar x = (Lx*i)/M;
77c4762a1bSJed Brown       PetscScalar y = (Ly*j)/N;
78c4762a1bSJed Brown       va[j][i] = PetscPowScalarInt(x-0.5*Lx,2)+PetscPowScalarInt(y-0.5*Ly,2);
79c4762a1bSJed Brown     }
80c4762a1bSJed Brown   }
815f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,v,&va));
825f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(v,view));
845f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&view));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&v));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
87c4762a1bSJed Brown   return 0;
88c4762a1bSJed Brown }
89c4762a1bSJed Brown 
90c4762a1bSJed Brown /*
91c4762a1bSJed Brown   Write 2D DMDA vector without coordinates in VTK VTR format
92c4762a1bSJed Brown 
93c4762a1bSJed Brown */
94c4762a1bSJed Brown PetscErrorCode test_2d_nocoord(const char filename[])
95c4762a1bSJed Brown {
96c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
97c4762a1bSJed Brown   const PetscInt    M=10, N=20, dof=1, sw=1;
98c4762a1bSJed Brown   const PetscScalar Lx=1.0, Ly=1.0;
99c4762a1bSJed Brown   DM                da;
100c4762a1bSJed Brown   Vec               v;
101c4762a1bSJed Brown   PetscViewer       view;
102c4762a1bSJed Brown   DMDALocalInfo     info;
103c4762a1bSJed Brown   PetscScalar       **va;
104c4762a1bSJed Brown   PetscInt          i,j;
105c4762a1bSJed Brown 
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(comm,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,M,N,PETSC_DECIDE,PETSC_DECIDE,dof,sw,NULL,NULL,&da));
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(da,&info));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&v));
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,v,&va));
112c4762a1bSJed Brown   for (j=info.ys; j<info.ys+info.ym; j++) {
113c4762a1bSJed Brown     for (i=info.xs; i<info.xs+info.xm; i++) {
114c4762a1bSJed Brown       PetscScalar x = (Lx*i)/M;
115c4762a1bSJed Brown       PetscScalar y = (Ly*j)/N;
116c4762a1bSJed Brown       va[j][i] = PetscPowScalarInt(x-0.5*Lx,2)+PetscPowScalarInt(y-0.5*Ly,2);
117c4762a1bSJed Brown     }
118c4762a1bSJed Brown   }
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,v,&va));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view));
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(v,view));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&view));
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&v));
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
125c4762a1bSJed Brown   return 0;
126c4762a1bSJed Brown }
127c4762a1bSJed Brown 
128c4762a1bSJed Brown /*
129c4762a1bSJed Brown   Write 3D DMDA vector without coordinates in VTK VTR format
130c4762a1bSJed Brown 
131c4762a1bSJed Brown */
132c4762a1bSJed Brown PetscErrorCode test_3d_nocoord(const char filename[])
133c4762a1bSJed Brown {
134c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
135c4762a1bSJed Brown   const PetscInt    M=10, N=20, P=30, dof=1, sw=1;
136c4762a1bSJed Brown   const PetscScalar Lx=1.0, Ly=1.0, Lz=1.0;
137c4762a1bSJed Brown   DM                da;
138c4762a1bSJed Brown   Vec               v;
139c4762a1bSJed Brown   PetscViewer       view;
140c4762a1bSJed Brown   DMDALocalInfo     info;
141c4762a1bSJed Brown   PetscScalar       ***va;
142c4762a1bSJed Brown   PetscInt          i,j,k;
143c4762a1bSJed Brown 
1445f80ce2aSJacob 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));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
147c4762a1bSJed Brown 
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(da,&info));
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&v));
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,v,&va));
151c4762a1bSJed Brown   for (k=info.zs; k<info.zs+info.zm; k++) {
152c4762a1bSJed Brown     for (j=info.ys; j<info.ys+info.ym; j++) {
153c4762a1bSJed Brown       for (i=info.xs; i<info.xs+info.xm; i++) {
154c4762a1bSJed Brown         PetscScalar x = (Lx*i)/M;
155c4762a1bSJed Brown         PetscScalar y = (Ly*j)/N;
156c4762a1bSJed Brown         PetscScalar z = (Lz*k)/P;
157c4762a1bSJed Brown         va[k][j][i] = PetscPowScalarInt(x-0.5*Lx,2)+PetscPowScalarInt(y-0.5*Ly,2)+PetscPowScalarInt(z-0.5*Lz,2);
158c4762a1bSJed Brown       }
159c4762a1bSJed Brown     }
160c4762a1bSJed Brown   }
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,v,&va));
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view));
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(v,view));
1645f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&view));
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&v));
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
167c4762a1bSJed Brown   return 0;
168c4762a1bSJed Brown }
169c4762a1bSJed Brown 
170c4762a1bSJed Brown int main(int argc, char *argv[])
171c4762a1bSJed Brown {
172c4762a1bSJed Brown 
173*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,0,help));
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(test_3d("3d.vtr"));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(test_2d("2d.vtr"));
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(test_2d_nocoord("2d_nocoord.vtr"));
1775f80ce2aSJacob Faibussowitsch   CHKERRQ(test_3d_nocoord("3d_nocoord.vtr"));
178*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
179*b122ec5aSJacob Faibussowitsch   return 0;
180c4762a1bSJed Brown }
181c4762a1bSJed Brown 
182c4762a1bSJed Brown /*TEST
183c4762a1bSJed Brown 
184c4762a1bSJed Brown    build:
185c4762a1bSJed Brown       requires: !complex
186c4762a1bSJed Brown 
187c4762a1bSJed Brown    test:
188c4762a1bSJed Brown       nsize: 2
189c4762a1bSJed Brown 
190c4762a1bSJed Brown TEST*/
191