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