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*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 165f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,NULL,"-detect_bug",&detect_bug)); 175f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 185f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 19c4762a1bSJed Brown n = 2*size; 20c4762a1bSJed Brown 215f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 235f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 245f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(C)); 25c4762a1bSJed Brown 26c4762a1bSJed Brown /* 27c4762a1bSJed Brown This is JUST to generate a nice test matrix, all processors fill up 28c4762a1bSJed Brown the entire matrix. This is not something one would ever do in practice. 29c4762a1bSJed Brown */ 305f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(C,&rstart,&rend)); 31c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 32c4762a1bSJed Brown for (j=0; j<m*n; j++) { 33c4762a1bSJed Brown v = i + j + 1; 345f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES)); 35c4762a1bSJed Brown } 36c4762a1bSJed Brown } 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON)); 405f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* 44c4762a1bSJed Brown Generate a new matrix consisting every row and column of the original matrix 45c4762a1bSJed Brown */ 465f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(C,&rstart,&rend)); 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* Create parallel IS with the rows we want on THIS processor */ 49dd400576SPatrick Sanan if (detect_bug && rank == 0) { 505f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_WORLD,1,rstart,1,&isrow)); 51c4762a1bSJed Brown } else { 525f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_WORLD,rend-rstart,rstart,1,&isrow)); 53c4762a1bSJed Brown } 545f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrix(C,isrow,NULL,MAT_INITIAL_MATRIX,&A)); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* Change C to test the case MAT_REUSE_MATRIX */ 57dd400576SPatrick Sanan if (rank == 0) { 58c4762a1bSJed Brown i = 0; j = 0; v = 100; 595f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES)); 60c4762a1bSJed Brown } 615f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 625f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 63c4762a1bSJed Brown 645f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrix(C,isrow,NULL,MAT_REUSE_MATRIX,&A)); 655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON)); 665f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* Test MatDuplicate */ 705f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(A,MAT_COPY_VALUES,&Adup)); 715f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Adup)); 72c4762a1bSJed Brown 735f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrow)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 76*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 77*b122ec5aSJacob Faibussowitsch return 0; 78c4762a1bSJed Brown } 79c4762a1bSJed Brown 80c4762a1bSJed Brown /*TEST 81c4762a1bSJed Brown 82c4762a1bSJed Brown test: 83c4762a1bSJed Brown nsize: 2 84c4762a1bSJed Brown filter: grep -v "Mat Object" 85c4762a1bSJed Brown requires: !complex 86c4762a1bSJed Brown 87c4762a1bSJed Brown test: 88c4762a1bSJed Brown suffix: 2 89c4762a1bSJed Brown nsize: 3 90c4762a1bSJed Brown args: -detect_bug 91c4762a1bSJed Brown filter: grep -v "Mat Object" 92c4762a1bSJed Brown requires: !complex 93c4762a1bSJed Brown 94c4762a1bSJed Brown TEST*/ 95