xref: /petsc/src/mat/tests/ex159.c (revision ebead697dbf761eb322f829370bbe90b3bd93fa3)
1 static const char help[] = "Test MatGetLocalSubMatrix() with multiple levels of nesting.\n\n";
2 
3 #include <petscmat.h>
4 
5 int main(int argc, char *argv[])
6 {
7   IS             is0a,is0b,is0,is1,isl0a,isl0b,isl0,isl1;
8   Mat            A,Aexplicit;
9   PetscBool      usenest;
10   PetscMPIInt    rank,size;
11   PetscInt       i,j;
12 
13   PetscFunctionBeginUser;
14   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
15   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
16   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
17 
18   {
19     PetscInt ix0a[1],ix0b[1],ix0[2],ix1[1];
20 
21     ix0a[0] = rank*2+0;
22     ix0b[0] = rank*2+1;
23     ix0[0]  = rank*3+0; ix0[1] = rank*3+1;
24     ix1[0]  = rank*3+2;
25     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,1,ix0a,PETSC_COPY_VALUES,&is0a));
26     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,1,ix0b,PETSC_COPY_VALUES,&is0b));
27     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,2,ix0,PETSC_COPY_VALUES,&is0));
28     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,1,ix1,PETSC_COPY_VALUES,&is1));
29   }
30   {
31     PetscCall(ISCreateStride(PETSC_COMM_SELF,6,0,1,&isl0));
32     PetscCall(ISCreateStride(PETSC_COMM_SELF,3,0,1,&isl0a));
33     PetscCall(ISCreateStride(PETSC_COMM_SELF,3,3,1,&isl0b));
34     PetscCall(ISCreateStride(PETSC_COMM_SELF,3,6,1,&isl1));
35   }
36 
37   usenest = PETSC_FALSE;
38   PetscCall(PetscOptionsGetBool(NULL,NULL,"-nest",&usenest,NULL));
39   if (usenest) {
40     ISLocalToGlobalMapping l2g;
41     PetscInt               l2gind[3];
42     Mat                    B[9];
43 
44     l2gind[0] = (rank-1+size)%size; l2gind[1] = rank; l2gind[2] = (rank+1)%size;
45     PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,3,l2gind,PETSC_COPY_VALUES,&l2g));
46     for (i=0; i<9; i++) {
47       PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,1,1,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,NULL,PETSC_DECIDE,NULL,&B[i]));
48       PetscCall(MatSetUp(B[i]));
49       PetscCall(MatSetLocalToGlobalMapping(B[i],l2g,l2g));
50     }
51     {
52       IS  isx[2];
53       Mat Bx00[4],Bx01[2],Bx10[2];
54       Mat B00,B01,B10;
55 
56       isx[0]  = is0a; isx[1] = is0b;
57       Bx00[0] = B[0]; Bx00[1] = B[1]; Bx00[2] = B[3]; Bx00[3] = B[4];
58       Bx01[0] = B[2]; Bx01[1] = B[5];
59       Bx10[0] = B[6]; Bx10[1] = B[7];
60 
61       PetscCall(MatCreateNest(PETSC_COMM_WORLD,2,isx,2,isx,Bx00,&B00));
62       PetscCall(MatSetUp(B00));
63       PetscCall(MatCreateNest(PETSC_COMM_WORLD,2,isx,1,NULL,Bx01,&B01));
64       PetscCall(MatSetUp(B01));
65       PetscCall(MatCreateNest(PETSC_COMM_WORLD,1,NULL,2,isx,Bx10,&B10));
66       PetscCall(MatSetUp(B10));
67       {
68         Mat By[4];
69         IS  isy[2];
70 
71         By[0]  = B00; By[1] = B01; By[2] = B10; By[3] = B[8];
72         isy[0] = is0; isy[1] = is1;
73 
74         PetscCall(MatCreateNest(PETSC_COMM_WORLD,2,isy,2,isy,By,&A));
75         PetscCall(MatSetUp(A));
76       }
77       PetscCall(MatDestroy(&B00));
78       PetscCall(MatDestroy(&B01));
79       PetscCall(MatDestroy(&B10));
80     }
81     for (i=0; i<9; i++) PetscCall(MatDestroy(&B[i]));
82     PetscCall(ISLocalToGlobalMappingDestroy(&l2g));
83   } else {
84     ISLocalToGlobalMapping l2g;
85     PetscInt               l2gind[9];
86     for (i=0; i<3; i++) for (j=0; j<3; j++) l2gind[3*i+j] = ((rank-1+j+size) % size)*3 + i;
87     PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,9,l2gind,PETSC_COPY_VALUES,&l2g));
88     PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,3,3,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,NULL,PETSC_DECIDE,NULL,&A));
89     PetscCall(MatSetLocalToGlobalMapping(A,l2g,l2g));
90     PetscCall(ISLocalToGlobalMappingDestroy(&l2g));
91   }
92 
93   {
94     Mat A00,A11,A0a0a,A0a0b;
95     PetscCall(MatGetLocalSubMatrix(A,isl0,isl0,&A00));
96     PetscCall(MatGetLocalSubMatrix(A,isl1,isl1,&A11));
97     PetscCall(MatGetLocalSubMatrix(A00,isl0a,isl0a,&A0a0a));
98     PetscCall(MatGetLocalSubMatrix(A00,isl0a,isl0b,&A0a0b));
99 
100     PetscCall(MatSetValueLocal(A0a0a,0,0,100*rank+1,ADD_VALUES));
101     PetscCall(MatSetValueLocal(A0a0a,0,1,100*rank+2,ADD_VALUES));
102     PetscCall(MatSetValueLocal(A0a0a,2,2,100*rank+9,ADD_VALUES));
103 
104     PetscCall(MatSetValueLocal(A0a0b,1,1,100*rank+50+5,ADD_VALUES));
105 
106     PetscCall(MatSetValueLocal(A11,0,0,1000*(rank+1)+1,ADD_VALUES));
107     PetscCall(MatSetValueLocal(A11,1,2,1000*(rank+1)+6,ADD_VALUES));
108 
109     PetscCall(MatRestoreLocalSubMatrix(A00,isl0a,isl0a,&A0a0a));
110     PetscCall(MatRestoreLocalSubMatrix(A00,isl0a,isl0b,&A0a0b));
111     PetscCall(MatRestoreLocalSubMatrix(A,isl0,isl0,&A00));
112     PetscCall(MatRestoreLocalSubMatrix(A,isl1,isl1,&A11));
113   }
114   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
115   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
116 
117   PetscCall(MatComputeOperator(A,MATAIJ,&Aexplicit));
118   PetscCall(MatView(Aexplicit,PETSC_VIEWER_STDOUT_WORLD));
119 
120   PetscCall(MatDestroy(&A));
121   PetscCall(MatDestroy(&Aexplicit));
122   PetscCall(ISDestroy(&is0a));
123   PetscCall(ISDestroy(&is0b));
124   PetscCall(ISDestroy(&is0));
125   PetscCall(ISDestroy(&is1));
126   PetscCall(ISDestroy(&isl0a));
127   PetscCall(ISDestroy(&isl0b));
128   PetscCall(ISDestroy(&isl0));
129   PetscCall(ISDestroy(&isl1));
130   PetscCall(PetscFinalize());
131   return 0;
132 }
133 
134 /*TEST
135 
136    test:
137       nsize: 3
138 
139    test:
140       suffix: nest
141       nsize: 3
142       args: -nest
143 
144 TEST*/
145