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 PetscErrorCode ierr; 11c4762a1bSJed Brown IS irow,icol; 12c4762a1bSJed Brown PetscScalar value = 1.0; 13c4762a1bSJed Brown PetscViewer sviewer; 14c4762a1bSJed Brown PetscBool allA = PETSC_FALSE; 15c4762a1bSJed Brown 16c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 17*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_COMMON)); 18*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF,PETSC_VIEWER_ASCII_COMMON)); 19c4762a1bSJed Brown 20*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&mat)); 21*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(mat,PETSC_DECIDE,PETSC_DECIDE,m,n)); 22*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(mat)); 23*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(mat)); 24*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(mat,&rstart,&rend)); 25c4762a1bSJed Brown for (i=rstart; i<rend; i++) { 26c4762a1bSJed Brown value = (PetscReal)i+1; tmp = i % 5; 27*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(mat,1,&tmp,1,&i,&value,INSERT_VALUES)); 28c4762a1bSJed Brown } 29*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(mat,MAT_FINAL_ASSEMBLY)); 30*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(mat,MAT_FINAL_ASSEMBLY)); 31*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"Original matrix\n")); 32*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(mat,PETSC_VIEWER_STDOUT_WORLD)); 33c4762a1bSJed Brown 34c4762a1bSJed Brown /* Test MatCreateSubMatrix_XXX_All(), i.e., submatrix = A */ 35*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-test_all",&allA,NULL)); 36c4762a1bSJed Brown if (allA) { 37*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_SELF,m,0,1,&irow)); 38*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_SELF,n,0,1,&icol)); 39*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrices(mat,1,&irow,&icol,MAT_INITIAL_MATRIX,&submatrices)); 40*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 44*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"\nSubmatrices with all\n")); 45*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"--------------------\n")); 46*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 47*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(submat,sviewer)); 48*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 49*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 50c4762a1bSJed Brown 51*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&irow)); 52*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&icol)); 53c4762a1bSJed Brown 54c4762a1bSJed Brown /* test getting a reference on a submat */ 55*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)submat)); 56*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroySubMatrices(1,&submatrices)); 57*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&submat)); 58c4762a1bSJed Brown } 59c4762a1bSJed Brown 60c4762a1bSJed Brown /* Form submatrix with rows 2-4 and columns 4-8 */ 61*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_SELF,3,2,1,&irow)); 62*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_SELF,5,4,1,&icol)); 63*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrices(mat,1,&irow,&icol,MAT_INITIAL_MATRIX,&submatrices)); 64c4762a1bSJed Brown submat = *submatrices; 65c4762a1bSJed Brown 66c4762a1bSJed Brown /* Test reuse submatrices */ 67*5f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 70*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"\nSubmatrices\n")); 71*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 72*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(submat,sviewer)); 73*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 74*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 75*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)submat)); 76*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroySubMatrices(1,&submatrices)); 77*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&submat)); 78c4762a1bSJed Brown 79c4762a1bSJed Brown /* Form submatrix with rows 2-4 and all columns */ 80*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&icol)); 81*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_SELF,10,0,1,&icol)); 82*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrices(mat,1,&irow,&icol,MAT_INITIAL_MATRIX,&submatrices)); 83*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrices(mat,1,&irow,&icol,MAT_REUSE_MATRIX,&submatrices)); 84c4762a1bSJed Brown submat = *submatrices; 85c4762a1bSJed Brown 86*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"\nSubmatrices with allcolumns\n")); 87*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 88*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(submat,sviewer)); 89*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD,PETSC_COMM_SELF,&sviewer)); 90*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(PETSC_VIEWER_STDOUT_WORLD)); 91c4762a1bSJed Brown 92c4762a1bSJed Brown /* Test MatDuplicate */ 93*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDuplicate(submat,MAT_COPY_VALUES,&submat1)); 94*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&submat1)); 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* Zero the original matrix */ 97*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD,"Original zeroed matrix\n")); 98*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroEntries(mat)); 99*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatView(mat,PETSC_VIEWER_STDOUT_WORLD)); 100c4762a1bSJed Brown 101*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&irow)); 102*5f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&icol)); 103*5f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)submat)); 104*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroySubMatrices(1,&submatrices)); 105*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&submat)); 106*5f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&mat)); 107c4762a1bSJed Brown ierr = PetscFinalize(); 108c4762a1bSJed Brown return ierr; 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