1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatCreateSubmatrix() with certain entire rows of matrix, modified from ex181.c."; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **args) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat C,A; 9c4762a1bSJed Brown PetscInt i,j,m = 3,n = 2,rstart,rend,cstart,cend; 10c4762a1bSJed Brown PetscMPIInt size,rank; 11c4762a1bSJed Brown PetscScalar v; 12c4762a1bSJed Brown IS isrow,iscol; 13c4762a1bSJed Brown 14*327415f7SBarry Smith PetscFunctionBeginUser; 159566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 18c4762a1bSJed Brown n = 2*size; 19c4762a1bSJed Brown 209566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 219566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 229566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 239566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 24c4762a1bSJed Brown 25c4762a1bSJed Brown /* 26c4762a1bSJed Brown This is JUST to generate a nice test matrix, all processors fill up 27c4762a1bSJed Brown the entire matrix. This is not something one would ever do in practice. 28c4762a1bSJed Brown */ 299566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C,&rstart,&rend)); 30c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 31c4762a1bSJed Brown for (j=0; j<m*n; j++) { 32c4762a1bSJed Brown v = i + j + 1; 339566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES)); 34c4762a1bSJed Brown } 35c4762a1bSJed Brown } 369566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 379566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 389566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON)); 399566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 409566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 41c4762a1bSJed Brown 42c4762a1bSJed Brown /* 43c4762a1bSJed Brown Generate a new matrix consisting every column of the original matrix 44c4762a1bSJed Brown */ 459566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C,&rstart,&rend)); 469566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(C,&cstart,&cend)); 47c4762a1bSJed Brown 489566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,rend-rstart > 0 ? rend-rstart-1 : 0,rstart,1,&isrow)); 499566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,cend-cstart,cstart,1,&iscol)); 509566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(C,isrow,NULL,MAT_INITIAL_MATRIX,&A)); 51c4762a1bSJed Brown 52c4762a1bSJed Brown /* Change C to test the case MAT_REUSE_MATRIX */ 53dd400576SPatrick Sanan if (rank == 0) { 54c4762a1bSJed Brown i = 0; j = 0; v = 100; 559566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES)); 56c4762a1bSJed Brown } 579566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 589566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 59c4762a1bSJed Brown 609566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(C,isrow,NULL,MAT_REUSE_MATRIX,&A)); 619566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON)); 629566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 639566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 64c4762a1bSJed Brown 659566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 669566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol)); 679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 699566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 70b122ec5aSJacob Faibussowitsch return 0; 71c4762a1bSJed Brown } 72c4762a1bSJed Brown 73c4762a1bSJed Brown /*TEST 74c4762a1bSJed Brown 75c4762a1bSJed Brown test: 76c4762a1bSJed Brown nsize: 2 77c4762a1bSJed Brown 78c4762a1bSJed Brown TEST*/ 79