xref: /petsc/src/mat/tests/ex16.c (revision a5225ed31aa9ee76d5b3e7690f8edc53c82a15a1)
1c4762a1bSJed Brown static char help[] = "Tests MatDenseGetArray() and MatView()/MatLoad() with binary viewers.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown #include <petscviewer.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown static PetscErrorCode CheckValues(Mat A,PetscBool one)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   const PetscScalar *array;
9c4762a1bSJed Brown   PetscInt          M,N,rstart,rend,lda,i,j;
10c4762a1bSJed Brown   PetscErrorCode    ierr;
11c4762a1bSJed Brown 
12c4762a1bSJed Brown   PetscFunctionBegin;
13c4762a1bSJed Brown   ierr = MatDenseGetArrayRead(A,&array);CHKERRQ(ierr);
14c4762a1bSJed Brown   ierr = MatDenseGetLDA(A,&lda);CHKERRQ(ierr);
15c4762a1bSJed Brown   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
16c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
17c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
18c4762a1bSJed Brown     for (j=0; j<N; j++) {
19c4762a1bSJed Brown       PetscInt ii = i - rstart, jj = j;
20c4762a1bSJed Brown       PetscReal v = (PetscReal)(one ? 1 : (1 + i + j*M));
21c4762a1bSJed Brown       PetscReal w = PetscRealPart(array[ii + jj*lda]);
22c4762a1bSJed Brown       if (PetscAbsReal(v-w) > 0) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Matrix entry (%D,%D) should be %g, got %g",i,j,(double)v,(double)w);
23c4762a1bSJed Brown     }
24c4762a1bSJed Brown   }
25c4762a1bSJed Brown   ierr = MatDenseRestoreArrayRead(A,&array);CHKERRQ(ierr);
26c4762a1bSJed Brown   PetscFunctionReturn(0);
27c4762a1bSJed Brown }
28c4762a1bSJed Brown 
29c4762a1bSJed Brown #define CheckValuesIJ(A)  CheckValues(A,PETSC_FALSE)
30c4762a1bSJed Brown #define CheckValuesOne(A) CheckValues(A,PETSC_TRUE)
31c4762a1bSJed Brown 
32c4762a1bSJed Brown int main(int argc,char **args)
33c4762a1bSJed Brown {
34c4762a1bSJed Brown   Mat            A;
35c4762a1bSJed Brown   PetscInt       i,j,M = 4,N = 3,rstart,rend;
36c4762a1bSJed Brown   PetscErrorCode ierr;
37c4762a1bSJed Brown   PetscScalar    *array;
38*a5225ed3SStefano Zampini   char           mattype[256];
39c4762a1bSJed Brown   PetscViewer    view;
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,NULL,help);if (ierr) return ierr;
42*a5225ed3SStefano Zampini   ierr = PetscStrcpy(mattype,MATMPIDENSE);CHKERRQ(ierr);
43*a5225ed3SStefano Zampini   ierr = PetscOptionsGetString(NULL,NULL,"-mat_type",mattype,sizeof(mattype),NULL);CHKERRQ(ierr);
44c4762a1bSJed Brown   /*
45c4762a1bSJed Brown       Create a parallel dense matrix shared by all processors
46c4762a1bSJed Brown   */
47c4762a1bSJed Brown   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,M,N,NULL,&A);CHKERRQ(ierr);
48*a5225ed3SStefano Zampini   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
49*a5225ed3SStefano Zampini   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
50*a5225ed3SStefano Zampini   ierr = MatConvert(A,mattype,MAT_INPLACE_MATRIX,&A);CHKERRQ(ierr);
51c4762a1bSJed Brown   /*
52c4762a1bSJed Brown      Set values into the matrix
53c4762a1bSJed Brown   */
54c4762a1bSJed Brown   for (i=0; i<M; i++) {
55c4762a1bSJed Brown     for (j=0; j<N; j++) {
56c4762a1bSJed Brown       PetscScalar v = (PetscReal)(1 + i + j*M);
57c4762a1bSJed Brown       ierr = MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);
58c4762a1bSJed Brown     }
59c4762a1bSJed Brown   }
60c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
61c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
62*a5225ed3SStefano Zampini   ierr = MatScale(A,2.0);CHKERRQ(ierr);
63*a5225ed3SStefano Zampini   ierr = MatScale(A,1.0/2.0);CHKERRQ(ierr);
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /*
66c4762a1bSJed Brown       Store the binary matrix to a file
67c4762a1bSJed Brown   */
68c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_WRITE, &view);CHKERRQ(ierr);
69c4762a1bSJed Brown   for (i=0; i<2; i++) {
70c4762a1bSJed Brown     ierr = MatView(A,view);CHKERRQ(ierr);
71c4762a1bSJed Brown     ierr = PetscViewerPushFormat(view,PETSC_VIEWER_NATIVE);CHKERRQ(ierr);
72c4762a1bSJed Brown     ierr = MatView(A,view);CHKERRQ(ierr);
73c4762a1bSJed Brown     ierr = PetscViewerPopFormat(view);CHKERRQ(ierr);
74c4762a1bSJed Brown   }
75c4762a1bSJed Brown   ierr = PetscViewerDestroy(&view);CHKERRQ(ierr);
76c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   /*
79c4762a1bSJed Brown       Now reload the matrix and check its values
80c4762a1bSJed Brown   */
81c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view);CHKERRQ(ierr);
82c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
83*a5225ed3SStefano Zampini   ierr = MatSetType(A,mattype);CHKERRQ(ierr);
84c4762a1bSJed Brown   for (i=0; i<4; i++) {
85c4762a1bSJed Brown     if (i > 0) {ierr = MatZeroEntries(A);CHKERRQ(ierr);}
86c4762a1bSJed Brown     ierr = MatLoad(A,view);CHKERRQ(ierr);
87c4762a1bSJed Brown     ierr = CheckValuesIJ(A);CHKERRQ(ierr);
88c4762a1bSJed Brown   }
89c4762a1bSJed Brown   ierr = PetscViewerDestroy(&view);CHKERRQ(ierr);
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr = PetscMalloc1((rend-rstart)*N,&array);CHKERRQ(ierr);
93c4762a1bSJed Brown   for (i=0; i<(rend-rstart)*N; i++) array[i] = (PetscReal)1;
94c4762a1bSJed Brown   ierr = MatDensePlaceArray(A,array);CHKERRQ(ierr);
95*a5225ed3SStefano Zampini   ierr = MatScale(A,2.0);CHKERRQ(ierr);
96*a5225ed3SStefano Zampini   ierr = MatScale(A,1.0/2.0);CHKERRQ(ierr);
97c4762a1bSJed Brown   ierr = CheckValuesOne(A);CHKERRQ(ierr);
98c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_WRITE,&view);CHKERRQ(ierr);
99c4762a1bSJed Brown   ierr = MatView(A,view);CHKERRQ(ierr);
100c4762a1bSJed Brown   ierr = MatDenseResetArray(A);CHKERRQ(ierr);
101c4762a1bSJed Brown   ierr = PetscFree(array);CHKERRQ(ierr);
102c4762a1bSJed Brown   ierr = CheckValuesIJ(A);CHKERRQ(ierr);
103c4762a1bSJed Brown   ierr = PetscViewerBinarySetSkipHeader(view,PETSC_TRUE);CHKERRQ(ierr);
104c4762a1bSJed Brown   ierr = MatView(A,view);CHKERRQ(ierr);
105c4762a1bSJed Brown   ierr = PetscViewerBinarySetSkipHeader(view,PETSC_FALSE);CHKERRQ(ierr);
106c4762a1bSJed Brown   ierr = PetscViewerDestroy(&view);CHKERRQ(ierr);
107c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
110*a5225ed3SStefano Zampini   ierr = MatSetType(A,mattype);CHKERRQ(ierr);
111c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view);CHKERRQ(ierr);
112c4762a1bSJed Brown   ierr = MatLoad(A,view);CHKERRQ(ierr);
113c4762a1bSJed Brown   ierr = CheckValuesOne(A);CHKERRQ(ierr);
114c4762a1bSJed Brown   ierr = MatZeroEntries(A);CHKERRQ(ierr);
115c4762a1bSJed Brown   ierr = PetscViewerBinarySetSkipHeader(view,PETSC_TRUE);CHKERRQ(ierr);
116c4762a1bSJed Brown   ierr = MatLoad(A,view);CHKERRQ(ierr);
117c4762a1bSJed Brown   ierr = PetscViewerBinarySetSkipHeader(view,PETSC_FALSE);CHKERRQ(ierr);
118c4762a1bSJed Brown   ierr = CheckValuesIJ(A);CHKERRQ(ierr);
119c4762a1bSJed Brown   ierr = PetscViewerDestroy(&view);CHKERRQ(ierr);
120c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   {
123c4762a1bSJed Brown     PetscInt m = PETSC_DECIDE, n = PETSC_DECIDE;
124c4762a1bSJed Brown     ierr = PetscSplitOwnership(PETSC_COMM_WORLD,&m,&M);CHKERRQ(ierr);
125c4762a1bSJed Brown     ierr = PetscSplitOwnership(PETSC_COMM_WORLD,&n,&N);CHKERRQ(ierr);
126c4762a1bSJed Brown     /* TODO: MatCreateDense requires data!=NULL at all processes! */
127c4762a1bSJed Brown     ierr = PetscMalloc1(m*N+1,&array);CHKERRQ(ierr);
128c4762a1bSJed Brown 
129c4762a1bSJed Brown     ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view);CHKERRQ(ierr);
130c4762a1bSJed Brown     ierr = MatCreateDense(PETSC_COMM_WORLD,m,n,M,N,array,&A);CHKERRQ(ierr);
131c4762a1bSJed Brown     ierr = MatLoad(A,view);CHKERRQ(ierr);
132c4762a1bSJed Brown     ierr = CheckValuesOne(A);CHKERRQ(ierr);
133c4762a1bSJed Brown     ierr = PetscViewerBinarySetSkipHeader(view,PETSC_TRUE);CHKERRQ(ierr);
134c4762a1bSJed Brown     ierr = MatLoad(A,view);CHKERRQ(ierr);
135c4762a1bSJed Brown     ierr = PetscViewerBinarySetSkipHeader(view,PETSC_FALSE);CHKERRQ(ierr);
136c4762a1bSJed Brown     ierr = CheckValuesIJ(A);CHKERRQ(ierr);
137c4762a1bSJed Brown     ierr = MatDestroy(&A);CHKERRQ(ierr);
138c4762a1bSJed Brown     ierr = PetscViewerDestroy(&view);CHKERRQ(ierr);
139c4762a1bSJed Brown 
140c4762a1bSJed Brown     ierr = MatCreateDense(PETSC_COMM_WORLD,m,n,M,N,array,&A);CHKERRQ(ierr);
141c4762a1bSJed Brown     ierr = CheckValuesIJ(A);CHKERRQ(ierr);
142c4762a1bSJed Brown     ierr = MatDestroy(&A);CHKERRQ(ierr);
143c4762a1bSJed Brown 
144c4762a1bSJed Brown     ierr = PetscFree(array);CHKERRQ(ierr);
145c4762a1bSJed Brown   }
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   ierr = PetscFinalize();
148c4762a1bSJed Brown   return ierr;
149c4762a1bSJed Brown }
150c4762a1bSJed Brown 
151c4762a1bSJed Brown 
152c4762a1bSJed Brown /*TEST
153c4762a1bSJed Brown 
154c4762a1bSJed Brown    testset:
155c4762a1bSJed Brown       args: -viewer_binary_mpiio 0
156c4762a1bSJed Brown       output_file: output/ex16.out
157c4762a1bSJed Brown       test:
158c4762a1bSJed Brown         suffix: stdio_1
159c4762a1bSJed Brown         nsize: 1
160*a5225ed3SStefano Zampini         args: -mat_type seqdense
161c4762a1bSJed Brown       test:
162c4762a1bSJed Brown         suffix: stdio_2
163c4762a1bSJed Brown         nsize: 2
164c4762a1bSJed Brown       test:
165c4762a1bSJed Brown         suffix: stdio_3
166c4762a1bSJed Brown         nsize: 3
167c4762a1bSJed Brown       test:
168c4762a1bSJed Brown         suffix: stdio_4
169c4762a1bSJed Brown         nsize: 4
170c4762a1bSJed Brown       test:
171c4762a1bSJed Brown         suffix: stdio_5
172c4762a1bSJed Brown         nsize: 5
173*a5225ed3SStefano Zampini       test:
174*a5225ed3SStefano Zampini         requires: cuda
175*a5225ed3SStefano Zampini         args: -mat_type seqdensecuda
176*a5225ed3SStefano Zampini         suffix: stdio_cuda_1
177*a5225ed3SStefano Zampini         nsize: 1
178*a5225ed3SStefano Zampini       test:
179*a5225ed3SStefano Zampini         requires: cuda
180*a5225ed3SStefano Zampini         args: -mat_type mpidensecuda
181*a5225ed3SStefano Zampini         suffix: stdio_cuda_2
182*a5225ed3SStefano Zampini         nsize: 2
183*a5225ed3SStefano Zampini       test:
184*a5225ed3SStefano Zampini         requires: cuda
185*a5225ed3SStefano Zampini         args: -mat_type mpidensecuda
186*a5225ed3SStefano Zampini         suffix: stdio_cuda_3
187*a5225ed3SStefano Zampini         nsize: 3
188*a5225ed3SStefano Zampini       test:
189*a5225ed3SStefano Zampini         requires: cuda
190*a5225ed3SStefano Zampini         args: -mat_type mpidensecuda
191*a5225ed3SStefano Zampini         suffix: stdio_cuda_4
192*a5225ed3SStefano Zampini         nsize: 4
193*a5225ed3SStefano Zampini       test:
194*a5225ed3SStefano Zampini         requires: cuda
195*a5225ed3SStefano Zampini         args: -mat_type mpidensecuda
196*a5225ed3SStefano Zampini         suffix: stdio_cuda_5
197*a5225ed3SStefano Zampini         nsize: 5
198c4762a1bSJed Brown 
199c4762a1bSJed Brown    testset:
200c4762a1bSJed Brown       requires: mpiio
201c4762a1bSJed Brown       args: -viewer_binary_mpiio 1
202c4762a1bSJed Brown       output_file: output/ex16.out
203c4762a1bSJed Brown       test:
204c4762a1bSJed Brown         suffix: mpiio_1
205c4762a1bSJed Brown         nsize: 1
206c4762a1bSJed Brown       test:
207c4762a1bSJed Brown         suffix: mpiio_2
208c4762a1bSJed Brown         nsize: 2
209c4762a1bSJed Brown       test:
210c4762a1bSJed Brown         suffix: mpiio_3
211c4762a1bSJed Brown         nsize: 3
212c4762a1bSJed Brown       test:
213c4762a1bSJed Brown         suffix: mpiio_4
214c4762a1bSJed Brown         nsize: 4
215c4762a1bSJed Brown       test:
216c4762a1bSJed Brown         suffix: mpiio_5
217c4762a1bSJed Brown         nsize: 5
218*a5225ed3SStefano Zampini       test:
219*a5225ed3SStefano Zampini         requires: cuda
220*a5225ed3SStefano Zampini         args: -mat_type mpidensecuda
221*a5225ed3SStefano Zampini         suffix: mpiio_cuda_1
222*a5225ed3SStefano Zampini         nsize: 1
223*a5225ed3SStefano Zampini       test:
224*a5225ed3SStefano Zampini         requires: cuda
225*a5225ed3SStefano Zampini         args: -mat_type mpidensecuda
226*a5225ed3SStefano Zampini         suffix: mpiio_cuda_2
227*a5225ed3SStefano Zampini         nsize: 2
228*a5225ed3SStefano Zampini       test:
229*a5225ed3SStefano Zampini         requires: cuda
230*a5225ed3SStefano Zampini         args: -mat_type mpidensecuda
231*a5225ed3SStefano Zampini         suffix: mpiio_cuda_3
232*a5225ed3SStefano Zampini         nsize: 3
233*a5225ed3SStefano Zampini       test:
234*a5225ed3SStefano Zampini         requires: cuda
235*a5225ed3SStefano Zampini         args: -mat_type mpidensecuda
236*a5225ed3SStefano Zampini         suffix: mpiio_cuda_4
237*a5225ed3SStefano Zampini         nsize: 4
238*a5225ed3SStefano Zampini       test:
239*a5225ed3SStefano Zampini         requires: cuda
240*a5225ed3SStefano Zampini         args: -mat_type mpidensecuda
241*a5225ed3SStefano Zampini         suffix: mpiio_cuda_5
242*a5225ed3SStefano Zampini         nsize: 5
243c4762a1bSJed Brown 
244c4762a1bSJed Brown TEST*/
245