xref: /petsc/src/mat/tests/ex35.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests MatCreateSubMatrices().\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown #include <petscmat.h>
5c4762a1bSJed Brown 
6c4762a1bSJed Brown int main(int argc,char **args)
7c4762a1bSJed Brown {
8c4762a1bSJed Brown   Mat            A,B,*Bsub;
9c4762a1bSJed Brown   PetscInt       i,j,m = 6,n = 6,N = 36,Ii,J;
10c4762a1bSJed Brown   PetscScalar    v;
11c4762a1bSJed Brown   IS             isrow;
12c4762a1bSJed Brown 
13*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help));
145f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,5,NULL,&A));
15c4762a1bSJed Brown   for (i=0; i<m; i++) {
16c4762a1bSJed Brown     for (j=0; j<n; j++) {
17c4762a1bSJed Brown       v = -1.0;  Ii = j + n*i;
185f80ce2aSJacob Faibussowitsch       if (i>0)   {J = Ii - n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
195f80ce2aSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
205f80ce2aSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
215f80ce2aSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; CHKERRQ(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
225f80ce2aSJacob Faibussowitsch       v = 4.0; CHKERRQ(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
23c4762a1bSJed Brown     }
24c4762a1bSJed Brown   }
255f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_SELF));
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   /* take the first diagonal block */
305f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateStride(PETSC_COMM_WORLD,m,0,1,&isrow));
315f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub));
32c4762a1bSJed Brown   B    = *Bsub;
335f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isrow));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_SELF));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroySubMatrices(1,&Bsub));
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   /* take a strided block */
385f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateStride(PETSC_COMM_WORLD,m,0,2,&isrow));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub));
40c4762a1bSJed Brown   B    = *Bsub;
415f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isrow));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_SELF));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroySubMatrices(1,&Bsub));
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /* take the last block */
465f80ce2aSJacob Faibussowitsch   CHKERRQ(ISCreateStride(PETSC_COMM_WORLD,m,N-m-1,1,&isrow));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub));
48c4762a1bSJed Brown   B    = *Bsub;
495f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isrow));
505f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(B,PETSC_VIEWER_STDOUT_SELF));
51c4762a1bSJed Brown 
525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroySubMatrices(1,&Bsub));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
54c4762a1bSJed Brown 
55*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
56*b122ec5aSJacob Faibussowitsch   return 0;
57c4762a1bSJed Brown }
58c4762a1bSJed Brown 
59c4762a1bSJed Brown /*TEST
60c4762a1bSJed Brown 
61c4762a1bSJed Brown    test:
62c4762a1bSJed Brown 
63c4762a1bSJed Brown TEST*/
64