xref: /petsc/src/dm/tutorials/ex15.c (revision 732aec7a18f2199fb53bb9a2f3aef439a834ce31)
19371c9d4SSatish Balay static char help[] = "Tests VecView() functionality with DMDA objects when using:"
2c4762a1bSJed Brown                      "(i) a PetscViewer binary with MPI-IO support; and (ii) when the binary header is skipped.\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscdm.h>
5c4762a1bSJed Brown #include <petscdmda.h>
6c4762a1bSJed Brown 
7c4762a1bSJed Brown #define DMDA_I 5
8c4762a1bSJed Brown #define DMDA_J 4
9c4762a1bSJed Brown #define DMDA_K 6
10c4762a1bSJed Brown 
11c4762a1bSJed Brown const PetscReal dmda_i_val[] = {1.10, 2.3006, 2.32444, 3.44006, 66.9009};
12c4762a1bSJed Brown const PetscReal dmda_j_val[] = {0.0, 0.25, 0.5, 0.75};
13c4762a1bSJed Brown const PetscReal dmda_k_val[] = {0.0, 1.1, 2.2, 3.3, 4.4, 5.5};
14c4762a1bSJed Brown 
MyVecDump(const char fname[],PetscBool skippheader,PetscBool usempiio,Vec x)15d71ae5a4SJacob Faibussowitsch PetscErrorCode MyVecDump(const char fname[], PetscBool skippheader, PetscBool usempiio, Vec x)
16d71ae5a4SJacob Faibussowitsch {
17c4762a1bSJed Brown   MPI_Comm    comm;
18c4762a1bSJed Brown   PetscViewer viewer;
19c4762a1bSJed Brown   PetscBool   ismpiio, isskip;
20c4762a1bSJed Brown 
21c4762a1bSJed Brown   PetscFunctionBeginUser;
229566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)x, &comm));
23c4762a1bSJed Brown 
249566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, &viewer));
259566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(viewer, PETSCVIEWERBINARY));
269566063dSJacob Faibussowitsch   if (skippheader) PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
279566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_WRITE));
289566063dSJacob Faibussowitsch   if (usempiio) PetscCall(PetscViewerBinarySetUseMPIIO(viewer, PETSC_TRUE));
299566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(viewer, fname));
30c4762a1bSJed Brown 
319566063dSJacob Faibussowitsch   PetscCall(VecView(x, viewer));
32c4762a1bSJed Brown 
339566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &ismpiio));
349566063dSJacob Faibussowitsch   if (ismpiio) PetscCall(PetscPrintf(comm, "*** PetscViewer[write] using MPI-IO ***\n"));
359566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &isskip));
369566063dSJacob Faibussowitsch   if (isskip) PetscCall(PetscPrintf(comm, "*** PetscViewer[write] skipping header ***\n"));
37c4762a1bSJed Brown 
389566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
40c4762a1bSJed Brown }
41c4762a1bSJed Brown 
MyVecLoad(const char fname[],PetscBool skippheader,PetscBool usempiio,Vec x)42d71ae5a4SJacob Faibussowitsch PetscErrorCode MyVecLoad(const char fname[], PetscBool skippheader, PetscBool usempiio, Vec x)
43d71ae5a4SJacob Faibussowitsch {
44c4762a1bSJed Brown   MPI_Comm    comm;
45c4762a1bSJed Brown   PetscViewer viewer;
46c4762a1bSJed Brown   PetscBool   ismpiio, isskip;
47c4762a1bSJed Brown 
48c4762a1bSJed Brown   PetscFunctionBeginUser;
499566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)x, &comm));
50c4762a1bSJed Brown 
519566063dSJacob Faibussowitsch   PetscCall(PetscViewerCreate(comm, &viewer));
529566063dSJacob Faibussowitsch   PetscCall(PetscViewerSetType(viewer, PETSCVIEWERBINARY));
539566063dSJacob Faibussowitsch   if (skippheader) PetscCall(PetscViewerBinarySetSkipHeader(viewer, PETSC_TRUE));
549566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetMode(viewer, FILE_MODE_READ));
559566063dSJacob Faibussowitsch   if (usempiio) PetscCall(PetscViewerBinarySetUseMPIIO(viewer, PETSC_TRUE));
569566063dSJacob Faibussowitsch   PetscCall(PetscViewerFileSetName(viewer, fname));
57c4762a1bSJed Brown 
589566063dSJacob Faibussowitsch   PetscCall(VecLoad(x, viewer));
59c4762a1bSJed Brown 
609566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetSkipHeader(viewer, &isskip));
619566063dSJacob Faibussowitsch   if (isskip) PetscCall(PetscPrintf(comm, "*** PetscViewer[load] skipping header ***\n"));
629566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryGetUseMPIIO(viewer, &ismpiio));
639566063dSJacob Faibussowitsch   if (ismpiio) PetscCall(PetscPrintf(comm, "*** PetscViewer[load] using MPI-IO ***\n"));
64c4762a1bSJed Brown 
659566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
67c4762a1bSJed Brown }
68c4762a1bSJed Brown 
DMDAVecGenerateEntries(DM dm,Vec a)69d71ae5a4SJacob Faibussowitsch PetscErrorCode DMDAVecGenerateEntries(DM dm, Vec a)
70d71ae5a4SJacob Faibussowitsch {
71c4762a1bSJed Brown   PetscScalar ****LA_v;
72c4762a1bSJed Brown   PetscInt        i, j, k, l, si, sj, sk, ni, nj, nk, M, N, dof;
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   PetscFunctionBeginUser;
759566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, NULL, &M, &N, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
769566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(dm, &si, &sj, &sk, &ni, &nj, &nk));
779566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayDOF(dm, a, &LA_v));
78c4762a1bSJed Brown   for (k = sk; k < sk + nk; k++) {
79c4762a1bSJed Brown     for (j = sj; j < sj + nj; j++) {
80c4762a1bSJed Brown       for (i = si; i < si + ni; i++) {
81c4762a1bSJed Brown         PetscScalar test_value_s;
82c4762a1bSJed Brown 
83c4762a1bSJed Brown         test_value_s = dmda_i_val[i] * ((PetscScalar)i) + dmda_j_val[j] * ((PetscScalar)(i + j * M)) + dmda_k_val[k] * ((PetscScalar)(i + j * M + k * M * N));
84ad540459SPierre Jolivet         for (l = 0; l < dof; l++) LA_v[k][j][i][l] = (PetscScalar)dof * test_value_s + (PetscScalar)l;
85c4762a1bSJed Brown       }
86c4762a1bSJed Brown     }
87c4762a1bSJed Brown   }
889566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOF(dm, a, &LA_v));
893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
90c4762a1bSJed Brown }
91c4762a1bSJed Brown 
HeaderlessBinaryReadCheck(DM dm,const char name[])92d71ae5a4SJacob Faibussowitsch PetscErrorCode HeaderlessBinaryReadCheck(DM dm, const char name[])
93d71ae5a4SJacob Faibussowitsch {
94c4762a1bSJed Brown   int         fdes;
95c4762a1bSJed Brown   PetscScalar buffer[DMDA_I * DMDA_J * DMDA_K * 10];
96c4762a1bSJed Brown   PetscInt    len, d, i, j, k, M, N, dof;
97c4762a1bSJed Brown   PetscMPIInt rank;
98c4762a1bSJed Brown   PetscBool   dataverified = PETSC_TRUE;
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   PetscFunctionBeginUser;
1019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1029566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(dm, NULL, &M, &N, NULL, NULL, NULL, NULL, &dof, NULL, NULL, NULL, NULL, NULL));
103c4762a1bSJed Brown   len = DMDA_I * DMDA_J * DMDA_K * dof;
104dd400576SPatrick Sanan   if (rank == 0) {
1059566063dSJacob Faibussowitsch     PetscCall(PetscBinaryOpen(name, FILE_MODE_READ, &fdes));
1069566063dSJacob Faibussowitsch     PetscCall(PetscBinaryRead(fdes, buffer, len, NULL, PETSC_SCALAR));
1079566063dSJacob Faibussowitsch     PetscCall(PetscBinaryClose(fdes));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown     for (k = 0; k < DMDA_K; k++) {
110c4762a1bSJed Brown       for (j = 0; j < DMDA_J; j++) {
111c4762a1bSJed Brown         for (i = 0; i < DMDA_I; i++) {
112c4762a1bSJed Brown           for (d = 0; d < dof; d++) {
113c4762a1bSJed Brown             PetscScalar v, test_value_s, test_value;
114c4762a1bSJed Brown             PetscInt    index;
115c4762a1bSJed Brown 
116c4762a1bSJed Brown             test_value_s = dmda_i_val[i] * ((PetscScalar)i) + dmda_j_val[j] * ((PetscScalar)(i + j * M)) + dmda_k_val[k] * ((PetscScalar)(i + j * M + k * M * N));
117c4762a1bSJed Brown             test_value   = (PetscScalar)dof * test_value_s + (PetscScalar)d;
118c4762a1bSJed Brown 
119c4762a1bSJed Brown             index = dof * (i + j * M + k * M * N) + d;
120c4762a1bSJed Brown             v     = PetscAbsScalar(test_value - buffer[index]);
121c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
122c4762a1bSJed Brown             if ((PetscRealPart(v) > 1.0e-10) || (PetscImaginaryPart(v) > 1.0e-10)) {
12363a3b9bcSJacob Faibussowitsch               PetscCall(PetscPrintf(PETSC_COMM_SELF, "ERROR: Difference > 1.0e-10 occurred (delta = (%+1.12e,%+1.12e) [loc %" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "(%" PetscInt_FMT ")])\n", (double)PetscRealPart(test_value), (double)PetscImaginaryPart(test_value), i, j, k, d));
124c4762a1bSJed Brown               dataverified = PETSC_FALSE;
125c4762a1bSJed Brown             }
126c4762a1bSJed Brown #else
127c4762a1bSJed Brown             if (PetscRealPart(v) > 1.0e-10) {
12863a3b9bcSJacob Faibussowitsch               PetscCall(PetscPrintf(PETSC_COMM_SELF, "ERROR: Difference > 1.0e-10 occurred (delta = %+1.12e [loc %" PetscInt_FMT ",%" PetscInt_FMT ",%" PetscInt_FMT "(%" PetscInt_FMT ")])\n", (double)PetscRealPart(test_value), i, j, k, d));
129c4762a1bSJed Brown               dataverified = PETSC_FALSE;
130c4762a1bSJed Brown             }
131c4762a1bSJed Brown #endif
132c4762a1bSJed Brown           }
133c4762a1bSJed Brown         }
134c4762a1bSJed Brown       }
135c4762a1bSJed Brown     }
13648a46eb9SPierre Jolivet     if (dataverified) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Headerless read of data verified for: %s\n", name));
137c4762a1bSJed Brown   }
1383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
139c4762a1bSJed Brown }
140c4762a1bSJed Brown 
VecCompare(Vec a,Vec b)141d71ae5a4SJacob Faibussowitsch PetscErrorCode VecCompare(Vec a, Vec b)
142d71ae5a4SJacob Faibussowitsch {
143c4762a1bSJed Brown   PetscInt  locmin[2], locmax[2];
144c4762a1bSJed Brown   PetscReal min[2], max[2];
145c4762a1bSJed Brown   Vec       ref;
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   PetscFunctionBeginUser;
1489566063dSJacob Faibussowitsch   PetscCall(VecMin(a, &locmin[0], &min[0]));
1499566063dSJacob Faibussowitsch   PetscCall(VecMax(a, &locmax[0], &max[0]));
150c4762a1bSJed Brown 
1519566063dSJacob Faibussowitsch   PetscCall(VecMin(b, &locmin[1], &min[1]));
1529566063dSJacob Faibussowitsch   PetscCall(VecMax(b, &locmax[1], &max[1]));
153c4762a1bSJed Brown 
1549566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "VecCompare\n"));
15563a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  min(a)   = %+1.2e [loc %" PetscInt_FMT "]\n", (double)min[0], locmin[0]));
15663a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  max(a)   = %+1.2e [loc %" PetscInt_FMT "]\n", (double)max[0], locmax[0]));
157c4762a1bSJed Brown 
15863a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  min(b)   = %+1.2e [loc %" PetscInt_FMT "]\n", (double)min[1], locmin[1]));
15963a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  max(b)   = %+1.2e [loc %" PetscInt_FMT "]\n", (double)max[1], locmax[1]));
160c4762a1bSJed Brown 
1619566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(a, &ref));
1629566063dSJacob Faibussowitsch   PetscCall(VecCopy(a, ref));
1639566063dSJacob Faibussowitsch   PetscCall(VecAXPY(ref, -1.0, b));
1649566063dSJacob Faibussowitsch   PetscCall(VecMin(ref, &locmin[0], &min[0]));
165c4762a1bSJed Brown   if (PetscAbsReal(min[0]) > 1.0e-10) {
1669566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  ERROR: min(a-b) > 1.0e-10\n"));
1679566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  min(a-b) = %+1.10e\n", (double)PetscAbsReal(min[0])));
168c4762a1bSJed Brown   } else {
1699566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "  min(a-b) < 1.0e-10\n"));
170c4762a1bSJed Brown   }
1719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ref));
1723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
173c4762a1bSJed Brown }
174c4762a1bSJed Brown 
TestDMDAVec(PetscBool usempiio)175d71ae5a4SJacob Faibussowitsch PetscErrorCode TestDMDAVec(PetscBool usempiio)
176d71ae5a4SJacob Faibussowitsch {
177c4762a1bSJed Brown   DM        dm;
178c4762a1bSJed Brown   Vec       x_ref, x_test;
179c4762a1bSJed Brown   PetscBool skipheader = PETSC_TRUE;
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   PetscFunctionBeginUser;
1829566063dSJacob Faibussowitsch   if (!usempiio) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s\n", PETSC_FUNCTION_NAME));
1839566063dSJacob Faibussowitsch   else PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s [using mpi-io]\n", PETSC_FUNCTION_NAME));
1849371c9d4SSatish Balay   PetscCall(DMDACreate3d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_BOX, DMDA_I, DMDA_J, DMDA_K, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, 3, 2, NULL, NULL, NULL, &dm));
1859566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(dm));
1869566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
187c4762a1bSJed Brown 
1889566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &x_ref));
1899566063dSJacob Faibussowitsch   PetscCall(DMDAVecGenerateEntries(dm, x_ref));
190c4762a1bSJed Brown 
191c4762a1bSJed Brown   if (!usempiio) {
1929566063dSJacob Faibussowitsch     PetscCall(MyVecDump("dmda.pbvec", skipheader, PETSC_FALSE, x_ref));
193c4762a1bSJed Brown   } else {
1949566063dSJacob Faibussowitsch     PetscCall(MyVecDump("dmda-mpiio.pbvec", skipheader, PETSC_TRUE, x_ref));
195c4762a1bSJed Brown   }
196c4762a1bSJed Brown 
1979566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(dm, &x_test));
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   if (!usempiio) {
2009566063dSJacob Faibussowitsch     PetscCall(MyVecLoad("dmda.pbvec", skipheader, usempiio, x_test));
201c4762a1bSJed Brown   } else {
2029566063dSJacob Faibussowitsch     PetscCall(MyVecLoad("dmda-mpiio.pbvec", skipheader, usempiio, x_test));
203c4762a1bSJed Brown   }
204c4762a1bSJed Brown 
2059566063dSJacob Faibussowitsch   PetscCall(VecCompare(x_ref, x_test));
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   if (!usempiio) {
2089566063dSJacob Faibussowitsch     PetscCall(HeaderlessBinaryReadCheck(dm, "dmda.pbvec"));
209c4762a1bSJed Brown   } else {
2109566063dSJacob Faibussowitsch     PetscCall(HeaderlessBinaryReadCheck(dm, "dmda-mpiio.pbvec"));
211c4762a1bSJed Brown   }
2129566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x_ref));
2139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x_test));
2149566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&dm));
2153ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
216c4762a1bSJed Brown }
217c4762a1bSJed Brown 
main(int argc,char ** args)218d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
219d71ae5a4SJacob Faibussowitsch {
220c4762a1bSJed Brown   PetscBool usempiio = PETSC_FALSE;
221c4762a1bSJed Brown 
222327415f7SBarry Smith   PetscFunctionBeginUser;
223*c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &args, NULL, help));
2249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-usempiio", &usempiio, NULL));
225c4762a1bSJed Brown   if (!usempiio) {
2269566063dSJacob Faibussowitsch     PetscCall(TestDMDAVec(PETSC_FALSE));
227c4762a1bSJed Brown   } else {
228c4762a1bSJed Brown #if defined(PETSC_HAVE_MPIIO)
2299566063dSJacob Faibussowitsch     PetscCall(TestDMDAVec(PETSC_TRUE));
230c4762a1bSJed Brown #else
2319566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Warning: Executing TestDMDAVec(PETSC_TRUE) requires a working MPI-2 implementation\n"));
232c4762a1bSJed Brown #endif
233c4762a1bSJed Brown   }
2349566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
235b122ec5aSJacob Faibussowitsch   return 0;
236c4762a1bSJed Brown }
237c4762a1bSJed Brown 
238c4762a1bSJed Brown /*TEST
239c4762a1bSJed Brown 
240c4762a1bSJed Brown    test:
241c4762a1bSJed Brown 
242c4762a1bSJed Brown    test:
243c4762a1bSJed Brown       suffix: 2
244c4762a1bSJed Brown       nsize: 12
245c4762a1bSJed Brown 
246c4762a1bSJed Brown    test:
247c4762a1bSJed Brown       suffix: 3
248c4762a1bSJed Brown       nsize: 12
249dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPIIO)
250c4762a1bSJed Brown       args: -usempiio
251c4762a1bSJed Brown 
252c4762a1bSJed Brown TEST*/
253