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