1c4762a1bSJed Brown static char help[] = "Tests MatCreateSubmatrix() with entire matrix, modified from ex59.c.";
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown
main(int argc,char ** args)5d71ae5a4SJacob Faibussowitsch int main(int argc, char **args)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown Mat C, A, Adup;
8c4762a1bSJed Brown PetscInt i, j, m = 3, n = 2, rstart, rend;
9c4762a1bSJed Brown PetscMPIInt size, rank;
10c4762a1bSJed Brown PetscScalar v;
11c4762a1bSJed Brown IS isrow;
12c4762a1bSJed Brown PetscBool detect_bug = PETSC_FALSE;
13c4762a1bSJed Brown
14327415f7SBarry Smith PetscFunctionBeginUser;
15*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
169566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-detect_bug", &detect_bug));
179566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
189566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
19c4762a1bSJed Brown n = 2 * size;
20c4762a1bSJed Brown
219566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
229566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
239566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C));
249566063dSJacob Faibussowitsch PetscCall(MatSetUp(C));
25c4762a1bSJed Brown
26c4762a1bSJed Brown /*
27c4762a1bSJed Brown This is JUST to generate a nice test matrix, all processors fill up
28c4762a1bSJed Brown the entire matrix. This is not something one would ever do in practice.
29c4762a1bSJed Brown */
309566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C, &rstart, &rend));
31c4762a1bSJed Brown for (i = rstart; i < rend; i++) {
32c4762a1bSJed Brown for (j = 0; j < m * n; j++) {
33c4762a1bSJed Brown v = i + j + 1;
349566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
35c4762a1bSJed Brown }
36c4762a1bSJed Brown }
379566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
389566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
399566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_COMMON));
409566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
419566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
42c4762a1bSJed Brown
43c4762a1bSJed Brown /*
44c4762a1bSJed Brown Generate a new matrix consisting every row and column of the original matrix
45c4762a1bSJed Brown */
469566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C, &rstart, &rend));
47c4762a1bSJed Brown
48c4762a1bSJed Brown /* Create parallel IS with the rows we want on THIS processor */
49dd400576SPatrick Sanan if (detect_bug && rank == 0) {
509566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, 1, rstart, 1, &isrow));
51c4762a1bSJed Brown } else {
529566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, rend - rstart, rstart, 1, &isrow));
53c4762a1bSJed Brown }
549566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(C, isrow, NULL, MAT_INITIAL_MATRIX, &A));
55c4762a1bSJed Brown
56c4762a1bSJed Brown /* Change C to test the case MAT_REUSE_MATRIX */
57dd400576SPatrick Sanan if (rank == 0) {
589371c9d4SSatish Balay i = 0;
599371c9d4SSatish Balay j = 0;
609371c9d4SSatish Balay v = 100;
619566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
62c4762a1bSJed Brown }
639566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
649566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
65c4762a1bSJed Brown
669566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(C, isrow, NULL, MAT_REUSE_MATRIX, &A));
679566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_COMMON));
689566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
699566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
70c4762a1bSJed Brown
71c4762a1bSJed Brown /* Test MatDuplicate */
729566063dSJacob Faibussowitsch PetscCall(MatDuplicate(A, MAT_COPY_VALUES, &Adup));
739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Adup));
74c4762a1bSJed Brown
759566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow));
769566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
789566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
79b122ec5aSJacob Faibussowitsch return 0;
80c4762a1bSJed Brown }
81c4762a1bSJed Brown
82c4762a1bSJed Brown /*TEST
83c4762a1bSJed Brown
84c4762a1bSJed Brown test:
85c4762a1bSJed Brown nsize: 2
86c4762a1bSJed Brown filter: grep -v "Mat Object"
87c4762a1bSJed Brown requires: !complex
88c4762a1bSJed Brown
89c4762a1bSJed Brown test:
90c4762a1bSJed Brown suffix: 2
91c4762a1bSJed Brown nsize: 3
92c4762a1bSJed Brown args: -detect_bug
93c4762a1bSJed Brown filter: grep -v "Mat Object"
94c4762a1bSJed Brown requires: !complex
95c4762a1bSJed Brown
96c4762a1bSJed Brown TEST*/
97