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*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 14*9566063dSJacob Faibussowitsch PetscCall(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; 18*9566063dSJacob Faibussowitsch if (i>0) {J = Ii - n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 19*9566063dSJacob Faibussowitsch if (i<m-1) {J = Ii + n; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 20*9566063dSJacob Faibussowitsch if (j>0) {J = Ii - 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 21*9566063dSJacob Faibussowitsch if (j<n-1) {J = Ii + 1; PetscCall(MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES));} 22*9566063dSJacob Faibussowitsch v = 4.0; PetscCall(MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES)); 23c4762a1bSJed Brown } 24c4762a1bSJed Brown } 25*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 26*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 27*9566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_SELF)); 28c4762a1bSJed Brown 29c4762a1bSJed Brown /* take the first diagonal block */ 30*9566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,m,0,1,&isrow)); 31*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub)); 32c4762a1bSJed Brown B = *Bsub; 33*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 34*9566063dSJacob Faibussowitsch PetscCall(MatView(B,PETSC_VIEWER_STDOUT_SELF)); 35*9566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(1,&Bsub)); 36c4762a1bSJed Brown 37c4762a1bSJed Brown /* take a strided block */ 38*9566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,m,0,2,&isrow)); 39*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub)); 40c4762a1bSJed Brown B = *Bsub; 41*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 42*9566063dSJacob Faibussowitsch PetscCall(MatView(B,PETSC_VIEWER_STDOUT_SELF)); 43*9566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(1,&Bsub)); 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* take the last block */ 46*9566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,m,N-m-1,1,&isrow)); 47*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub)); 48c4762a1bSJed Brown B = *Bsub; 49*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 50*9566063dSJacob Faibussowitsch PetscCall(MatView(B,PETSC_VIEWER_STDOUT_SELF)); 51c4762a1bSJed Brown 52*9566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(1,&Bsub)); 53*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 54c4762a1bSJed Brown 55*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 56b122ec5aSJacob Faibussowitsch return 0; 57c4762a1bSJed Brown } 58c4762a1bSJed Brown 59c4762a1bSJed Brown /*TEST 60c4762a1bSJed Brown 61c4762a1bSJed Brown test: 62c4762a1bSJed Brown 63c4762a1bSJed Brown TEST*/ 64