xref: /petsc/src/mat/tests/ex159.c (revision 39b6f3f95800db5f8164a53bdb9c6f06ef9d762f)
1c4762a1bSJed Brown static const char help[] = "Test MatGetLocalSubMatrix() with multiple levels of nesting.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5d71ae5a4SJacob Faibussowitsch int main(int argc, char *argv[])
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown   IS          is0a, is0b, is0, is1, isl0a, isl0b, isl0, isl1;
8c4762a1bSJed Brown   Mat         A, Aexplicit;
9*39b6f3f9SStefano Zampini   PetscBool   usenest, test_mat_is;
10c4762a1bSJed Brown   PetscMPIInt rank, size;
11c4762a1bSJed Brown   PetscInt    i, j;
12c4762a1bSJed Brown 
13327415f7SBarry Smith   PetscFunctionBeginUser;
149566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
159566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
17c4762a1bSJed Brown 
18*39b6f3f9SStefano Zampini   usenest     = PETSC_FALSE;
19*39b6f3f9SStefano Zampini   test_mat_is = PETSC_FALSE;
20*39b6f3f9SStefano Zampini   PetscCall(PetscOptionsGetBool(NULL, NULL, "-nest", &usenest, NULL));
21*39b6f3f9SStefano Zampini   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test_mat_is", &test_mat_is, NULL));
22*39b6f3f9SStefano Zampini 
23c4762a1bSJed Brown   {
24c4762a1bSJed Brown     PetscInt ix0a[1], ix0b[1], ix0[2], ix1[1];
25c4762a1bSJed Brown 
26c4762a1bSJed Brown     ix0a[0] = rank * 2 + 0;
27c4762a1bSJed Brown     ix0b[0] = rank * 2 + 1;
289371c9d4SSatish Balay     ix0[0]  = rank * 3 + 0;
299371c9d4SSatish Balay     ix0[1]  = rank * 3 + 1;
30c4762a1bSJed Brown     ix1[0]  = rank * 3 + 2;
319566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 1, ix0a, PETSC_COPY_VALUES, &is0a));
329566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 1, ix0b, PETSC_COPY_VALUES, &is0b));
339566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 2, ix0, PETSC_COPY_VALUES, &is0));
349566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 1, ix1, PETSC_COPY_VALUES, &is1));
35c4762a1bSJed Brown   }
36c4762a1bSJed Brown   {
379566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, 6, 0, 1, &isl0));
389566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, 3, 0, 1, &isl0a));
399566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, 3, 3, 1, &isl0b));
409566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, 3, 6, 1, &isl1));
41c4762a1bSJed Brown   }
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   if (usenest) {
44c4762a1bSJed Brown     ISLocalToGlobalMapping l2g;
45c4762a1bSJed Brown     PetscInt               l2gind[3];
46c4762a1bSJed Brown     Mat                    B[9];
47c4762a1bSJed Brown 
489371c9d4SSatish Balay     l2gind[0] = (rank - 1 + size) % size;
499371c9d4SSatish Balay     l2gind[1] = rank;
509371c9d4SSatish Balay     l2gind[2] = (rank + 1) % size;
519566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 1, 3, l2gind, PETSC_COPY_VALUES, &l2g));
52c4762a1bSJed Brown     for (i = 0; i < 9; i++) {
53*39b6f3f9SStefano Zampini       if (test_mat_is) {
54*39b6f3f9SStefano Zampini         PetscCall(MatCreateIS(PETSC_COMM_WORLD, 1, 1, 1, PETSC_DECIDE, PETSC_DECIDE, l2g, l2g, &B[i]));
55*39b6f3f9SStefano Zampini       } else {
569566063dSJacob Faibussowitsch         PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, 1, 1, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, NULL, PETSC_DECIDE, NULL, &B[i]));
579566063dSJacob Faibussowitsch         PetscCall(MatSetUp(B[i]));
589566063dSJacob Faibussowitsch         PetscCall(MatSetLocalToGlobalMapping(B[i], l2g, l2g));
59c4762a1bSJed Brown       }
60*39b6f3f9SStefano Zampini     }
61c4762a1bSJed Brown     {
62c4762a1bSJed Brown       IS  isx[2];
63c4762a1bSJed Brown       Mat Bx00[4], Bx01[2], Bx10[2];
64c4762a1bSJed Brown       Mat B00, B01, B10;
65c4762a1bSJed Brown 
669371c9d4SSatish Balay       isx[0]  = is0a;
679371c9d4SSatish Balay       isx[1]  = is0b;
689371c9d4SSatish Balay       Bx00[0] = B[0];
699371c9d4SSatish Balay       Bx00[1] = B[1];
709371c9d4SSatish Balay       Bx00[2] = B[3];
719371c9d4SSatish Balay       Bx00[3] = B[4];
729371c9d4SSatish Balay       Bx01[0] = B[2];
739371c9d4SSatish Balay       Bx01[1] = B[5];
749371c9d4SSatish Balay       Bx10[0] = B[6];
759371c9d4SSatish Balay       Bx10[1] = B[7];
76c4762a1bSJed Brown 
779566063dSJacob Faibussowitsch       PetscCall(MatCreateNest(PETSC_COMM_WORLD, 2, isx, 2, isx, Bx00, &B00));
789566063dSJacob Faibussowitsch       PetscCall(MatSetUp(B00));
799566063dSJacob Faibussowitsch       PetscCall(MatCreateNest(PETSC_COMM_WORLD, 2, isx, 1, NULL, Bx01, &B01));
809566063dSJacob Faibussowitsch       PetscCall(MatSetUp(B01));
819566063dSJacob Faibussowitsch       PetscCall(MatCreateNest(PETSC_COMM_WORLD, 1, NULL, 2, isx, Bx10, &B10));
829566063dSJacob Faibussowitsch       PetscCall(MatSetUp(B10));
83c4762a1bSJed Brown       {
84c4762a1bSJed Brown         Mat By[4];
85c4762a1bSJed Brown         IS  isy[2];
86c4762a1bSJed Brown 
879371c9d4SSatish Balay         By[0]  = B00;
889371c9d4SSatish Balay         By[1]  = B01;
899371c9d4SSatish Balay         By[2]  = B10;
909371c9d4SSatish Balay         By[3]  = B[8];
919371c9d4SSatish Balay         isy[0] = is0;
929371c9d4SSatish Balay         isy[1] = is1;
93c4762a1bSJed Brown 
949566063dSJacob Faibussowitsch         PetscCall(MatCreateNest(PETSC_COMM_WORLD, 2, isy, 2, isy, By, &A));
959566063dSJacob Faibussowitsch         PetscCall(MatSetUp(A));
96c4762a1bSJed Brown       }
979566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&B00));
989566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&B01));
999566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&B10));
100c4762a1bSJed Brown     }
1019566063dSJacob Faibussowitsch     for (i = 0; i < 9; i++) PetscCall(MatDestroy(&B[i]));
1029566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingDestroy(&l2g));
103c4762a1bSJed Brown   } else {
104c4762a1bSJed Brown     ISLocalToGlobalMapping l2g;
105c4762a1bSJed Brown     PetscInt               l2gind[9];
1069371c9d4SSatish Balay     for (i = 0; i < 3; i++)
1079371c9d4SSatish Balay       for (j = 0; j < 3; j++) l2gind[3 * i + j] = ((rank - 1 + j + size) % size) * 3 + i;
1089566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD, 1, 9, l2gind, PETSC_COPY_VALUES, &l2g));
109*39b6f3f9SStefano Zampini     if (test_mat_is) {
110*39b6f3f9SStefano Zampini       PetscCall(MatCreateIS(PETSC_COMM_WORLD, 1, 3, 3, PETSC_DECIDE, PETSC_DECIDE, l2g, l2g, &A));
111*39b6f3f9SStefano Zampini     } else {
1129566063dSJacob Faibussowitsch       PetscCall(MatCreateAIJ(PETSC_COMM_WORLD, 3, 3, PETSC_DECIDE, PETSC_DECIDE, PETSC_DECIDE, NULL, PETSC_DECIDE, NULL, &A));
1139566063dSJacob Faibussowitsch       PetscCall(MatSetLocalToGlobalMapping(A, l2g, l2g));
114*39b6f3f9SStefano Zampini     }
1159566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingDestroy(&l2g));
116c4762a1bSJed Brown   }
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   {
119c4762a1bSJed Brown     Mat A00, A11, A0a0a, A0a0b;
1209566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSubMatrix(A, isl0, isl0, &A00));
1219566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSubMatrix(A, isl1, isl1, &A11));
1229566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSubMatrix(A00, isl0a, isl0a, &A0a0a));
1239566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSubMatrix(A00, isl0a, isl0b, &A0a0b));
124c4762a1bSJed Brown 
1259566063dSJacob Faibussowitsch     PetscCall(MatSetValueLocal(A0a0a, 0, 0, 100 * rank + 1, ADD_VALUES));
1269566063dSJacob Faibussowitsch     PetscCall(MatSetValueLocal(A0a0a, 0, 1, 100 * rank + 2, ADD_VALUES));
1279566063dSJacob Faibussowitsch     PetscCall(MatSetValueLocal(A0a0a, 2, 2, 100 * rank + 9, ADD_VALUES));
128c4762a1bSJed Brown 
1299566063dSJacob Faibussowitsch     PetscCall(MatSetValueLocal(A0a0b, 1, 1, 100 * rank + 50 + 5, ADD_VALUES));
130c4762a1bSJed Brown 
1319566063dSJacob Faibussowitsch     PetscCall(MatSetValueLocal(A11, 0, 0, 1000 * (rank + 1) + 1, ADD_VALUES));
1329566063dSJacob Faibussowitsch     PetscCall(MatSetValueLocal(A11, 1, 2, 1000 * (rank + 1) + 6, ADD_VALUES));
133c4762a1bSJed Brown 
1349566063dSJacob Faibussowitsch     PetscCall(MatRestoreLocalSubMatrix(A00, isl0a, isl0a, &A0a0a));
1359566063dSJacob Faibussowitsch     PetscCall(MatRestoreLocalSubMatrix(A00, isl0a, isl0b, &A0a0b));
1369566063dSJacob Faibussowitsch     PetscCall(MatRestoreLocalSubMatrix(A, isl0, isl0, &A00));
1379566063dSJacob Faibussowitsch     PetscCall(MatRestoreLocalSubMatrix(A, isl1, isl1, &A11));
138c4762a1bSJed Brown   }
1399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
141c4762a1bSJed Brown 
1429566063dSJacob Faibussowitsch   PetscCall(MatComputeOperator(A, MATAIJ, &Aexplicit));
1439566063dSJacob Faibussowitsch   PetscCall(MatView(Aexplicit, PETSC_VIEWER_STDOUT_WORLD));
144c4762a1bSJed Brown 
1459566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1469566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Aexplicit));
1479566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is0a));
1489566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is0b));
1499566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is0));
1509566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is1));
1519566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isl0a));
1529566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isl0b));
1539566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isl0));
1549566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isl1));
1559566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
156b122ec5aSJacob Faibussowitsch   return 0;
157c4762a1bSJed Brown }
158c4762a1bSJed Brown 
159c4762a1bSJed Brown /*TEST
160c4762a1bSJed Brown 
161c4762a1bSJed Brown    test:
162c4762a1bSJed Brown       nsize: 3
163*39b6f3f9SStefano Zampini       args: -test_mat_is {{0 1}}
164*39b6f3f9SStefano Zampini       diff_args: -j
165c4762a1bSJed Brown 
166c4762a1bSJed Brown    test:
167c4762a1bSJed Brown       suffix: nest
168c4762a1bSJed Brown       nsize: 3
169*39b6f3f9SStefano Zampini       args: -nest -test_mat_is {{0 1}}
170*39b6f3f9SStefano Zampini       diff_args: -j
171c4762a1bSJed Brown 
172c4762a1bSJed Brown TEST*/
173