xref: /petsc/src/mat/tests/ex35.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 
2 static char help[] = "Tests MatCreateSubMatrices().\n\n";
3 
4 #include <petscmat.h>
5 
6 int main(int argc,char **args)
7 {
8   Mat            A,B,*Bsub;
9   PetscInt       i,j,m = 6,n = 6,N = 36,Ii,J;
10   PetscErrorCode ierr;
11   PetscScalar    v;
12   IS             isrow;
13 
14   ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr;
15   CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,5,NULL,&A));
16   for (i=0; i<m; i++) {
17     for (j=0; j<n; j++) {
18       v = -1.0;  Ii = j + n*i;
19       if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
20       if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
21       if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
22       if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
23       v = 4.0; CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
24     }
25   }
26   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
27   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
28   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_SELF));
29 
30   /* take the first diagonal block */
31   CHKERRQ(ISCreateStride(PETSC_COMM_WORLD,m,0,1,&isrow));
32   CHKERRQ(MatCreateSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub));
33   B    = *Bsub;
34   CHKERRQ(ISDestroy(&isrow));
35   CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_SELF));
36   CHKERRQ(MatDestroySubMatrices(1,&Bsub));
37 
38   /* take a strided block */
39   CHKERRQ(ISCreateStride(PETSC_COMM_WORLD,m,0,2,&isrow));
40   CHKERRQ(MatCreateSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub));
41   B    = *Bsub;
42   CHKERRQ(ISDestroy(&isrow));
43   CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_SELF));
44   CHKERRQ(MatDestroySubMatrices(1,&Bsub));
45 
46   /* take the last block */
47   CHKERRQ(ISCreateStride(PETSC_COMM_WORLD,m,N-m-1,1,&isrow));
48   CHKERRQ(MatCreateSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub));
49   B    = *Bsub;
50   CHKERRQ(ISDestroy(&isrow));
51   CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_SELF));
52 
53   CHKERRQ(MatDestroySubMatrices(1,&Bsub));
54   CHKERRQ(MatDestroy(&A));
55 
56   ierr = PetscFinalize();
57   return ierr;
58 }
59 
60 /*TEST
61 
62    test:
63 
64 TEST*/
65