1c4762a1bSJed Brown static char help[] = "Test VTK structured (.vts) and rectilinear (.vtr) viewer support with multi-dof DMDAs.\n\ 2c4762a1bSJed Brown Supply the -namefields flag to test with field names.\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscdm.h> 5c4762a1bSJed Brown #include <petscdmda.h> 6c4762a1bSJed Brown 7c4762a1bSJed Brown /* Helper function to name DMDA fields */ 8d71ae5a4SJacob Faibussowitsch PetscErrorCode NameFields(DM da, PetscInt dof) 9d71ae5a4SJacob Faibussowitsch { 10c4762a1bSJed Brown PetscInt c; 11c4762a1bSJed Brown 12c4762a1bSJed Brown PetscFunctionBeginUser; 13c4762a1bSJed Brown for (c = 0; c < dof; ++c) { 14c4762a1bSJed Brown char fieldname[256]; 1563a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(fieldname, sizeof(fieldname), "field_%" PetscInt_FMT, c)); 169566063dSJacob Faibussowitsch PetscCall(DMDASetFieldName(da, c, fieldname)); 17c4762a1bSJed Brown } 18*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 19c4762a1bSJed Brown } 20c4762a1bSJed Brown 21c4762a1bSJed Brown /* 22c4762a1bSJed Brown Write 3D DMDA vector with coordinates in VTK format 23c4762a1bSJed Brown */ 24d71ae5a4SJacob Faibussowitsch PetscErrorCode test_3d(const char filename[], PetscInt dof, PetscBool namefields) 25d71ae5a4SJacob Faibussowitsch { 26c4762a1bSJed Brown MPI_Comm comm = MPI_COMM_WORLD; 27c4762a1bSJed Brown const PetscInt M = 10, N = 15, P = 30, sw = 1; 28c4762a1bSJed Brown const PetscScalar Lx = 1.0, Ly = 1.0, Lz = 1.0; 29c4762a1bSJed Brown DM da; 30c4762a1bSJed Brown Vec v; 31c4762a1bSJed Brown PetscViewer view; 32c4762a1bSJed Brown DMDALocalInfo info; 33c4762a1bSJed Brown PetscScalar ****va; 34c4762a1bSJed Brown PetscInt i, j, k, c; 35c4762a1bSJed Brown 369566063dSJacob 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)); 379566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 389566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 399566063dSJacob Faibussowitsch if (namefields) PetscCall(NameFields(da, dof)); 40c4762a1bSJed Brown 419566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, Lx, 0.0, Ly, 0.0, Lz)); 429566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 439566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &v)); 449566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOF(da, v, &va)); 45c4762a1bSJed Brown for (k = info.zs; k < info.zs + info.zm; k++) { 46c4762a1bSJed Brown for (j = info.ys; j < info.ys + info.ym; j++) { 47c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; i++) { 48c4762a1bSJed Brown const PetscScalar x = (Lx * i) / M; 49c4762a1bSJed Brown const PetscScalar y = (Ly * j) / N; 50c4762a1bSJed Brown const PetscScalar z = (Lz * k) / P; 51ad540459SPierre Jolivet for (c = 0; c < dof; ++c) 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; 52c4762a1bSJed Brown } 53c4762a1bSJed Brown } 54c4762a1bSJed Brown } 559566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOF(da, v, &va)); 569566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKOpen(comm, filename, FILE_MODE_WRITE, &view)); 579566063dSJacob Faibussowitsch PetscCall(VecView(v, view)); 589566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 599566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 609566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 61*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* 65c4762a1bSJed Brown Write 2D DMDA vector with coordinates in VTK format 66c4762a1bSJed Brown */ 67d71ae5a4SJacob Faibussowitsch PetscErrorCode test_2d(const char filename[], PetscInt dof, PetscBool namefields) 68d71ae5a4SJacob Faibussowitsch { 69c4762a1bSJed Brown MPI_Comm comm = MPI_COMM_WORLD; 70c4762a1bSJed Brown const PetscInt M = 10, N = 20, sw = 1; 71c4762a1bSJed Brown const PetscScalar Lx = 1.0, Ly = 1.0, Lz = 1.0; 72c4762a1bSJed Brown DM da; 73c4762a1bSJed Brown Vec v; 74c4762a1bSJed Brown PetscViewer view; 75c4762a1bSJed Brown DMDALocalInfo info; 76c4762a1bSJed Brown PetscScalar ***va; 77c4762a1bSJed Brown PetscInt i, j, c; 78c4762a1bSJed Brown 799566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, M, N, PETSC_DECIDE, PETSC_DECIDE, dof, sw, NULL, NULL, &da)); 809566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 819566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 829566063dSJacob Faibussowitsch if (namefields) PetscCall(NameFields(da, dof)); 839566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, Lx, 0.0, Ly, 0.0, Lz)); 849566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 859566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &v)); 869566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOF(da, v, &va)); 87c4762a1bSJed Brown for (j = info.ys; j < info.ys + info.ym; j++) { 88c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; i++) { 89c4762a1bSJed Brown const PetscScalar x = (Lx * i) / M; 90c4762a1bSJed Brown const PetscScalar y = (Ly * j) / N; 91ad540459SPierre Jolivet for (c = 0; c < dof; ++c) va[j][i][c] = PetscPowScalarInt(x - 0.5 * Lx, 2) + PetscPowScalarInt(y - 0.5 * Ly, 2) + 10.0 * c; 92c4762a1bSJed Brown } 93c4762a1bSJed Brown } 949566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOF(da, v, &va)); 959566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKOpen(comm, filename, FILE_MODE_WRITE, &view)); 969566063dSJacob Faibussowitsch PetscCall(VecView(v, view)); 979566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 989566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 999566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 100*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 101c4762a1bSJed Brown } 102c4762a1bSJed Brown 103c4762a1bSJed Brown /* 104c4762a1bSJed Brown Write a scalar and a vector field from two compatible 3d DMDAs 105c4762a1bSJed Brown */ 106d71ae5a4SJacob Faibussowitsch PetscErrorCode test_3d_compat(const char filename[], PetscInt dof, PetscBool namefields) 107d71ae5a4SJacob Faibussowitsch { 108c4762a1bSJed Brown MPI_Comm comm = MPI_COMM_WORLD; 109c4762a1bSJed Brown const PetscInt M = 10, N = 15, P = 30, sw = 1; 110c4762a1bSJed Brown const PetscScalar Lx = 1.0, Ly = 1.0, Lz = 1.0; 111c4762a1bSJed Brown DM da, daVector; 112c4762a1bSJed Brown Vec v, vVector; 113c4762a1bSJed Brown PetscViewer view; 114c4762a1bSJed Brown DMDALocalInfo info; 115c4762a1bSJed Brown PetscScalar ***va, ****vVectora; 116c4762a1bSJed Brown PetscInt i, j, k, c; 117c4762a1bSJed Brown 1189566063dSJacob 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:*/ 1, sw, NULL, NULL, NULL, &da)); 1199566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 1209566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1219566063dSJacob Faibussowitsch if (namefields) PetscCall(NameFields(da, 1)); 122c4762a1bSJed Brown 1239566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, Lx, 0.0, Ly, 0.0, Lz)); 1249566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 1259566063dSJacob Faibussowitsch PetscCall(DMDACreateCompatibleDMDA(da, dof, &daVector)); 1269566063dSJacob Faibussowitsch if (namefields) PetscCall(NameFields(daVector, dof)); 1279566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &v)); 1289566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daVector, &vVector)); 1299566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, v, &va)); 1309566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOF(daVector, vVector, &vVectora)); 131c4762a1bSJed Brown for (k = info.zs; k < info.zs + info.zm; k++) { 132c4762a1bSJed Brown for (j = info.ys; j < info.ys + info.ym; j++) { 133c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; i++) { 134c4762a1bSJed Brown const PetscScalar x = (Lx * i) / M; 135c4762a1bSJed Brown const PetscScalar y = (Ly * j) / N; 136c4762a1bSJed Brown const PetscScalar z = (Lz * k) / P; 137c4762a1bSJed Brown va[k][j][i] = PetscPowScalarInt(x - 0.5 * Lx, 2) + PetscPowScalarInt(y - 0.5 * Ly, 2) + PetscPowScalarInt(z - 0.5 * Lz, 2); 138ad540459SPierre Jolivet for (c = 0; c < dof; ++c) 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; 139c4762a1bSJed Brown } 140c4762a1bSJed Brown } 141c4762a1bSJed Brown } 1429566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, v, &va)); 1439566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOF(da, v, &vVectora)); 1449566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKOpen(comm, filename, FILE_MODE_WRITE, &view)); 1459566063dSJacob Faibussowitsch PetscCall(VecView(v, view)); 1469566063dSJacob Faibussowitsch PetscCall(VecView(vVector, view)); 1479566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 1489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 1499566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vVector)); 1509566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1519566063dSJacob Faibussowitsch PetscCall(DMDestroy(&daVector)); 152*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 153c4762a1bSJed Brown } 154c4762a1bSJed Brown 155c4762a1bSJed Brown /* 156c4762a1bSJed Brown Write a scalar and a vector field from two compatible 2d DMDAs 157c4762a1bSJed Brown */ 158d71ae5a4SJacob Faibussowitsch PetscErrorCode test_2d_compat(const char filename[], PetscInt dof, PetscBool namefields) 159d71ae5a4SJacob Faibussowitsch { 160c4762a1bSJed Brown MPI_Comm comm = MPI_COMM_WORLD; 161c4762a1bSJed Brown const PetscInt M = 10, N = 20, sw = 1; 162c4762a1bSJed Brown const PetscScalar Lx = 1.0, Ly = 1.0, Lz = 1.0; 163c4762a1bSJed Brown DM da, daVector; 164c4762a1bSJed Brown Vec v, vVector; 165c4762a1bSJed Brown PetscViewer view; 166c4762a1bSJed Brown DMDALocalInfo info; 167c4762a1bSJed Brown PetscScalar **va, ***vVectora; 168c4762a1bSJed Brown PetscInt i, j, c; 169c4762a1bSJed Brown 1709566063dSJacob Faibussowitsch PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, M, N, PETSC_DECIDE, PETSC_DECIDE, /* dof:*/ 1, sw, NULL, NULL, &da)); 1719566063dSJacob Faibussowitsch PetscCall(DMSetFromOptions(da)); 1729566063dSJacob Faibussowitsch PetscCall(DMSetUp(da)); 1739566063dSJacob Faibussowitsch if (namefields) PetscCall(NameFields(da, 1)); 1749566063dSJacob Faibussowitsch PetscCall(DMDASetUniformCoordinates(da, 0.0, Lx, 0.0, Ly, 0.0, Lz)); 1759566063dSJacob Faibussowitsch PetscCall(DMDACreateCompatibleDMDA(da, dof, &daVector)); 1769566063dSJacob Faibussowitsch if (namefields) PetscCall(NameFields(daVector, dof)); 1779566063dSJacob Faibussowitsch PetscCall(DMDAGetLocalInfo(da, &info)); 1789566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(da, &v)); 1799566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(daVector, &vVector)); 1809566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArray(da, v, &va)); 1819566063dSJacob Faibussowitsch PetscCall(DMDAVecGetArrayDOF(daVector, vVector, &vVectora)); 182c4762a1bSJed Brown for (j = info.ys; j < info.ys + info.ym; j++) { 183c4762a1bSJed Brown for (i = info.xs; i < info.xs + info.xm; i++) { 184c4762a1bSJed Brown const PetscScalar x = (Lx * i) / M; 185c4762a1bSJed Brown const PetscScalar y = (Ly * j) / N; 186c4762a1bSJed Brown va[j][i] = PetscPowScalarInt(x - 0.5 * Lx, 2) + PetscPowScalarInt(y - 0.5 * Ly, 2); 187ad540459SPierre Jolivet for (c = 0; c < dof; ++c) vVectora[j][i][c] = PetscPowScalarInt(x - 0.5 * Lx, 2) + PetscPowScalarInt(y - 0.5 * Ly, 2) + 10.0 * c; 188c4762a1bSJed Brown } 189c4762a1bSJed Brown } 1909566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArray(da, v, &va)); 1919566063dSJacob Faibussowitsch PetscCall(DMDAVecRestoreArrayDOF(daVector, vVector, &vVectora)); 1929566063dSJacob Faibussowitsch PetscCall(PetscViewerVTKOpen(comm, filename, FILE_MODE_WRITE, &view)); 1939566063dSJacob Faibussowitsch PetscCall(VecView(v, view)); 1949566063dSJacob Faibussowitsch PetscCall(VecView(vVector, view)); 1959566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&view)); 1969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&v)); 1979566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vVector)); 1989566063dSJacob Faibussowitsch PetscCall(DMDestroy(&da)); 1999566063dSJacob Faibussowitsch PetscCall(DMDestroy(&daVector)); 200*3ba16761SJacob Faibussowitsch return PETSC_SUCCESS; 201c4762a1bSJed Brown } 202c4762a1bSJed Brown 203d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[]) 204d71ae5a4SJacob Faibussowitsch { 205c4762a1bSJed Brown PetscInt dof; 206c4762a1bSJed Brown PetscBool namefields; 207c4762a1bSJed Brown 208327415f7SBarry Smith PetscFunctionBeginUser; 2099566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, 0, help)); 210c4762a1bSJed Brown dof = 2; 2119566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-dof", &dof, NULL)); 212c4762a1bSJed Brown namefields = PETSC_FALSE; 2139566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-namefields", &namefields, NULL)); 2149566063dSJacob Faibussowitsch PetscCall(test_3d("3d.vtr", dof, namefields)); 2159566063dSJacob Faibussowitsch PetscCall(test_2d("2d.vtr", dof, namefields)); 2169566063dSJacob Faibussowitsch PetscCall(test_3d_compat("3d_compat.vtr", dof, namefields)); 2179566063dSJacob Faibussowitsch PetscCall(test_2d_compat("2d_compat.vtr", dof, namefields)); 2189566063dSJacob Faibussowitsch PetscCall(test_3d("3d.vts", dof, namefields)); 2199566063dSJacob Faibussowitsch PetscCall(test_2d("2d.vts", dof, namefields)); 2209566063dSJacob Faibussowitsch PetscCall(test_3d_compat("3d_compat.vts", dof, namefields)); 2219566063dSJacob Faibussowitsch PetscCall(test_2d_compat("2d_compat.vts", dof, namefields)); 2229566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 223b122ec5aSJacob Faibussowitsch return 0; 224c4762a1bSJed Brown } 225c4762a1bSJed Brown 226c4762a1bSJed Brown /*TEST 227c4762a1bSJed Brown 228c4762a1bSJed Brown build: 229c4762a1bSJed Brown requires: !complex 230c4762a1bSJed Brown 231c4762a1bSJed Brown test: 232c4762a1bSJed Brown suffix: 1 233c4762a1bSJed Brown nsize: 2 234c4762a1bSJed Brown args: -dof 2 235c4762a1bSJed Brown 236c4762a1bSJed Brown test: 237c4762a1bSJed Brown suffix: 2 238c4762a1bSJed Brown nsize: 2 239c4762a1bSJed Brown args: -dof 2 240c4762a1bSJed Brown 241c4762a1bSJed Brown test: 242c4762a1bSJed Brown suffix: 3 243c4762a1bSJed Brown nsize: 2 244c4762a1bSJed Brown args: -dof 3 245c4762a1bSJed Brown 246c4762a1bSJed Brown test: 247c4762a1bSJed Brown suffix: 4 248c4762a1bSJed Brown nsize: 1 249c4762a1bSJed Brown args: -dof 2 -namefields 250c4762a1bSJed Brown 251c4762a1bSJed Brown test: 252c4762a1bSJed Brown suffix: 5 253c4762a1bSJed Brown nsize: 2 254c4762a1bSJed Brown args: -dof 4 -namefields 255c4762a1bSJed Brown 256c4762a1bSJed Brown TEST*/ 257