xref: /petsc/src/mat/tests/ex47.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscInt          i,j,row,m,n,ncols1,ncols2,ct,m2,n2;
14c4762a1bSJed Brown   const PetscInt    *cols1,*cols2;
15c4762a1bSJed Brown   char              file[PETSC_MAX_PATH_LEN];
16c4762a1bSJed Brown   PetscBool         tflg;
17c4762a1bSJed Brown   PetscScalar       rval;
18c4762a1bSJed Brown   const PetscScalar *vals1,*vals2;
19c4762a1bSJed Brown   PetscReal         norm1,norm2,rnorm;
20c4762a1bSJed Brown   PetscRandom       r;
21c4762a1bSJed Brown 
22*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),NULL));
24c4762a1bSJed Brown 
25c4762a1bSJed Brown   /* Load the matrix as AIJ format */
265f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&va));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATSEQAIJ));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(A,va));
305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&va));
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   /* Load the matrix as BAIJ format */
335f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&vb));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(B,MATSEQBAIJ));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(B,vb));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&vb));
38c4762a1bSJed Brown 
39c4762a1bSJed Brown   /* Load the matrix as BAIJ format */
405f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&vc));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(C,MATSEQBAIJ));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(C));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatLoad(C,vc));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerDestroy(&vc));
46c4762a1bSJed Brown 
475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(A,&m,&n));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetSize(B,&m2,&n2));
492c71b3e2SJacob Faibussowitsch   PetscCheckFalse(m!=m2,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices are of different size. Cannot run this example");
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* Test MatEqual() */
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatEqual(B,C,&tflg));
5328b400f6SJacob Faibussowitsch   PetscCheck(tflg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatEqual() failed");
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   /* Test MatGetDiagonal() */
565f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m,&x));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m,&y));
58c4762a1bSJed Brown 
595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetDiagonal(A,x));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatGetDiagonal(B,y));
61c4762a1bSJed Brown 
625f80ce2aSJacob Faibussowitsch   CHKERRQ(VecEqual(x,y,&tflg));
6328b400f6SJacob Faibussowitsch   PetscCheck(tflg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetDiagonal() failed");
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* Test MatDiagonalScale() */
665f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&r));
675f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomSetFromOptions(r));
685f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(x,r));
695f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSetRandom(y,r));
70c4762a1bSJed Brown 
715f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDiagonalScale(A,x,y));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDiagonalScale(B,x,y));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(A,x,y));
745f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(y,NORM_2,&norm1));
755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(B,x,y));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(y,NORM_2,&norm2));
77c4762a1bSJed Brown   rnorm = ((norm1-norm2)*100)/norm1;
782c71b3e2SJacob Faibussowitsch   PetscCheckFalse(rnorm<-0.1 || rnorm>0.01,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatDiagonalScale() failed Norm1 %g Norm2 %g",(double)norm1,(double)norm2);
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* Test MatGetRow()/ MatRestoreRow() */
81c4762a1bSJed Brown   for (ct=0; ct<100; ct++) {
825f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscRandomGetValue(r,&rval));
83c4762a1bSJed Brown     row  = (int)(rval*m);
845f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetRow(A,row,&ncols1,&cols1,&vals1));
855f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetRow(B,row,&ncols2,&cols2,&vals2));
86c4762a1bSJed Brown 
87c4762a1bSJed Brown     for (i=0,j=0; i<ncols1 && j<ncols2; i++) {
88c4762a1bSJed Brown       while (cols2[j] != cols1[i]) j++;
892c71b3e2SJacob Faibussowitsch       PetscCheckFalse(vals1[i] != vals2[j],PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetRow() failed - vals incorrect.");
90c4762a1bSJed Brown     }
912c71b3e2SJacob Faibussowitsch     PetscCheckFalse(i<ncols1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetRow() failed - cols incorrect");
92c4762a1bSJed Brown 
935f80ce2aSJacob Faibussowitsch     CHKERRQ(MatRestoreRow(A,row,&ncols1,&cols1,&vals1));
945f80ce2aSJacob Faibussowitsch     CHKERRQ(MatRestoreRow(B,row,&ncols2,&cols2,&vals2));
95c4762a1bSJed Brown   }
96c4762a1bSJed Brown 
975f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
985f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&B));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&C));
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&y));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscRandomDestroy(&r));
103*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
104*b122ec5aSJacob Faibussowitsch   return 0;
105c4762a1bSJed Brown }
106c4762a1bSJed Brown 
107c4762a1bSJed Brown /*TEST
108c4762a1bSJed Brown 
109c4762a1bSJed Brown    build:
110c4762a1bSJed Brown       requires: !complex
111c4762a1bSJed Brown 
112c4762a1bSJed Brown    test:
113c4762a1bSJed Brown       args: -f ${DATAFILESPATH}/matrices/cfd.1.10 -mat_block_size 5
114dfd57a17SPierre Jolivet       requires: !complex double datafilespath !defined(PETSC_USE_64BIT_INDICES)
115c4762a1bSJed Brown 
116c4762a1bSJed Brown TEST*/
117