xref: /petsc/src/mat/tests/ex35.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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*327415f7SBarry Smith   PetscFunctionBeginUser;
149566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&args,(char*)0,help));
159566063dSJacob Faibussowitsch   PetscCall(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;
199566063dSJacob Faibussowitsch       if (i>0)   {J = Ii - n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
209566063dSJacob Faibussowitsch       if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
219566063dSJacob Faibussowitsch       if (j>0)   {J = Ii - 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
229566063dSJacob Faibussowitsch       if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));}
239566063dSJacob Faibussowitsch       v = 4.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES));
24c4762a1bSJed Brown     }
25c4762a1bSJed Brown   }
269566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
279566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
289566063dSJacob Faibussowitsch   PetscCall(MatView(A,PETSC_VIEWER_STDOUT_SELF));
29c4762a1bSJed Brown 
30c4762a1bSJed Brown   /* take the first diagonal block */
319566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_WORLD,m,0,1,&isrow));
329566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub));
33c4762a1bSJed Brown   B    = *Bsub;
349566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrow));
359566063dSJacob Faibussowitsch   PetscCall(MatView(B,PETSC_VIEWER_STDOUT_SELF));
369566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(1,&Bsub));
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   /* take a strided block */
399566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_WORLD,m,0,2,&isrow));
409566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub));
41c4762a1bSJed Brown   B    = *Bsub;
429566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrow));
439566063dSJacob Faibussowitsch   PetscCall(MatView(B,PETSC_VIEWER_STDOUT_SELF));
449566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(1,&Bsub));
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   /* take the last block */
479566063dSJacob Faibussowitsch   PetscCall(ISCreateStride(PETSC_COMM_WORLD,m,N-m-1,1,&isrow));
489566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub));
49c4762a1bSJed Brown   B    = *Bsub;
509566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isrow));
519566063dSJacob Faibussowitsch   PetscCall(MatView(B,PETSC_VIEWER_STDOUT_SELF));
52c4762a1bSJed Brown 
539566063dSJacob Faibussowitsch   PetscCall(MatDestroySubMatrices(1,&Bsub));
549566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
55c4762a1bSJed Brown 
569566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
57b122ec5aSJacob Faibussowitsch   return 0;
58c4762a1bSJed Brown }
59c4762a1bSJed Brown 
60c4762a1bSJed Brown /*TEST
61c4762a1bSJed Brown 
62c4762a1bSJed Brown    test:
63c4762a1bSJed Brown 
64c4762a1bSJed Brown TEST*/
65