1c4762a1bSJed Brown static char help[] = "Creates a matrix, inserts some values, and tests MatCreateSubMatrices() and MatZeroEntries().\n\n";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown
main(int argc,char ** argv)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown Mat mat, submat, submat1, *submatrices;
8c4762a1bSJed Brown PetscInt m = 10, n = 10, i = 4, tmp, rstart, rend;
9c4762a1bSJed Brown IS irow, icol;
10c4762a1bSJed Brown PetscScalar value = 1.0;
11c4762a1bSJed Brown PetscViewer sviewer;
12c4762a1bSJed Brown PetscBool allA = PETSC_FALSE;
13c4762a1bSJed Brown
14327415f7SBarry Smith PetscFunctionBeginUser;
15*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &argv, NULL, help));
169566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_COMMON));
179566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_SELF, PETSC_VIEWER_ASCII_COMMON));
18c4762a1bSJed Brown
199566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &mat));
209566063dSJacob Faibussowitsch PetscCall(MatSetSizes(mat, PETSC_DECIDE, PETSC_DECIDE, m, n));
219566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(mat));
229566063dSJacob Faibussowitsch PetscCall(MatSetUp(mat));
239566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(mat, &rstart, &rend));
24c4762a1bSJed Brown for (i = rstart; i < rend; i++) {
259371c9d4SSatish Balay value = (PetscReal)i + 1;
269371c9d4SSatish Balay 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));
49c4762a1bSJed Brown
509566063dSJacob Faibussowitsch PetscCall(ISDestroy(&irow));
519566063dSJacob Faibussowitsch PetscCall(ISDestroy(&icol));
52c4762a1bSJed Brown
53c4762a1bSJed Brown /* test getting a reference on a submat */
549566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)submat));
559566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(1, &submatrices));
569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&submat));
57c4762a1bSJed Brown }
58c4762a1bSJed Brown
59c4762a1bSJed Brown /* Form submatrix with rows 2-4 and columns 4-8 */
609566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, 3, 2, 1, &irow));
619566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, 5, 4, 1, &icol));
629566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(mat, 1, &irow, &icol, MAT_INITIAL_MATRIX, &submatrices));
63c4762a1bSJed Brown submat = *submatrices;
64c4762a1bSJed Brown
65c4762a1bSJed Brown /* Test reuse submatrices */
669566063dSJacob 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 */
699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "\nSubmatrices\n"));
709566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &sviewer));
719566063dSJacob Faibussowitsch PetscCall(MatView(submat, sviewer));
729566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &sviewer));
739566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)submat));
749566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(1, &submatrices));
759566063dSJacob Faibussowitsch PetscCall(MatDestroy(&submat));
76c4762a1bSJed Brown
77c4762a1bSJed Brown /* Form submatrix with rows 2-4 and all columns */
789566063dSJacob Faibussowitsch PetscCall(ISDestroy(&icol));
799566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, 10, 0, 1, &icol));
809566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(mat, 1, &irow, &icol, MAT_INITIAL_MATRIX, &submatrices));
819566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(mat, 1, &irow, &icol, MAT_REUSE_MATRIX, &submatrices));
82c4762a1bSJed Brown submat = *submatrices;
83c4762a1bSJed Brown
849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "\nSubmatrices with allcolumns\n"));
859566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &sviewer));
869566063dSJacob Faibussowitsch PetscCall(MatView(submat, sviewer));
879566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(PETSC_VIEWER_STDOUT_WORLD, PETSC_COMM_SELF, &sviewer));
88c4762a1bSJed Brown
89c4762a1bSJed Brown /* Test MatDuplicate */
909566063dSJacob Faibussowitsch PetscCall(MatDuplicate(submat, MAT_COPY_VALUES, &submat1));
919566063dSJacob Faibussowitsch PetscCall(MatDestroy(&submat1));
92c4762a1bSJed Brown
93c4762a1bSJed Brown /* Zero the original matrix */
949566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(PETSC_VIEWER_STDOUT_WORLD, "Original zeroed matrix\n"));
959566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(mat));
969566063dSJacob Faibussowitsch PetscCall(MatView(mat, PETSC_VIEWER_STDOUT_WORLD));
97c4762a1bSJed Brown
989566063dSJacob Faibussowitsch PetscCall(ISDestroy(&irow));
999566063dSJacob Faibussowitsch PetscCall(ISDestroy(&icol));
1009566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)submat));
1019566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(1, &submatrices));
1029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&submat));
1039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&mat));
1049566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
105b122ec5aSJacob Faibussowitsch return 0;
106c4762a1bSJed Brown }
107c4762a1bSJed Brown
108c4762a1bSJed Brown /*TEST
109c4762a1bSJed Brown
110c4762a1bSJed Brown test:
111c4762a1bSJed Brown args: -mat_type aij
112c4762a1bSJed Brown
113c4762a1bSJed Brown test:
114c4762a1bSJed Brown suffix: 2
115c4762a1bSJed Brown args: -mat_type dense
116c4762a1bSJed Brown
117c4762a1bSJed Brown test:
118c4762a1bSJed Brown suffix: 3
119c4762a1bSJed Brown nsize: 3
120c4762a1bSJed Brown args: -mat_type aij
121c4762a1bSJed Brown
122c4762a1bSJed Brown test:
123c4762a1bSJed Brown suffix: 4
124c4762a1bSJed Brown nsize: 3
125c4762a1bSJed Brown args: -mat_type dense
126c4762a1bSJed Brown
127c4762a1bSJed Brown test:
128c4762a1bSJed Brown suffix: 5
129c4762a1bSJed Brown nsize: 3
130c4762a1bSJed Brown args: -mat_type aij -test_all
131c4762a1bSJed Brown
132c4762a1bSJed Brown TEST*/
133