xref: /petsc/src/mat/tests/ex16.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Tests MatDenseGetArray() and MatView()/MatLoad() with binary viewers.\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petscmat.h>
4*c4762a1bSJed Brown #include <petscviewer.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown static PetscErrorCode CheckValues(Mat A,PetscBool one)
7*c4762a1bSJed Brown {
8*c4762a1bSJed Brown   const PetscScalar *array;
9*c4762a1bSJed Brown   PetscInt          M,N,rstart,rend,lda,i,j;
10*c4762a1bSJed Brown   PetscErrorCode    ierr;
11*c4762a1bSJed Brown 
12*c4762a1bSJed Brown   PetscFunctionBegin;
13*c4762a1bSJed Brown   ierr = MatDenseGetArrayRead(A,&array);CHKERRQ(ierr);
14*c4762a1bSJed Brown   ierr = MatDenseGetLDA(A,&lda);CHKERRQ(ierr);
15*c4762a1bSJed Brown   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
16*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
17*c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
18*c4762a1bSJed Brown     for (j=0; j<N; j++) {
19*c4762a1bSJed Brown       PetscInt ii = i - rstart, jj = j;
20*c4762a1bSJed Brown       PetscReal v = (PetscReal)(one ? 1 : (1 + i + j*M));
21*c4762a1bSJed Brown       PetscReal w = PetscRealPart(array[ii + jj*lda]);
22*c4762a1bSJed 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);
23*c4762a1bSJed Brown     }
24*c4762a1bSJed Brown   }
25*c4762a1bSJed Brown   ierr = MatDenseRestoreArrayRead(A,&array);CHKERRQ(ierr);
26*c4762a1bSJed Brown   PetscFunctionReturn(0);
27*c4762a1bSJed Brown }
28*c4762a1bSJed Brown 
29*c4762a1bSJed Brown #define CheckValuesIJ(A)  CheckValues(A,PETSC_FALSE)
30*c4762a1bSJed Brown #define CheckValuesOne(A) CheckValues(A,PETSC_TRUE)
31*c4762a1bSJed Brown 
32*c4762a1bSJed Brown int main(int argc,char **args)
33*c4762a1bSJed Brown {
34*c4762a1bSJed Brown   Mat            A;
35*c4762a1bSJed Brown   PetscInt       i,j,M = 4,N = 3,rstart,rend;
36*c4762a1bSJed Brown   PetscErrorCode ierr;
37*c4762a1bSJed Brown   PetscScalar    *array;
38*c4762a1bSJed Brown   PetscViewer    view;
39*c4762a1bSJed Brown 
40*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,NULL,help);if (ierr) return ierr;
41*c4762a1bSJed Brown   /*
42*c4762a1bSJed Brown       Create a parallel dense matrix shared by all processors
43*c4762a1bSJed Brown   */
44*c4762a1bSJed Brown   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,M,N,NULL,&A);CHKERRQ(ierr);
45*c4762a1bSJed Brown 
46*c4762a1bSJed Brown   /*
47*c4762a1bSJed Brown      Set values into the matrix
48*c4762a1bSJed Brown   */
49*c4762a1bSJed Brown   for (i=0; i<M; i++) {
50*c4762a1bSJed Brown     for (j=0; j<N; j++) {
51*c4762a1bSJed Brown       PetscScalar v = (PetscReal)(1 + i + j*M);
52*c4762a1bSJed Brown       ierr = MatSetValues(A,1,&i,1,&j,&v,INSERT_VALUES);CHKERRQ(ierr);
53*c4762a1bSJed Brown     }
54*c4762a1bSJed Brown   }
55*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
56*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
57*c4762a1bSJed Brown 
58*c4762a1bSJed Brown   /*
59*c4762a1bSJed Brown       Store the binary matrix to a file
60*c4762a1bSJed Brown   */
61*c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_WRITE, &view);CHKERRQ(ierr);
62*c4762a1bSJed Brown   for (i=0; i<2; i++) {
63*c4762a1bSJed Brown     ierr = MatView(A,view);CHKERRQ(ierr);
64*c4762a1bSJed Brown     ierr = PetscViewerPushFormat(view,PETSC_VIEWER_NATIVE);CHKERRQ(ierr);
65*c4762a1bSJed Brown     ierr = MatView(A,view);CHKERRQ(ierr);
66*c4762a1bSJed Brown     ierr = PetscViewerPopFormat(view);CHKERRQ(ierr);
67*c4762a1bSJed Brown   }
68*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&view);CHKERRQ(ierr);
69*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
70*c4762a1bSJed Brown 
71*c4762a1bSJed Brown   /*
72*c4762a1bSJed Brown       Now reload the matrix and check its values
73*c4762a1bSJed Brown   */
74*c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view);CHKERRQ(ierr);
75*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
76*c4762a1bSJed Brown   ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
77*c4762a1bSJed Brown   for (i=0; i<4; i++) {
78*c4762a1bSJed Brown     if (i > 0) {ierr = MatZeroEntries(A);CHKERRQ(ierr);}
79*c4762a1bSJed Brown     ierr = MatLoad(A,view);CHKERRQ(ierr);
80*c4762a1bSJed Brown     ierr = CheckValuesIJ(A);CHKERRQ(ierr);
81*c4762a1bSJed Brown   }
82*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&view);CHKERRQ(ierr);
83*c4762a1bSJed Brown 
84*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
85*c4762a1bSJed Brown   ierr = PetscMalloc1((rend-rstart)*N,&array);CHKERRQ(ierr);
86*c4762a1bSJed Brown   for (i=0; i<(rend-rstart)*N; i++) array[i] = (PetscReal)1;
87*c4762a1bSJed Brown   ierr = MatDensePlaceArray(A,array);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = CheckValuesOne(A);CHKERRQ(ierr);
89*c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_WRITE,&view);CHKERRQ(ierr);
90*c4762a1bSJed Brown   ierr = MatView(A,view);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = MatDenseResetArray(A);CHKERRQ(ierr);
92*c4762a1bSJed Brown   ierr = PetscFree(array);CHKERRQ(ierr);
93*c4762a1bSJed Brown   ierr = CheckValuesIJ(A);CHKERRQ(ierr);
94*c4762a1bSJed Brown   ierr = PetscViewerBinarySetSkipHeader(view,PETSC_TRUE);CHKERRQ(ierr);
95*c4762a1bSJed Brown   ierr = MatView(A,view);CHKERRQ(ierr);
96*c4762a1bSJed Brown   ierr = PetscViewerBinarySetSkipHeader(view,PETSC_FALSE);CHKERRQ(ierr);
97*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&view);CHKERRQ(ierr);
98*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
99*c4762a1bSJed Brown 
100*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
101*c4762a1bSJed Brown   ierr = MatSetType(A,MATMPIDENSE);CHKERRQ(ierr);
102*c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view);CHKERRQ(ierr);
103*c4762a1bSJed Brown   ierr = MatLoad(A,view);CHKERRQ(ierr);
104*c4762a1bSJed Brown   ierr = CheckValuesOne(A);CHKERRQ(ierr);
105*c4762a1bSJed Brown   ierr = MatZeroEntries(A);CHKERRQ(ierr);
106*c4762a1bSJed Brown   ierr = PetscViewerBinarySetSkipHeader(view,PETSC_TRUE);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = MatLoad(A,view);CHKERRQ(ierr);
108*c4762a1bSJed Brown   ierr = PetscViewerBinarySetSkipHeader(view,PETSC_FALSE);CHKERRQ(ierr);
109*c4762a1bSJed Brown   ierr = CheckValuesIJ(A);CHKERRQ(ierr);
110*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&view);CHKERRQ(ierr);
111*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
112*c4762a1bSJed Brown 
113*c4762a1bSJed Brown   {
114*c4762a1bSJed Brown     PetscInt m = PETSC_DECIDE, n = PETSC_DECIDE;
115*c4762a1bSJed Brown     ierr = PetscSplitOwnership(PETSC_COMM_WORLD,&m,&M);CHKERRQ(ierr);
116*c4762a1bSJed Brown     ierr = PetscSplitOwnership(PETSC_COMM_WORLD,&n,&N);CHKERRQ(ierr);
117*c4762a1bSJed Brown     /* TODO: MatCreateDense requires data!=NULL at all processes! */
118*c4762a1bSJed Brown     ierr = PetscMalloc1(m*N+1,&array);CHKERRQ(ierr);
119*c4762a1bSJed Brown 
120*c4762a1bSJed Brown     ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view);CHKERRQ(ierr);
121*c4762a1bSJed Brown     ierr = MatCreateDense(PETSC_COMM_WORLD,m,n,M,N,array,&A);CHKERRQ(ierr);
122*c4762a1bSJed Brown     ierr = MatLoad(A,view);CHKERRQ(ierr);
123*c4762a1bSJed Brown     ierr = CheckValuesOne(A);CHKERRQ(ierr);
124*c4762a1bSJed Brown     ierr = PetscViewerBinarySetSkipHeader(view,PETSC_TRUE);CHKERRQ(ierr);
125*c4762a1bSJed Brown     ierr = MatLoad(A,view);CHKERRQ(ierr);
126*c4762a1bSJed Brown     ierr = PetscViewerBinarySetSkipHeader(view,PETSC_FALSE);CHKERRQ(ierr);
127*c4762a1bSJed Brown     ierr = CheckValuesIJ(A);CHKERRQ(ierr);
128*c4762a1bSJed Brown     ierr = MatDestroy(&A);CHKERRQ(ierr);
129*c4762a1bSJed Brown     ierr = PetscViewerDestroy(&view);CHKERRQ(ierr);
130*c4762a1bSJed Brown 
131*c4762a1bSJed Brown     ierr = MatCreateDense(PETSC_COMM_WORLD,m,n,M,N,array,&A);CHKERRQ(ierr);
132*c4762a1bSJed Brown     ierr = CheckValuesIJ(A);CHKERRQ(ierr);
133*c4762a1bSJed Brown     ierr = MatDestroy(&A);CHKERRQ(ierr);
134*c4762a1bSJed Brown 
135*c4762a1bSJed Brown     ierr = PetscFree(array);CHKERRQ(ierr);
136*c4762a1bSJed Brown   }
137*c4762a1bSJed Brown 
138*c4762a1bSJed Brown   ierr = PetscFinalize();
139*c4762a1bSJed Brown   return ierr;
140*c4762a1bSJed Brown }
141*c4762a1bSJed Brown 
142*c4762a1bSJed Brown 
143*c4762a1bSJed Brown /*TEST
144*c4762a1bSJed Brown 
145*c4762a1bSJed Brown    testset:
146*c4762a1bSJed Brown       args: -viewer_binary_mpiio 0
147*c4762a1bSJed Brown       output_file: output/ex16.out
148*c4762a1bSJed Brown       test:
149*c4762a1bSJed Brown         suffix: stdio_1
150*c4762a1bSJed Brown         nsize: 1
151*c4762a1bSJed Brown       test:
152*c4762a1bSJed Brown         suffix: stdio_2
153*c4762a1bSJed Brown         nsize: 2
154*c4762a1bSJed Brown       test:
155*c4762a1bSJed Brown         suffix: stdio_3
156*c4762a1bSJed Brown         nsize: 3
157*c4762a1bSJed Brown       test:
158*c4762a1bSJed Brown         suffix: stdio_4
159*c4762a1bSJed Brown         nsize: 4
160*c4762a1bSJed Brown       test:
161*c4762a1bSJed Brown         suffix: stdio_5
162*c4762a1bSJed Brown         nsize: 5
163*c4762a1bSJed Brown 
164*c4762a1bSJed Brown    testset:
165*c4762a1bSJed Brown       requires: mpiio
166*c4762a1bSJed Brown       args: -viewer_binary_mpiio 1
167*c4762a1bSJed Brown       output_file: output/ex16.out
168*c4762a1bSJed Brown       test:
169*c4762a1bSJed Brown         suffix: mpiio_1
170*c4762a1bSJed Brown         nsize: 1
171*c4762a1bSJed Brown       test:
172*c4762a1bSJed Brown         suffix: mpiio_2
173*c4762a1bSJed Brown         nsize: 2
174*c4762a1bSJed Brown       test:
175*c4762a1bSJed Brown         suffix: mpiio_3
176*c4762a1bSJed Brown         nsize: 3
177*c4762a1bSJed Brown       test:
178*c4762a1bSJed Brown         suffix: mpiio_4
179*c4762a1bSJed Brown         nsize: 4
180*c4762a1bSJed Brown       test:
181*c4762a1bSJed Brown         suffix: mpiio_5
182*c4762a1bSJed Brown         nsize: 5
183*c4762a1bSJed Brown 
184*c4762a1bSJed Brown 
185*c4762a1bSJed Brown TEST*/
186