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