1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests MatCreateSubMatrices() for SBAIJ matrices\n\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown #include <petscmat.h> 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown int main(int argc,char **args) 7*c4762a1bSJed Brown { 8*c4762a1bSJed Brown Mat BAIJ,SBAIJ,*subBAIJ,*subSBAIJ; 9*c4762a1bSJed Brown PetscViewer viewer; 10*c4762a1bSJed Brown char file[PETSC_MAX_PATH_LEN]; 11*c4762a1bSJed Brown PetscBool flg; 12*c4762a1bSJed Brown PetscErrorCode ierr; 13*c4762a1bSJed Brown PetscInt n = 2,issize,M,N; 14*c4762a1bSJed Brown PetscMPIInt rank; 15*c4762a1bSJed Brown IS isrow,iscol,irow[n],icol[n]; 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 18*c4762a1bSJed Brown ierr = PetscOptionsGetString(NULL,NULL,"-f",file,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 19*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&viewer);CHKERRQ(ierr); 20*c4762a1bSJed Brown 21*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&BAIJ);CHKERRQ(ierr); 22*c4762a1bSJed Brown ierr = MatSetType(BAIJ,MATMPIBAIJ);CHKERRQ(ierr); 23*c4762a1bSJed Brown ierr = MatLoad(BAIJ,viewer);CHKERRQ(ierr); 24*c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 25*c4762a1bSJed Brown 26*c4762a1bSJed Brown ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,file,FILE_MODE_READ,&viewer);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&SBAIJ);CHKERRQ(ierr); 28*c4762a1bSJed Brown ierr = MatSetType(SBAIJ,MATMPISBAIJ);CHKERRQ(ierr); 29*c4762a1bSJed Brown ierr = MatLoad(SBAIJ,viewer);CHKERRQ(ierr); 30*c4762a1bSJed Brown ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown ierr = MatGetSize(BAIJ,&M,&N);CHKERRQ(ierr); 33*c4762a1bSJed Brown issize = M/4; 34*c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,issize,0,1,&isrow);CHKERRQ(ierr); 35*c4762a1bSJed Brown irow[0] = irow[1] = isrow; 36*c4762a1bSJed Brown issize = N/2; 37*c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,issize,0,1,&iscol);CHKERRQ(ierr); 38*c4762a1bSJed Brown icol[0] = icol[1] = iscol; 39*c4762a1bSJed Brown ierr = MatCreateSubMatrices(BAIJ,n,irow,icol,MAT_INITIAL_MATRIX,&subBAIJ);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = MatCreateSubMatrices(BAIJ,n,irow,icol,MAT_REUSE_MATRIX,&subBAIJ);CHKERRQ(ierr); 41*c4762a1bSJed Brown 42*c4762a1bSJed Brown /* irow and icol must be same for SBAIJ matrices! */ 43*c4762a1bSJed Brown icol[0] = icol[1] = isrow; 44*c4762a1bSJed Brown ierr = MatCreateSubMatrices(SBAIJ,n,irow,icol,MAT_INITIAL_MATRIX,&subSBAIJ);CHKERRQ(ierr); 45*c4762a1bSJed Brown ierr = MatCreateSubMatrices(SBAIJ,n,irow,icol,MAT_REUSE_MATRIX,&subSBAIJ);CHKERRQ(ierr); 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr); 48*c4762a1bSJed Brown if (!rank) { 49*c4762a1bSJed Brown ierr = MatView(subBAIJ[0],PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 50*c4762a1bSJed Brown ierr = MatView(subSBAIJ[0],PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 51*c4762a1bSJed Brown } 52*c4762a1bSJed Brown 53*c4762a1bSJed Brown /* Free data structures */ 54*c4762a1bSJed Brown ierr = ISDestroy(&isrow);CHKERRQ(ierr); 55*c4762a1bSJed Brown ierr = ISDestroy(&iscol);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = MatDestroySubMatrices(n,&subBAIJ);CHKERRQ(ierr); 57*c4762a1bSJed Brown ierr = MatDestroySubMatrices(n,&subSBAIJ);CHKERRQ(ierr); 58*c4762a1bSJed Brown ierr = MatDestroy(&BAIJ);CHKERRQ(ierr); 59*c4762a1bSJed Brown ierr = MatDestroy(&SBAIJ);CHKERRQ(ierr); 60*c4762a1bSJed Brown 61*c4762a1bSJed Brown ierr = PetscFinalize(); 62*c4762a1bSJed Brown return ierr; 63*c4762a1bSJed Brown } 64*c4762a1bSJed Brown 65*c4762a1bSJed Brown 66