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