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 24*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)); 25*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 26*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 27c4762a1bSJed Brown 28*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(da,0.0,Lx,0.0,Ly,0.0,Lz)); 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&v)); 31*5f80ce2aSJacob 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 } 42*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,v,&va)); 43*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view)); 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(v,view)); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v)); 47*5f80ce2aSJacob 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 67*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)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 69*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDASetUniformCoordinates(da,0.0,Lx,0.0,Ly,0.0,Lz)); 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&v)); 73*5f80ce2aSJacob 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 } 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,v,&va)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(v,view)); 84*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 85*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v)); 86*5f80ce2aSJacob 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 106*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)); 107*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 108*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 109*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 110*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&v)); 111*5f80ce2aSJacob 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 } 119*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,v,&va)); 120*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view)); 121*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(v,view)); 122*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 123*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v)); 124*5f80ce2aSJacob 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 144*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)); 145*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetFromOptions(da)); 146*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMSetUp(da)); 147c4762a1bSJed Brown 148*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAGetLocalInfo(da,&info)); 149*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateGlobalVector(da,&v)); 150*5f80ce2aSJacob 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 } 161*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDAVecRestoreArray(da,v,&va)); 162*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerVTKOpen(comm,filename,FILE_MODE_WRITE,&view)); 163*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(v,view)); 164*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&view)); 165*5f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&v)); 166*5f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(&da)); 167c4762a1bSJed Brown return 0; 168c4762a1bSJed Brown } 169c4762a1bSJed Brown 170c4762a1bSJed Brown int main(int argc, char *argv[]) 171c4762a1bSJed Brown { 172c4762a1bSJed Brown PetscErrorCode ierr; 173c4762a1bSJed Brown 174c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,0,help);if (ierr) return ierr; 175*5f80ce2aSJacob Faibussowitsch CHKERRQ(test_3d("3d.vtr")); 176*5f80ce2aSJacob Faibussowitsch CHKERRQ(test_2d("2d.vtr")); 177*5f80ce2aSJacob Faibussowitsch CHKERRQ(test_2d_nocoord("2d_nocoord.vtr")); 178*5f80ce2aSJacob Faibussowitsch CHKERRQ(test_3d_nocoord("3d_nocoord.vtr")); 179c4762a1bSJed Brown ierr = PetscFinalize(); 180c4762a1bSJed Brown return ierr; 181c4762a1bSJed Brown } 182c4762a1bSJed Brown 183c4762a1bSJed Brown /*TEST 184c4762a1bSJed Brown 185c4762a1bSJed Brown build: 186c4762a1bSJed Brown requires: !complex 187c4762a1bSJed Brown 188c4762a1bSJed Brown test: 189c4762a1bSJed Brown nsize: 2 190c4762a1bSJed Brown 191c4762a1bSJed Brown TEST*/ 192