1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests MatCreateSubmatrix() in parallel."; 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; 10c4762a1bSJed Brown PetscMPIInt size,rank; 11c4762a1bSJed Brown PetscScalar v; 12c4762a1bSJed Brown IS isrow,iscol; 13c4762a1bSJed Brown 14*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&args,(char*)0,help)); 155f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 165f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 17c4762a1bSJed Brown n = 2*size; 18c4762a1bSJed Brown 195f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&C)); 205f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,m*n,m*n)); 215f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(C)); 225f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(C)); 23c4762a1bSJed Brown 24c4762a1bSJed Brown /* 25c4762a1bSJed Brown This is JUST to generate a nice test matrix, all processors fill up 26c4762a1bSJed Brown the entire matrix. This is not something one would ever do in practice. 27c4762a1bSJed Brown */ 285f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(C,&rstart,&rend)); 29c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 30c4762a1bSJed Brown for (j=0; j<m*n; j++) { 31c4762a1bSJed Brown v = i + j + 1; 325f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(C,1,&i,1,&j,&v,INSERT_VALUES)); 33c4762a1bSJed Brown } 34c4762a1bSJed Brown } 35c4762a1bSJed Brown 365f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY)); 375f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY)); 385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON)); 395f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(C,PETSC_VIEWER_STDOUT_WORLD)); 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* 42c4762a1bSJed Brown Generate a new matrix consisting of every second row and column of 43c4762a1bSJed Brown the original matrix 44c4762a1bSJed Brown */ 455f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(C,&rstart,&rend)); 46c4762a1bSJed Brown /* Create parallel IS with the rows we want on THIS processor */ 475f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_WORLD,(rend-rstart)/2,rstart,2,&isrow)); 48c4762a1bSJed Brown /* Create parallel IS with the rows we want on THIS processor (same as rows for now) */ 495f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_WORLD,(rend-rstart)/2,rstart,2,&iscol)); 50c4762a1bSJed Brown 515f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrix(C,isrow,iscol,MAT_INITIAL_MATRIX,&A)); 525f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrix(C,isrow,iscol,MAT_REUSE_MATRIX,&A)); 535f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(A,PETSC_VIEWER_STDOUT_WORLD)); 54c4762a1bSJed Brown 555f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isrow)); 565f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&iscol)); 575f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&C)); 59*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 60*b122ec5aSJacob Faibussowitsch return 0; 61c4762a1bSJed Brown } 62c4762a1bSJed Brown 63c4762a1bSJed Brown /*TEST 64c4762a1bSJed Brown 65c4762a1bSJed Brown test: 66c4762a1bSJed Brown 67c4762a1bSJed Brown test: 68c4762a1bSJed Brown suffix: 2 69c4762a1bSJed Brown nsize: 3 70c4762a1bSJed Brown 71c4762a1bSJed Brown test: 72c4762a1bSJed Brown suffix: 2_baij 73c4762a1bSJed Brown nsize: 3 74c4762a1bSJed Brown args: -mat_type baij 75c4762a1bSJed Brown 76c4762a1bSJed Brown test: 77c4762a1bSJed Brown suffix: 2_sbaij 78c4762a1bSJed Brown nsize: 3 79c4762a1bSJed Brown args: -mat_type sbaij 80c4762a1bSJed Brown 81c4762a1bSJed Brown test: 82c4762a1bSJed Brown suffix: baij 83c4762a1bSJed Brown args: -mat_type baij 84c4762a1bSJed Brown output_file: output/ex59_1_baij.out 85c4762a1bSJed Brown 86c4762a1bSJed Brown test: 87c4762a1bSJed Brown suffix: sbaij 88c4762a1bSJed Brown args: -mat_type sbaij 89c4762a1bSJed Brown output_file: output/ex59_1_sbaij.out 90c4762a1bSJed Brown 91c4762a1bSJed Brown TEST*/ 92