xref: /petsc/src/mat/tests/ex45.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Tests MatView()/MatLoad() with binary viewers for BAIJ matrices.\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include <petscmat.h>
4*c4762a1bSJed Brown #include <petscviewer.h>
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown #include <petsc/private/hashtable.h>
7*c4762a1bSJed Brown static PetscReal MakeValue(PetscInt i,PetscInt j,PetscInt M)
8*c4762a1bSJed Brown {
9*c4762a1bSJed Brown   PetscHash_t h = PetscHashCombine(PetscHashInt(i),PetscHashInt(j));
10*c4762a1bSJed Brown   return (PetscReal) ((h % 5 == 0) ? (1 + i + j*M) : 0);
11*c4762a1bSJed Brown }
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown static PetscErrorCode CheckValuesAIJ(Mat A)
14*c4762a1bSJed Brown {
15*c4762a1bSJed Brown   PetscInt        M,N,rstart,rend,i,j;
16*c4762a1bSJed Brown   PetscReal       v,w;
17*c4762a1bSJed Brown   PetscScalar     val;
18*c4762a1bSJed Brown   PetscErrorCode  ierr;
19*c4762a1bSJed Brown 
20*c4762a1bSJed Brown   PetscFunctionBegin;
21*c4762a1bSJed Brown   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
22*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
23*c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
24*c4762a1bSJed Brown     for (j=0; j<N; j++) {
25*c4762a1bSJed Brown       ierr = MatGetValue(A,i,j,&val);CHKERRQ(ierr);
26*c4762a1bSJed Brown       v = MakeValue(i,j,M); w = PetscRealPart(val);
27*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);
28*c4762a1bSJed Brown     }
29*c4762a1bSJed Brown   }
30*c4762a1bSJed Brown   PetscFunctionReturn(0);
31*c4762a1bSJed Brown }
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown int main(int argc,char **args)
34*c4762a1bSJed Brown {
35*c4762a1bSJed Brown   Mat            A;
36*c4762a1bSJed Brown   PetscInt       M = 24,N = 48,bs = 2;
37*c4762a1bSJed Brown   PetscInt       rstart,rend,i,j;
38*c4762a1bSJed Brown   PetscErrorCode ierr;
39*c4762a1bSJed Brown   PetscViewer    view;
40*c4762a1bSJed Brown 
41*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,NULL,help);if (ierr) return ierr;
42*c4762a1bSJed Brown   /*
43*c4762a1bSJed Brown       Create a parallel BAIJ matrix shared by all processors
44*c4762a1bSJed Brown   */
45*c4762a1bSJed Brown   ierr = MatCreateBAIJ(PETSC_COMM_WORLD,
46*c4762a1bSJed Brown                        bs,
47*c4762a1bSJed Brown                        PETSC_DECIDE,PETSC_DECIDE,
48*c4762a1bSJed Brown                        M,N,
49*c4762a1bSJed Brown                        PETSC_DECIDE,NULL,
50*c4762a1bSJed Brown                        PETSC_DECIDE,NULL,
51*c4762a1bSJed Brown                        &A);CHKERRQ(ierr);
52*c4762a1bSJed Brown 
53*c4762a1bSJed Brown   /*
54*c4762a1bSJed Brown       Set values into the matrix
55*c4762a1bSJed Brown   */
56*c4762a1bSJed Brown   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
57*c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
58*c4762a1bSJed Brown   for (i=rstart; i<rend; i++) {
59*c4762a1bSJed Brown     for (j=0; j<N; j++) {
60*c4762a1bSJed Brown       PetscReal v = MakeValue(i,j,M);
61*c4762a1bSJed Brown       if (PetscAbsReal(v) > 0) {
62*c4762a1bSJed Brown         ierr = MatSetValue(A,i,j,v,INSERT_VALUES);CHKERRQ(ierr);
63*c4762a1bSJed Brown       }
64*c4762a1bSJed Brown     }
65*c4762a1bSJed Brown   }
66*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
67*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68*c4762a1bSJed Brown   ierr = MatViewFromOptions(A,NULL,"-mat_base_view");CHKERRQ(ierr);
69*c4762a1bSJed Brown 
70*c4762a1bSJed Brown   /*
71*c4762a1bSJed Brown       Store the binary matrix to a file
72*c4762a1bSJed Brown   */
73*c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD, "matrix.dat", FILE_MODE_WRITE, &view);CHKERRQ(ierr);
74*c4762a1bSJed Brown   for (i=0; i<3; i++) {
75*c4762a1bSJed Brown     ierr = MatView(A,view);CHKERRQ(ierr);
76*c4762a1bSJed Brown   }
77*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&view);CHKERRQ(ierr);
78*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
79*c4762a1bSJed Brown 
80*c4762a1bSJed Brown   /*
81*c4762a1bSJed Brown       Now reload the matrix and check its values
82*c4762a1bSJed Brown   */
83*c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view);CHKERRQ(ierr);
84*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
85*c4762a1bSJed Brown   ierr = MatSetType(A,MATBAIJ);CHKERRQ(ierr);
86*c4762a1bSJed Brown   for (i=0; i<3; i++) {
87*c4762a1bSJed Brown     if (i > 0) {ierr = MatZeroEntries(A);CHKERRQ(ierr);}
88*c4762a1bSJed Brown     ierr = MatLoad(A,view);CHKERRQ(ierr);
89*c4762a1bSJed Brown     ierr = CheckValuesAIJ(A);CHKERRQ(ierr);
90*c4762a1bSJed Brown   }
91*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&view);CHKERRQ(ierr);
92*c4762a1bSJed Brown   ierr = MatViewFromOptions(A,NULL,"-mat_load_view");CHKERRQ(ierr);
93*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown   /*
96*c4762a1bSJed Brown       Reload in SEQBAIJ matrix and check its values
97*c4762a1bSJed Brown   */
98*c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_SELF,"matrix.dat",FILE_MODE_READ,&view);CHKERRQ(ierr);
99*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_SELF,&A);CHKERRQ(ierr);
100*c4762a1bSJed Brown   ierr = MatSetType(A,MATSEQBAIJ);CHKERRQ(ierr);
101*c4762a1bSJed Brown   for (i=0; i<3; i++) {
102*c4762a1bSJed Brown     if (i > 0) {ierr = MatZeroEntries(A);CHKERRQ(ierr);}
103*c4762a1bSJed Brown     ierr = MatLoad(A,view);CHKERRQ(ierr);
104*c4762a1bSJed Brown     ierr = CheckValuesAIJ(A);CHKERRQ(ierr);
105*c4762a1bSJed Brown   }
106*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&view);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
108*c4762a1bSJed Brown 
109*c4762a1bSJed Brown   /*
110*c4762a1bSJed Brown      Reload in MPIBAIJ matrix and check its values
111*c4762a1bSJed Brown   */
112*c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,"matrix.dat",FILE_MODE_READ,&view);CHKERRQ(ierr);
113*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
114*c4762a1bSJed Brown   ierr = MatSetType(A,MATMPIBAIJ);CHKERRQ(ierr);
115*c4762a1bSJed Brown   for (i=0; i<3; i++) {
116*c4762a1bSJed Brown     if (i > 0) {ierr = MatZeroEntries(A);CHKERRQ(ierr);}
117*c4762a1bSJed Brown     ierr = MatLoad(A,view);CHKERRQ(ierr);
118*c4762a1bSJed Brown     ierr = CheckValuesAIJ(A);CHKERRQ(ierr);
119*c4762a1bSJed Brown   }
120*c4762a1bSJed Brown   ierr = PetscViewerDestroy(&view);CHKERRQ(ierr);
121*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
122*c4762a1bSJed Brown 
123*c4762a1bSJed Brown   ierr = PetscFinalize();
124*c4762a1bSJed Brown   return ierr;
125*c4762a1bSJed Brown }
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown /*TEST
128*c4762a1bSJed Brown 
129*c4762a1bSJed Brown    testset:
130*c4762a1bSJed Brown       args: -viewer_binary_mpiio 0
131*c4762a1bSJed Brown       output_file: output/ex45.out
132*c4762a1bSJed Brown       test:
133*c4762a1bSJed Brown         suffix: stdio_1
134*c4762a1bSJed Brown         nsize: 1
135*c4762a1bSJed Brown       test:
136*c4762a1bSJed Brown         suffix: stdio_2
137*c4762a1bSJed Brown         nsize: 2
138*c4762a1bSJed Brown       test:
139*c4762a1bSJed Brown         suffix: stdio_3
140*c4762a1bSJed Brown         nsize: 3
141*c4762a1bSJed Brown       test:
142*c4762a1bSJed Brown         suffix: stdio_4
143*c4762a1bSJed Brown         nsize: 4
144*c4762a1bSJed Brown       test:
145*c4762a1bSJed Brown         suffix: stdio_5
146*c4762a1bSJed Brown         nsize: 4
147*c4762a1bSJed Brown 
148*c4762a1bSJed Brown    testset:
149*c4762a1bSJed Brown       requires: mpiio
150*c4762a1bSJed Brown       args: -viewer_binary_mpiio 1
151*c4762a1bSJed Brown       output_file: output/ex45.out
152*c4762a1bSJed Brown       test:
153*c4762a1bSJed Brown         suffix: mpiio_1
154*c4762a1bSJed Brown         nsize: 1
155*c4762a1bSJed Brown       test:
156*c4762a1bSJed Brown         suffix: mpiio_2
157*c4762a1bSJed Brown         nsize: 2
158*c4762a1bSJed Brown       test:
159*c4762a1bSJed Brown         suffix: mpiio_3
160*c4762a1bSJed Brown         nsize: 3
161*c4762a1bSJed Brown       test:
162*c4762a1bSJed Brown         suffix: mpiio_4
163*c4762a1bSJed Brown         nsize: 4
164*c4762a1bSJed Brown       test:
165*c4762a1bSJed Brown         suffix: mpiio_5
166*c4762a1bSJed Brown         nsize: 5
167*c4762a1bSJed Brown 
168*c4762a1bSJed Brown TEST*/
169