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