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