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