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 */ 12d71ae5a4SJacob Faibussowitsch PetscErrorCode test_3d(const char filename[]) 13d71ae5a4SJacob Faibussowitsch { 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 249566063dSJacob Faibussowitsch PetscCall(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)); 259566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 269566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 27c4762a1bSJed Brown 289566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, Lx, 0.0, Ly, 0.0, Lz)); 299566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 309566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &v)); 319566063dSJacob Faibussowitsch PetscCall(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 } 429566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, v, &va)); 439566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKOpen(comm, filename, FILE_MODE_WRITE, &view)); 449566063dSJacob Faibussowitsch PetscCall(VecView(v, view)); 459566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 479566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 483ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 49c4762a1bSJed Brown } 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* 52c4762a1bSJed Brown Write 2D DMDA vector with coordinates in VTK VTR format 53c4762a1bSJed Brown 54c4762a1bSJed Brown */ 55d71ae5a4SJacob Faibussowitsch PetscErrorCode test_2d(const char filename[]) 56d71ae5a4SJacob Faibussowitsch { 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 679566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, M, N, PETSC_DECIDE, PETSC_DECIDE, dof, sw, NULL, NULL, &da)); 689566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 699566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 709566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, Lx, 0.0, Ly, 0.0, Lz)); 719566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 729566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &v)); 739566063dSJacob Faibussowitsch PetscCall(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 } 819566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, v, &va)); 829566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKOpen(comm, filename, FILE_MODE_WRITE, &view)); 839566063dSJacob Faibussowitsch PetscCall(VecView(v, view)); 849566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 859566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 869566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 873ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 88c4762a1bSJed Brown } 89c4762a1bSJed Brown 90c4762a1bSJed Brown /* 91c4762a1bSJed Brown Write 2D DMDA vector without coordinates in VTK VTR format 92c4762a1bSJed Brown 93c4762a1bSJed Brown */ 94d71ae5a4SJacob Faibussowitsch PetscErrorCode test_2d_nocoord(const char filename[]) 95d71ae5a4SJacob Faibussowitsch { 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 1069566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, M, N, PETSC_DECIDE, PETSC_DECIDE, dof, sw, NULL, NULL, &da)); 1079566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 1089566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1099566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 1109566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &v)); 1119566063dSJacob Faibussowitsch PetscCall(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 } 1199566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, v, &va)); 1209566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKOpen(comm, filename, FILE_MODE_WRITE, &view)); 1219566063dSJacob Faibussowitsch PetscCall(VecView(v, view)); 1229566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 1239566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 1249566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1253ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 126c4762a1bSJed Brown } 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* 129c4762a1bSJed Brown Write 3D DMDA vector without coordinates in VTK VTR format 130c4762a1bSJed Brown 131c4762a1bSJed Brown */ 132d71ae5a4SJacob Faibussowitsch PetscErrorCode test_3d_nocoord(const char filename[]) 133d71ae5a4SJacob Faibussowitsch { 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 1449566063dSJacob Faibussowitsch PetscCall(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)); 1459566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 1469566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 147c4762a1bSJed Brown 1489566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 1499566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &v)); 1509566063dSJacob Faibussowitsch PetscCall(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 } 1619566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, v, &va)); 1629566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKOpen(comm, filename, FILE_MODE_WRITE, &view)); 1639566063dSJacob Faibussowitsch PetscCall(VecView(v, view)); 1649566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 1659566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 1669566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1673ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 168c4762a1bSJed Brown } 169c4762a1bSJed Brown 170d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[]) 171d71ae5a4SJacob Faibussowitsch { 172327415f7SBarry Smith PetscFunctionBeginUser; 1739566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help)); 1749566063dSJacob Faibussowitsch PetscCall(test_3d("3d.vtr")); 1759566063dSJacob Faibussowitsch PetscCall(test_2d("2d.vtr")); 1769566063dSJacob Faibussowitsch PetscCall(test_2d_nocoord("2d_nocoord.vtr")); 1779566063dSJacob Faibussowitsch PetscCall(test_3d_nocoord("3d_nocoord.vtr")); 1789566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 179b122ec5aSJacob 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 189*3886731fSPierre Jolivet output_file: output/empty.out 190c4762a1bSJed Brown 191c4762a1bSJed Brown TEST*/ 192