xref: /petsc/src/mat/impls/nest/matnest.c (revision 7874fa86d92623c024d734b6b3975594867d60ee)
1d8588912SDave May 
2ec9811eeSJed Brown #include "matnestimpl.h" /*I   "petscmat.h"   I*/
3d8588912SDave May 
4c8883902SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]);
5*7874fa86SDave May static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left);
6c8883902SJed Brown 
7d8588912SDave May /* private functions */
8d8588912SDave May #undef __FUNCT__
9c8883902SJed Brown #define __FUNCT__ "MatNestGetSize_Private"
10c8883902SJed Brown static PetscErrorCode MatNestGetSize_Private(Mat A,PetscBool globalsize,PetscInt *M,PetscInt *N)
11d8588912SDave May {
12d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
13b9dfa011SHong Zhang   PetscInt       sizeM,sizeN,i,j,index = -1;
14d8588912SDave May   PetscErrorCode ierr;
15d8588912SDave May 
16d8588912SDave May   PetscFunctionBegin;
17d8588912SDave May   /* rows */
18d8588912SDave May   for (i=0; i<bA->nr; i++) {
19d8588912SDave May     for (j=0; j<bA->nc; j++) {
20d8588912SDave May       if (bA->m[i][j]) { index = j; break; }
21d8588912SDave May     }
22d8588912SDave May     if (bA->m[i][index]) {
23d8588912SDave May       if (globalsize) { ierr = MatGetSize(bA->m[i][index],&sizeM,&sizeN);CHKERRQ(ierr); }
24d8588912SDave May       else {            ierr = MatGetLocalSize(bA->m[i][index],&sizeM,&sizeN);CHKERRQ(ierr); }
25d8588912SDave May       *M = *M + sizeM;
26d8588912SDave May     }
27d8588912SDave May   }
28d8588912SDave May 
29d8588912SDave May   /* cols */
30d8588912SDave May   /* now descend recursively */
31d8588912SDave May   for (j=0; j<bA->nc; j++) {
32d8588912SDave May     for (i=0; i<bA->nr; i++) {
33d8588912SDave May       if (bA->m[i][j]) { index = i; break; }
34d8588912SDave May     }
35d8588912SDave May     if (bA->m[index][j]) {
36d8588912SDave May       if (globalsize) { ierr = MatGetSize(bA->m[index][j],&sizeM,&sizeN);CHKERRQ(ierr); }
37d8588912SDave May       else {            ierr = MatGetLocalSize(bA->m[index][j],&sizeM,&sizeN);CHKERRQ(ierr); }
38d8588912SDave May       *N = *N + sizeN;
39d8588912SDave May     }
40d8588912SDave May   }
41d8588912SDave May   PetscFunctionReturn(0);
42d8588912SDave May }
43d8588912SDave May 
44d8588912SDave May #undef __FUNCT__
45d8588912SDave May #define __FUNCT__ "PETSc_MatNest_CheckMatVecCompatibility2"
46207556f9SJed Brown PETSC_UNUSED
47d8588912SDave May static PetscErrorCode PETSc_MatNest_CheckMatVecCompatibility2(Mat A,Vec x,Vec y)
48d8588912SDave May {
49d8588912SDave May   PetscBool      isANest, isxNest, isyNest;
50d8588912SDave May   PetscErrorCode ierr;
51d8588912SDave May 
52d8588912SDave May   PetscFunctionBegin;
53d8588912SDave May   isANest = isxNest = isyNest = PETSC_FALSE;
54d8588912SDave May   ierr = PetscTypeCompare( (PetscObject)A, MATNEST, &isANest );CHKERRQ(ierr);
55d8588912SDave May   ierr = PetscTypeCompare( (PetscObject)x, VECNEST, &isxNest );CHKERRQ(ierr);
56d8588912SDave May   ierr = PetscTypeCompare( (PetscObject)y, VECNEST, &isyNest );CHKERRQ(ierr);
57d8588912SDave May 
58d8588912SDave May   if (isANest && isxNest && isyNest) {
59d8588912SDave May     PetscFunctionReturn(0);
60d8588912SDave May   } else {
61d8588912SDave May     SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_SUP, "Operation only valid for Mat(Nest)-Vec(Nest) combinations");
62d8588912SDave May     PetscFunctionReturn(0);
63d8588912SDave May   }
64d8588912SDave May   PetscFunctionReturn(0);
65d8588912SDave May }
66d8588912SDave May 
67d8588912SDave May /*
68d8588912SDave May  This function is called every time we insert a sub matrix the Nest.
69d8588912SDave May  It ensures that every Nest along row r, and col c has the same dimensions
70d8588912SDave May  as the submat being inserted.
71d8588912SDave May */
72d8588912SDave May #undef __FUNCT__
73d8588912SDave May #define __FUNCT__ "PETSc_MatNest_CheckConsistency"
74063a28d4SJed Brown PETSC_UNUSED
75d8588912SDave May static PetscErrorCode PETSc_MatNest_CheckConsistency(Mat A,Mat submat,PetscInt r,PetscInt c)
76d8588912SDave May {
77d8588912SDave May   Mat_Nest       *b = (Mat_Nest*)A->data;
78d8588912SDave May   PetscInt       i,j;
79d8588912SDave May   PetscInt       nr,nc;
80d8588912SDave May   PetscInt       sM,sN,M,N;
81d8588912SDave May   Mat            mat;
82d8588912SDave May   PetscErrorCode ierr;
83d8588912SDave May 
84d8588912SDave May   PetscFunctionBegin;
85d8588912SDave May   nr = b->nr;
86d8588912SDave May   nc = b->nc;
87d8588912SDave May   ierr = MatGetSize(submat,&sM,&sN);CHKERRQ(ierr);
88d8588912SDave May   for (i=0; i<nr; i++) {
89d8588912SDave May     mat = b->m[i][c];
90d8588912SDave May     if (mat) {
91d8588912SDave May       ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
92d8588912SDave May       /* Check columns match */
93d8588912SDave May       if (sN != N) {
94d8588912SDave May         SETERRQ3(((PetscObject)A)->comm,PETSC_ERR_SUP,"Inserted incompatible submatrix into Nest at (%D,%D). Submatrix must have %D rows",r,c,N );
95d8588912SDave May       }
96d8588912SDave May     }
97d8588912SDave May   }
98d8588912SDave May 
99d8588912SDave May   for (j=0; j<nc; j++) {
100d8588912SDave May     mat = b->m[r][j];
101d8588912SDave May     if (mat) {
102d8588912SDave May       ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
103d8588912SDave May       /* Check rows match */
104d8588912SDave May       if (sM != M) {
105d8588912SDave May         SETERRQ3(((PetscObject)A)->comm,PETSC_ERR_SUP,"Inserted incompatible submatrix into Nest at (%D,%D). Submatrix must have %D cols",r,c,M );
106d8588912SDave May       }
107d8588912SDave May     }
108d8588912SDave May   }
109d8588912SDave May   PetscFunctionReturn(0);
110d8588912SDave May }
111d8588912SDave May 
112d8588912SDave May /*
113d8588912SDave May  Checks if the row/col's contain a complete set of IS's.
114d8588912SDave May  Once they are filled, the offsets are computed. This saves having to update
115d8588912SDave May  the index set entries each time we insert something new.
116d8588912SDave May  */
117d8588912SDave May #undef __FUNCT__
118d8588912SDave May #define __FUNCT__ "PETSc_MatNest_UpdateStructure"
119063a28d4SJed Brown PETSC_UNUSED
120d8588912SDave May static PetscErrorCode PETSc_MatNest_UpdateStructure(Mat A,PetscInt ridx,PetscInt cidx)
121d8588912SDave May {
122d8588912SDave May   Mat_Nest       *b = (Mat_Nest*)A->data;
123d8588912SDave May   PetscInt       m,n,i,j;
124d8588912SDave May   PetscBool      fullrow,fullcol;
125d8588912SDave May   PetscErrorCode ierr;
126d8588912SDave May 
127d8588912SDave May   PetscFunctionBegin;
128d8588912SDave May   ierr = MatGetLocalSize(b->m[ridx][cidx],&m,&n);CHKERRQ(ierr);
129d8588912SDave May   b->row_len[ridx] = m;
130d8588912SDave May   b->col_len[cidx] = n;
131d8588912SDave May 
132d8588912SDave May   fullrow = PETSC_TRUE;
133d8588912SDave May   for (i=0; i<b->nr; i++) {
134d8588912SDave May     if (b->row_len[i]==-1) { fullrow = PETSC_FALSE; break; }
135d8588912SDave May   }
136f349c1fdSJed Brown   if ( (fullrow) && (!b->isglobal.row[0]) ){
137d8588912SDave May     PetscInt cnt = 0;
138d8588912SDave May     for (i=0; i<b->nr; i++) {
139f349c1fdSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,b->row_len[i],cnt,1,&b->isglobal.row[i]);CHKERRQ(ierr);
140d8588912SDave May       cnt = cnt + b->row_len[i];
141d8588912SDave May     }
142d8588912SDave May   }
143d8588912SDave May 
144d8588912SDave May   fullcol = PETSC_TRUE;
145d8588912SDave May   for (j=0; j<b->nc; j++) {
146d8588912SDave May     if (b->col_len[j]==-1) { fullcol = PETSC_FALSE; break; }
147d8588912SDave May   }
148f349c1fdSJed Brown   if( (fullcol) && (!b->isglobal.col[0]) ){
149d8588912SDave May     PetscInt cnt = 0;
150d8588912SDave May     for (j=0; j<b->nc; j++) {
151f349c1fdSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,b->col_len[j],cnt,1,&b->isglobal.col[j]);CHKERRQ(ierr);
152d8588912SDave May       cnt = cnt + b->col_len[j];
153d8588912SDave May     }
154d8588912SDave May   }
155d8588912SDave May   PetscFunctionReturn(0);
156d8588912SDave May }
157d8588912SDave May 
158d8588912SDave May /* operations */
159d8588912SDave May #undef __FUNCT__
160d8588912SDave May #define __FUNCT__ "MatMult_Nest"
161207556f9SJed Brown static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y)
162d8588912SDave May {
163d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
164207556f9SJed Brown   Vec            *bx = bA->right,*by = bA->left;
165207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
166d8588912SDave May   PetscErrorCode ierr;
167d8588912SDave May 
168d8588912SDave May   PetscFunctionBegin;
169207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
170207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
171207556f9SJed Brown   for (i=0; i<nr; i++) {
172d8588912SDave May     ierr = VecZeroEntries(by[i]);CHKERRQ(ierr);
173207556f9SJed Brown     for (j=0; j<nc; j++) {
174207556f9SJed Brown       if (!bA->m[i][j]) continue;
175d8588912SDave May       /* y[i] <- y[i] + A[i][j] * x[j] */
176d8588912SDave May       ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr);
177d8588912SDave May     }
178d8588912SDave May   }
179207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
180207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
181d8588912SDave May   PetscFunctionReturn(0);
182d8588912SDave May }
183d8588912SDave May 
184d8588912SDave May #undef __FUNCT__
185d8588912SDave May #define __FUNCT__ "MatMultTranspose_Nest"
186207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y)
187d8588912SDave May {
188d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
189207556f9SJed Brown   Vec            *bx = bA->left,*by = bA->right;
190207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
191d8588912SDave May   PetscErrorCode ierr;
192d8588912SDave May 
193d8588912SDave May   PetscFunctionBegin;
194d8588912SDave May   if (A->symmetric) {
195d8588912SDave May     ierr = MatMult_Nest(A,x,y);CHKERRQ(ierr);
196d8588912SDave May     PetscFunctionReturn(0);
197d8588912SDave May   }
198207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
199207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
200207556f9SJed Brown   for (i=0; i<nr; i++) {
201d8588912SDave May     ierr = VecZeroEntries(by[i]);CHKERRQ(ierr);
202207556f9SJed Brown     for (j=0; j<nc; j++) {
203207556f9SJed Brown       if (!bA->m[j][i]) continue;
204d8588912SDave May       /* y[i] <- y[i] + A^T[i][j] * x[j], so we swap i,j in mat[][] */
205d8588912SDave May       ierr = MatMultTransposeAdd(bA->m[j][i],bx[j],by[i],by[i]);CHKERRQ(ierr);
206d8588912SDave May     }
207d8588912SDave May   }
208207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
209207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
210d8588912SDave May   PetscFunctionReturn(0);
211d8588912SDave May }
212d8588912SDave May 
213d8588912SDave May #undef __FUNCT__
214e2d7f03fSJed Brown #define __FUNCT__ "MatNestDestroyISList"
215e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list)
216e2d7f03fSJed Brown {
217e2d7f03fSJed Brown   PetscErrorCode ierr;
218e2d7f03fSJed Brown   IS             *lst = *list;
219e2d7f03fSJed Brown   PetscInt       i;
220e2d7f03fSJed Brown 
221e2d7f03fSJed Brown   PetscFunctionBegin;
222e2d7f03fSJed Brown   if (!lst) PetscFunctionReturn(0);
223e2d7f03fSJed Brown   for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(lst[i]);CHKERRQ(ierr);}
224e2d7f03fSJed Brown   ierr = PetscFree(lst);CHKERRQ(ierr);
225e2d7f03fSJed Brown   *list = PETSC_NULL;
226e2d7f03fSJed Brown   PetscFunctionReturn(0);
227e2d7f03fSJed Brown }
228e2d7f03fSJed Brown 
229e2d7f03fSJed Brown #undef __FUNCT__
230d8588912SDave May #define __FUNCT__ "MatDestroy_Nest"
231207556f9SJed Brown static PetscErrorCode MatDestroy_Nest(Mat A)
232d8588912SDave May {
233d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
234d8588912SDave May   PetscInt       i,j;
235d8588912SDave May   PetscErrorCode ierr;
236d8588912SDave May 
237d8588912SDave May   PetscFunctionBegin;
238d8588912SDave May   /* release the matrices and the place holders */
239e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr);
240e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr);
241e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
242e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
243d8588912SDave May 
244d8588912SDave May   ierr = PetscFree(vs->row_len);CHKERRQ(ierr);
245d8588912SDave May   ierr = PetscFree(vs->col_len);CHKERRQ(ierr);
246d8588912SDave May 
247207556f9SJed Brown   ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr);
248207556f9SJed Brown 
249d8588912SDave May   /* release the matrices and the place holders */
250d8588912SDave May   if (vs->m) {
251d8588912SDave May     for (i=0; i<vs->nr; i++) {
252d8588912SDave May       for (j=0; j<vs->nc; j++) {
253d8588912SDave May 
254d8588912SDave May         if (vs->m[i][j]) {
255d8588912SDave May           ierr = MatDestroy(vs->m[i][j]);CHKERRQ(ierr);
256d8588912SDave May           vs->m[i][j] = PETSC_NULL;
257d8588912SDave May         }
258d8588912SDave May       }
259d8588912SDave May       ierr = PetscFree( vs->m[i] );CHKERRQ(ierr);
260d8588912SDave May     }
261d8588912SDave May     ierr = PetscFree(vs->m);CHKERRQ(ierr);
262d8588912SDave May   }
263d8588912SDave May   ierr = PetscFree(vs);CHKERRQ(ierr);
264d8588912SDave May 
265c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C",   "",0);CHKERRQ(ierr);
266c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C",  "",0);CHKERRQ(ierr);
267c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C",     "",0);CHKERRQ(ierr);
268c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C",  "",0);CHKERRQ(ierr);
269c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C",   "",0);CHKERRQ(ierr);
270d8588912SDave May   PetscFunctionReturn(0);
271d8588912SDave May }
272d8588912SDave May 
273d8588912SDave May #undef __FUNCT__
274d8588912SDave May #define __FUNCT__ "MatAssemblyBegin_Nest"
275207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
276d8588912SDave May {
277d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
278d8588912SDave May   PetscInt       i,j;
279d8588912SDave May   PetscErrorCode ierr;
280d8588912SDave May 
281d8588912SDave May   PetscFunctionBegin;
282d8588912SDave May   for (i=0; i<vs->nr; i++) {
283d8588912SDave May     for (j=0; j<vs->nc; j++) {
284d8588912SDave May       if (vs->m[i][j]) { ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); }
285d8588912SDave May     }
286d8588912SDave May   }
287d8588912SDave May   PetscFunctionReturn(0);
288d8588912SDave May }
289d8588912SDave May 
290d8588912SDave May #undef __FUNCT__
291d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest"
292207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
293d8588912SDave May {
294d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
295d8588912SDave May   PetscInt       i,j;
296d8588912SDave May   PetscErrorCode ierr;
297d8588912SDave May 
298d8588912SDave May   PetscFunctionBegin;
299d8588912SDave May   for (i=0; i<vs->nr; i++) {
300d8588912SDave May     for (j=0; j<vs->nc; j++) {
301d8588912SDave May       if (vs->m[i][j]) { ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); }
302d8588912SDave May     }
303d8588912SDave May   }
304d8588912SDave May   PetscFunctionReturn(0);
305d8588912SDave May }
306d8588912SDave May 
307d8588912SDave May #undef __FUNCT__
308f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatRow"
309f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
310d8588912SDave May {
311207556f9SJed Brown   PetscErrorCode ierr;
312f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
313f349c1fdSJed Brown   PetscInt       j;
314f349c1fdSJed Brown   Mat            sub;
315d8588912SDave May 
316d8588912SDave May   PetscFunctionBegin;
317f349c1fdSJed Brown   sub = vs->m[row][row];        /* Prefer to find on the diagonal */
318f349c1fdSJed Brown   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
319f349c1fdSJed Brown   if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",row);
320207556f9SJed Brown   ierr = MatPreallocated(sub);CHKERRQ(ierr); /* Ensure that the sizes are available */
321f349c1fdSJed Brown   *B = sub;
322f349c1fdSJed Brown   PetscFunctionReturn(0);
323d8588912SDave May }
324d8588912SDave May 
325f349c1fdSJed Brown #undef __FUNCT__
326f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatCol"
327f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
328f349c1fdSJed Brown {
329207556f9SJed Brown   PetscErrorCode ierr;
330f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
331f349c1fdSJed Brown   PetscInt       i;
332f349c1fdSJed Brown   Mat            sub;
333f349c1fdSJed Brown 
334f349c1fdSJed Brown   PetscFunctionBegin;
335f349c1fdSJed Brown   sub = vs->m[col][col];        /* Prefer to find on the diagonal */
336f349c1fdSJed Brown   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
337f349c1fdSJed Brown   if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",col);
338207556f9SJed Brown   ierr = MatPreallocated(sub);CHKERRQ(ierr); /* Ensure that the sizes are available */
339f349c1fdSJed Brown   *B = sub;
340f349c1fdSJed Brown   PetscFunctionReturn(0);
341d8588912SDave May }
342d8588912SDave May 
343f349c1fdSJed Brown #undef __FUNCT__
344f349c1fdSJed Brown #define __FUNCT__ "MatNestFindIS"
345f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
346f349c1fdSJed Brown {
347f349c1fdSJed Brown   PetscErrorCode ierr;
348f349c1fdSJed Brown   PetscInt       i;
349f349c1fdSJed Brown   PetscBool      flg;
350f349c1fdSJed Brown 
351f349c1fdSJed Brown   PetscFunctionBegin;
352f349c1fdSJed Brown   PetscValidPointer(list,3);
353f349c1fdSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,4);
354f349c1fdSJed Brown   PetscValidIntPointer(found,5);
355f349c1fdSJed Brown   *found = -1;
356f349c1fdSJed Brown   for (i=0; i<n; i++) {
357207556f9SJed Brown     if (!list[i]) continue;
358f349c1fdSJed Brown     ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr);
359f349c1fdSJed Brown     if (flg) {
360f349c1fdSJed Brown       *found = i;
361f349c1fdSJed Brown       PetscFunctionReturn(0);
362f349c1fdSJed Brown     }
363f349c1fdSJed Brown   }
364f349c1fdSJed Brown   SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Could not find index set");
365f349c1fdSJed Brown   PetscFunctionReturn(0);
366f349c1fdSJed Brown }
367f349c1fdSJed Brown 
368f349c1fdSJed Brown #undef __FUNCT__
369f349c1fdSJed Brown #define __FUNCT__ "MatNestFindSubMat"
370f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
371f349c1fdSJed Brown {
372f349c1fdSJed Brown   PetscErrorCode ierr;
373f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
374f349c1fdSJed Brown   PetscInt       row,col;
375f349c1fdSJed Brown 
376f349c1fdSJed Brown   PetscFunctionBegin;
377f349c1fdSJed Brown   ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
378f349c1fdSJed Brown   ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr);
379f349c1fdSJed Brown   *B = vs->m[row][col];
380f349c1fdSJed Brown   PetscFunctionReturn(0);
381f349c1fdSJed Brown }
382f349c1fdSJed Brown 
383f349c1fdSJed Brown #undef __FUNCT__
384f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest"
385f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
386f349c1fdSJed Brown {
387f349c1fdSJed Brown   PetscErrorCode ierr;
388f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
389f349c1fdSJed Brown   Mat            sub;
390f349c1fdSJed Brown 
391f349c1fdSJed Brown   PetscFunctionBegin;
392f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr);
393f349c1fdSJed Brown   switch (reuse) {
394f349c1fdSJed Brown   case MAT_INITIAL_MATRIX:
395*7874fa86SDave May     if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); }
396f349c1fdSJed Brown     *B = sub;
397f349c1fdSJed Brown     break;
398f349c1fdSJed Brown   case MAT_REUSE_MATRIX:
399f349c1fdSJed Brown     if (sub != *B) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
400f349c1fdSJed Brown     break;
401f349c1fdSJed Brown   case MAT_IGNORE_MATRIX:       /* Nothing to do */
402f349c1fdSJed Brown     break;
403f349c1fdSJed Brown   }
404f349c1fdSJed Brown   PetscFunctionReturn(0);
405f349c1fdSJed Brown }
406f349c1fdSJed Brown 
407f349c1fdSJed Brown #undef __FUNCT__
408f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest"
409f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
410f349c1fdSJed Brown {
411f349c1fdSJed Brown   PetscErrorCode ierr;
412f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
413f349c1fdSJed Brown   Mat            sub;
414f349c1fdSJed Brown 
415f349c1fdSJed Brown   PetscFunctionBegin;
416f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
417f349c1fdSJed Brown   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
418f349c1fdSJed Brown   if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);}
419f349c1fdSJed Brown   *B = sub;
420d8588912SDave May   PetscFunctionReturn(0);
421d8588912SDave May }
422d8588912SDave May 
423d8588912SDave May #undef __FUNCT__
424d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest"
425207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
426d8588912SDave May {
427d8588912SDave May   PetscErrorCode ierr;
428f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
429f349c1fdSJed Brown   Mat            sub;
430d8588912SDave May 
431d8588912SDave May   PetscFunctionBegin;
432f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
433f349c1fdSJed Brown   if (*B != sub) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
434f349c1fdSJed Brown   if (sub) {
435f349c1fdSJed Brown     if (((PetscObject)sub)->refct <= 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
436f349c1fdSJed Brown     ierr = MatDestroy(*B);CHKERRQ(ierr);
437d8588912SDave May   }
438f349c1fdSJed Brown   *B = PETSC_NULL;
439d8588912SDave May   PetscFunctionReturn(0);
440d8588912SDave May }
441d8588912SDave May 
442d8588912SDave May #undef __FUNCT__
443*7874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest"
444*7874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
445*7874fa86SDave May {
446*7874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
447*7874fa86SDave May   Vec            *bdiag;
448*7874fa86SDave May   PetscInt       i;
449*7874fa86SDave May   PetscErrorCode ierr;
450*7874fa86SDave May 
451*7874fa86SDave May   PetscFunctionBegin;
452*7874fa86SDave May //  ierr = MatGetVecs_Nest(A,&diag,PETSC_NULL);CHKERRQ(ierr);
453*7874fa86SDave May //  ierr = VecNestGetSubVecs(diag,PETSC_NULL,&bdiag);CHKERRQ(ierr);
454*7874fa86SDave May   ierr = VecNestGetSubVecs(v,PETSC_NULL,&bdiag);CHKERRQ(ierr);
455*7874fa86SDave May   for (i=0; i<bA->nr; i++) {
456*7874fa86SDave May     if (bA->m[i][i]) {
457*7874fa86SDave May       ierr = MatGetDiagonal(bA->m[i][i],bdiag[i]);CHKERRQ(ierr);
458*7874fa86SDave May     } else {
459*7874fa86SDave May       ierr = VecSet(bdiag[i],1.0);CHKERRQ(ierr);
460*7874fa86SDave May     }
461*7874fa86SDave May   }
462*7874fa86SDave May   PetscFunctionReturn(0);
463*7874fa86SDave May }
464*7874fa86SDave May 
465*7874fa86SDave May #undef __FUNCT__
466*7874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest2"
467*7874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest2(Mat A,Vec v)
468*7874fa86SDave May {
469*7874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
470*7874fa86SDave May   Vec            diag,*bdiag;
471*7874fa86SDave May   VecScatter     *vscat;
472*7874fa86SDave May   PetscInt       i;
473*7874fa86SDave May   PetscErrorCode ierr;
474*7874fa86SDave May 
475*7874fa86SDave May   PetscFunctionBegin;
476*7874fa86SDave May   ierr = MatGetVecs_Nest(A,&diag,PETSC_NULL);CHKERRQ(ierr);
477*7874fa86SDave May   ierr = MatGetDiagonal_Nest(A,diag);CHKERRQ(ierr);
478*7874fa86SDave May 
479*7874fa86SDave May   /* scatter diag into v */
480*7874fa86SDave May   ierr = VecNestGetSubVecs(diag,PETSC_NULL,&bdiag);CHKERRQ(ierr);
481*7874fa86SDave May   ierr = PetscMalloc( sizeof(VecScatter) * bA->nr, &vscat );CHKERRQ(ierr);
482*7874fa86SDave May   for (i=0; i<bA->nr; i++) {
483*7874fa86SDave May     ierr = VecScatterCreate(v,bA->isglobal.row[i], bdiag[i],PETSC_NULL,&vscat[i]);CHKERRQ(ierr);
484*7874fa86SDave May     ierr = VecScatterBegin(vscat[i],bdiag[i],v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
485*7874fa86SDave May   }
486*7874fa86SDave May   for (i=0; i<bA->nr; i++) {
487*7874fa86SDave May     ierr = VecScatterEnd(vscat[i],bdiag[i],v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
488*7874fa86SDave May   }
489*7874fa86SDave May   for (i=0; i<bA->nr; i++) {
490*7874fa86SDave May     ierr = VecScatterDestroy(vscat[i]);CHKERRQ(ierr);
491*7874fa86SDave May   }
492*7874fa86SDave May   ierr = PetscFree(vscat);CHKERRQ(ierr);
493*7874fa86SDave May   ierr = VecDestroy(diag);CHKERRQ(ierr);
494*7874fa86SDave May   PetscFunctionReturn(0);
495*7874fa86SDave May }
496*7874fa86SDave May 
497*7874fa86SDave May #undef __FUNCT__
498*7874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest"
499*7874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
500*7874fa86SDave May {
501*7874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
502*7874fa86SDave May   Vec            *bl,*br;
503*7874fa86SDave May   PetscInt       i,j;
504*7874fa86SDave May   PetscErrorCode ierr;
505*7874fa86SDave May 
506*7874fa86SDave May   PetscFunctionBegin;
507*7874fa86SDave May   ierr = VecNestGetSubVecs(l,PETSC_NULL,&bl);CHKERRQ(ierr);
508*7874fa86SDave May   ierr = VecNestGetSubVecs(r,PETSC_NULL,&br);CHKERRQ(ierr);
509*7874fa86SDave May   for (i=0; i<bA->nr; i++) {
510*7874fa86SDave May     for (j=0; j<bA->nc; j++) {
511*7874fa86SDave May       if (bA->m[i][j]) {
512*7874fa86SDave May         ierr = MatDiagonalScale(bA->m[i][j],bl[i],br[j]);CHKERRQ(ierr);
513*7874fa86SDave May       }
514*7874fa86SDave May     }
515*7874fa86SDave May   }
516*7874fa86SDave May   PetscFunctionReturn(0);
517*7874fa86SDave May }
518*7874fa86SDave May 
519*7874fa86SDave May #undef __FUNCT__
520*7874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest2"
521*7874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest2(Mat A,Vec l,Vec r)
522*7874fa86SDave May {
523*7874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
524*7874fa86SDave May   Vec            bl,br,*ble,*bre;
525*7874fa86SDave May   VecScatter     *vscatl,*vscatr;
526*7874fa86SDave May   PetscInt       i;
527*7874fa86SDave May   PetscErrorCode ierr;
528*7874fa86SDave May 
529*7874fa86SDave May   PetscFunctionBegin;
530*7874fa86SDave May   /* scatter l,r into bl,br */
531*7874fa86SDave May   ierr = MatGetVecs_Nest(A,&bl,&br);CHKERRQ(ierr);
532*7874fa86SDave May   ierr = VecNestGetSubVecs(bl,PETSC_NULL,&ble);CHKERRQ(ierr);
533*7874fa86SDave May   ierr = VecNestGetSubVecs(br,PETSC_NULL,&bre);CHKERRQ(ierr);
534*7874fa86SDave May 
535*7874fa86SDave May   /* row */
536*7874fa86SDave May   ierr = PetscMalloc( sizeof(VecScatter) * bA->nr, &vscatl );CHKERRQ(ierr);
537*7874fa86SDave May   for (i=0; i<bA->nr; i++) {
538*7874fa86SDave May     ierr = VecScatterCreate(l,bA->isglobal.row[i], ble[i],PETSC_NULL,&vscatl[i]);CHKERRQ(ierr);
539*7874fa86SDave May     ierr = VecScatterBegin(vscatl[i],l,ble[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
540*7874fa86SDave May   }
541*7874fa86SDave May   for (i=0; i<bA->nr; i++) {
542*7874fa86SDave May     ierr = VecScatterEnd(vscatl[i],l,ble[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
543*7874fa86SDave May   }
544*7874fa86SDave May   /* col */
545*7874fa86SDave May   ierr = PetscMalloc( sizeof(VecScatter) * bA->nc, &vscatr );CHKERRQ(ierr);
546*7874fa86SDave May   for (i=0; i<bA->nc; i++) {
547*7874fa86SDave May     ierr = VecScatterCreate(r,bA->isglobal.col[i], bre[i],PETSC_NULL,&vscatr[i]);CHKERRQ(ierr);
548*7874fa86SDave May     ierr = VecScatterBegin(vscatr[i],l,bre[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
549*7874fa86SDave May   }
550*7874fa86SDave May   for (i=0; i<bA->nc; i++) {
551*7874fa86SDave May     ierr = VecScatterEnd(vscatr[i],r,bre[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
552*7874fa86SDave May   }
553*7874fa86SDave May 
554*7874fa86SDave May   ierr = MatDiagonalScale_Nest(A,bl,br);CHKERRQ(ierr);
555*7874fa86SDave May 
556*7874fa86SDave May   for (i=0; i<bA->nr; i++) {
557*7874fa86SDave May     ierr = VecScatterDestroy(vscatl[i]);CHKERRQ(ierr);
558*7874fa86SDave May   }
559*7874fa86SDave May   for (i=0; i<bA->nc; i++) {
560*7874fa86SDave May     ierr = VecScatterDestroy(vscatr[i]);CHKERRQ(ierr);
561*7874fa86SDave May   }
562*7874fa86SDave May   ierr = PetscFree(vscatl);CHKERRQ(ierr);
563*7874fa86SDave May   ierr = PetscFree(vscatr);CHKERRQ(ierr);
564*7874fa86SDave May   ierr = VecDestroy(bl);CHKERRQ(ierr);
565*7874fa86SDave May   ierr = VecDestroy(br);CHKERRQ(ierr);
566*7874fa86SDave May   PetscFunctionReturn(0);
567*7874fa86SDave May }
568*7874fa86SDave May 
569*7874fa86SDave May #undef __FUNCT__
570d8588912SDave May #define __FUNCT__ "MatGetVecs_Nest"
571207556f9SJed Brown static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left)
572d8588912SDave May {
573d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
574d8588912SDave May   Vec            *L,*R;
575d8588912SDave May   MPI_Comm       comm;
576d8588912SDave May   PetscInt       i,j;
577d8588912SDave May   PetscErrorCode ierr;
578d8588912SDave May 
579d8588912SDave May   PetscFunctionBegin;
580d8588912SDave May   comm = ((PetscObject)A)->comm;
581d8588912SDave May   if (right) {
582d8588912SDave May     /* allocate R */
583d8588912SDave May     ierr = PetscMalloc( sizeof(Vec) * bA->nc, &R );CHKERRQ(ierr);
584d8588912SDave May     /* Create the right vectors */
585d8588912SDave May     for (j=0; j<bA->nc; j++) {
586d8588912SDave May       for (i=0; i<bA->nr; i++) {
587d8588912SDave May         if (bA->m[i][j]) {
588d8588912SDave May           ierr = MatGetVecs(bA->m[i][j],&R[j],PETSC_NULL);CHKERRQ(ierr);
589d8588912SDave May           break;
590d8588912SDave May         }
591d8588912SDave May       }
592d8588912SDave May       if (i==bA->nr) {
593d8588912SDave May         /* have an empty column */
594d8588912SDave May         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
595d8588912SDave May       }
596d8588912SDave May     }
597f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
598d8588912SDave May     /* hand back control to the nest vector */
599d8588912SDave May     for (j=0; j<bA->nc; j++) {
600d8588912SDave May       ierr = VecDestroy(R[j]);CHKERRQ(ierr);
601d8588912SDave May     }
602d8588912SDave May     ierr = PetscFree(R);CHKERRQ(ierr);
603d8588912SDave May   }
604d8588912SDave May 
605d8588912SDave May   if (left) {
606d8588912SDave May     /* allocate L */
607d8588912SDave May     ierr = PetscMalloc( sizeof(Vec) * bA->nr, &L );CHKERRQ(ierr);
608d8588912SDave May     /* Create the left vectors */
609d8588912SDave May     for (i=0; i<bA->nr; i++) {
610d8588912SDave May       for (j=0; j<bA->nc; j++) {
611d8588912SDave May         if (bA->m[i][j]) {
612d8588912SDave May           ierr = MatGetVecs(bA->m[i][j],PETSC_NULL,&L[i]);CHKERRQ(ierr);
613d8588912SDave May           break;
614d8588912SDave May         }
615d8588912SDave May       }
616d8588912SDave May       if (j==bA->nc) {
617d8588912SDave May         /* have an empty row */
618d8588912SDave May         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
619d8588912SDave May       }
620d8588912SDave May     }
621d8588912SDave May 
622f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
623d8588912SDave May     for (i=0; i<bA->nr; i++) {
624d8588912SDave May       ierr = VecDestroy(L[i]);CHKERRQ(ierr);
625d8588912SDave May     }
626d8588912SDave May 
627d8588912SDave May     ierr = PetscFree(L);CHKERRQ(ierr);
628d8588912SDave May   }
629d8588912SDave May   PetscFunctionReturn(0);
630d8588912SDave May }
631d8588912SDave May 
632d8588912SDave May #undef __FUNCT__
633d8588912SDave May #define __FUNCT__ "MatView_Nest"
634207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
635d8588912SDave May {
636d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
637d8588912SDave May   PetscBool      isascii;
638d8588912SDave May   PetscInt       i,j;
639d8588912SDave May   PetscErrorCode ierr;
640d8588912SDave May 
641d8588912SDave May   PetscFunctionBegin;
642d8588912SDave May   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
643d8588912SDave May   if (isascii) {
644d8588912SDave May 
645d8588912SDave May     PetscViewerASCIIPrintf(viewer,"Matrix object: \n" );
646d8588912SDave May     PetscViewerASCIIPushTab( viewer );    /* push0 */
647d8588912SDave May     PetscViewerASCIIPrintf( viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc);
648d8588912SDave May 
649d8588912SDave May     PetscViewerASCIIPrintf(viewer,"MatNest structure: \n" );
650d8588912SDave May     for (i=0; i<bA->nr; i++) {
651d8588912SDave May       for (j=0; j<bA->nc; j++) {
652d8588912SDave May         const MatType type;
653d8588912SDave May         const char *name;
654d8588912SDave May         PetscInt NR,NC;
655d8588912SDave May         PetscBool isNest = PETSC_FALSE;
656d8588912SDave May 
657d8588912SDave May         if (!bA->m[i][j]) {
658d8588912SDave May           PetscViewerASCIIPrintf( viewer, "(%D,%D) : PETSC_NULL \n",i,j);
659d8588912SDave May           continue;
660d8588912SDave May         }
661d8588912SDave May         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
662d8588912SDave May         ierr = MatGetType( bA->m[i][j], &type );CHKERRQ(ierr);
663d8588912SDave May         name = ((PetscObject)bA->m[i][j])->prefix;
664d8588912SDave May         ierr = PetscTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
665d8588912SDave May 
666d8588912SDave May         PetscViewerASCIIPrintf( viewer, "(%D,%D) : name=\"%s\", type=%s, rows=%D, cols=%D \n",i,j,name,type,NR,NC);
667d8588912SDave May 
668d8588912SDave May         if (isNest) {
669d8588912SDave May           PetscViewerASCIIPushTab( viewer );  /* push1 */
670d8588912SDave May           ierr = MatView( bA->m[i][j], viewer );CHKERRQ(ierr);
671d8588912SDave May           PetscViewerASCIIPopTab(viewer);    /* pop1 */
672d8588912SDave May         }
673d8588912SDave May       }
674d8588912SDave May     }
675d8588912SDave May     PetscViewerASCIIPopTab(viewer);    /* pop0 */
676d8588912SDave May   }
677d8588912SDave May   PetscFunctionReturn(0);
678d8588912SDave May }
679d8588912SDave May 
680d8588912SDave May #undef __FUNCT__
681d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest"
682207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A)
683d8588912SDave May {
684d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
685d8588912SDave May   PetscInt       i,j;
686d8588912SDave May   PetscErrorCode ierr;
687d8588912SDave May 
688d8588912SDave May   PetscFunctionBegin;
689d8588912SDave May   for (i=0; i<bA->nr; i++) {
690d8588912SDave May     for (j=0; j<bA->nc; j++) {
691d8588912SDave May       if (!bA->m[i][j]) continue;
692d8588912SDave May       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
693d8588912SDave May     }
694d8588912SDave May   }
695d8588912SDave May   PetscFunctionReturn(0);
696d8588912SDave May }
697d8588912SDave May 
698d8588912SDave May #undef __FUNCT__
699d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest"
700207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
701d8588912SDave May {
702d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
703841e96a3SJed Brown   Mat            *b;
704841e96a3SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
705d8588912SDave May   PetscErrorCode ierr;
706d8588912SDave May 
707d8588912SDave May   PetscFunctionBegin;
708841e96a3SJed Brown   ierr = PetscMalloc(nr*nc*sizeof(Mat),&b);CHKERRQ(ierr);
709841e96a3SJed Brown   for (i=0; i<nr; i++) {
710841e96a3SJed Brown     for (j=0; j<nc; j++) {
711841e96a3SJed Brown       if (bA->m[i][j]) {
712841e96a3SJed Brown         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
713841e96a3SJed Brown       } else {
714841e96a3SJed Brown         b[i*nc+j] = PETSC_NULL;
715d8588912SDave May       }
716d8588912SDave May     }
717d8588912SDave May   }
718f349c1fdSJed Brown   ierr = MatCreateNest(((PetscObject)A)->comm,nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
719841e96a3SJed Brown   /* Give the new MatNest exclusive ownership */
720841e96a3SJed Brown   for (i=0; i<nr*nc; i++) {
721841e96a3SJed Brown     if (b[i]) {ierr = MatDestroy(b[i]);CHKERRQ(ierr);}
722d8588912SDave May   }
723d8588912SDave May   ierr = PetscFree(b);CHKERRQ(ierr);
724d8588912SDave May 
725841e96a3SJed Brown   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
726841e96a3SJed Brown   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
727d8588912SDave May   PetscFunctionReturn(0);
728d8588912SDave May }
729d8588912SDave May 
730d8588912SDave May /* nest api */
731d8588912SDave May EXTERN_C_BEGIN
732d8588912SDave May #undef __FUNCT__
733d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest"
734d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
735d8588912SDave May {
736d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
737d8588912SDave May   PetscFunctionBegin;
738d8588912SDave May   if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
739d8588912SDave May   if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
740d8588912SDave May   *mat = bA->m[idxm][jdxm];
741d8588912SDave May   PetscFunctionReturn(0);
742d8588912SDave May }
743d8588912SDave May EXTERN_C_END
744d8588912SDave May 
745d8588912SDave May #undef __FUNCT__
746d8588912SDave May #define __FUNCT__ "MatNestGetSubMat"
747d8588912SDave May /*@C
748d8588912SDave May  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
749d8588912SDave May 
750d8588912SDave May  Not collective
751d8588912SDave May 
752d8588912SDave May  Input Parameters:
753d8588912SDave May  .  A  - nest matrix
754d8588912SDave May  .  idxm - index of the matrix within the nest matrix
755d8588912SDave May  .  jdxm - index of the matrix within the nest matrix
756d8588912SDave May 
757d8588912SDave May  Output Parameter:
758d8588912SDave May  .  sub - matrix at index idxm,jdxm within the nest matrix
759d8588912SDave May 
760d8588912SDave May  Notes:
761d8588912SDave May 
762d8588912SDave May  Level: developer
763d8588912SDave May 
764d8588912SDave May  .seealso: MatNestGetSize(), MatNestGetSubMats()
765d8588912SDave May @*/
7667087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
767d8588912SDave May {
768699a902aSJed Brown   PetscErrorCode ierr;
769d8588912SDave May 
770d8588912SDave May   PetscFunctionBegin;
771699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
772d8588912SDave May   PetscFunctionReturn(0);
773d8588912SDave May }
774d8588912SDave May 
775d8588912SDave May EXTERN_C_BEGIN
776d8588912SDave May #undef __FUNCT__
777d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest"
778d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
779d8588912SDave May {
780d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
781d8588912SDave May   PetscFunctionBegin;
782d8588912SDave May   if (M)   { *M   = bA->nr; }
783d8588912SDave May   if (N)   { *N   = bA->nc; }
784d8588912SDave May   if (mat) { *mat = bA->m;  }
785d8588912SDave May   PetscFunctionReturn(0);
786d8588912SDave May }
787d8588912SDave May EXTERN_C_END
788d8588912SDave May 
789d8588912SDave May #undef __FUNCT__
790d8588912SDave May #define __FUNCT__ "MatNestGetSubMats"
791d8588912SDave May /*@C
792d8588912SDave May  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
793d8588912SDave May 
794d8588912SDave May  Not collective
795d8588912SDave May 
796d8588912SDave May  Input Parameters:
797d8588912SDave May  .  A  - nest matri
798d8588912SDave May 
799d8588912SDave May  Output Parameter:
800d8588912SDave May  .  M - number of rows in the nest matrix
801d8588912SDave May  .  N - number of cols in the nest matrix
802d8588912SDave May  .  mat - 2d array of matrices
803d8588912SDave May 
804d8588912SDave May  Notes:
805d8588912SDave May 
806d8588912SDave May  The user should not free the array mat.
807d8588912SDave May 
808d8588912SDave May  Level: developer
809d8588912SDave May 
810d8588912SDave May  .seealso: MatNestGetSize(), MatNestGetSubMat()
811d8588912SDave May @*/
8127087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
813d8588912SDave May {
814699a902aSJed Brown   PetscErrorCode ierr;
815d8588912SDave May 
816d8588912SDave May   PetscFunctionBegin;
817699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
818d8588912SDave May   PetscFunctionReturn(0);
819d8588912SDave May }
820d8588912SDave May 
821d8588912SDave May EXTERN_C_BEGIN
822d8588912SDave May #undef __FUNCT__
823d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest"
8247087cfbeSBarry Smith PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
825d8588912SDave May {
826d8588912SDave May   Mat_Nest  *bA = (Mat_Nest*)A->data;
827d8588912SDave May 
828d8588912SDave May   PetscFunctionBegin;
829d8588912SDave May   if (M) { *M  = bA->nr; }
830d8588912SDave May   if (N) { *N  = bA->nc; }
831d8588912SDave May   PetscFunctionReturn(0);
832d8588912SDave May }
833d8588912SDave May EXTERN_C_END
834d8588912SDave May 
835d8588912SDave May #undef __FUNCT__
836d8588912SDave May #define __FUNCT__ "MatNestGetSize"
837d8588912SDave May /*@C
838d8588912SDave May  MatNestGetSize - Returns the size of the nest matrix.
839d8588912SDave May 
840d8588912SDave May  Not collective
841d8588912SDave May 
842d8588912SDave May  Input Parameters:
843d8588912SDave May  .  A  - nest matrix
844d8588912SDave May 
845d8588912SDave May  Output Parameter:
846d8588912SDave May  .  M - number of rows in the nested mat
847d8588912SDave May  .  N - number of cols in the nested mat
848d8588912SDave May 
849d8588912SDave May  Notes:
850d8588912SDave May 
851d8588912SDave May  Level: developer
852d8588912SDave May 
853d8588912SDave May  .seealso: MatNestGetSubMat(), MatNestGetSubMats()
854d8588912SDave May @*/
8557087cfbeSBarry Smith PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
856d8588912SDave May {
857699a902aSJed Brown   PetscErrorCode ierr;
858d8588912SDave May 
859d8588912SDave May   PetscFunctionBegin;
860699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
861d8588912SDave May   PetscFunctionReturn(0);
862d8588912SDave May }
863d8588912SDave May 
864207556f9SJed Brown EXTERN_C_BEGIN
865207556f9SJed Brown #undef __FUNCT__
866207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest"
8677087cfbeSBarry Smith PetscErrorCode  MatNestSetVecType_Nest(Mat A,const VecType vtype)
868207556f9SJed Brown {
869207556f9SJed Brown   PetscErrorCode ierr;
870207556f9SJed Brown   PetscBool      flg;
871207556f9SJed Brown 
872207556f9SJed Brown   PetscFunctionBegin;
873207556f9SJed Brown   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
874207556f9SJed Brown   /* In reality, this only distinguishes VECNEST and "other" */
875207556f9SJed Brown   A->ops->getvecs       = flg ? MatGetVecs_Nest       : 0;
876*7874fa86SDave May   A->ops->getdiagonal   = flg ? MatGetDiagonal_Nest   : 0;
877*7874fa86SDave May   A->ops->diagonalscale = flg ? MatDiagonalScale_Nest : 0;
878*7874fa86SDave May 
879207556f9SJed Brown  PetscFunctionReturn(0);
880207556f9SJed Brown }
881207556f9SJed Brown EXTERN_C_END
882207556f9SJed Brown 
883207556f9SJed Brown #undef __FUNCT__
884207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType"
885207556f9SJed Brown /*@C
886207556f9SJed Brown  MatNestSetVecType - Sets the type of Vec returned by MatGetVecs()
887207556f9SJed Brown 
888207556f9SJed Brown  Not collective
889207556f9SJed Brown 
890207556f9SJed Brown  Input Parameters:
891207556f9SJed Brown +  A  - nest matrix
892207556f9SJed Brown -  vtype - type to use for creating vectors
893207556f9SJed Brown 
894207556f9SJed Brown  Notes:
895207556f9SJed Brown 
896207556f9SJed Brown  Level: developer
897207556f9SJed Brown 
898207556f9SJed Brown  .seealso: MatGetVecs()
899207556f9SJed Brown @*/
9007087cfbeSBarry Smith PetscErrorCode  MatNestSetVecType(Mat A,const VecType vtype)
901207556f9SJed Brown {
902207556f9SJed Brown   PetscErrorCode ierr;
903207556f9SJed Brown 
904207556f9SJed Brown   PetscFunctionBegin;
905207556f9SJed Brown   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,const VecType),(A,vtype));CHKERRQ(ierr);
906207556f9SJed Brown   PetscFunctionReturn(0);
907207556f9SJed Brown }
908207556f9SJed Brown 
909c8883902SJed Brown EXTERN_C_BEGIN
910d8588912SDave May #undef __FUNCT__
911c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats_Nest"
912c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
913d8588912SDave May {
914c8883902SJed Brown   Mat_Nest       *s = (Mat_Nest*)A->data;
915c8883902SJed Brown   PetscInt       i,j,m,n,M,N;
916d8588912SDave May   PetscErrorCode ierr;
917d8588912SDave May 
918d8588912SDave May   PetscFunctionBegin;
919c8883902SJed Brown   s->nr = nr;
920c8883902SJed Brown   s->nc = nc;
921d8588912SDave May 
922c8883902SJed Brown   /* Create space for submatrices */
923c8883902SJed Brown   ierr = PetscMalloc(sizeof(Mat*)*nr,&s->m);CHKERRQ(ierr);
924c8883902SJed Brown   for (i=0; i<nr; i++) {
925c8883902SJed Brown     ierr = PetscMalloc(sizeof(Mat)*nc,&s->m[i]);CHKERRQ(ierr);
926d8588912SDave May   }
927c8883902SJed Brown   for (i=0; i<nr; i++) {
928c8883902SJed Brown     for (j=0; j<nc; j++) {
929c8883902SJed Brown       s->m[i][j] = a[i*nc+j];
930c8883902SJed Brown       if (a[i*nc+j]) {
931c8883902SJed Brown         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
932d8588912SDave May       }
933d8588912SDave May     }
934d8588912SDave May   }
935d8588912SDave May 
936c8883902SJed Brown   ierr = PetscMalloc(sizeof(IS)*nr,&s->isglobal.row);CHKERRQ(ierr);
937c8883902SJed Brown   ierr = PetscMalloc(sizeof(IS)*nc,&s->isglobal.col);CHKERRQ(ierr);
938d8588912SDave May 
939c8883902SJed Brown   ierr = PetscMalloc(sizeof(PetscInt)*nr,&s->row_len);CHKERRQ(ierr);
940c8883902SJed Brown   ierr = PetscMalloc(sizeof(PetscInt)*nc,&s->col_len);CHKERRQ(ierr);
941c8883902SJed Brown   for (i=0; i<nr; i++) s->row_len[i]=-1;
942c8883902SJed Brown   for (j=0; j<nc; j++) s->col_len[j]=-1;
943d8588912SDave May 
944c8883902SJed Brown   m = n = 0;
945c8883902SJed Brown   M = N = 0;
946c8883902SJed Brown   ierr = MatNestGetSize_Private(A,PETSC_TRUE,&M,&N);CHKERRQ(ierr);
947c8883902SJed Brown   ierr = MatNestGetSize_Private(A,PETSC_FALSE,&m,&n);CHKERRQ(ierr);
948d8588912SDave May 
949c8883902SJed Brown   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
950c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
951c8883902SJed Brown   ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr);
952c8883902SJed Brown   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
953c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
954c8883902SJed Brown   ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr);
955c8883902SJed Brown 
956c8883902SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
957c8883902SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
958c8883902SJed Brown 
959c8883902SJed Brown   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
960c8883902SJed Brown   ierr = PetscMalloc2(nr,Vec,&s->left,nc,Vec,&s->right);CHKERRQ(ierr);
961c8883902SJed Brown   ierr = PetscMemzero(s->left,nr*sizeof(Vec));CHKERRQ(ierr);
962c8883902SJed Brown   ierr = PetscMemzero(s->right,nc*sizeof(Vec));CHKERRQ(ierr);
963d8588912SDave May   PetscFunctionReturn(0);
964d8588912SDave May }
965c8883902SJed Brown EXTERN_C_END
966d8588912SDave May 
967c8883902SJed Brown #undef __FUNCT__
968c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats"
969c8883902SJed Brown /*@
970c8883902SJed Brown    MatNestSetSubMats - Sets the nested submatrices
971c8883902SJed Brown 
972c8883902SJed Brown    Collective on Mat
973c8883902SJed Brown 
974c8883902SJed Brown    Input Parameter:
975c8883902SJed Brown +  N - nested matrix
976c8883902SJed Brown .  nr - number of nested row blocks
977c8883902SJed Brown .  is_row - index sets for each nested row block, or PETSC_NULL to make contiguous
978c8883902SJed Brown .  nc - number of nested column blocks
979c8883902SJed Brown .  is_col - index sets for each nested column block, or PETSC_NULL to make contiguous
980c8883902SJed Brown -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL
981c8883902SJed Brown 
982c8883902SJed Brown    Level: advanced
983c8883902SJed Brown 
984c8883902SJed Brown .seealso: MatCreateNest(), MATNEST
985c8883902SJed Brown @*/
986c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
987c8883902SJed Brown {
988c8883902SJed Brown   PetscErrorCode ierr;
989c8883902SJed Brown   PetscInt i;
990c8883902SJed Brown 
991c8883902SJed Brown   PetscFunctionBegin;
992c8883902SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
993c8883902SJed Brown   if (nr < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
994c8883902SJed Brown   if (nr && is_row) {
995c8883902SJed Brown     PetscValidPointer(is_row,3);
996c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
997c8883902SJed Brown   }
998c8883902SJed Brown   if (nc < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
999c8883902SJed Brown   if (nc && is_row) {
1000c8883902SJed Brown     PetscValidPointer(is_col,5);
1001c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
1002c8883902SJed Brown   }
1003c8883902SJed Brown   if (nr*nc) PetscValidPointer(a,6);
1004c8883902SJed Brown   ierr = PetscUseMethod(A,"MatNestSetSubMats_C",(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]),(A,nr,is_row,nc,is_col,a));CHKERRQ(ierr);
1005c8883902SJed Brown   PetscFunctionReturn(0);
1006c8883902SJed Brown }
1007d8588912SDave May 
1008d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1009d8588912SDave May /*
1010d8588912SDave May   nprocessors = NP
1011d8588912SDave May   Nest x^T = ( (g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1) )
1012d8588912SDave May        proc 0: => (g_0,h_0,)
1013d8588912SDave May        proc 1: => (g_1,h_1,)
1014d8588912SDave May        ...
1015d8588912SDave May        proc nprocs-1: => (g_NP-1,h_NP-1,)
1016d8588912SDave May 
1017d8588912SDave May             proc 0:                      proc 1:                    proc nprocs-1:
1018d8588912SDave May     is[0] = ( 0,1,2,...,nlocal(g_0)-1 )  ( 0,1,...,nlocal(g_1)-1 )  ( 0,1,...,nlocal(g_NP-1) )
1019d8588912SDave May 
1020d8588912SDave May             proc 0:
1021d8588912SDave May     is[1] = ( nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1 )
1022d8588912SDave May             proc 1:
1023d8588912SDave May     is[1] = ( nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1 )
1024d8588912SDave May 
1025d8588912SDave May             proc NP-1:
1026d8588912SDave May     is[1] = ( nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1 )
1027d8588912SDave May */
1028d8588912SDave May #undef __FUNCT__
1029d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private"
1030841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1031d8588912SDave May {
1032e2d7f03fSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
10332ae74bdbSJed Brown   PetscInt       i,j,offset,n,bs;
1034d8588912SDave May   PetscErrorCode ierr;
10352ae74bdbSJed Brown   Mat            sub;
1036d8588912SDave May 
1037d8588912SDave May   PetscFunctionBegin;
1038d8588912SDave May   if (is_row) { /* valid IS is passed in */
1039d8588912SDave May     /* refs on is[] are incremeneted */
1040e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1041d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
1042e2d7f03fSJed Brown       vs->isglobal.row[i] = is_row[i];
1043d8588912SDave May     }
10442ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
10452ae74bdbSJed Brown     offset = A->rmap->rstart;
1046e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1047f349c1fdSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
10482ae74bdbSJed Brown       ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr);
1049207556f9SJed Brown       if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
10502ae74bdbSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1051e2d7f03fSJed Brown       ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
1052e2d7f03fSJed Brown       ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
10532ae74bdbSJed Brown       offset += n;
1054d8588912SDave May     }
1055d8588912SDave May   }
1056d8588912SDave May 
1057d8588912SDave May   if (is_col) { /* valid IS is passed in */
1058d8588912SDave May     /* refs on is[] are incremeneted */
1059e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1060d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
1061e2d7f03fSJed Brown       vs->isglobal.col[j] = is_col[j];
1062d8588912SDave May     }
10632ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
10642ae74bdbSJed Brown     offset = A->cmap->rstart;
1065e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1066f349c1fdSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
10672ae74bdbSJed Brown       ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr);
1068207556f9SJed Brown       if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
10692ae74bdbSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1070e2d7f03fSJed Brown       ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
1071e2d7f03fSJed Brown       ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
10722ae74bdbSJed Brown       offset += n;
1073d8588912SDave May     }
1074d8588912SDave May   }
1075e2d7f03fSJed Brown 
1076e2d7f03fSJed Brown   /* Set up the local ISs */
1077e2d7f03fSJed Brown   ierr = PetscMalloc(vs->nr*sizeof(IS),&vs->islocal.row);CHKERRQ(ierr);
1078e2d7f03fSJed Brown   ierr = PetscMalloc(vs->nc*sizeof(IS),&vs->islocal.col);CHKERRQ(ierr);
1079e2d7f03fSJed Brown   for (i=0,offset=0; i<vs->nr; i++) {
1080e2d7f03fSJed Brown     IS                     isloc;
1081e2d7f03fSJed Brown     ISLocalToGlobalMapping rmap;
1082e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1083e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1084e2d7f03fSJed Brown     ierr = MatGetLocalToGlobalMapping(sub,&rmap,PETSC_NULL);CHKERRQ(ierr);
1085207556f9SJed Brown     if (rmap) {
1086e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1087e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
1088e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1089e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1090207556f9SJed Brown     } else {
1091207556f9SJed Brown       nlocal = 0;
1092207556f9SJed Brown       isloc  = PETSC_NULL;
1093207556f9SJed Brown     }
1094e2d7f03fSJed Brown     vs->islocal.row[i] = isloc;
1095e2d7f03fSJed Brown     offset += nlocal;
1096e2d7f03fSJed Brown   }
1097e2d7f03fSJed Brown   for (i=0,offset=0; i<vs->nr; i++) {
1098e2d7f03fSJed Brown     IS                     isloc;
1099e2d7f03fSJed Brown     ISLocalToGlobalMapping cmap;
1100e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1101e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
1102e2d7f03fSJed Brown     ierr = MatGetLocalToGlobalMapping(sub,PETSC_NULL,&cmap);CHKERRQ(ierr);
1103207556f9SJed Brown     if (cmap) {
1104e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1105e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr);
1106e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1107e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1108207556f9SJed Brown     } else {
1109207556f9SJed Brown       nlocal = 0;
1110207556f9SJed Brown       isloc  = PETSC_NULL;
1111207556f9SJed Brown     }
1112e2d7f03fSJed Brown     vs->islocal.col[i] = isloc;
1113e2d7f03fSJed Brown     offset += nlocal;
1114e2d7f03fSJed Brown   }
1115d8588912SDave May   PetscFunctionReturn(0);
1116d8588912SDave May }
1117d8588912SDave May 
1118d8588912SDave May #undef __FUNCT__
1119d8588912SDave May #define __FUNCT__ "MatCreateNest"
1120659c6bb0SJed Brown /*@
1121659c6bb0SJed Brown    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
1122659c6bb0SJed Brown 
1123659c6bb0SJed Brown    Collective on Mat
1124659c6bb0SJed Brown 
1125659c6bb0SJed Brown    Input Parameter:
1126659c6bb0SJed Brown +  comm - Communicator for the new Mat
1127659c6bb0SJed Brown .  nr - number of nested row blocks
1128659c6bb0SJed Brown .  is_row - index sets for each nested row block, or PETSC_NULL to make contiguous
1129659c6bb0SJed Brown .  nc - number of nested column blocks
1130659c6bb0SJed Brown .  is_col - index sets for each nested column block, or PETSC_NULL to make contiguous
1131659c6bb0SJed Brown -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL
1132659c6bb0SJed Brown 
1133659c6bb0SJed Brown    Output Parameter:
1134659c6bb0SJed Brown .  B - new matrix
1135659c6bb0SJed Brown 
1136659c6bb0SJed Brown    Level: advanced
1137659c6bb0SJed Brown 
1138659c6bb0SJed Brown .seealso: MatCreate(), VecCreateNest(), DMGetMatrix(), MATNEST
1139659c6bb0SJed Brown @*/
11407087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1141d8588912SDave May {
1142d8588912SDave May   Mat            A;
1143d8588912SDave May   PetscErrorCode ierr;
1144d8588912SDave May 
1145d8588912SDave May   PetscFunctionBegin;
1146c8883902SJed Brown   *B = 0;
1147d8588912SDave May   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1148c8883902SJed Brown   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
1149c8883902SJed Brown   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1150d8588912SDave May   *B = A;
1151d8588912SDave May   PetscFunctionReturn(0);
1152d8588912SDave May }
1153659c6bb0SJed Brown 
1154659c6bb0SJed Brown /*MC
1155659c6bb0SJed Brown   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
1156659c6bb0SJed Brown 
1157659c6bb0SJed Brown   Level: intermediate
1158659c6bb0SJed Brown 
1159659c6bb0SJed Brown   Notes:
1160659c6bb0SJed Brown   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1161659c6bb0SJed Brown   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1162659c6bb0SJed Brown   It is usually used with DMComposite and DMGetMatrix()
1163659c6bb0SJed Brown 
1164659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest()
1165659c6bb0SJed Brown M*/
1166c8883902SJed Brown EXTERN_C_BEGIN
1167c8883902SJed Brown #undef __FUNCT__
1168c8883902SJed Brown #define __FUNCT__ "MatCreate_Nest"
1169c8883902SJed Brown PetscErrorCode MatCreate_Nest(Mat A)
1170c8883902SJed Brown {
1171c8883902SJed Brown   Mat_Nest       *s;
1172c8883902SJed Brown   PetscErrorCode ierr;
1173c8883902SJed Brown 
1174c8883902SJed Brown   PetscFunctionBegin;
1175c8883902SJed Brown   ierr = PetscNewLog(A,Mat_Nest,&s);CHKERRQ(ierr);
1176c8883902SJed Brown   A->data         = (void*)s;
1177c8883902SJed Brown   s->nr = s->nc   = -1;
1178c8883902SJed Brown   s->m            = PETSC_NULL;
1179c8883902SJed Brown 
1180c8883902SJed Brown   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
1181c8883902SJed Brown   A->ops->mult                  = MatMult_Nest;
1182c8883902SJed Brown   A->ops->multadd               = 0;
1183c8883902SJed Brown   A->ops->multtranspose         = MatMultTranspose_Nest;
1184c8883902SJed Brown   A->ops->multtransposeadd      = 0;
1185c8883902SJed Brown   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1186c8883902SJed Brown   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1187c8883902SJed Brown   A->ops->zeroentries           = MatZeroEntries_Nest;
1188c8883902SJed Brown   A->ops->duplicate             = MatDuplicate_Nest;
1189c8883902SJed Brown   A->ops->getsubmatrix          = MatGetSubMatrix_Nest;
1190c8883902SJed Brown   A->ops->destroy               = MatDestroy_Nest;
1191c8883902SJed Brown   A->ops->view                  = MatView_Nest;
1192c8883902SJed Brown   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1193c8883902SJed Brown   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1194c8883902SJed Brown   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
1195*7874fa86SDave May   A->ops->getdiagonal           = MatGetDiagonal_Nest2; /* VECNEST version activated by calling MatNestSetVecType(A,VECNEST) */
1196*7874fa86SDave May   A->ops->diagonalscale         = MatDiagonalScale_Nest2; /* VECNEST version activated by calling MatNestSetVecType(A,VECNEST) */
1197c8883902SJed Brown 
1198c8883902SJed Brown   A->spptr        = 0;
1199c8883902SJed Brown   A->same_nonzero = PETSC_FALSE;
1200c8883902SJed Brown   A->assembled    = PETSC_FALSE;
1201c8883902SJed Brown 
1202c8883902SJed Brown   /* expose Nest api's */
1203c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C",   "MatNestGetSubMat_Nest",   MatNestGetSubMat_Nest);CHKERRQ(ierr);
1204c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C",  "MatNestGetSubMats_Nest",  MatNestGetSubMats_Nest);CHKERRQ(ierr);
1205c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C",     "MatNestGetSize_Nest",     MatNestGetSize_Nest);CHKERRQ(ierr);
1206c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C",  "MatNestSetVecType_Nest",  MatNestSetVecType_Nest);CHKERRQ(ierr);
1207c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C",  "MatNestSetSubMats_Nest",  MatNestSetSubMats_Nest);CHKERRQ(ierr);
1208c8883902SJed Brown 
1209c8883902SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
1210c8883902SJed Brown   PetscFunctionReturn(0);
1211c8883902SJed Brown }
121286e80854SHong Zhang EXTERN_C_END
1213