1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Creates a matrix, inserts some values, and tests MatCreateSubMatrices() and MatZeroEntries().\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscmat.h> 5c4762a1bSJed Brown 6c4762a1bSJed Brown int main(int argc,char **argv) 7c4762a1bSJed Brown { 8c4762a1bSJed Brown Mat mat,submat,submat1,*submatrices; 9c4762a1bSJed Brown PetscInt m = 10,n = 10,i = 4,tmp,rstart,rend; 10c4762a1bSJed Brown IS irow,icol; 11c4762a1bSJed Brown PetscScalar value = 1.0; 12c4762a1bSJed Brown PetscViewer sviewer; 13c4762a1bSJed Brown PetscBool allA = PETSC_FALSE; 14c4762a1bSJed Brown 15*327415f7SBarry Smith PetscFunctionBeginUser; 169566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 179566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON)); 189566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_COMMON)); 19c4762a1bSJed Brown 209566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&mat)); 219566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat,PETSC_DECIDE,PETSC_DECIDE,m,n)); 229566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(mat)); 239566063dSJacob Faibussowitsch PetscCall(MatSetUp(mat)); 249566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(mat,&rstart,&rend)); 25c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 26c4762a1bSJed Brown value = (PetscReal)i+1; tmp = i % 5; 279566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat,1,&tmp,1,&i,&value,INSERT_VALUES)); 28c4762a1bSJed Brown } 299566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY)); 309566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY)); 319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"Original matrix\n")); 329566063dSJacob Faibussowitsch PetscCall(MatView(mat,PETSC_VIEWER_STDOUT_WORLD)); 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* Test MatCreateSubMatrix_XXX_All(), i.e., submatrix = A */ 359566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_all",&allA,NULL)); 36c4762a1bSJed Brown if (allA) { 379566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,m,0,1,&irow)); 389566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,n,0,1,&icol)); 399566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(mat,1,&irow,&icol,MAT_INITIAL_MATRIX,&submatrices)); 409566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(mat,1,&irow,&icol,MAT_REUSE_MATRIX,&submatrices)); 41c4762a1bSJed Brown submat = *submatrices; 42c4762a1bSJed Brown 43c4762a1bSJed Brown /* sviewer will cause the submatrices (one per processor) to be printed in the correct order */ 449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"\nSubmatrices with all\n")); 459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"--------------------\n")); 469566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 479566063dSJacob Faibussowitsch PetscCall(MatView(submat,sviewer)); 489566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 499566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 50c4762a1bSJed Brown 519566063dSJacob Faibussowitsch PetscCall(ISDestroy(&irow)); 529566063dSJacob Faibussowitsch PetscCall(ISDestroy(&icol)); 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* test getting a reference on a submat */ 559566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)submat)); 569566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(1,&submatrices)); 579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&submat)); 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Form submatrix with rows 2-4 and columns 4-8 */ 619566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,3,2,1,&irow)); 629566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,5,4,1,&icol)); 639566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(mat,1,&irow,&icol,MAT_INITIAL_MATRIX,&submatrices)); 64c4762a1bSJed Brown submat = *submatrices; 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* Test reuse submatrices */ 679566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(mat,1,&irow,&icol,MAT_REUSE_MATRIX,&submatrices)); 68c4762a1bSJed Brown 69c4762a1bSJed Brown /* sviewer will cause the submatrices (one per processor) to be printed in the correct order */ 709566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"\nSubmatrices\n")); 719566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 729566063dSJacob Faibussowitsch PetscCall(MatView(submat,sviewer)); 739566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 749566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 759566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)submat)); 769566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(1,&submatrices)); 779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&submat)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* Form submatrix with rows 2-4 and all columns */ 809566063dSJacob Faibussowitsch PetscCall(ISDestroy(&icol)); 819566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,10,0,1,&icol)); 829566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(mat,1,&irow,&icol,MAT_INITIAL_MATRIX,&submatrices)); 839566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(mat,1,&irow,&icol,MAT_REUSE_MATRIX,&submatrices)); 84c4762a1bSJed Brown submat = *submatrices; 85c4762a1bSJed Brown 869566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"\nSubmatrices with allcolumns\n")); 879566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 889566063dSJacob Faibussowitsch PetscCall(MatView(submat,sviewer)); 899566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 909566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* Test MatDuplicate */ 939566063dSJacob Faibussowitsch PetscCall(MatDuplicate(submat,MAT_COPY_VALUES,&submat1)); 949566063dSJacob Faibussowitsch PetscCall(MatDestroy(&submat1)); 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* Zero the original matrix */ 979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"Original zeroed matrix\n")); 989566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(mat)); 999566063dSJacob Faibussowitsch PetscCall(MatView(mat,PETSC_VIEWER_STDOUT_WORLD)); 100c4762a1bSJed Brown 1019566063dSJacob Faibussowitsch PetscCall(ISDestroy(&irow)); 1029566063dSJacob Faibussowitsch PetscCall(ISDestroy(&icol)); 1039566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)submat)); 1049566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(1,&submatrices)); 1059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&submat)); 1069566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 1079566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 108b122ec5aSJacob Faibussowitsch return 0; 109c4762a1bSJed Brown } 110c4762a1bSJed Brown 111c4762a1bSJed Brown /*TEST 112c4762a1bSJed Brown 113c4762a1bSJed Brown test: 114c4762a1bSJed Brown args: -mat_type aij 115c4762a1bSJed Brown 116c4762a1bSJed Brown test: 117c4762a1bSJed Brown suffix: 2 118c4762a1bSJed Brown args: -mat_type dense 119c4762a1bSJed Brown 120c4762a1bSJed Brown test: 121c4762a1bSJed Brown suffix: 3 122c4762a1bSJed Brown nsize: 3 123c4762a1bSJed Brown args: -mat_type aij 124c4762a1bSJed Brown 125c4762a1bSJed Brown test: 126c4762a1bSJed Brown suffix: 4 127c4762a1bSJed Brown nsize: 3 128c4762a1bSJed Brown args: -mat_type dense 129c4762a1bSJed Brown 130c4762a1bSJed Brown test: 131c4762a1bSJed Brown suffix: 5 132c4762a1bSJed Brown nsize: 3 133c4762a1bSJed Brown args: -mat_type aij -test_all 134c4762a1bSJed Brown 135c4762a1bSJed Brown TEST*/ 136