static char help[] = "Tests converting a SBAIJ matrix to BAIJ format with MatConvert. Modified from ex55.c\n\n"; #include /* Usage: ./ex141 -mat_view */ int main(int argc,char **args) { Mat C,B; PetscInt i,bs=2,mbs,m,block,d_nz=6,col[3]; PetscMPIInt size; char file[PETSC_MAX_PATH_LEN]; PetscViewer fd; PetscBool equal,loadmat; PetscScalar value[4]; PetscInt ridx[2],cidx[2]; PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); PetscCall(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&loadmat)); PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); /* input matrix C */ if (loadmat) { /* Open binary file. Load a sbaij matrix, then destroy the viewer. */ PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&fd)); PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); PetscCall(MatSetType(C,MATSEQSBAIJ)); PetscCall(MatLoad(C,fd)); PetscCall(PetscViewerDestroy(&fd)); } else { /* Create a sbaij mat with bs>1 */ mbs =8; PetscCall(PetscOptionsGetInt(NULL,NULL,"-mbs",&mbs,NULL)); m = mbs*bs; PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m,m)); PetscCall(MatSetType(C,MATSBAIJ)); PetscCall(MatSetFromOptions(C)); PetscCall(MatSeqSBAIJSetPreallocation(C,bs,d_nz,NULL)); PetscCall(MatSetUp(C)); PetscCall(MatSetOption(C,MAT_IGNORE_LOWER_TRIANGULAR,PETSC_TRUE)); for (block=0; block