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; PetscFunctionBeginUser; PetscCall(PetscObjectBaseTypeCompare((PetscObject)A,MATSEQBAIJ,&isseqbaij)); PetscCall(PetscObjectBaseTypeCompare((PetscObject)A,MATSEQSBAIJ,&isseqsbaij)); if (isseqbaij || isseqsbaij) { PetscCall(MatGetBlockSize(A,&bs)); } PetscCall(MatGetType(A,&type)); PetscCall(MatGetRowIJ(A,shift,symmetric,compressed,&nr,&ia,&ja,&done)); PetscCall(PetscPrintf(PETSC_COMM_SELF,"===========================================================\n")); PetscCall(PetscPrintf(PETSC_COMM_SELF,"CSR for %s: shift %" PetscInt_FMT " symmetric %" PetscInt_FMT " compressed %" PetscInt_FMT "\n",type,shift,(PetscInt)symmetric,(PetscInt)compressed)); for (i=0;i=0) { PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES)); PetscCall(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES)); PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); } J = Ii + n; if (J=0) { PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES)); PetscCall(MatSetValues(B,1,&Ii,1,&J,&v,INSERT_VALUES)); PetscCall(MatSetValues(C,1,&Ii,1,&J,&v,INSERT_VALUES)); } J = Ii + 1; if (J0) ? PETSC_FALSE : PETSC_TRUE); for (k=0;k<2;k++) { PetscBool compressed = ((k>0) ? PETSC_FALSE : PETSC_TRUE); PetscCall(DumpCSR(A,shift,symmetric,compressed)); PetscCall(DumpCSR(B,shift,symmetric,compressed)); PetscCall(DumpCSR(C,shift,symmetric,compressed)); } } } PetscCall(MatDestroy(&A)); PetscCall(MatDestroy(&B)); PetscCall(MatDestroy(&C)); PetscCall(PetscFinalize()); return 0; } /*TEST test: test: suffix: 2 args: -bs 2 TEST*/