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; 245f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),NULL)); 25c4762a1bSJed Brown 26c4762a1bSJed Brown /* Load the matrix as AIJ format */ 275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&va)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATSEQAIJ)); 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(A,va)); 315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&va)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* Load the matrix as BAIJ format */ 345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&vb)); 355f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&B)); 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(B,MATSEQBAIJ)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(B,vb)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&vb)); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* Load the matrix as BAIJ format */ 415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&vc)); 425f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 435f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(C,MATSEQBAIJ)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(C,vc)); 465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&vc)); 47c4762a1bSJed Brown 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(A,&m,&n)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(B,&m2,&n2)); 502c71b3e2SJacob Faibussowitsch PetscCheckFalse(m!=m2,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices are of different size. Cannot run this example"); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Test MatEqual() */ 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatEqual(B,C,&tflg)); 54*28b400f6SJacob Faibussowitsch PetscCheck(tflg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatEqual() failed"); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* Test MatGetDiagonal() */ 575f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m,&x)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,m,&y)); 59c4762a1bSJed Brown 605f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(A,x)); 615f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonal(B,y)); 62c4762a1bSJed Brown 635f80ce2aSJacob Faibussowitsch CHKERRQ(VecEqual(x,y,&tflg)); 64*28b400f6SJacob Faibussowitsch PetscCheck(tflg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetDiagonal() failed"); 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* Test MatDiagonalScale() */ 675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomCreate(PETSC_COMM_SELF,&r)); 685f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomSetFromOptions(r)); 695f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(x,r)); 705f80ce2aSJacob Faibussowitsch CHKERRQ(VecSetRandom(y,r)); 71c4762a1bSJed Brown 725f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(A,x,y)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatDiagonalScale(B,x,y)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(A,x,y)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&norm1)); 765f80ce2aSJacob Faibussowitsch CHKERRQ(MatMult(B,x,y)); 775f80ce2aSJacob Faibussowitsch CHKERRQ(VecNorm(y,NORM_2,&norm2)); 78c4762a1bSJed Brown rnorm = ((norm1-norm2)*100)/norm1; 792c71b3e2SJacob Faibussowitsch PetscCheckFalse(rnorm<-0.1 || rnorm>0.01,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++) { 835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomGetValue(r,&rval)); 84c4762a1bSJed Brown row = (int)(rval*m); 855f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRow(A,row,&ncols1,&cols1,&vals1)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRow(B,row,&ncols2,&cols2,&vals2)); 87c4762a1bSJed Brown 88c4762a1bSJed Brown for (i=0,j=0; i<ncols1 && j<ncols2; i++) { 89c4762a1bSJed Brown while (cols2[j] != cols1[i]) j++; 902c71b3e2SJacob Faibussowitsch PetscCheckFalse(vals1[i] != vals2[j],PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetRow() failed - vals incorrect."); 91c4762a1bSJed Brown } 922c71b3e2SJacob Faibussowitsch PetscCheckFalse(i<ncols1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetRow() failed - cols incorrect"); 93c4762a1bSJed Brown 945f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRow(A,row,&ncols1,&cols1,&vals1)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRow(B,row,&ncols2,&cols2,&vals2)); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 985f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 995f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&B)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 1015f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 1025f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 1035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscRandomDestroy(&r)); 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 115dfd57a17SPierre Jolivet requires: !complex double datafilespath !defined(PETSC_USE_64BIT_INDICES) 116c4762a1bSJed Brown 117c4762a1bSJed Brown TEST*/ 118