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*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 16*9566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON)); 17*9566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_COMMON)); 18c4762a1bSJed Brown 19*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&mat)); 20*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat,PETSC_DECIDE,PETSC_DECIDE,m,n)); 21*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(mat)); 22*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(mat)); 23*9566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(mat,&rstart,&rend)); 24c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 25c4762a1bSJed Brown value = (PetscReal)i+1; tmp = i % 5; 26*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(mat,1,&tmp,1,&i,&value,INSERT_VALUES)); 27c4762a1bSJed Brown } 28*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY)); 29*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY)); 30*9566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"Original matrix\n")); 31*9566063dSJacob Faibussowitsch PetscCall(MatView(mat,PETSC_VIEWER_STDOUT_WORLD)); 32c4762a1bSJed Brown 33c4762a1bSJed Brown /* Test MatCreateSubMatrix_XXX_All(), i.e., submatrix = A */ 34*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,NULL,"-test_all",&allA,NULL)); 35c4762a1bSJed Brown if (allA) { 36*9566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,m,0,1,&irow)); 37*9566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,n,0,1,&icol)); 38*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(mat,1,&irow,&icol,MAT_INITIAL_MATRIX,&submatrices)); 39*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(mat,1,&irow,&icol,MAT_REUSE_MATRIX,&submatrices)); 40c4762a1bSJed Brown submat = *submatrices; 41c4762a1bSJed Brown 42c4762a1bSJed Brown /* sviewer will cause the submatrices (one per processor) to be printed in the correct order */ 43*9566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"\nSubmatrices with all\n")); 44*9566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"--------------------\n")); 45*9566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 46*9566063dSJacob Faibussowitsch PetscCall(MatView(submat,sviewer)); 47*9566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 48*9566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 49c4762a1bSJed Brown 50*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&irow)); 51*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&icol)); 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* test getting a reference on a submat */ 54*9566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)submat)); 55*9566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(1,&submatrices)); 56*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&submat)); 57c4762a1bSJed Brown } 58c4762a1bSJed Brown 59c4762a1bSJed Brown /* Form submatrix with rows 2-4 and columns 4-8 */ 60*9566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,3,2,1,&irow)); 61*9566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,5,4,1,&icol)); 62*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(mat,1,&irow,&icol,MAT_INITIAL_MATRIX,&submatrices)); 63c4762a1bSJed Brown submat = *submatrices; 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* Test reuse submatrices */ 66*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(mat,1,&irow,&icol,MAT_REUSE_MATRIX,&submatrices)); 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* sviewer will cause the submatrices (one per processor) to be printed in the correct order */ 69*9566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"\nSubmatrices\n")); 70*9566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 71*9566063dSJacob Faibussowitsch PetscCall(MatView(submat,sviewer)); 72*9566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 73*9566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 74*9566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)submat)); 75*9566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(1,&submatrices)); 76*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&submat)); 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* Form submatrix with rows 2-4 and all columns */ 79*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&icol)); 80*9566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,10,0,1,&icol)); 81*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(mat,1,&irow,&icol,MAT_INITIAL_MATRIX,&submatrices)); 82*9566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(mat,1,&irow,&icol,MAT_REUSE_MATRIX,&submatrices)); 83c4762a1bSJed Brown submat = *submatrices; 84c4762a1bSJed Brown 85*9566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"\nSubmatrices with allcolumns\n")); 86*9566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 87*9566063dSJacob Faibussowitsch PetscCall(MatView(submat,sviewer)); 88*9566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 89*9566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 90c4762a1bSJed Brown 91c4762a1bSJed Brown /* Test MatDuplicate */ 92*9566063dSJacob Faibussowitsch PetscCall(MatDuplicate(submat,MAT_COPY_VALUES,&submat1)); 93*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&submat1)); 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* Zero the original matrix */ 96*9566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"Original zeroed matrix\n")); 97*9566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(mat)); 98*9566063dSJacob Faibussowitsch PetscCall(MatView(mat,PETSC_VIEWER_STDOUT_WORLD)); 99c4762a1bSJed Brown 100*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&irow)); 101*9566063dSJacob Faibussowitsch PetscCall(ISDestroy(&icol)); 102*9566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)submat)); 103*9566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(1,&submatrices)); 104*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&submat)); 105*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat)); 106*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 107b122ec5aSJacob Faibussowitsch return 0; 108c4762a1bSJed Brown } 109c4762a1bSJed Brown 110c4762a1bSJed Brown /*TEST 111c4762a1bSJed Brown 112c4762a1bSJed Brown test: 113c4762a1bSJed Brown args: -mat_type aij 114c4762a1bSJed Brown 115c4762a1bSJed Brown test: 116c4762a1bSJed Brown suffix: 2 117c4762a1bSJed Brown args: -mat_type dense 118c4762a1bSJed Brown 119c4762a1bSJed Brown test: 120c4762a1bSJed Brown suffix: 3 121c4762a1bSJed Brown nsize: 3 122c4762a1bSJed Brown args: -mat_type aij 123c4762a1bSJed Brown 124c4762a1bSJed Brown test: 125c4762a1bSJed Brown suffix: 4 126c4762a1bSJed Brown nsize: 3 127c4762a1bSJed Brown args: -mat_type dense 128c4762a1bSJed Brown 129c4762a1bSJed Brown test: 130c4762a1bSJed Brown suffix: 5 131c4762a1bSJed Brown nsize: 3 132c4762a1bSJed Brown args: -mat_type aij -test_all 133c4762a1bSJed Brown 134c4762a1bSJed Brown TEST*/ 135