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