1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatCreateSubmatrix() with entire matrix, modified from ex59.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,Adup; 9c4762a1bSJed Brown PetscInt i,j,m = 3,n = 2,rstart,rend; 10c4762a1bSJed Brown PetscMPIInt size,rank; 11c4762a1bSJed Brown PetscScalar v; 12c4762a1bSJed Brown IS isrow; 13c4762a1bSJed Brown PetscBool detect_bug = PETSC_FALSE; 14c4762a1bSJed Brown 15*327415f7SBarry Smith PetscFunctionBeginUser; 169566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&args,(char*)0,help)); 179566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,NULL,"-detect_bug",&detect_bug)); 189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 199566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 20c4762a1bSJed Brown n = 2*size; 21c4762a1bSJed Brown 229566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&C)); 239566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 249566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C)); 259566063dSJacob Faibussowitsch PetscCall(MatSetUp(C)); 26c4762a1bSJed Brown 27c4762a1bSJed Brown /* 28c4762a1bSJed Brown This is JUST to generate a nice test matrix, all processors fill up 29c4762a1bSJed Brown the entire matrix. This is not something one would ever do in practice. 30c4762a1bSJed Brown */ 319566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C,&rstart,&rend)); 32c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 33c4762a1bSJed Brown for (j=0; j<m*n; j++) { 34c4762a1bSJed Brown v = i + j + 1; 359566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES)); 36c4762a1bSJed Brown } 37c4762a1bSJed Brown } 389566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 399566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 409566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON)); 419566063dSJacob Faibussowitsch PetscCall(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 429566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 43c4762a1bSJed Brown 44c4762a1bSJed Brown /* 45c4762a1bSJed Brown Generate a new matrix consisting every row and column of the original matrix 46c4762a1bSJed Brown */ 479566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C,&rstart,&rend)); 48c4762a1bSJed Brown 49c4762a1bSJed Brown /* Create parallel IS with the rows we want on THIS processor */ 50dd400576SPatrick Sanan if (detect_bug && rank == 0) { 519566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,1,rstart,1,&isrow)); 52c4762a1bSJed Brown } else { 539566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD,rend-rstart,rstart,1,&isrow)); 54c4762a1bSJed Brown } 559566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(C,isrow,NULL,MAT_INITIAL_MATRIX,&A)); 56c4762a1bSJed Brown 57c4762a1bSJed Brown /* Change C to test the case MAT_REUSE_MATRIX */ 58dd400576SPatrick Sanan if (rank == 0) { 59c4762a1bSJed Brown i = 0; j = 0; v = 100; 609566063dSJacob Faibussowitsch PetscCall(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES)); 61c4762a1bSJed Brown } 629566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 639566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 64c4762a1bSJed Brown 659566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(C,isrow,NULL,MAT_REUSE_MATRIX,&A)); 669566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON)); 679566063dSJacob Faibussowitsch PetscCall(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 689566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 69c4762a1bSJed Brown 70c4762a1bSJed Brown /* Test MatDuplicate */ 719566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A,MAT_COPY_VALUES,&Adup)); 729566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Adup)); 73c4762a1bSJed Brown 749566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow)); 759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C)); 779566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 78b122ec5aSJacob Faibussowitsch return 0; 79c4762a1bSJed Brown } 80c4762a1bSJed Brown 81c4762a1bSJed Brown /*TEST 82c4762a1bSJed Brown 83c4762a1bSJed Brown test: 84c4762a1bSJed Brown nsize: 2 85c4762a1bSJed Brown filter: grep -v "Mat Object" 86c4762a1bSJed Brown requires: !complex 87c4762a1bSJed Brown 88c4762a1bSJed Brown test: 89c4762a1bSJed Brown suffix: 2 90c4762a1bSJed Brown nsize: 3 91c4762a1bSJed Brown args: -detect_bug 92c4762a1bSJed Brown filter: grep -v "Mat Object" 93c4762a1bSJed Brown requires: !complex 94c4762a1bSJed Brown 95c4762a1bSJed Brown TEST*/ 96