xref: /petsc/src/mat/tests/ex159.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode ierr;
8c4762a1bSJed Brown   IS             is0a,is0b,is0,is1,isl0a,isl0b,isl0,isl1;
9c4762a1bSJed Brown   Mat            A,Aexplicit;
10c4762a1bSJed Brown   PetscBool      usenest;
11c4762a1bSJed Brown   PetscMPIInt    rank,size;
12c4762a1bSJed Brown   PetscInt       i,j;
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
15*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
16*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(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;
25*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,1,ix0a,PETSC_COPY_VALUES,&is0a));
26*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,1,ix0b,PETSC_COPY_VALUES,&is0b));
27*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,2,ix0,PETSC_COPY_VALUES,&is0));
28*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateGeneral(PETSC_COMM_WORLD,1,ix1,PETSC_COPY_VALUES,&is1));
29c4762a1bSJed Brown   }
30c4762a1bSJed Brown   {
31*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,6,0,1,&isl0));
32*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,3,0,1,&isl0a));
33*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,3,3,1,&isl0b));
34*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISCreateStride(PETSC_COMM_SELF,3,6,1,&isl1));
35c4762a1bSJed Brown   }
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   usenest = PETSC_FALSE;
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
45*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,3,l2gind,PETSC_COPY_VALUES,&l2g));
46c4762a1bSJed Brown     for (i=0; i<9; i++) {
47*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateAIJ(PETSC_COMM_WORLD,1,1,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,NULL,PETSC_DECIDE,NULL,&B[i]));
48*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetUp(B[i]));
49*5f80ce2aSJacob Faibussowitsch       CHKERRQ(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 
61*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateNest(PETSC_COMM_WORLD,2,isx,2,isx,Bx00,&B00));
62*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetUp(B00));
63*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateNest(PETSC_COMM_WORLD,2,isx,1,NULL,Bx01,&B01));
64*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetUp(B01));
65*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatCreateNest(PETSC_COMM_WORLD,1,NULL,2,isx,Bx10,&B10));
66*5f80ce2aSJacob Faibussowitsch       CHKERRQ(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 
74*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatCreateNest(PETSC_COMM_WORLD,2,isy,2,isy,By,&A));
75*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetUp(A));
76c4762a1bSJed Brown       }
77*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&B00));
78*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&B01));
79*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatDestroy(&B10));
80c4762a1bSJed Brown     }
81*5f80ce2aSJacob Faibussowitsch     for (i=0; i<9; i++) CHKERRQ(MatDestroy(&B[i]));
82*5f80ce2aSJacob Faibussowitsch     CHKERRQ(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;
87*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,9,l2gind,PETSC_COPY_VALUES,&l2g));
88*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateAIJ(PETSC_COMM_WORLD,3,3,PETSC_DECIDE,PETSC_DECIDE,PETSC_DECIDE,NULL,PETSC_DECIDE,NULL,&A));
89*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetLocalToGlobalMapping(A,l2g,l2g));
90*5f80ce2aSJacob Faibussowitsch     CHKERRQ(ISLocalToGlobalMappingDestroy(&l2g));
91c4762a1bSJed Brown   }
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   {
94c4762a1bSJed Brown     Mat A00,A11,A0a0a,A0a0b;
95*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetLocalSubMatrix(A,isl0,isl0,&A00));
96*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetLocalSubMatrix(A,isl1,isl1,&A11));
97*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetLocalSubMatrix(A00,isl0a,isl0a,&A0a0a));
98*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetLocalSubMatrix(A00,isl0a,isl0b,&A0a0b));
99c4762a1bSJed Brown 
100*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValueLocal(A0a0a,0,0,100*rank+1,ADD_VALUES));
101*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValueLocal(A0a0a,0,1,100*rank+2,ADD_VALUES));
102*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValueLocal(A0a0a,2,2,100*rank+9,ADD_VALUES));
103c4762a1bSJed Brown 
104*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValueLocal(A0a0b,1,1,100*rank+50+5,ADD_VALUES));
105c4762a1bSJed Brown 
106*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValueLocal(A11,0,0,1000*(rank+1)+1,ADD_VALUES));
107*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValueLocal(A11,1,2,1000*(rank+1)+6,ADD_VALUES));
108c4762a1bSJed Brown 
109*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatRestoreLocalSubMatrix(A00,isl0a,isl0a,&A0a0a));
110*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatRestoreLocalSubMatrix(A00,isl0a,isl0b,&A0a0b));
111*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatRestoreLocalSubMatrix(A,isl0,isl0,&A00));
112*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatRestoreLocalSubMatrix(A,isl1,isl1,&A11));
113c4762a1bSJed Brown   }
114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
115*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
116c4762a1bSJed Brown 
117*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatComputeOperator(A,MATAIJ,&Aexplicit));
118*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(Aexplicit,PETSC_VIEWER_STDOUT_WORLD));
119c4762a1bSJed Brown 
120*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
121*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Aexplicit));
122*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is0a));
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is0b));
124*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is0));
125*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&is1));
126*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isl0a));
127*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isl0b));
128*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isl0));
129*5f80ce2aSJacob Faibussowitsch   CHKERRQ(ISDestroy(&isl1));
130c4762a1bSJed Brown   ierr = PetscFinalize();
131c4762a1bSJed Brown   return ierr;
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