xref: /petsc/src/dm/tutorials/ex15.c (revision dfd57a172ac9fa6c7b5fe6de6ab5df85cefc2996)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests VecView() functionality with DMDA objects when using:"\
3c4762a1bSJed Brown "(i) a PetscViewer binary with MPI-IO support; and (ii) when the binary header is skipped.\n\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown #include <petscdm.h>
6c4762a1bSJed Brown #include <petscdmda.h>
7c4762a1bSJed Brown 
8c4762a1bSJed Brown #define DMDA_I 5
9c4762a1bSJed Brown #define DMDA_J 4
10c4762a1bSJed Brown #define DMDA_K 6
11c4762a1bSJed Brown 
12c4762a1bSJed Brown const PetscReal dmda_i_val[] = { 1.10, 2.3006, 2.32444, 3.44006, 66.9009 };
13c4762a1bSJed Brown const PetscReal dmda_j_val[] = { 0.0, 0.25, 0.5, 0.75 };
14c4762a1bSJed Brown const PetscReal dmda_k_val[] = { 0.0, 1.1, 2.2, 3.3, 4.4, 5.5 };
15c4762a1bSJed Brown 
16c4762a1bSJed Brown PetscErrorCode MyVecDump(const char fname[],PetscBool skippheader,PetscBool usempiio,Vec x)
17c4762a1bSJed Brown {
18c4762a1bSJed Brown   MPI_Comm       comm;
19c4762a1bSJed Brown   PetscViewer    viewer;
20c4762a1bSJed Brown   PetscBool      ismpiio,isskip;
21c4762a1bSJed Brown   PetscErrorCode ierr;
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   PetscFunctionBeginUser;
24c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr);
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   ierr = PetscViewerCreate(comm,&viewer);CHKERRQ(ierr);
27c4762a1bSJed Brown   ierr = PetscViewerSetType(viewer,PETSCVIEWERBINARY);CHKERRQ(ierr);
28c4762a1bSJed Brown   if (skippheader) { ierr = PetscViewerBinarySetSkipHeader(viewer,PETSC_TRUE);CHKERRQ(ierr); }
29c4762a1bSJed Brown   ierr = PetscViewerFileSetMode(viewer,FILE_MODE_WRITE);CHKERRQ(ierr);
30c4762a1bSJed Brown   if (usempiio) { ierr = PetscViewerBinarySetUseMPIIO(viewer,PETSC_TRUE);CHKERRQ(ierr); }
31c4762a1bSJed Brown   ierr = PetscViewerFileSetName(viewer,fname);CHKERRQ(ierr);
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   ierr = VecView(x,viewer);CHKERRQ(ierr);
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   ierr = PetscViewerBinaryGetUseMPIIO(viewer,&ismpiio);CHKERRQ(ierr);
36c4762a1bSJed Brown   if (ismpiio) { ierr = PetscPrintf(comm,"*** PetscViewer[write] using MPI-IO ***\n");CHKERRQ(ierr); }
37c4762a1bSJed Brown   ierr = PetscViewerBinaryGetSkipHeader(viewer,&isskip);CHKERRQ(ierr);
38c4762a1bSJed Brown   if (isskip) { ierr = PetscPrintf(comm,"*** PetscViewer[write] skipping header ***\n");CHKERRQ(ierr); }
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
41c4762a1bSJed Brown   PetscFunctionReturn(0);
42c4762a1bSJed Brown }
43c4762a1bSJed Brown 
44c4762a1bSJed Brown PetscErrorCode MyVecLoad(const char fname[],PetscBool skippheader,PetscBool usempiio,Vec x)
45c4762a1bSJed Brown {
46c4762a1bSJed Brown   MPI_Comm       comm;
47c4762a1bSJed Brown   PetscViewer    viewer;
48c4762a1bSJed Brown   PetscBool      ismpiio,isskip;
49c4762a1bSJed Brown   PetscErrorCode ierr;
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   PetscFunctionBeginUser;
52c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)x,&comm);CHKERRQ(ierr);
53c4762a1bSJed Brown 
54c4762a1bSJed Brown   ierr = PetscViewerCreate(comm,&viewer);CHKERRQ(ierr);
55c4762a1bSJed Brown   ierr = PetscViewerSetType(viewer,PETSCVIEWERBINARY);CHKERRQ(ierr);
56c4762a1bSJed Brown   if (skippheader) { ierr = PetscViewerBinarySetSkipHeader(viewer,PETSC_TRUE);CHKERRQ(ierr); }
57c4762a1bSJed Brown   ierr = PetscViewerFileSetMode(viewer,FILE_MODE_READ);CHKERRQ(ierr);
58c4762a1bSJed Brown   if (usempiio) { ierr = PetscViewerBinarySetUseMPIIO(viewer,PETSC_TRUE);CHKERRQ(ierr); }
59c4762a1bSJed Brown   ierr = PetscViewerFileSetName(viewer,fname);CHKERRQ(ierr);
60c4762a1bSJed Brown 
61c4762a1bSJed Brown   ierr = VecLoad(x,viewer);CHKERRQ(ierr);
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   ierr = PetscViewerBinaryGetSkipHeader(viewer,&isskip);CHKERRQ(ierr);
64c4762a1bSJed Brown   if (isskip) { ierr = PetscPrintf(comm,"*** PetscViewer[load] skipping header ***\n");CHKERRQ(ierr); }
65c4762a1bSJed Brown   ierr = PetscViewerBinaryGetUseMPIIO(viewer,&ismpiio);CHKERRQ(ierr);
66c4762a1bSJed Brown   if (ismpiio) { ierr = PetscPrintf(comm,"*** PetscViewer[load] using MPI-IO ***\n");CHKERRQ(ierr); }
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
69c4762a1bSJed Brown   PetscFunctionReturn(0);
70c4762a1bSJed Brown }
71c4762a1bSJed Brown 
72c4762a1bSJed Brown PetscErrorCode DMDAVecGenerateEntries(DM dm,Vec a)
73c4762a1bSJed Brown {
74c4762a1bSJed Brown   PetscScalar    ****LA_v;
75c4762a1bSJed Brown   PetscInt       i,j,k,l,si,sj,sk,ni,nj,nk,M,N,dof;
76c4762a1bSJed Brown   PetscErrorCode ierr;
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   PetscFunctionBeginUser;
79c4762a1bSJed Brown   ierr = DMDAGetInfo(dm,NULL,&M,&N,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
80c4762a1bSJed Brown   ierr = DMDAGetCorners(dm,&si,&sj,&sk,&ni,&nj,&nk);CHKERRQ(ierr);
81c4762a1bSJed Brown   ierr = DMDAVecGetArrayDOF(dm,a,&LA_v);CHKERRQ(ierr);
82c4762a1bSJed Brown   for (k=sk; k<sk+nk; k++) {
83c4762a1bSJed Brown     for (j=sj; j<sj+nj; j++) {
84c4762a1bSJed Brown       for (i=si; i<si+ni; i++) {
85c4762a1bSJed Brown         PetscScalar test_value_s;
86c4762a1bSJed Brown 
87c4762a1bSJed 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));
88c4762a1bSJed Brown         for (l=0; l<dof; l++) {
89c4762a1bSJed Brown           LA_v[k][j][i][l] = (PetscScalar)dof * test_value_s + (PetscScalar)l;
90c4762a1bSJed Brown         }
91c4762a1bSJed Brown       }
92c4762a1bSJed Brown     }
93c4762a1bSJed Brown   }
94c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayDOF(dm,a,&LA_v);CHKERRQ(ierr);
95c4762a1bSJed Brown   PetscFunctionReturn(0);
96c4762a1bSJed Brown }
97c4762a1bSJed Brown 
98c4762a1bSJed Brown PetscErrorCode HeaderlessBinaryReadCheck(DM dm,const char name[])
99c4762a1bSJed Brown {
100c4762a1bSJed Brown   PetscErrorCode ierr;
101c4762a1bSJed Brown   int            fdes;
102c4762a1bSJed Brown   PetscScalar    buffer[DMDA_I*DMDA_J*DMDA_K*10];
103c4762a1bSJed Brown   PetscInt       len,d,i,j,k,M,N,dof;
104c4762a1bSJed Brown   PetscMPIInt    rank;
105c4762a1bSJed Brown   PetscBool      dataverified = PETSC_TRUE;
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   PetscFunctionBeginUser;
108ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
109c4762a1bSJed Brown   ierr = DMDAGetInfo(dm,NULL,&M,&N,NULL,NULL,NULL,NULL,&dof,NULL,NULL,NULL,NULL,NULL);CHKERRQ(ierr);
110c4762a1bSJed Brown   len = DMDA_I*DMDA_J*DMDA_K*dof;
111c4762a1bSJed Brown   if (!rank) {
112c4762a1bSJed Brown     ierr = PetscBinaryOpen(name,FILE_MODE_READ,&fdes);CHKERRQ(ierr);
113c4762a1bSJed Brown     ierr = PetscBinaryRead(fdes,buffer,len,NULL,PETSC_SCALAR);CHKERRQ(ierr);
114c4762a1bSJed Brown     ierr = PetscBinaryClose(fdes);CHKERRQ(ierr);
115c4762a1bSJed Brown 
116c4762a1bSJed Brown     for (k=0; k<DMDA_K; k++) {
117c4762a1bSJed Brown       for (j=0; j<DMDA_J; j++) {
118c4762a1bSJed Brown         for (i=0; i<DMDA_I; i++) {
119c4762a1bSJed Brown           for (d=0; d<dof; d++) {
120c4762a1bSJed Brown             PetscScalar v,test_value_s,test_value;
121c4762a1bSJed Brown             PetscInt    index;
122c4762a1bSJed Brown 
123c4762a1bSJed 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));
124c4762a1bSJed Brown             test_value = (PetscScalar)dof * test_value_s + (PetscScalar)d;
125c4762a1bSJed Brown 
126c4762a1bSJed Brown             index = dof*(i + j*M + k*M*N) + d;
127c4762a1bSJed Brown             v = PetscAbsScalar(test_value-buffer[index]);
128c4762a1bSJed Brown #if defined(PETSC_USE_COMPLEX)
129c4762a1bSJed Brown             if ((PetscRealPart(v) > 1.0e-10) || (PetscImaginaryPart(v) > 1.0e-10)) {
130c4762a1bSJed Brown               ierr = PetscPrintf(PETSC_COMM_SELF,"ERROR: Difference > 1.0e-10 occurred (delta = (%+1.12e,%+1.12e) [loc %D,%D,%D(%D)])\n",(double)PetscRealPart(test_value),(double)PetscImaginaryPart(test_value),i,j,k,d);CHKERRQ(ierr);
131c4762a1bSJed Brown               dataverified = PETSC_FALSE;
132c4762a1bSJed Brown             }
133c4762a1bSJed Brown #else
134c4762a1bSJed Brown             if (PetscRealPart(v) > 1.0e-10) {
135c4762a1bSJed Brown               ierr = PetscPrintf(PETSC_COMM_SELF,"ERROR: Difference > 1.0e-10 occurred (delta = %+1.12e [loc %D,%D,%D(%D)])\n",(double)PetscRealPart(test_value),i,j,k,d);CHKERRQ(ierr);
136c4762a1bSJed Brown               dataverified = PETSC_FALSE;
137c4762a1bSJed Brown             }
138c4762a1bSJed Brown #endif
139c4762a1bSJed Brown           }
140c4762a1bSJed Brown         }
141c4762a1bSJed Brown       }
142c4762a1bSJed Brown     }
143c4762a1bSJed Brown     if (dataverified) {
144c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"Headerless read of data verified for: %s\n",name);CHKERRQ(ierr);
145c4762a1bSJed Brown     }
146c4762a1bSJed Brown   }
147c4762a1bSJed Brown   PetscFunctionReturn(0);
148c4762a1bSJed Brown }
149c4762a1bSJed Brown 
150c4762a1bSJed Brown PetscErrorCode VecCompare(Vec a,Vec b)
151c4762a1bSJed Brown {
152c4762a1bSJed Brown   PetscInt       locmin[2],locmax[2];
153c4762a1bSJed Brown   PetscReal      min[2],max[2];
154c4762a1bSJed Brown   Vec            ref;
155c4762a1bSJed Brown   PetscErrorCode ierr;
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   PetscFunctionBeginUser;
158c4762a1bSJed Brown   ierr = VecMin(a,&locmin[0],&min[0]);CHKERRQ(ierr);
159c4762a1bSJed Brown   ierr = VecMax(a,&locmax[0],&max[0]);CHKERRQ(ierr);
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   ierr = VecMin(b,&locmin[1],&min[1]);CHKERRQ(ierr);
162c4762a1bSJed Brown   ierr = VecMax(b,&locmax[1],&max[1]);CHKERRQ(ierr);
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"VecCompare\n");CHKERRQ(ierr);
165c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"  min(a)   = %+1.2e [loc %D]\n",(double)min[0],locmin[0]);CHKERRQ(ierr);
166c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"  max(a)   = %+1.2e [loc %D]\n",(double)max[0],locmax[0]);CHKERRQ(ierr);
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"  min(b)   = %+1.2e [loc %D]\n",(double)min[1],locmin[1]);CHKERRQ(ierr);
169c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"  max(b)   = %+1.2e [loc %D]\n",(double)max[1],locmax[1]);CHKERRQ(ierr);
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   ierr = VecDuplicate(a,&ref);CHKERRQ(ierr);
172c4762a1bSJed Brown   ierr = VecCopy(a,ref);CHKERRQ(ierr);
173c4762a1bSJed Brown   ierr = VecAXPY(ref,-1.0,b);CHKERRQ(ierr);
174c4762a1bSJed Brown   ierr = VecMin(ref,&locmin[0],&min[0]);CHKERRQ(ierr);
175c4762a1bSJed Brown   if (PetscAbsReal(min[0]) > 1.0e-10) {
176c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"  ERROR: min(a-b) > 1.0e-10\n");CHKERRQ(ierr);
177c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"  min(a-b) = %+1.10e\n",(double)PetscAbsReal(min[0]));CHKERRQ(ierr);
178c4762a1bSJed Brown   } else {
179c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"  min(a-b) < 1.0e-10\n");CHKERRQ(ierr);
180c4762a1bSJed Brown   }
181c4762a1bSJed Brown   ierr = VecDestroy(&ref);CHKERRQ(ierr);
182c4762a1bSJed Brown   PetscFunctionReturn(0);
183c4762a1bSJed Brown }
184c4762a1bSJed Brown 
185c4762a1bSJed Brown PetscErrorCode TestDMDAVec(PetscBool usempiio)
186c4762a1bSJed Brown {
187c4762a1bSJed Brown   DM             dm;
188c4762a1bSJed Brown   Vec            x_ref,x_test;
189c4762a1bSJed Brown   PetscBool      skipheader = PETSC_TRUE;
190c4762a1bSJed Brown   PetscErrorCode ierr;
191c4762a1bSJed Brown 
192c4762a1bSJed Brown   PetscFunctionBeginUser;
193c4762a1bSJed Brown   if (!usempiio) { ierr = PetscPrintf(PETSC_COMM_WORLD,"%s\n",PETSC_FUNCTION_NAME);CHKERRQ(ierr); }
194c4762a1bSJed Brown   else { ierr = PetscPrintf(PETSC_COMM_WORLD,"%s [using mpi-io]\n",PETSC_FUNCTION_NAME);CHKERRQ(ierr); }
195c4762a1bSJed Brown   ierr = 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,
196c4762a1bSJed Brown                         3,2,NULL,NULL,NULL,&dm);CHKERRQ(ierr);
197c4762a1bSJed Brown   ierr = DMSetFromOptions(dm);CHKERRQ(ierr);
198c4762a1bSJed Brown   ierr = DMSetUp(dm);CHKERRQ(ierr);
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dm,&x_ref);CHKERRQ(ierr);
201c4762a1bSJed Brown   ierr = DMDAVecGenerateEntries(dm,x_ref);CHKERRQ(ierr);
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   if (!usempiio) {
204c4762a1bSJed Brown     ierr = MyVecDump("dmda.pbvec",skipheader,PETSC_FALSE,x_ref);CHKERRQ(ierr);
205c4762a1bSJed Brown   } else {
206c4762a1bSJed Brown     ierr = MyVecDump("dmda-mpiio.pbvec",skipheader,PETSC_TRUE,x_ref);CHKERRQ(ierr);
207c4762a1bSJed Brown   }
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   ierr = DMCreateGlobalVector(dm,&x_test);CHKERRQ(ierr);
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   if (!usempiio) {
212c4762a1bSJed Brown     ierr = MyVecLoad("dmda.pbvec",skipheader,usempiio,x_test);CHKERRQ(ierr);
213c4762a1bSJed Brown   } else {
214c4762a1bSJed Brown     ierr = MyVecLoad("dmda-mpiio.pbvec",skipheader,usempiio,x_test);CHKERRQ(ierr);
215c4762a1bSJed Brown   }
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   ierr = VecCompare(x_ref,x_test);CHKERRQ(ierr);
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   if (!usempiio) {
220c4762a1bSJed Brown     ierr = HeaderlessBinaryReadCheck(dm,"dmda.pbvec");CHKERRQ(ierr);
221c4762a1bSJed Brown   } else {
222c4762a1bSJed Brown     ierr = HeaderlessBinaryReadCheck(dm,"dmda-mpiio.pbvec");CHKERRQ(ierr);
223c4762a1bSJed Brown   }
224c4762a1bSJed Brown   ierr = VecDestroy(&x_ref);CHKERRQ(ierr);
225c4762a1bSJed Brown   ierr = VecDestroy(&x_test);CHKERRQ(ierr);
226c4762a1bSJed Brown   ierr = DMDestroy(&dm);CHKERRQ(ierr);
227c4762a1bSJed Brown   PetscFunctionReturn(0);
228c4762a1bSJed Brown }
229c4762a1bSJed Brown 
230c4762a1bSJed Brown int main(int argc,char **args)
231c4762a1bSJed Brown {
232c4762a1bSJed Brown   PetscErrorCode ierr;
233c4762a1bSJed Brown   PetscBool      usempiio = PETSC_FALSE;
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
236c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-usempiio",&usempiio,NULL);CHKERRQ(ierr);
237c4762a1bSJed Brown   if (!usempiio) {
238c4762a1bSJed Brown     ierr = TestDMDAVec(PETSC_FALSE);CHKERRQ(ierr);
239c4762a1bSJed Brown   } else {
240c4762a1bSJed Brown #if defined(PETSC_HAVE_MPIIO)
241c4762a1bSJed Brown     ierr = TestDMDAVec(PETSC_TRUE);CHKERRQ(ierr);
242c4762a1bSJed Brown #else
243c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Warning: Executing TestDMDAVec(PETSC_TRUE) requires a working MPI-2 implementation\n");CHKERRQ(ierr);
244c4762a1bSJed Brown #endif
245c4762a1bSJed Brown   }
246c4762a1bSJed Brown   ierr = PetscFinalize();
247c4762a1bSJed Brown   return ierr;
248c4762a1bSJed Brown }
249c4762a1bSJed Brown 
250c4762a1bSJed Brown /*TEST
251c4762a1bSJed Brown 
252c4762a1bSJed Brown    test:
253c4762a1bSJed Brown 
254c4762a1bSJed Brown    test:
255c4762a1bSJed Brown       suffix: 2
256c4762a1bSJed Brown       nsize: 12
257c4762a1bSJed Brown 
258c4762a1bSJed Brown    test:
259c4762a1bSJed Brown       suffix: 3
260c4762a1bSJed Brown       nsize: 12
261*dfd57a17SPierre Jolivet       requires: defined(PETSC_HAVE_MPIIO)
262c4762a1bSJed Brown       args: -usempiio
263c4762a1bSJed Brown 
264c4762a1bSJed Brown TEST*/
265c4762a1bSJed Brown 
266