xref: /petsc/src/mat/tests/ex47.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
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   PetscInt          i,j,row,m,n,ncols1,ncols2,ct,m2,n2;
14   const PetscInt    *cols1,*cols2;
15   char              file[PETSC_MAX_PATH_LEN];
16   PetscBool         tflg;
17   PetscScalar       rval;
18   const PetscScalar *vals1,*vals2;
19   PetscReal         norm1,norm2,rnorm;
20   PetscRandom       r;
21 
22   PetscFunctionBeginUser;
23   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
24   PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),NULL));
25 
26   /* Load the matrix as AIJ format */
27   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&va));
28   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
29   PetscCall(MatSetType(A,MATSEQAIJ));
30   PetscCall(MatLoad(A,va));
31   PetscCall(PetscViewerDestroy(&va));
32 
33   /* Load the matrix as BAIJ format */
34   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&vb));
35   PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
36   PetscCall(MatSetType(B,MATSEQBAIJ));
37   PetscCall(MatLoad(B,vb));
38   PetscCall(PetscViewerDestroy(&vb));
39 
40   /* Load the matrix as BAIJ format */
41   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&vc));
42   PetscCall(MatCreate(PETSC_COMM_WORLD,&C));
43   PetscCall(MatSetType(C,MATSEQBAIJ));
44   PetscCall(MatSetFromOptions(C));
45   PetscCall(MatLoad(C,vc));
46   PetscCall(PetscViewerDestroy(&vc));
47 
48   PetscCall(MatGetSize(A,&m,&n));
49   PetscCall(MatGetSize(B,&m2,&n2));
50   PetscCheck(m==m2,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices are of different size. Cannot run this example");
51 
52   /* Test MatEqual() */
53   PetscCall(MatEqual(B,C,&tflg));
54   PetscCheck(tflg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatEqual() failed");
55 
56   /* Test MatGetDiagonal() */
57   PetscCall(VecCreateSeq(PETSC_COMM_SELF,m,&x));
58   PetscCall(VecCreateSeq(PETSC_COMM_SELF,m,&y));
59 
60   PetscCall(MatGetDiagonal(A,x));
61   PetscCall(MatGetDiagonal(B,y));
62 
63   PetscCall(VecEqual(x,y,&tflg));
64   PetscCheck(tflg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetDiagonal() failed");
65 
66   /* Test MatDiagonalScale() */
67   PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&r));
68   PetscCall(PetscRandomSetFromOptions(r));
69   PetscCall(VecSetRandom(x,r));
70   PetscCall(VecSetRandom(y,r));
71 
72   PetscCall(MatDiagonalScale(A,x,y));
73   PetscCall(MatDiagonalScale(B,x,y));
74   PetscCall(MatMult(A,x,y));
75   PetscCall(VecNorm(y,NORM_2,&norm1));
76   PetscCall(MatMult(B,x,y));
77   PetscCall(VecNorm(y,NORM_2,&norm2));
78   rnorm = ((norm1-norm2)*100)/norm1;
79   PetscCheck(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     PetscCall(PetscRandomGetValue(r,&rval));
84     row  = (int)(rval*m);
85     PetscCall(MatGetRow(A,row,&ncols1,&cols1,&vals1));
86     PetscCall(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       PetscCheck(vals1[i] == vals2[j],PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetRow() failed - vals incorrect.");
91     }
92     PetscCheck(i>=ncols1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetRow() failed - cols incorrect");
93 
94     PetscCall(MatRestoreRow(A,row,&ncols1,&cols1,&vals1));
95     PetscCall(MatRestoreRow(B,row,&ncols2,&cols2,&vals2));
96   }
97 
98   PetscCall(MatDestroy(&A));
99   PetscCall(MatDestroy(&B));
100   PetscCall(MatDestroy(&C));
101   PetscCall(VecDestroy(&x));
102   PetscCall(VecDestroy(&y));
103   PetscCall(PetscRandomDestroy(&r));
104   PetscCall(PetscFinalize());
105   return 0;
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