xref: /petsc/src/dm/tests/ex48.c (revision ad540459ab38c4a232710a68077e487eb6627fe2)
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 */
89371c9d4SSatish Balay PetscErrorCode NameFields(DM da, PetscInt dof) {
9c4762a1bSJed Brown   PetscInt c;
10c4762a1bSJed Brown 
11c4762a1bSJed Brown   PetscFunctionBeginUser;
12c4762a1bSJed Brown   for (c = 0; c < dof; ++c) {
13c4762a1bSJed Brown     char fieldname[256];
1463a3b9bcSJacob Faibussowitsch     PetscCall(PetscSNPrintf(fieldname, sizeof(fieldname), "field_%" PetscInt_FMT, c));
159566063dSJacob Faibussowitsch     PetscCall(DMDASetFieldName(da, c, fieldname));
16c4762a1bSJed Brown   }
17c4762a1bSJed Brown   PetscFunctionReturn(0);
18c4762a1bSJed Brown }
19c4762a1bSJed Brown 
20c4762a1bSJed Brown /*
21c4762a1bSJed Brown   Write 3D DMDA vector with coordinates in VTK format
22c4762a1bSJed Brown */
239371c9d4SSatish Balay PetscErrorCode test_3d(const char filename[], PetscInt dof, PetscBool namefields) {
24c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
25c4762a1bSJed Brown   const PetscInt    M = 10, N = 15, P = 30, sw = 1;
26c4762a1bSJed Brown   const PetscScalar Lx = 1.0, Ly = 1.0, Lz = 1.0;
27c4762a1bSJed Brown   DM                da;
28c4762a1bSJed Brown   Vec               v;
29c4762a1bSJed Brown   PetscViewer       view;
30c4762a1bSJed Brown   DMDALocalInfo     info;
31c4762a1bSJed Brown   PetscScalar   ****va;
32c4762a1bSJed Brown   PetscInt          i, j, k, c;
33c4762a1bSJed Brown 
349566063dSJacob 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));
359566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
369566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
379566063dSJacob Faibussowitsch   if (namefields) PetscCall(NameFields(da, dof));
38c4762a1bSJed Brown 
399566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, 0.0, Lx, 0.0, Ly, 0.0, Lz));
409566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
419566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &v));
429566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayDOF(da, v, &va));
43c4762a1bSJed Brown   for (k = info.zs; k < info.zs + info.zm; k++) {
44c4762a1bSJed Brown     for (j = info.ys; j < info.ys + info.ym; j++) {
45c4762a1bSJed Brown       for (i = info.xs; i < info.xs + info.xm; i++) {
46c4762a1bSJed Brown         const PetscScalar x = (Lx * i) / M;
47c4762a1bSJed Brown         const PetscScalar y = (Ly * j) / N;
48c4762a1bSJed Brown         const PetscScalar z = (Lz * k) / P;
49*ad540459SPierre 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;
50c4762a1bSJed Brown       }
51c4762a1bSJed Brown     }
52c4762a1bSJed Brown   }
539566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOF(da, v, &va));
549566063dSJacob Faibussowitsch   PetscCall(PetscViewerVTKOpen(comm, filename, FILE_MODE_WRITE, &view));
559566063dSJacob Faibussowitsch   PetscCall(VecView(v, view));
569566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&view));
579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v));
589566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
59c4762a1bSJed Brown   return 0;
60c4762a1bSJed Brown }
61c4762a1bSJed Brown 
62c4762a1bSJed Brown /*
63c4762a1bSJed Brown   Write 2D DMDA vector with coordinates in VTK format
64c4762a1bSJed Brown */
659371c9d4SSatish Balay PetscErrorCode test_2d(const char filename[], PetscInt dof, PetscBool namefields) {
66c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
67c4762a1bSJed Brown   const PetscInt    M = 10, N = 20, sw = 1;
68c4762a1bSJed Brown   const PetscScalar Lx = 1.0, Ly = 1.0, Lz = 1.0;
69c4762a1bSJed Brown   DM                da;
70c4762a1bSJed Brown   Vec               v;
71c4762a1bSJed Brown   PetscViewer       view;
72c4762a1bSJed Brown   DMDALocalInfo     info;
73c4762a1bSJed Brown   PetscScalar    ***va;
74c4762a1bSJed Brown   PetscInt          i, j, c;
75c4762a1bSJed Brown 
769566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(comm, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, M, N, PETSC_DECIDE, PETSC_DECIDE, dof, sw, NULL, NULL, &da));
779566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
789566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
799566063dSJacob Faibussowitsch   if (namefields) PetscCall(NameFields(da, dof));
809566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, 0.0, Lx, 0.0, Ly, 0.0, Lz));
819566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
829566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &v));
839566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayDOF(da, v, &va));
84c4762a1bSJed Brown   for (j = info.ys; j < info.ys + info.ym; j++) {
85c4762a1bSJed Brown     for (i = info.xs; i < info.xs + info.xm; i++) {
86c4762a1bSJed Brown       const PetscScalar x = (Lx * i) / M;
87c4762a1bSJed Brown       const PetscScalar y = (Ly * j) / N;
88*ad540459SPierre 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;
89c4762a1bSJed Brown     }
90c4762a1bSJed Brown   }
919566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOF(da, v, &va));
929566063dSJacob Faibussowitsch   PetscCall(PetscViewerVTKOpen(comm, filename, FILE_MODE_WRITE, &view));
939566063dSJacob Faibussowitsch   PetscCall(VecView(v, view));
949566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&view));
959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v));
969566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
97c4762a1bSJed Brown   return 0;
98c4762a1bSJed Brown }
99c4762a1bSJed Brown 
100c4762a1bSJed Brown /*
101c4762a1bSJed Brown   Write a scalar and a vector field from two compatible 3d DMDAs
102c4762a1bSJed Brown */
1039371c9d4SSatish Balay PetscErrorCode test_3d_compat(const char filename[], PetscInt dof, PetscBool namefields) {
104c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
105c4762a1bSJed Brown   const PetscInt    M = 10, N = 15, P = 30, sw = 1;
106c4762a1bSJed Brown   const PetscScalar Lx = 1.0, Ly = 1.0, Lz = 1.0;
107c4762a1bSJed Brown   DM                da, daVector;
108c4762a1bSJed Brown   Vec               v, vVector;
109c4762a1bSJed Brown   PetscViewer       view;
110c4762a1bSJed Brown   DMDALocalInfo     info;
111c4762a1bSJed Brown   PetscScalar    ***va, ****vVectora;
112c4762a1bSJed Brown   PetscInt          i, j, k, c;
113c4762a1bSJed Brown 
1149566063dSJacob 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));
1159566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
1169566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
1179566063dSJacob Faibussowitsch   if (namefields) PetscCall(NameFields(da, 1));
118c4762a1bSJed Brown 
1199566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, 0.0, Lx, 0.0, Ly, 0.0, Lz));
1209566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
1219566063dSJacob Faibussowitsch   PetscCall(DMDACreateCompatibleDMDA(da, dof, &daVector));
1229566063dSJacob Faibussowitsch   if (namefields) PetscCall(NameFields(daVector, dof));
1239566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &v));
1249566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(daVector, &vVector));
1259566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, v, &va));
1269566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayDOF(daVector, vVector, &vVectora));
127c4762a1bSJed Brown   for (k = info.zs; k < info.zs + info.zm; k++) {
128c4762a1bSJed Brown     for (j = info.ys; j < info.ys + info.ym; j++) {
129c4762a1bSJed Brown       for (i = info.xs; i < info.xs + info.xm; i++) {
130c4762a1bSJed Brown         const PetscScalar x = (Lx * i) / M;
131c4762a1bSJed Brown         const PetscScalar y = (Ly * j) / N;
132c4762a1bSJed Brown         const PetscScalar z = (Lz * k) / P;
133c4762a1bSJed Brown         va[k][j][i]         = PetscPowScalarInt(x - 0.5 * Lx, 2) + PetscPowScalarInt(y - 0.5 * Ly, 2) + PetscPowScalarInt(z - 0.5 * Lz, 2);
134*ad540459SPierre 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;
135c4762a1bSJed Brown       }
136c4762a1bSJed Brown     }
137c4762a1bSJed Brown   }
1389566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, v, &va));
1399566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOF(da, v, &vVectora));
1409566063dSJacob Faibussowitsch   PetscCall(PetscViewerVTKOpen(comm, filename, FILE_MODE_WRITE, &view));
1419566063dSJacob Faibussowitsch   PetscCall(VecView(v, view));
1429566063dSJacob Faibussowitsch   PetscCall(VecView(vVector, view));
1439566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&view));
1449566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v));
1459566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vVector));
1469566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1479566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&daVector));
148c4762a1bSJed Brown   return 0;
149c4762a1bSJed Brown }
150c4762a1bSJed Brown 
151c4762a1bSJed Brown /*
152c4762a1bSJed Brown   Write a scalar and a vector field from two compatible 2d DMDAs
153c4762a1bSJed Brown */
1549371c9d4SSatish Balay PetscErrorCode test_2d_compat(const char filename[], PetscInt dof, PetscBool namefields) {
155c4762a1bSJed Brown   MPI_Comm          comm = MPI_COMM_WORLD;
156c4762a1bSJed Brown   const PetscInt    M = 10, N = 20, sw = 1;
157c4762a1bSJed Brown   const PetscScalar Lx = 1.0, Ly = 1.0, Lz = 1.0;
158c4762a1bSJed Brown   DM                da, daVector;
159c4762a1bSJed Brown   Vec               v, vVector;
160c4762a1bSJed Brown   PetscViewer       view;
161c4762a1bSJed Brown   DMDALocalInfo     info;
162c4762a1bSJed Brown   PetscScalar     **va, ***vVectora;
163c4762a1bSJed Brown   PetscInt          i, j, c;
164c4762a1bSJed Brown 
1659566063dSJacob 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));
1669566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
1679566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
1689566063dSJacob Faibussowitsch   if (namefields) PetscCall(NameFields(da, 1));
1699566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(da, 0.0, Lx, 0.0, Ly, 0.0, Lz));
1709566063dSJacob Faibussowitsch   PetscCall(DMDACreateCompatibleDMDA(da, dof, &daVector));
1719566063dSJacob Faibussowitsch   if (namefields) PetscCall(NameFields(daVector, dof));
1729566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
1739566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &v));
1749566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(daVector, &vVector));
1759566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, v, &va));
1769566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayDOF(daVector, vVector, &vVectora));
177c4762a1bSJed Brown   for (j = info.ys; j < info.ys + info.ym; j++) {
178c4762a1bSJed Brown     for (i = info.xs; i < info.xs + info.xm; i++) {
179c4762a1bSJed Brown       const PetscScalar x = (Lx * i) / M;
180c4762a1bSJed Brown       const PetscScalar y = (Ly * j) / N;
181c4762a1bSJed Brown       va[j][i]            = PetscPowScalarInt(x - 0.5 * Lx, 2) + PetscPowScalarInt(y - 0.5 * Ly, 2);
182*ad540459SPierre 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;
183c4762a1bSJed Brown     }
184c4762a1bSJed Brown   }
1859566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, v, &va));
1869566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOF(daVector, vVector, &vVectora));
1879566063dSJacob Faibussowitsch   PetscCall(PetscViewerVTKOpen(comm, filename, FILE_MODE_WRITE, &view));
1889566063dSJacob Faibussowitsch   PetscCall(VecView(v, view));
1899566063dSJacob Faibussowitsch   PetscCall(VecView(vVector, view));
1909566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&view));
1919566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&v));
1929566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&vVector));
1939566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1949566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&daVector));
195c4762a1bSJed Brown   return 0;
196c4762a1bSJed Brown }
197c4762a1bSJed Brown 
1989371c9d4SSatish Balay int main(int argc, char *argv[]) {
199c4762a1bSJed Brown   PetscInt  dof;
200c4762a1bSJed Brown   PetscBool namefields;
201c4762a1bSJed Brown 
202327415f7SBarry Smith   PetscFunctionBeginUser;
2039566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, 0, help));
204c4762a1bSJed Brown   dof = 2;
2059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-dof", &dof, NULL));
206c4762a1bSJed Brown   namefields = PETSC_FALSE;
2079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-namefields", &namefields, NULL));
2089566063dSJacob Faibussowitsch   PetscCall(test_3d("3d.vtr", dof, namefields));
2099566063dSJacob Faibussowitsch   PetscCall(test_2d("2d.vtr", dof, namefields));
2109566063dSJacob Faibussowitsch   PetscCall(test_3d_compat("3d_compat.vtr", dof, namefields));
2119566063dSJacob Faibussowitsch   PetscCall(test_2d_compat("2d_compat.vtr", dof, namefields));
2129566063dSJacob Faibussowitsch   PetscCall(test_3d("3d.vts", dof, namefields));
2139566063dSJacob Faibussowitsch   PetscCall(test_2d("2d.vts", dof, namefields));
2149566063dSJacob Faibussowitsch   PetscCall(test_3d_compat("3d_compat.vts", dof, namefields));
2159566063dSJacob Faibussowitsch   PetscCall(test_2d_compat("2d_compat.vts", dof, namefields));
2169566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
217b122ec5aSJacob Faibussowitsch   return 0;
218c4762a1bSJed Brown }
219c4762a1bSJed Brown 
220c4762a1bSJed Brown /*TEST
221c4762a1bSJed Brown 
222c4762a1bSJed Brown    build:
223c4762a1bSJed Brown       requires: !complex
224c4762a1bSJed Brown 
225c4762a1bSJed Brown    test:
226c4762a1bSJed Brown       suffix: 1
227c4762a1bSJed Brown       nsize: 2
228c4762a1bSJed Brown       args: -dof 2
229c4762a1bSJed Brown 
230c4762a1bSJed Brown    test:
231c4762a1bSJed Brown       suffix: 2
232c4762a1bSJed Brown       nsize: 2
233c4762a1bSJed Brown       args: -dof 2
234c4762a1bSJed Brown 
235c4762a1bSJed Brown    test:
236c4762a1bSJed Brown       suffix: 3
237c4762a1bSJed Brown       nsize: 2
238c4762a1bSJed Brown       args: -dof 3
239c4762a1bSJed Brown 
240c4762a1bSJed Brown    test:
241c4762a1bSJed Brown       suffix: 4
242c4762a1bSJed Brown       nsize: 1
243c4762a1bSJed Brown       args: -dof 2 -namefields
244c4762a1bSJed Brown 
245c4762a1bSJed Brown    test:
246c4762a1bSJed Brown       suffix: 5
247c4762a1bSJed Brown       nsize: 2
248c4762a1bSJed Brown       args: -dof 4 -namefields
249c4762a1bSJed Brown 
250c4762a1bSJed Brown TEST*/
251