xref: /petsc/src/mat/tests/ex159.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown static const char help[] = "Test MatGetLocalSubMatrix() with multiple levels of nesting.\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown #include <petscmat.h>
4c4762a1bSJed Brown 
5c4762a1bSJed Brown int main(int argc, char *argv[])
6c4762a1bSJed Brown {
7c4762a1bSJed Brown   IS             is0a,is0b,is0,is1,isl0a,isl0b,isl0,isl1;
8c4762a1bSJed Brown   Mat            A,Aexplicit;
9c4762a1bSJed Brown   PetscBool      usenest;
10c4762a1bSJed Brown   PetscMPIInt    rank,size;
11c4762a1bSJed Brown   PetscInt       i,j;
12c4762a1bSJed Brown 
13*327415f7SBarry 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 
18c4762a1bSJed Brown   {
19c4762a1bSJed Brown     PetscInt ix0a[1],ix0b[1],ix0[2],ix1[1];
20c4762a1bSJed Brown 
21c4762a1bSJed Brown     ix0a[0] = rank*2+0;
22c4762a1bSJed Brown     ix0b[0] = rank*2+1;
23c4762a1bSJed Brown     ix0[0]  = rank*3+0; ix0[1] = rank*3+1;
24c4762a1bSJed Brown     ix1[0]  = rank*3+2;
259566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,1,ix0a,PETSC_COPY_VALUES,&is0a));
269566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,1,ix0b,PETSC_COPY_VALUES,&is0b));
279566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,2,ix0,PETSC_COPY_VALUES,&is0));
289566063dSJacob Faibussowitsch     PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,1,ix1,PETSC_COPY_VALUES,&is1));
29c4762a1bSJed Brown   }
30c4762a1bSJed Brown   {
319566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,6,0,1,&isl0));
329566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,3,0,1,&isl0a));
339566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,3,3,1,&isl0b));
349566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF,3,6,1,&isl1));
35c4762a1bSJed Brown   }
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   usenest = PETSC_FALSE;
389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-nest",&usenest,NULL));
39c4762a1bSJed Brown   if (usenest) {
40c4762a1bSJed Brown     ISLocalToGlobalMapping l2g;
41c4762a1bSJed Brown     PetscInt               l2gind[3];
42c4762a1bSJed Brown     Mat                    B[9];
43c4762a1bSJed Brown 
44c4762a1bSJed Brown     l2gind[0] = (rank-1+size)%size; l2gind[1] = rank; l2gind[2] = (rank+1)%size;
459566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,3,l2gind,PETSC_COPY_VALUES,&l2g));
46c4762a1bSJed Brown     for (i=0; i<9; i++) {
479566063dSJacob Faibussowitsch       PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,1,1,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,NULL,PETSC_DECIDE,NULL,&B[i]));
489566063dSJacob Faibussowitsch       PetscCall(MatSetUp(B[i]));
499566063dSJacob Faibussowitsch       PetscCall(MatSetLocalToGlobalMapping(B[i],l2g,l2g));
50c4762a1bSJed Brown     }
51c4762a1bSJed Brown     {
52c4762a1bSJed Brown       IS  isx[2];
53c4762a1bSJed Brown       Mat Bx00[4],Bx01[2],Bx10[2];
54c4762a1bSJed Brown       Mat B00,B01,B10;
55c4762a1bSJed Brown 
56c4762a1bSJed Brown       isx[0]  = is0a; isx[1] = is0b;
57c4762a1bSJed Brown       Bx00[0] = B[0]; Bx00[1] = B[1]; Bx00[2] = B[3]; Bx00[3] = B[4];
58c4762a1bSJed Brown       Bx01[0] = B[2]; Bx01[1] = B[5];
59c4762a1bSJed Brown       Bx10[0] = B[6]; Bx10[1] = B[7];
60c4762a1bSJed Brown 
619566063dSJacob Faibussowitsch       PetscCall(MatCreateNest(PETSC_COMM_WORLD,2,isx,2,isx,Bx00,&B00));
629566063dSJacob Faibussowitsch       PetscCall(MatSetUp(B00));
639566063dSJacob Faibussowitsch       PetscCall(MatCreateNest(PETSC_COMM_WORLD,2,isx,1,NULL,Bx01,&B01));
649566063dSJacob Faibussowitsch       PetscCall(MatSetUp(B01));
659566063dSJacob Faibussowitsch       PetscCall(MatCreateNest(PETSC_COMM_WORLD,1,NULL,2,isx,Bx10,&B10));
669566063dSJacob Faibussowitsch       PetscCall(MatSetUp(B10));
67c4762a1bSJed Brown       {
68c4762a1bSJed Brown         Mat By[4];
69c4762a1bSJed Brown         IS  isy[2];
70c4762a1bSJed Brown 
71c4762a1bSJed Brown         By[0]  = B00; By[1] = B01; By[2] = B10; By[3] = B[8];
72c4762a1bSJed Brown         isy[0] = is0; isy[1] = is1;
73c4762a1bSJed Brown 
749566063dSJacob Faibussowitsch         PetscCall(MatCreateNest(PETSC_COMM_WORLD,2,isy,2,isy,By,&A));
759566063dSJacob Faibussowitsch         PetscCall(MatSetUp(A));
76c4762a1bSJed Brown       }
779566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&B00));
789566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&B01));
799566063dSJacob Faibussowitsch       PetscCall(MatDestroy(&B10));
80c4762a1bSJed Brown     }
819566063dSJacob Faibussowitsch     for (i=0; i<9; i++) PetscCall(MatDestroy(&B[i]));
829566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingDestroy(&l2g));
83c4762a1bSJed Brown   } else {
84c4762a1bSJed Brown     ISLocalToGlobalMapping l2g;
85c4762a1bSJed Brown     PetscInt               l2gind[9];
86c4762a1bSJed Brown     for (i=0; i<3; i++) for (j=0; j<3; j++) l2gind[3*i+j] = ((rank-1+j+size) % size)*3 + i;
879566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,9,l2gind,PETSC_COPY_VALUES,&l2g));
889566063dSJacob Faibussowitsch     PetscCall(MatCreateAIJ(PETSC_COMM_WORLD,3,3,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,NULL,PETSC_DECIDE,NULL,&A));
899566063dSJacob Faibussowitsch     PetscCall(MatSetLocalToGlobalMapping(A,l2g,l2g));
909566063dSJacob Faibussowitsch     PetscCall(ISLocalToGlobalMappingDestroy(&l2g));
91c4762a1bSJed Brown   }
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   {
94c4762a1bSJed Brown     Mat A00,A11,A0a0a,A0a0b;
959566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSubMatrix(A,isl0,isl0,&A00));
969566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSubMatrix(A,isl1,isl1,&A11));
979566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSubMatrix(A00,isl0a,isl0a,&A0a0a));
989566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSubMatrix(A00,isl0a,isl0b,&A0a0b));
99c4762a1bSJed Brown 
1009566063dSJacob Faibussowitsch     PetscCall(MatSetValueLocal(A0a0a,0,0,100*rank+1,ADD_VALUES));
1019566063dSJacob Faibussowitsch     PetscCall(MatSetValueLocal(A0a0a,0,1,100*rank+2,ADD_VALUES));
1029566063dSJacob Faibussowitsch     PetscCall(MatSetValueLocal(A0a0a,2,2,100*rank+9,ADD_VALUES));
103c4762a1bSJed Brown 
1049566063dSJacob Faibussowitsch     PetscCall(MatSetValueLocal(A0a0b,1,1,100*rank+50+5,ADD_VALUES));
105c4762a1bSJed Brown 
1069566063dSJacob Faibussowitsch     PetscCall(MatSetValueLocal(A11,0,0,1000*(rank+1)+1,ADD_VALUES));
1079566063dSJacob Faibussowitsch     PetscCall(MatSetValueLocal(A11,1,2,1000*(rank+1)+6,ADD_VALUES));
108c4762a1bSJed Brown 
1099566063dSJacob Faibussowitsch     PetscCall(MatRestoreLocalSubMatrix(A00,isl0a,isl0a,&A0a0a));
1109566063dSJacob Faibussowitsch     PetscCall(MatRestoreLocalSubMatrix(A00,isl0a,isl0b,&A0a0b));
1119566063dSJacob Faibussowitsch     PetscCall(MatRestoreLocalSubMatrix(A,isl0,isl0,&A00));
1129566063dSJacob Faibussowitsch     PetscCall(MatRestoreLocalSubMatrix(A,isl1,isl1,&A11));
113c4762a1bSJed Brown   }
1149566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1159566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
116c4762a1bSJed Brown 
1179566063dSJacob Faibussowitsch   PetscCall(MatComputeOperator(A,MATAIJ,&Aexplicit));
1189566063dSJacob Faibussowitsch   PetscCall(MatView(Aexplicit,PETSC_VIEWER_STDOUT_WORLD));
119c4762a1bSJed Brown 
1209566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1219566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Aexplicit));
1229566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is0a));
1239566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is0b));
1249566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is0));
1259566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&is1));
1269566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isl0a));
1279566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isl0b));
1289566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isl0));
1299566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&isl1));
1309566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
131b122ec5aSJacob Faibussowitsch   return 0;
132c4762a1bSJed Brown }
133c4762a1bSJed Brown 
134c4762a1bSJed Brown /*TEST
135c4762a1bSJed Brown 
136c4762a1bSJed Brown    test:
137c4762a1bSJed Brown       nsize: 3
138c4762a1bSJed Brown 
139c4762a1bSJed Brown    test:
140c4762a1bSJed Brown       suffix: nest
141c4762a1bSJed Brown       nsize: 3
142c4762a1bSJed Brown       args: -nest
143c4762a1bSJed Brown 
144c4762a1bSJed Brown TEST*/
145