1c4762a1bSJed Brown static char help[] = "Tests MatCreateSubmatrix() with certain entire rows of matrix, modified from ex181.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;
8c4762a1bSJed Brown PetscInt i, j, m = 3, n = 2, rstart, rend, cstart, cend;
9c4762a1bSJed Brown PetscMPIInt size, rank;
10c4762a1bSJed Brown PetscScalar v;
11c4762a1bSJed Brown IS isrow, iscol;
12c4762a1bSJed Brown
13327415f7SBarry Smith PetscFunctionBeginUser;
14*c8025a54SPierre Jolivet PetscCall(PetscInitialize(&argc, &args, NULL, help));
159566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
169566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
17c4762a1bSJed Brown n = 2 * size;
18c4762a1bSJed Brown
199566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &C));
209566063dSJacob Faibussowitsch PetscCall(MatSetSizes(C, PETSC_DECIDE, PETSC_DECIDE, m * n, m * n));
219566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(C));
229566063dSJacob Faibussowitsch PetscCall(MatSetUp(C));
23c4762a1bSJed Brown
24c4762a1bSJed Brown /*
25c4762a1bSJed Brown This is JUST to generate a nice test matrix, all processors fill up
26c4762a1bSJed Brown the entire matrix. This is not something one would ever do in practice.
27c4762a1bSJed Brown */
289566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C, &rstart, &rend));
29c4762a1bSJed Brown for (i = rstart; i < rend; i++) {
30c4762a1bSJed Brown for (j = 0; j < m * n; j++) {
31c4762a1bSJed Brown v = i + j + 1;
329566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
33c4762a1bSJed Brown }
34c4762a1bSJed Brown }
359566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
369566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
379566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_COMMON));
389566063dSJacob Faibussowitsch PetscCall(MatView(C, PETSC_VIEWER_STDOUT_WORLD));
399566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
40c4762a1bSJed Brown
41c4762a1bSJed Brown /*
42c4762a1bSJed Brown Generate a new matrix consisting every column of the original matrix
43c4762a1bSJed Brown */
449566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(C, &rstart, &rend));
459566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRangeColumn(C, &cstart, &cend));
46c4762a1bSJed Brown
479566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, rend - rstart > 0 ? rend - rstart - 1 : 0, rstart, 1, &isrow));
489566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_WORLD, cend - cstart, cstart, 1, &iscol));
499566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(C, isrow, NULL, MAT_INITIAL_MATRIX, &A));
50c4762a1bSJed Brown
51c4762a1bSJed Brown /* Change C to test the case MAT_REUSE_MATRIX */
52dd400576SPatrick Sanan if (rank == 0) {
539371c9d4SSatish Balay i = 0;
549371c9d4SSatish Balay j = 0;
559371c9d4SSatish Balay v = 100;
569566063dSJacob Faibussowitsch PetscCall(MatSetValues(C, 1, &i, 1, &j, &v, INSERT_VALUES));
57c4762a1bSJed Brown }
589566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(C, MAT_FINAL_ASSEMBLY));
599566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(C, MAT_FINAL_ASSEMBLY));
60c4762a1bSJed Brown
619566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(C, isrow, NULL, MAT_REUSE_MATRIX, &A));
629566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_COMMON));
639566063dSJacob Faibussowitsch PetscCall(MatView(A, PETSC_VIEWER_STDOUT_WORLD));
649566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
65c4762a1bSJed Brown
669566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isrow));
679566063dSJacob Faibussowitsch PetscCall(ISDestroy(&iscol));
689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&C));
709566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
71b122ec5aSJacob Faibussowitsch return 0;
72c4762a1bSJed Brown }
73c4762a1bSJed Brown
74c4762a1bSJed Brown /*TEST
75c4762a1bSJed Brown
76c4762a1bSJed Brown test:
77c4762a1bSJed Brown nsize: 2
78c4762a1bSJed Brown
79c4762a1bSJed Brown TEST*/
80