1*c4762a1bSJed Brown 2*c4762a1bSJed Brown static char help[] = "Tests MatCreateSubMatrices().\n\n"; 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown #include <petscmat.h> 5*c4762a1bSJed Brown 6*c4762a1bSJed Brown int main(int argc,char **args) 7*c4762a1bSJed Brown { 8*c4762a1bSJed Brown Mat A,B,*Bsub; 9*c4762a1bSJed Brown PetscInt i,j,m = 6,n = 6,N = 36,Ii,J; 10*c4762a1bSJed Brown PetscErrorCode ierr; 11*c4762a1bSJed Brown PetscScalar v; 12*c4762a1bSJed Brown IS isrow; 13*c4762a1bSJed Brown 14*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 15*c4762a1bSJed Brown ierr = MatCreateSeqAIJ(PETSC_COMM_WORLD,N,N,5,NULL,&A);CHKERRQ(ierr); 16*c4762a1bSJed Brown for (i=0; i<m; i++) { 17*c4762a1bSJed Brown for (j=0; j<n; j++) { 18*c4762a1bSJed Brown v = -1.0; Ii = j + n*i; 19*c4762a1bSJed Brown if (i>0) {J = Ii - n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 20*c4762a1bSJed Brown if (i<m-1) {J = Ii + n; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 21*c4762a1bSJed Brown if (j>0) {J = Ii - 1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 22*c4762a1bSJed Brown if (j<n-1) {J = Ii + 1; ierr = MatSetValues(A,1,&Ii,1,&J,&v,INSERT_VALUES);CHKERRQ(ierr);} 23*c4762a1bSJed Brown v = 4.0; ierr = MatSetValues(A,1,&Ii,1,&Ii,&v,INSERT_VALUES);CHKERRQ(ierr); 24*c4762a1bSJed Brown } 25*c4762a1bSJed Brown } 26*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 27*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 28*c4762a1bSJed Brown ierr = MatView(A,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 29*c4762a1bSJed Brown 30*c4762a1bSJed Brown /* take the first diagonal block */ 31*c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,m,0,1,&isrow);CHKERRQ(ierr); 32*c4762a1bSJed Brown ierr = MatCreateSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub);CHKERRQ(ierr); 33*c4762a1bSJed Brown B = *Bsub; 34*c4762a1bSJed Brown ierr = ISDestroy(&isrow);CHKERRQ(ierr); 35*c4762a1bSJed Brown ierr = MatView(B,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = MatDestroySubMatrices(1,&Bsub);CHKERRQ(ierr); 37*c4762a1bSJed Brown 38*c4762a1bSJed Brown /* take a strided block */ 39*c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,m,0,2,&isrow);CHKERRQ(ierr); 40*c4762a1bSJed Brown ierr = MatCreateSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub);CHKERRQ(ierr); 41*c4762a1bSJed Brown B = *Bsub; 42*c4762a1bSJed Brown ierr = ISDestroy(&isrow);CHKERRQ(ierr); 43*c4762a1bSJed Brown ierr = MatView(B,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 44*c4762a1bSJed Brown ierr = MatDestroySubMatrices(1,&Bsub);CHKERRQ(ierr); 45*c4762a1bSJed Brown 46*c4762a1bSJed Brown /* take the last block */ 47*c4762a1bSJed Brown ierr = ISCreateStride(PETSC_COMM_WORLD,m,N-m-1,1,&isrow);CHKERRQ(ierr); 48*c4762a1bSJed Brown ierr = MatCreateSubMatrices(A,1,&isrow,&isrow,MAT_INITIAL_MATRIX,&Bsub);CHKERRQ(ierr); 49*c4762a1bSJed Brown B = *Bsub; 50*c4762a1bSJed Brown ierr = ISDestroy(&isrow);CHKERRQ(ierr); 51*c4762a1bSJed Brown ierr = MatView(B,PETSC_VIEWER_STDOUT_SELF);CHKERRQ(ierr); 52*c4762a1bSJed Brown 53*c4762a1bSJed Brown ierr = MatDestroySubMatrices(1,&Bsub);CHKERRQ(ierr); 54*c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 55*c4762a1bSJed Brown 56*c4762a1bSJed Brown ierr = PetscFinalize(); 57*c4762a1bSJed Brown return ierr; 58*c4762a1bSJed Brown } 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown 61*c4762a1bSJed Brown 62*c4762a1bSJed Brown /*TEST 63*c4762a1bSJed Brown 64*c4762a1bSJed Brown test: 65*c4762a1bSJed Brown 66*c4762a1bSJed Brown TEST*/ 67