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 PetscErrorCode ierr; 12c4762a1bSJed Brown PetscScalar v; 13c4762a1bSJed Brown IS isrow,iscol; 14c4762a1bSJed Brown 15c4762a1bSJed Brown ierr = PetscInitialize(&argc,&args,(char*)0,help);if (ierr) return ierr; 16*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 17*5f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 18c4762a1bSJed Brown n = 2*size; 19c4762a1bSJed Brown 20*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 21*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 22*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 23*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(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; 33*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES)); 34c4762a1bSJed Brown } 35c4762a1bSJed Brown } 36*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 38*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON)); 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 41c4762a1bSJed Brown 42c4762a1bSJed Brown /* 43c4762a1bSJed Brown Generate a new matrix consisting every column of the original matrix 44c4762a1bSJed Brown */ 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(C,&rstart,&rend)); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRangeColumn(C,&cstart,&cend)); 47c4762a1bSJed Brown 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_WORLD,rend-rstart > 0 ? rend-rstart-1 : 0,rstart,1,&isrow)); 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_WORLD,cend-cstart,cstart,1,&iscol)); 50*5f80ce2aSJacob Faibussowitsch CHKERRQ(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; 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES)); 56c4762a1bSJed Brown } 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 58*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 59c4762a1bSJed Brown 60*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrix(C,isrow,NULL,MAT_REUSE_MATRIX,&A)); 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON)); 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 64c4762a1bSJed Brown 65*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrow)); 66*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iscol)); 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 68*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 69c4762a1bSJed Brown ierr = PetscFinalize(); 70c4762a1bSJed Brown return ierr; 71c4762a1bSJed Brown } 72c4762a1bSJed Brown 73c4762a1bSJed Brown /*TEST 74c4762a1bSJed Brown 75c4762a1bSJed Brown test: 76c4762a1bSJed Brown nsize: 2 77c4762a1bSJed Brown 78c4762a1bSJed Brown TEST*/ 79