xref: /petsc/src/mat/tests/ex53.c (revision dfd57a172ac9fa6c7b5fe6de6ab5df85cefc2996)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests various routines in MatMPIBAIJ format.\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown #define IMAX 15
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat               A,B,C,At,Bt;
9c4762a1bSJed Brown   PetscViewer       fd;
10c4762a1bSJed Brown   char              file[PETSC_MAX_PATH_LEN];
11c4762a1bSJed Brown   PetscRandom       rand;
12c4762a1bSJed Brown   Vec               xx,yy,s1,s2;
13c4762a1bSJed Brown   PetscReal         s1norm,s2norm,rnorm,tol=1.e-10;
14c4762a1bSJed Brown   PetscInt          rstart,rend,rows[2],cols[2],m,n,i,j,M,N,ct,row,ncols1,ncols2,bs;
15c4762a1bSJed Brown   PetscMPIInt       rank,size;
16c4762a1bSJed Brown   PetscErrorCode    ierr = 0;
17c4762a1bSJed Brown   const PetscInt    *cols1,*cols2;
18c4762a1bSJed Brown   PetscScalar       vals1[4],vals2[4],v;
19c4762a1bSJed Brown   const PetscScalar *v1,*v2;
20c4762a1bSJed Brown   PetscBool         flg;
21c4762a1bSJed Brown 
22c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
23ffc4695bSBarry Smith   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRMPI(ierr);
24ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
25c4762a1bSJed Brown 
26c4762a1bSJed Brown   /* Check out if MatLoad() works */
27589a23caSBarry Smith   ierr = PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg);CHKERRQ(ierr);
28d8185827SBarry Smith   if (!flg) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_USER_INPUT,"Input file not specified");
29c4762a1bSJed Brown   ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd);CHKERRQ(ierr);
30c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
31c4762a1bSJed Brown   ierr = MatSetType(A,MATBAIJ);CHKERRQ(ierr);
32c4762a1bSJed Brown   ierr = MatLoad(A,fd);CHKERRQ(ierr);
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   ierr = MatConvert(A,MATAIJ,MAT_INITIAL_MATRIX,&B);CHKERRQ(ierr);
35c4762a1bSJed Brown   ierr = PetscViewerDestroy(&fd);CHKERRQ(ierr);
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_WORLD,&rand);CHKERRQ(ierr);
38c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(rand);CHKERRQ(ierr);
39c4762a1bSJed Brown   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
40c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&xx);CHKERRQ(ierr);
41c4762a1bSJed Brown   ierr = VecSetSizes(xx,m,PETSC_DECIDE);CHKERRQ(ierr);
42c4762a1bSJed Brown   ierr = VecSetFromOptions(xx);CHKERRQ(ierr);
43c4762a1bSJed Brown   ierr = VecDuplicate(xx,&s1);CHKERRQ(ierr);
44c4762a1bSJed Brown   ierr = VecDuplicate(xx,&s2);CHKERRQ(ierr);
45c4762a1bSJed Brown   ierr = VecDuplicate(xx,&yy);CHKERRQ(ierr);
46c4762a1bSJed Brown   ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr);
47c4762a1bSJed Brown 
48c4762a1bSJed Brown   /* Test MatNorm() */
49c4762a1bSJed Brown   ierr  = MatNorm(A,NORM_FROBENIUS,&s1norm);CHKERRQ(ierr);
50c4762a1bSJed Brown   ierr  = MatNorm(B,NORM_FROBENIUS,&s2norm);CHKERRQ(ierr);
51c4762a1bSJed Brown   rnorm = PetscAbsScalar(s2norm-s1norm)/s2norm;
52c4762a1bSJed Brown   if (rnorm>tol) {
53c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatNorm_FROBENIUS()- NormA=%16.14e NormB=%16.14e bs = %D\n",rank,s1norm,s2norm,bs);CHKERRQ(ierr);
54c4762a1bSJed Brown   }
55c4762a1bSJed Brown   ierr  = MatNorm(A,NORM_INFINITY,&s1norm);CHKERRQ(ierr);
56c4762a1bSJed Brown   ierr  = MatNorm(B,NORM_INFINITY,&s2norm);CHKERRQ(ierr);
57c4762a1bSJed Brown   rnorm = PetscAbsScalar(s2norm-s1norm)/s2norm;
58c4762a1bSJed Brown   if (rnorm>tol) {
59c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatNorm_INFINITY()- NormA=%16.14e NormB=%16.14e bs = %D\n",rank,s1norm,s2norm,bs);CHKERRQ(ierr);
60c4762a1bSJed Brown   }
61c4762a1bSJed Brown   ierr  = MatNorm(A,NORM_1,&s1norm);CHKERRQ(ierr);
62c4762a1bSJed Brown   ierr  = MatNorm(B,NORM_1,&s2norm);CHKERRQ(ierr);
63c4762a1bSJed Brown   rnorm = PetscAbsScalar(s2norm-s1norm)/s2norm;
64c4762a1bSJed Brown   if (rnorm>tol) {
65c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatNorm_NORM_1()- NormA=%16.14e NormB=%16.14e bs = %D\n",rank,s1norm,s2norm,bs);CHKERRQ(ierr);
66c4762a1bSJed Brown   }
67c4762a1bSJed Brown 
68c4762a1bSJed Brown   /* Test MatMult() */
69c4762a1bSJed Brown   for (i=0; i<IMAX; i++) {
70c4762a1bSJed Brown     ierr = VecSetRandom(xx,rand);CHKERRQ(ierr);
71c4762a1bSJed Brown     ierr = MatMult(A,xx,s1);CHKERRQ(ierr);
72c4762a1bSJed Brown     ierr = MatMult(B,xx,s2);CHKERRQ(ierr);
73c4762a1bSJed Brown     ierr = VecAXPY(s2,-1.0,s1);CHKERRQ(ierr);
74c4762a1bSJed Brown     ierr = VecNorm(s2,NORM_2,&rnorm);CHKERRQ(ierr);
75c4762a1bSJed Brown     if (rnorm>tol) {
76c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatMult - Norm2=%16.14e bs = %D\n",rank,rnorm,bs);CHKERRQ(ierr);
77c4762a1bSJed Brown     }
78c4762a1bSJed Brown   }
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* test MatMultAdd() */
81c4762a1bSJed Brown   for (i=0; i<IMAX; i++) {
82c4762a1bSJed Brown     ierr = VecSetRandom(xx,rand);CHKERRQ(ierr);
83c4762a1bSJed Brown     ierr = VecSetRandom(yy,rand);CHKERRQ(ierr);
84c4762a1bSJed Brown     ierr = MatMultAdd(A,xx,yy,s1);CHKERRQ(ierr);
85c4762a1bSJed Brown     ierr = MatMultAdd(B,xx,yy,s2);CHKERRQ(ierr);
86c4762a1bSJed Brown     ierr = VecAXPY(s2,-1.0,s1);CHKERRQ(ierr);
87c4762a1bSJed Brown     ierr = VecNorm(s2,NORM_2,&rnorm);CHKERRQ(ierr);
88c4762a1bSJed Brown     if (rnorm>tol) {
89c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatMultAdd - Norm2=%16.14e bs = %D\n",rank,rnorm,bs);CHKERRQ(ierr);
90c4762a1bSJed Brown     }
91c4762a1bSJed Brown   }
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   /* Test MatMultTranspose() */
94c4762a1bSJed Brown   for (i=0; i<IMAX; i++) {
95c4762a1bSJed Brown     ierr  = VecSetRandom(xx,rand);CHKERRQ(ierr);
96c4762a1bSJed Brown     ierr  = MatMultTranspose(A,xx,s1);CHKERRQ(ierr);
97c4762a1bSJed Brown     ierr  = MatMultTranspose(B,xx,s2);CHKERRQ(ierr);
98c4762a1bSJed Brown     ierr  = VecNorm(s1,NORM_2,&s1norm);CHKERRQ(ierr);
99c4762a1bSJed Brown     ierr  = VecNorm(s2,NORM_2,&s2norm);CHKERRQ(ierr);
100c4762a1bSJed Brown     rnorm = s2norm-s1norm;
101c4762a1bSJed Brown     if (rnorm<-tol || rnorm>tol) {
102c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatMultTranspose - Norm1=%16.14e Norm2=%16.14e bs = %D\n",rank,s1norm,s2norm,bs);CHKERRQ(ierr);
103c4762a1bSJed Brown     }
104c4762a1bSJed Brown   }
105c4762a1bSJed Brown   /* Test MatMultTransposeAdd() */
106c4762a1bSJed Brown   for (i=0; i<IMAX; i++) {
107c4762a1bSJed Brown     ierr  = VecSetRandom(xx,rand);CHKERRQ(ierr);
108c4762a1bSJed Brown     ierr  = VecSetRandom(yy,rand);CHKERRQ(ierr);
109c4762a1bSJed Brown     ierr  = MatMultTransposeAdd(A,xx,yy,s1);CHKERRQ(ierr);
110c4762a1bSJed Brown     ierr  = MatMultTransposeAdd(B,xx,yy,s2);CHKERRQ(ierr);
111c4762a1bSJed Brown     ierr  = VecNorm(s1,NORM_2,&s1norm);CHKERRQ(ierr);
112c4762a1bSJed Brown     ierr  = VecNorm(s2,NORM_2,&s2norm);CHKERRQ(ierr);
113c4762a1bSJed Brown     rnorm = s2norm-s1norm;
114c4762a1bSJed Brown     if (rnorm<-tol || rnorm>tol) {
115c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Error: MatMultTransposeAdd - Norm1=%16.14e Norm2=%16.14e bs = %D\n",rank,s1norm,s2norm,bs);CHKERRQ(ierr);
116c4762a1bSJed Brown     }
117c4762a1bSJed Brown   }
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /* Check MatGetValues() */
120c4762a1bSJed Brown   ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr);
121c4762a1bSJed Brown   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   for (i=0; i<IMAX; i++) {
124c4762a1bSJed Brown     /* Create random row numbers ad col numbers */
125c4762a1bSJed Brown     ierr    = PetscRandomGetValue(rand,&v);CHKERRQ(ierr);
126c4762a1bSJed Brown     cols[0] = (int)(PetscRealPart(v)*N);
127c4762a1bSJed Brown     ierr    = PetscRandomGetValue(rand,&v);CHKERRQ(ierr);
128c4762a1bSJed Brown     cols[1] = (int)(PetscRealPart(v)*N);
129c4762a1bSJed Brown     ierr    = PetscRandomGetValue(rand,&v);CHKERRQ(ierr);
130c4762a1bSJed Brown     rows[0] = rstart + (int)(PetscRealPart(v)*m);
131c4762a1bSJed Brown     ierr    = PetscRandomGetValue(rand,&v);CHKERRQ(ierr);
132c4762a1bSJed Brown     rows[1] = rstart + (int)(PetscRealPart(v)*m);
133c4762a1bSJed Brown 
134c4762a1bSJed Brown     ierr = MatGetValues(A,2,rows,2,cols,vals1);CHKERRQ(ierr);
135c4762a1bSJed Brown     ierr = MatGetValues(B,2,rows,2,cols,vals2);CHKERRQ(ierr);
136c4762a1bSJed Brown 
137c4762a1bSJed Brown     for (j=0; j<4; j++) {
138c4762a1bSJed Brown       if (vals1[j] != vals2[j]) {
139c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF,"[%d]: Error: MatGetValues rstart = %2d  row = %2d col = %2d val1 = %e val2 = %e bs = %D\n",rank,rstart,rows[j/2],cols[j%2],PetscRealPart(vals1[j]),PetscRealPart(vals2[j]),bs);CHKERRQ(ierr);
140c4762a1bSJed Brown       }
141c4762a1bSJed Brown     }
142c4762a1bSJed Brown   }
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* Test MatGetRow()/ MatRestoreRow() */
145c4762a1bSJed Brown   for (ct=0; ct<100; ct++) {
146c4762a1bSJed Brown     ierr = PetscRandomGetValue(rand,&v);CHKERRQ(ierr);
147c4762a1bSJed Brown     row  = rstart + (int)(PetscRealPart(v)*m);
148c4762a1bSJed Brown     ierr = MatGetRow(A,row,&ncols1,&cols1,&v1);CHKERRQ(ierr);
149c4762a1bSJed Brown     ierr = MatGetRow(B,row,&ncols2,&cols2,&v2);CHKERRQ(ierr);
150c4762a1bSJed Brown 
151c4762a1bSJed Brown     for (i=0,j=0; i<ncols1 && j<ncols2; j++) {
152c4762a1bSJed Brown       while (cols2[j] != cols1[i]) i++;
153c4762a1bSJed Brown       if (v1[i] != v2[j]) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetRow() failed - vals incorrect.");
154c4762a1bSJed Brown     }
155c4762a1bSJed Brown     if (j<ncols2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetRow() failed - cols incorrect");
156c4762a1bSJed Brown 
157c4762a1bSJed Brown     ierr = MatRestoreRow(A,row,&ncols1,&cols1,&v1);CHKERRQ(ierr);
158c4762a1bSJed Brown     ierr = MatRestoreRow(B,row,&ncols2,&cols2,&v2);CHKERRQ(ierr);
159c4762a1bSJed Brown   }
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /* Test MatConvert() */
162c4762a1bSJed Brown   ierr = MatConvert(A,MATSAME,MAT_INITIAL_MATRIX,&C);CHKERRQ(ierr);
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   /* See if MatMult Says both are same */
165c4762a1bSJed Brown   for (i=0; i<IMAX; i++) {
166c4762a1bSJed Brown     ierr  = VecSetRandom(xx,rand);CHKERRQ(ierr);
167c4762a1bSJed Brown     ierr  = MatMult(A,xx,s1);CHKERRQ(ierr);
168c4762a1bSJed Brown     ierr  = MatMult(C,xx,s2);CHKERRQ(ierr);
169c4762a1bSJed Brown     ierr  = VecNorm(s1,NORM_2,&s1norm);CHKERRQ(ierr);
170c4762a1bSJed Brown     ierr  = VecNorm(s2,NORM_2,&s2norm);CHKERRQ(ierr);
171c4762a1bSJed Brown     rnorm = s2norm-s1norm;
172c4762a1bSJed Brown     if (rnorm<-tol || rnorm>tol) {
173c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Error in MatConvert: MatMult - Norm1=%16.14e Norm2=%16.14e bs = %D\n",rank,s1norm,s2norm,bs);CHKERRQ(ierr);
174c4762a1bSJed Brown     }
175c4762a1bSJed Brown   }
176c4762a1bSJed Brown   ierr = MatDestroy(&C);CHKERRQ(ierr);
177c4762a1bSJed Brown 
178c4762a1bSJed Brown   /* Test MatTranspose() */
179c4762a1bSJed Brown   ierr = MatTranspose(A,MAT_INITIAL_MATRIX,&At);CHKERRQ(ierr);
180c4762a1bSJed Brown   ierr = MatTranspose(B,MAT_INITIAL_MATRIX,&Bt);CHKERRQ(ierr);
181c4762a1bSJed Brown   for (i=0; i<IMAX; i++) {
182c4762a1bSJed Brown     ierr  = VecSetRandom(xx,rand);CHKERRQ(ierr);
183c4762a1bSJed Brown     ierr  = MatMult(At,xx,s1);CHKERRQ(ierr);
184c4762a1bSJed Brown     ierr  = MatMult(Bt,xx,s2);CHKERRQ(ierr);
185c4762a1bSJed Brown     ierr  = VecNorm(s1,NORM_2,&s1norm);CHKERRQ(ierr);
186c4762a1bSJed Brown     ierr  = VecNorm(s2,NORM_2,&s2norm);CHKERRQ(ierr);
187c4762a1bSJed Brown     rnorm = s2norm-s1norm;
188c4762a1bSJed Brown     if (rnorm<-tol || rnorm>tol) {
189c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Error in MatConvert:MatMult - Norm1=%16.14e Norm2=%16.14e bs = %D\n",rank,s1norm,s2norm,bs);CHKERRQ(ierr);
190c4762a1bSJed Brown     }
191c4762a1bSJed Brown   }
192c4762a1bSJed Brown   ierr = MatDestroy(&At);CHKERRQ(ierr);
193c4762a1bSJed Brown   ierr = MatDestroy(&Bt);CHKERRQ(ierr);
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
196c4762a1bSJed Brown   ierr = MatDestroy(&B);CHKERRQ(ierr);
197c4762a1bSJed Brown   ierr = VecDestroy(&xx);CHKERRQ(ierr);
198c4762a1bSJed Brown   ierr = VecDestroy(&yy);CHKERRQ(ierr);
199c4762a1bSJed Brown   ierr = VecDestroy(&s1);CHKERRQ(ierr);
200c4762a1bSJed Brown   ierr = VecDestroy(&s2);CHKERRQ(ierr);
201c4762a1bSJed Brown   ierr = PetscRandomDestroy(&rand);CHKERRQ(ierr);
202c4762a1bSJed Brown   ierr = PetscFinalize();
203c4762a1bSJed Brown   return ierr;
204c4762a1bSJed Brown }
205c4762a1bSJed Brown 
206c4762a1bSJed Brown /*TEST
207c4762a1bSJed Brown 
208c4762a1bSJed Brown    build:
209c4762a1bSJed Brown       requires: !complex
210c4762a1bSJed Brown 
211c4762a1bSJed Brown    test:
212*dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
213c4762a1bSJed Brown       nsize: 3
214c4762a1bSJed Brown       args: -matload_block_size 1 -f ${DATAFILESPATH}/matrices/small
215c4762a1bSJed Brown 
216c4762a1bSJed Brown    test:
217c4762a1bSJed Brown       suffix: 2
218*dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
219c4762a1bSJed Brown       nsize: 3
220c4762a1bSJed Brown       args: -matload_block_size 2 -f ${DATAFILESPATH}/matrices/small
221c4762a1bSJed Brown       output_file: output/ex53_1.out
222c4762a1bSJed Brown 
223c4762a1bSJed Brown    test:
224c4762a1bSJed Brown       suffix: 3
225*dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
226c4762a1bSJed Brown       nsize: 3
227c4762a1bSJed Brown       args: -matload_block_size 4 -f ${DATAFILESPATH}/matrices/small
228c4762a1bSJed Brown       output_file: output/ex53_1.out
229c4762a1bSJed Brown 
230c4762a1bSJed Brown    test:
231c4762a1bSJed Brown       TODO: Matrix row/column sizes are not compatible with block size
232c4762a1bSJed Brown       suffix: 4
233*dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
234c4762a1bSJed Brown       nsize: 3
235c4762a1bSJed Brown       args: -matload_block_size 5 -f ${DATAFILESPATH}/matrices/small
236c4762a1bSJed Brown       output_file: output/ex53_1.out
237c4762a1bSJed Brown 
238c4762a1bSJed Brown    test:
239c4762a1bSJed Brown       suffix: 5
240*dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
241c4762a1bSJed Brown       nsize: 3
242c4762a1bSJed Brown       args: -matload_block_size 6 -f ${DATAFILESPATH}/matrices/small
243c4762a1bSJed Brown       output_file: output/ex53_1.out
244c4762a1bSJed Brown 
245c4762a1bSJed Brown    test:
246c4762a1bSJed Brown       TODO: Matrix row/column sizes are not compatible with block size
247c4762a1bSJed Brown       suffix: 6
248*dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
249c4762a1bSJed Brown       nsize: 3
250c4762a1bSJed Brown       args: -matload_block_size 7 -f ${DATAFILESPATH}/matrices/small
251c4762a1bSJed Brown       output_file: output/ex53_1.out
252c4762a1bSJed Brown 
253c4762a1bSJed Brown    test:
254c4762a1bSJed Brown       TODO: Matrix row/column sizes are not compatible with block size
255c4762a1bSJed Brown       suffix: 7
256*dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
257c4762a1bSJed Brown       nsize: 3
258c4762a1bSJed Brown       args: -matload_block_size 8 -f ${DATAFILESPATH}/matrices/small
259c4762a1bSJed Brown       output_file: output/ex53_1.out
260c4762a1bSJed Brown 
261c4762a1bSJed Brown    test:
262c4762a1bSJed Brown       suffix: 8
263*dfd57a17SPierre Jolivet       requires: datafilespath !complex double !defined(PETSC_USE_64BIT_INDICES)
264c4762a1bSJed Brown       nsize: 4
265c4762a1bSJed Brown       args: -matload_block_size 3 -f ${DATAFILESPATH}/matrices/small
266c4762a1bSJed Brown       output_file: output/ex53_1.out
267c4762a1bSJed Brown 
268c4762a1bSJed Brown TEST*/
269