static char help[] = "Tests MatGetRowIJ for SeqAIJ, SeqBAIJ and SeqSBAIJ\n\n"; #include PetscErrorCode DumpCSR(Mat A,PetscInt shift,PetscBool symmetric,PetscBool compressed) { MatType type; PetscInt i,j,nr,bs = 1; const PetscInt *ia,*ja; PetscBool done,isseqbaij,isseqsbaij; PetscErrorCode ierr; PetscFunctionBeginUser; ierr = PetscObjectBaseTypeCompare((PetscObject)A,MATSEQBAIJ,&isseqbaij);CHKERRQ(ierr); ierr = PetscObjectBaseTypeCompare((PetscObject)A,MATSEQSBAIJ,&isseqsbaij);CHKERRQ(ierr); if (isseqbaij || isseqsbaij) { ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); } ierr = MatGetType(A,&type);CHKERRQ(ierr); ierr = MatGetRowIJ(A,shift,symmetric,compressed,&nr,&ia,&ja,&done);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_SELF,"===========================================================\n");CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_SELF,"CSR for %s: shift %D symmetric %D compressed %D\n",type,shift,(PetscInt)symmetric,(PetscInt)compressed);CHKERRQ(ierr); for (i=0;i=0) { ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr); ierr = MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr); ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr); } J = Ii + n; if (J=0) { ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr); ierr = MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr); ierr = MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr); } J = Ii + 1; if (J0) ? PETSC_FALSE : PETSC_TRUE); for (k=0;k<2;k++) { PetscBool compressed = ((k>0) ? PETSC_FALSE : PETSC_TRUE); ierr = DumpCSR(A,shift,symmetric,compressed);CHKERRQ(ierr); ierr = DumpCSR(B,shift,symmetric,compressed);CHKERRQ(ierr); ierr = DumpCSR(C,shift,symmetric,compressed);CHKERRQ(ierr); } } } ierr = MatDestroy(&A);CHKERRQ(ierr); ierr = MatDestroy(&B);CHKERRQ(ierr); ierr = MatDestroy(&C);CHKERRQ(ierr); ierr = PetscFinalize(); return ierr; } /*TEST test: test: suffix: 2 args: -bs 2 TEST*/