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*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 23*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),NULL)); 24c4762a1bSJed Brown 25c4762a1bSJed Brown /* Load the matrix as AIJ format */ 26*9566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&va)); 27*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 28*9566063dSJacob Faibussowitsch PetscCall(MatSetType(A,MATSEQAIJ)); 29*9566063dSJacob Faibussowitsch PetscCall(MatLoad(A,va)); 30*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&va)); 31c4762a1bSJed Brown 32c4762a1bSJed Brown /* Load the matrix as BAIJ format */ 33*9566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&vb)); 34*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&B)); 35*9566063dSJacob Faibussowitsch PetscCall(MatSetType(B,MATSEQBAIJ)); 36*9566063dSJacob Faibussowitsch PetscCall(MatLoad(B,vb)); 37*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vb)); 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* Load the matrix as BAIJ format */ 40*9566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&vc)); 41*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 42*9566063dSJacob Faibussowitsch PetscCall(MatSetType(C,MATSEQBAIJ)); 43*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 44*9566063dSJacob Faibussowitsch PetscCall(MatLoad(C,vc)); 45*9566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&vc)); 46c4762a1bSJed Brown 47*9566063dSJacob Faibussowitsch PetscCall(MatGetSize(A,&m,&n)); 48*9566063dSJacob Faibussowitsch PetscCall(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() */ 52*9566063dSJacob Faibussowitsch PetscCall(MatEqual(B,C,&tflg)); 5328b400f6SJacob Faibussowitsch PetscCheck(tflg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatEqual() failed"); 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* Test MatGetDiagonal() */ 56*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,m,&x)); 57*9566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,m,&y)); 58c4762a1bSJed Brown 59*9566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(A,x)); 60*9566063dSJacob Faibussowitsch PetscCall(MatGetDiagonal(B,y)); 61c4762a1bSJed Brown 62*9566063dSJacob Faibussowitsch PetscCall(VecEqual(x,y,&tflg)); 6328b400f6SJacob Faibussowitsch PetscCheck(tflg,PETSC_COMM_SELF,PETSC_ERR_PLIB,"MatGetDiagonal() failed"); 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* Test MatDiagonalScale() */ 66*9566063dSJacob Faibussowitsch PetscCall(PetscRandomCreate(PETSC_COMM_SELF,&r)); 67*9566063dSJacob Faibussowitsch PetscCall(PetscRandomSetFromOptions(r)); 68*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(x,r)); 69*9566063dSJacob Faibussowitsch PetscCall(VecSetRandom(y,r)); 70c4762a1bSJed Brown 71*9566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(A,x,y)); 72*9566063dSJacob Faibussowitsch PetscCall(MatDiagonalScale(B,x,y)); 73*9566063dSJacob Faibussowitsch PetscCall(MatMult(A,x,y)); 74*9566063dSJacob Faibussowitsch PetscCall(VecNorm(y,NORM_2,&norm1)); 75*9566063dSJacob Faibussowitsch PetscCall(MatMult(B,x,y)); 76*9566063dSJacob Faibussowitsch PetscCall(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++) { 82*9566063dSJacob Faibussowitsch PetscCall(PetscRandomGetValue(r,&rval)); 83c4762a1bSJed Brown row = (int)(rval*m); 84*9566063dSJacob Faibussowitsch PetscCall(MatGetRow(A,row,&ncols1,&cols1,&vals1)); 85*9566063dSJacob Faibussowitsch PetscCall(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 93*9566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(A,row,&ncols1,&cols1,&vals1)); 94*9566063dSJacob Faibussowitsch PetscCall(MatRestoreRow(B,row,&ncols2,&cols2,&vals2)); 95c4762a1bSJed Brown } 96c4762a1bSJed Brown 97*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 98*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&B)); 99*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 100*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 101*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 102*9566063dSJacob Faibussowitsch PetscCall(PetscRandomDestroy(&r)); 103*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 104b122ec5aSJacob 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