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*327415f7SBarry Smith PetscFunctionBeginUser; 239566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 249566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),NULL)); 25c4762a1bSJed Brown 26c4762a1bSJed Brown /* Load the matrix as AIJ format */ 279566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&va)); 289566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 299566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATSEQAIJ)); 309566063dSJacob Faibussowitsch PetscCall(MatLoad(A,va)); 319566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&va)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* Load the matrix as BAIJ format */ 349566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&vb)); 359566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); 369566063dSJacob Faibussowitsch PetscCall(MatSetType(B,MATSEQBAIJ)); 379566063dSJacob Faibussowitsch PetscCall(MatLoad(B,vb)); 389566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vb)); 39c4762a1bSJed Brown 40c4762a1bSJed Brown /* Load the matrix as BAIJ format */ 419566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&vc)); 429566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 439566063dSJacob Faibussowitsch PetscCall(MatSetType(C,MATSEQBAIJ)); 449566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 459566063dSJacob Faibussowitsch PetscCall(MatLoad(C,vc)); 469566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vc)); 47c4762a1bSJed Brown 489566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&m,&n)); 499566063dSJacob Faibussowitsch PetscCall(MatGetSize(B,&m2,&n2)); 5008401ef6SPierre Jolivet PetscCheck(m==m2,PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrices are of different size. Cannot run this example"); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Test MatEqual() */ 539566063dSJacob Faibussowitsch PetscCall(MatEqual(B,C,&tflg)); 5428b400f6SJacob Faibussowitsch PetscCheck(tflg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatEqual() failed"); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* Test MatGetDiagonal() */ 579566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,m,&x)); 589566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,m,&y)); 59c4762a1bSJed Brown 609566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A,x)); 619566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(B,y)); 62c4762a1bSJed Brown 639566063dSJacob Faibussowitsch PetscCall(VecEqual(x,y,&tflg)); 6428b400f6SJacob Faibussowitsch PetscCheck(tflg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetDiagonal() failed"); 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* Test MatDiagonalScale() */ 679566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&r)); 689566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(r)); 699566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,r)); 709566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y,r)); 71c4762a1bSJed Brown 729566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A,x,y)); 739566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(B,x,y)); 749566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,y)); 759566063dSJacob Faibussowitsch PetscCall(VecNorm(y,NORM_2,&norm1)); 769566063dSJacob Faibussowitsch PetscCall(MatMult(B,x,y)); 779566063dSJacob Faibussowitsch PetscCall(VecNorm(y,NORM_2,&norm2)); 78c4762a1bSJed Brown rnorm = ((norm1-norm2)*100)/norm1; 79bc30f867SBarry Smith PetscCheck(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++) { 839566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r,&rval)); 84c4762a1bSJed Brown row = (int)(rval*m); 859566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,row,&ncols1,&cols1,&vals1)); 869566063dSJacob Faibussowitsch PetscCall(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++; 9008401ef6SPierre Jolivet PetscCheck(vals1[i] == vals2[j],PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetRow() failed - vals incorrect."); 91c4762a1bSJed Brown } 9208401ef6SPierre Jolivet PetscCheck(i>=ncols1,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetRow() failed - cols incorrect"); 93c4762a1bSJed Brown 949566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,row,&ncols1,&cols1,&vals1)); 959566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(B,row,&ncols2,&cols2,&vals2)); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 1009566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 1019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 1029566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 1039566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&r)); 1049566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 105b122ec5aSJacob Faibussowitsch return 0; 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