1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatCreateSubMatrices() for SBAIJ matrices\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat BAIJ,SBAIJ,*subBAIJ,*subSBAIJ; 9c4762a1bSJed Brown PetscViewer viewer; 10c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; 11c4762a1bSJed Brown PetscBool flg; 12c4762a1bSJed Brown PetscInt n = 2,issize,M,N; 13c4762a1bSJed Brown PetscMPIInt rank; 14c4762a1bSJed Brown IS isrow,iscol,irow[n],icol[n]; 15c4762a1bSJed Brown 16*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 175f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,NULL,"-f",file,sizeof(file),&flg)); 185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&viewer)); 19c4762a1bSJed Brown 205f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&BAIJ)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(BAIJ,MATMPIBAIJ)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(BAIJ,viewer)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 24c4762a1bSJed Brown 255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&viewer)); 265f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&SBAIJ)); 275f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(SBAIJ,MATMPISBAIJ)); 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatLoad(SBAIJ,viewer)); 295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&viewer)); 30c4762a1bSJed Brown 315f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(BAIJ,&M,&N)); 32c4762a1bSJed Brown issize = M/4; 335f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_SELF,issize,0,1,&isrow)); 34c4762a1bSJed Brown irow[0] = irow[1] = isrow; 35c4762a1bSJed Brown issize = N/2; 365f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_SELF,issize,0,1,&iscol)); 37c4762a1bSJed Brown icol[0] = icol[1] = iscol; 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrices(BAIJ,n,irow,icol,MAT_INITIAL_MATRIX,&subBAIJ)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrices(BAIJ,n,irow,icol,MAT_REUSE_MATRIX,&subBAIJ)); 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* irow and icol must be same for SBAIJ matrices! */ 42c4762a1bSJed Brown icol[0] = icol[1] = isrow; 435f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrices(SBAIJ,n,irow,icol,MAT_INITIAL_MATRIX,&subSBAIJ)); 445f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrices(SBAIJ,n,irow,icol,MAT_REUSE_MATRIX,&subSBAIJ)); 45c4762a1bSJed Brown 465f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 47dd400576SPatrick Sanan if (rank == 0) { 485f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(subBAIJ[0],PETSC_VIEWER_STDOUT_SELF)); 495f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(subSBAIJ[0],PETSC_VIEWER_STDOUT_SELF)); 50c4762a1bSJed Brown } 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Free data structures */ 535f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrow)); 545f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iscol)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroySubMatrices(n,&subBAIJ)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroySubMatrices(n,&subSBAIJ)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&BAIJ)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&SBAIJ)); 59c4762a1bSJed Brown 60*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 61*b122ec5aSJacob Faibussowitsch return 0; 62c4762a1bSJed Brown } 63