xref: /petsc/src/mat/tests/ex47.c (revision 28b400f66ebc7ae0049166a2294dfcd3df27e64b)
1 
2 static char help[] = "Tests the various routines in MatBAIJ format.\n\
3 Input arguments are:\n\
4   -f <input_file> : file to load.  For example see $PETSC_DIR/share/petsc/datafiles/matrices\n\n";
5 
6 #include <petscmat.h>
7 
8 int main(int argc,char **args)
9 {
10   Mat               A,B,C;
11   PetscViewer       va,vb,vc;
12   Vec               x,y;
13   PetscErrorCode    ierr;
14   PetscInt          i,j,row,m,n,ncols1,ncols2,ct,m2,n2;
15   const PetscInt    *cols1,*cols2;
16   char              file[PETSC_MAX_PATH_LEN];
17   PetscBool         tflg;
18   PetscScalar       rval;
19   const PetscScalar *vals1,*vals2;
20   PetscReal         norm1,norm2,rnorm;
21   PetscRandom       r;
22 
23   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
24   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),NULL));
25 
26   /* Load the matrix as AIJ format */
27   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&va));
28   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
29   CHKERRQ(MatSetType(A,MATSEQAIJ));
30   CHKERRQ(MatLoad(A,va));
31   CHKERRQ(PetscViewerDestroy(&va));
32 
33   /* Load the matrix as BAIJ format */
34   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&vb));
35   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B));
36   CHKERRQ(MatSetType(B,MATSEQBAIJ));
37   CHKERRQ(MatLoad(B,vb));
38   CHKERRQ(PetscViewerDestroy(&vb));
39 
40   /* Load the matrix as BAIJ format */
41   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&vc));
42   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
43   CHKERRQ(MatSetType(C,MATSEQBAIJ));
44   CHKERRQ(MatSetFromOptions(C));
45   CHKERRQ(MatLoad(C,vc));
46   CHKERRQ(PetscViewerDestroy(&vc));
47 
48   CHKERRQ(MatGetSize(A,&m,&n));
49   CHKERRQ(MatGetSize(B,&m2,&n2));
50   PetscCheckFalse(m!=m2,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices are of different size. Cannot run this example");
51 
52   /* Test MatEqual() */
53   CHKERRQ(MatEqual(B,C,&tflg));
54   PetscCheck(tflg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatEqual() failed");
55 
56   /* Test MatGetDiagonal() */
57   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m,&x));
58   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m,&y));
59 
60   CHKERRQ(MatGetDiagonal(A,x));
61   CHKERRQ(MatGetDiagonal(B,y));
62 
63   CHKERRQ(VecEqual(x,y,&tflg));
64   PetscCheck(tflg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetDiagonal() failed");
65 
66   /* Test MatDiagonalScale() */
67   CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&r));
68   CHKERRQ(PetscRandomSetFromOptions(r));
69   CHKERRQ(VecSetRandom(x,r));
70   CHKERRQ(VecSetRandom(y,r));
71 
72   CHKERRQ(MatDiagonalScale(A,x,y));
73   CHKERRQ(MatDiagonalScale(B,x,y));
74   CHKERRQ(MatMult(A,x,y));
75   CHKERRQ(VecNorm(y,NORM_2,&norm1));
76   CHKERRQ(MatMult(B,x,y));
77   CHKERRQ(VecNorm(y,NORM_2,&norm2));
78   rnorm = ((norm1-norm2)*100)/norm1;
79   PetscCheckFalse(rnorm<-0.1 || rnorm>0.01,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatDiagonalScale() failed Norm1 %g Norm2 %g",(double)norm1,(double)norm2);
80 
81   /* Test MatGetRow()/ MatRestoreRow() */
82   for (ct=0; ct<100; ct++) {
83     CHKERRQ(PetscRandomGetValue(r,&rval));
84     row  = (int)(rval*m);
85     CHKERRQ(MatGetRow(A,row,&ncols1,&cols1,&vals1));
86     CHKERRQ(MatGetRow(B,row,&ncols2,&cols2,&vals2));
87 
88     for (i=0,j=0; i<ncols1 && j<ncols2; i++) {
89       while (cols2[j] != cols1[i]) j++;
90       PetscCheckFalse(vals1[i] != vals2[j],PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetRow() failed - vals incorrect.");
91     }
92     PetscCheckFalse(i<ncols1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetRow() failed - cols incorrect");
93 
94     CHKERRQ(MatRestoreRow(A,row,&ncols1,&cols1,&vals1));
95     CHKERRQ(MatRestoreRow(B,row,&ncols2,&cols2,&vals2));
96   }
97 
98   CHKERRQ(MatDestroy(&A));
99   CHKERRQ(MatDestroy(&B));
100   CHKERRQ(MatDestroy(&C));
101   CHKERRQ(VecDestroy(&x));
102   CHKERRQ(VecDestroy(&y));
103   CHKERRQ(PetscRandomDestroy(&r));
104   ierr = PetscFinalize();
105   return ierr;
106 }
107 
108 /*TEST
109 
110    build:
111       requires: !complex
112 
113    test:
114       args: -f ${DATAFILESPATH}/matrices/cfd.1.10 -mat_block_size 5
115       requires: !complex double datafilespath !defined(PETSC_USE_64BIT_INDICES)
116 
117 TEST*/
118