xref: /petsc/src/mat/impls/nest/matnest.c (revision 8188e55adec81f61b4511e7984c970fbc87f0a5c)
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[]);
57874fa86SDave May static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left);
6c8883902SJed Brown 
7d8588912SDave May /* private functions */
8d8588912SDave May #undef __FUNCT__
9*8188e55aSJed Brown #define __FUNCT__ "MatNestGetSizes_Private"
10*8188e55aSJed Brown static PetscErrorCode MatNestGetSizes_Private(Mat A,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N)
11d8588912SDave May {
12d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
13*8188e55aSJed Brown   PetscInt       i,j;
14d8588912SDave May   PetscErrorCode ierr;
15d8588912SDave May 
16d8588912SDave May   PetscFunctionBegin;
17*8188e55aSJed Brown   *m = *n = *M = *N = 0;
18*8188e55aSJed Brown   for (i=0; i<bA->nr; i++) {  /* rows */
19*8188e55aSJed Brown     PetscInt sm,sM;
20*8188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.row[i],&sm);CHKERRQ(ierr);
21*8188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.row[i],&sM);CHKERRQ(ierr);
22*8188e55aSJed Brown     *m += sm;
23*8188e55aSJed Brown     *M += sM;
24d8588912SDave May   }
25*8188e55aSJed Brown   for (j=0; j<bA->nc; j++) {  /* cols */
26*8188e55aSJed Brown     PetscInt sn,sN;
27*8188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.col[j],&sn);CHKERRQ(ierr);
28*8188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.col[j],&sN);CHKERRQ(ierr);
29*8188e55aSJed Brown     *n += sn;
30*8188e55aSJed Brown     *N += sN;
31d8588912SDave May   }
32d8588912SDave May   PetscFunctionReturn(0);
33d8588912SDave May }
34d8588912SDave May 
35d8588912SDave May #undef __FUNCT__
36d8588912SDave May #define __FUNCT__ "PETSc_MatNest_CheckMatVecCompatibility2"
37207556f9SJed Brown PETSC_UNUSED
38d8588912SDave May static PetscErrorCode PETSc_MatNest_CheckMatVecCompatibility2(Mat A,Vec x,Vec y)
39d8588912SDave May {
40d8588912SDave May   PetscBool      isANest, isxNest, isyNest;
41d8588912SDave May   PetscErrorCode ierr;
42d8588912SDave May 
43d8588912SDave May   PetscFunctionBegin;
44d8588912SDave May   isANest = isxNest = isyNest = PETSC_FALSE;
45d8588912SDave May   ierr = PetscTypeCompare( (PetscObject)A, MATNEST, &isANest );CHKERRQ(ierr);
46d8588912SDave May   ierr = PetscTypeCompare( (PetscObject)x, VECNEST, &isxNest );CHKERRQ(ierr);
47d8588912SDave May   ierr = PetscTypeCompare( (PetscObject)y, VECNEST, &isyNest );CHKERRQ(ierr);
48d8588912SDave May 
49d8588912SDave May   if (isANest && isxNest && isyNest) {
50d8588912SDave May     PetscFunctionReturn(0);
51d8588912SDave May   } else {
52d8588912SDave May     SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_SUP, "Operation only valid for Mat(Nest)-Vec(Nest) combinations");
53d8588912SDave May     PetscFunctionReturn(0);
54d8588912SDave May   }
55d8588912SDave May   PetscFunctionReturn(0);
56d8588912SDave May }
57d8588912SDave May 
58d8588912SDave May /*
59d8588912SDave May  This function is called every time we insert a sub matrix the Nest.
60d8588912SDave May  It ensures that every Nest along row r, and col c has the same dimensions
61d8588912SDave May  as the submat being inserted.
62d8588912SDave May */
63d8588912SDave May #undef __FUNCT__
64d8588912SDave May #define __FUNCT__ "PETSc_MatNest_CheckConsistency"
65063a28d4SJed Brown PETSC_UNUSED
66d8588912SDave May static PetscErrorCode PETSc_MatNest_CheckConsistency(Mat A,Mat submat,PetscInt r,PetscInt c)
67d8588912SDave May {
68d8588912SDave May   Mat_Nest       *b = (Mat_Nest*)A->data;
69d8588912SDave May   PetscInt       i,j;
70d8588912SDave May   PetscInt       nr,nc;
71d8588912SDave May   PetscInt       sM,sN,M,N;
72d8588912SDave May   Mat            mat;
73d8588912SDave May   PetscErrorCode ierr;
74d8588912SDave May 
75d8588912SDave May   PetscFunctionBegin;
76d8588912SDave May   nr = b->nr;
77d8588912SDave May   nc = b->nc;
78d8588912SDave May   ierr = MatGetSize(submat,&sM,&sN);CHKERRQ(ierr);
79d8588912SDave May   for (i=0; i<nr; i++) {
80d8588912SDave May     mat = b->m[i][c];
81d8588912SDave May     if (mat) {
82d8588912SDave May       ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
83d8588912SDave May       /* Check columns match */
84d8588912SDave May       if (sN != N) {
85d8588912SDave May         SETERRQ3(((PetscObject)A)->comm,PETSC_ERR_SUP,"Inserted incompatible submatrix into Nest at (%D,%D). Submatrix must have %D rows",r,c,N );
86d8588912SDave May       }
87d8588912SDave May     }
88d8588912SDave May   }
89d8588912SDave May 
90d8588912SDave May   for (j=0; j<nc; j++) {
91d8588912SDave May     mat = b->m[r][j];
92d8588912SDave May     if (mat) {
93d8588912SDave May       ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
94d8588912SDave May       /* Check rows match */
95d8588912SDave May       if (sM != M) {
96d8588912SDave May         SETERRQ3(((PetscObject)A)->comm,PETSC_ERR_SUP,"Inserted incompatible submatrix into Nest at (%D,%D). Submatrix must have %D cols",r,c,M );
97d8588912SDave May       }
98d8588912SDave May     }
99d8588912SDave May   }
100d8588912SDave May   PetscFunctionReturn(0);
101d8588912SDave May }
102d8588912SDave May 
103d8588912SDave May /*
104d8588912SDave May  Checks if the row/col's contain a complete set of IS's.
105d8588912SDave May  Once they are filled, the offsets are computed. This saves having to update
106d8588912SDave May  the index set entries each time we insert something new.
107d8588912SDave May  */
108d8588912SDave May #undef __FUNCT__
109d8588912SDave May #define __FUNCT__ "PETSc_MatNest_UpdateStructure"
110063a28d4SJed Brown PETSC_UNUSED
111d8588912SDave May static PetscErrorCode PETSc_MatNest_UpdateStructure(Mat A,PetscInt ridx,PetscInt cidx)
112d8588912SDave May {
113d8588912SDave May   Mat_Nest       *b = (Mat_Nest*)A->data;
114d8588912SDave May   PetscInt       m,n,i,j;
115d8588912SDave May   PetscBool      fullrow,fullcol;
116d8588912SDave May   PetscErrorCode ierr;
117d8588912SDave May 
118d8588912SDave May   PetscFunctionBegin;
119d8588912SDave May   ierr = MatGetLocalSize(b->m[ridx][cidx],&m,&n);CHKERRQ(ierr);
120d8588912SDave May   b->row_len[ridx] = m;
121d8588912SDave May   b->col_len[cidx] = n;
122d8588912SDave May 
123d8588912SDave May   fullrow = PETSC_TRUE;
124d8588912SDave May   for (i=0; i<b->nr; i++) {
125d8588912SDave May     if (b->row_len[i]==-1) { fullrow = PETSC_FALSE; break; }
126d8588912SDave May   }
127f349c1fdSJed Brown   if ( (fullrow) && (!b->isglobal.row[0]) ){
128d8588912SDave May     PetscInt cnt = 0;
129d8588912SDave May     for (i=0; i<b->nr; i++) {
130f349c1fdSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,b->row_len[i],cnt,1,&b->isglobal.row[i]);CHKERRQ(ierr);
131d8588912SDave May       cnt = cnt + b->row_len[i];
132d8588912SDave May     }
133d8588912SDave May   }
134d8588912SDave May 
135d8588912SDave May   fullcol = PETSC_TRUE;
136d8588912SDave May   for (j=0; j<b->nc; j++) {
137d8588912SDave May     if (b->col_len[j]==-1) { fullcol = PETSC_FALSE; break; }
138d8588912SDave May   }
139f349c1fdSJed Brown   if( (fullcol) && (!b->isglobal.col[0]) ){
140d8588912SDave May     PetscInt cnt = 0;
141d8588912SDave May     for (j=0; j<b->nc; j++) {
142f349c1fdSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,b->col_len[j],cnt,1,&b->isglobal.col[j]);CHKERRQ(ierr);
143d8588912SDave May       cnt = cnt + b->col_len[j];
144d8588912SDave May     }
145d8588912SDave May   }
146d8588912SDave May   PetscFunctionReturn(0);
147d8588912SDave May }
148d8588912SDave May 
149d8588912SDave May /* operations */
150d8588912SDave May #undef __FUNCT__
151d8588912SDave May #define __FUNCT__ "MatMult_Nest"
152207556f9SJed Brown static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y)
153d8588912SDave May {
154d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
155207556f9SJed Brown   Vec            *bx = bA->right,*by = bA->left;
156207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
157d8588912SDave May   PetscErrorCode ierr;
158d8588912SDave May 
159d8588912SDave May   PetscFunctionBegin;
160207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
161207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
162207556f9SJed Brown   for (i=0; i<nr; i++) {
163d8588912SDave May     ierr = VecZeroEntries(by[i]);CHKERRQ(ierr);
164207556f9SJed Brown     for (j=0; j<nc; j++) {
165207556f9SJed Brown       if (!bA->m[i][j]) continue;
166d8588912SDave May       /* y[i] <- y[i] + A[i][j] * x[j] */
167d8588912SDave May       ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr);
168d8588912SDave May     }
169d8588912SDave May   }
170207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
171207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
172d8588912SDave May   PetscFunctionReturn(0);
173d8588912SDave May }
174d8588912SDave May 
175d8588912SDave May #undef __FUNCT__
176d8588912SDave May #define __FUNCT__ "MatMultTranspose_Nest"
177207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y)
178d8588912SDave May {
179d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
180207556f9SJed Brown   Vec            *bx = bA->left,*by = bA->right;
181207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
182d8588912SDave May   PetscErrorCode ierr;
183d8588912SDave May 
184d8588912SDave May   PetscFunctionBegin;
185d8588912SDave May   if (A->symmetric) {
186d8588912SDave May     ierr = MatMult_Nest(A,x,y);CHKERRQ(ierr);
187d8588912SDave May     PetscFunctionReturn(0);
188d8588912SDave May   }
189207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
190207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
191207556f9SJed Brown   for (i=0; i<nr; i++) {
192d8588912SDave May     ierr = VecZeroEntries(by[i]);CHKERRQ(ierr);
193207556f9SJed Brown     for (j=0; j<nc; j++) {
194207556f9SJed Brown       if (!bA->m[j][i]) continue;
195d8588912SDave May       /* y[i] <- y[i] + A^T[i][j] * x[j], so we swap i,j in mat[][] */
196d8588912SDave May       ierr = MatMultTransposeAdd(bA->m[j][i],bx[j],by[i],by[i]);CHKERRQ(ierr);
197d8588912SDave May     }
198d8588912SDave May   }
199207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
200207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
201d8588912SDave May   PetscFunctionReturn(0);
202d8588912SDave May }
203d8588912SDave May 
204d8588912SDave May #undef __FUNCT__
205e2d7f03fSJed Brown #define __FUNCT__ "MatNestDestroyISList"
206e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list)
207e2d7f03fSJed Brown {
208e2d7f03fSJed Brown   PetscErrorCode ierr;
209e2d7f03fSJed Brown   IS             *lst = *list;
210e2d7f03fSJed Brown   PetscInt       i;
211e2d7f03fSJed Brown 
212e2d7f03fSJed Brown   PetscFunctionBegin;
213e2d7f03fSJed Brown   if (!lst) PetscFunctionReturn(0);
2146bf464f9SBarry Smith   for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);}
215e2d7f03fSJed Brown   ierr = PetscFree(lst);CHKERRQ(ierr);
216e2d7f03fSJed Brown   *list = PETSC_NULL;
217e2d7f03fSJed Brown   PetscFunctionReturn(0);
218e2d7f03fSJed Brown }
219e2d7f03fSJed Brown 
220e2d7f03fSJed Brown #undef __FUNCT__
221d8588912SDave May #define __FUNCT__ "MatDestroy_Nest"
222207556f9SJed Brown static PetscErrorCode MatDestroy_Nest(Mat A)
223d8588912SDave May {
224d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
225d8588912SDave May   PetscInt       i,j;
226d8588912SDave May   PetscErrorCode ierr;
227d8588912SDave May 
228d8588912SDave May   PetscFunctionBegin;
229d8588912SDave May   /* release the matrices and the place holders */
230e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr);
231e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr);
232e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
233e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
234d8588912SDave May 
235d8588912SDave May   ierr = PetscFree(vs->row_len);CHKERRQ(ierr);
236d8588912SDave May   ierr = PetscFree(vs->col_len);CHKERRQ(ierr);
237d8588912SDave May 
238207556f9SJed Brown   ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr);
239207556f9SJed Brown 
240d8588912SDave May   /* release the matrices and the place holders */
241d8588912SDave May   if (vs->m) {
242d8588912SDave May     for (i=0; i<vs->nr; i++) {
243d8588912SDave May       for (j=0; j<vs->nc; j++) {
2446bf464f9SBarry Smith         ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr);
245d8588912SDave May       }
246d8588912SDave May       ierr = PetscFree( vs->m[i] );CHKERRQ(ierr);
247d8588912SDave May     }
248d8588912SDave May     ierr = PetscFree(vs->m);CHKERRQ(ierr);
249d8588912SDave May   }
250bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
251d8588912SDave May 
252c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C",   "",0);CHKERRQ(ierr);
253c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C",  "",0);CHKERRQ(ierr);
254c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C",     "",0);CHKERRQ(ierr);
255c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C",  "",0);CHKERRQ(ierr);
256c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C",   "",0);CHKERRQ(ierr);
257d8588912SDave May   PetscFunctionReturn(0);
258d8588912SDave May }
259d8588912SDave May 
260d8588912SDave May #undef __FUNCT__
261d8588912SDave May #define __FUNCT__ "MatAssemblyBegin_Nest"
262207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
263d8588912SDave May {
264d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
265d8588912SDave May   PetscInt       i,j;
266d8588912SDave May   PetscErrorCode ierr;
267d8588912SDave May 
268d8588912SDave May   PetscFunctionBegin;
269d8588912SDave May   for (i=0; i<vs->nr; i++) {
270d8588912SDave May     for (j=0; j<vs->nc; j++) {
271d8588912SDave May       if (vs->m[i][j]) { ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); }
272d8588912SDave May     }
273d8588912SDave May   }
274d8588912SDave May   PetscFunctionReturn(0);
275d8588912SDave May }
276d8588912SDave May 
277d8588912SDave May #undef __FUNCT__
278d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest"
279207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
280d8588912SDave May {
281d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
282d8588912SDave May   PetscInt       i,j;
283d8588912SDave May   PetscErrorCode ierr;
284d8588912SDave May 
285d8588912SDave May   PetscFunctionBegin;
286d8588912SDave May   for (i=0; i<vs->nr; i++) {
287d8588912SDave May     for (j=0; j<vs->nc; j++) {
288d8588912SDave May       if (vs->m[i][j]) { ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); }
289d8588912SDave May     }
290d8588912SDave May   }
291d8588912SDave May   PetscFunctionReturn(0);
292d8588912SDave May }
293d8588912SDave May 
294d8588912SDave May #undef __FUNCT__
295f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatRow"
296f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
297d8588912SDave May {
298207556f9SJed Brown   PetscErrorCode ierr;
299f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
300f349c1fdSJed Brown   PetscInt       j;
301f349c1fdSJed Brown   Mat            sub;
302d8588912SDave May 
303d8588912SDave May   PetscFunctionBegin;
304*8188e55aSJed Brown   sub = (row < vs->nc) ? vs->m[row][row] : PETSC_NULL; /* Prefer to find on the diagonal */
305f349c1fdSJed Brown   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
306*8188e55aSJed Brown   if (sub) {ierr = MatPreallocated(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */
307f349c1fdSJed Brown   *B = sub;
308f349c1fdSJed Brown   PetscFunctionReturn(0);
309d8588912SDave May }
310d8588912SDave May 
311f349c1fdSJed Brown #undef __FUNCT__
312f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatCol"
313f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
314f349c1fdSJed Brown {
315207556f9SJed Brown   PetscErrorCode ierr;
316f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
317f349c1fdSJed Brown   PetscInt       i;
318f349c1fdSJed Brown   Mat            sub;
319f349c1fdSJed Brown 
320f349c1fdSJed Brown   PetscFunctionBegin;
321*8188e55aSJed Brown   sub = (col < vs->nr) ? vs->m[col][col] : PETSC_NULL; /* Prefer to find on the diagonal */
322f349c1fdSJed Brown   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
323*8188e55aSJed Brown   if (sub) {ierr = MatPreallocated(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */
324f349c1fdSJed Brown   *B = sub;
325f349c1fdSJed Brown   PetscFunctionReturn(0);
326d8588912SDave May }
327d8588912SDave May 
328f349c1fdSJed Brown #undef __FUNCT__
329f349c1fdSJed Brown #define __FUNCT__ "MatNestFindIS"
330f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
331f349c1fdSJed Brown {
332f349c1fdSJed Brown   PetscErrorCode ierr;
333f349c1fdSJed Brown   PetscInt       i;
334f349c1fdSJed Brown   PetscBool      flg;
335f349c1fdSJed Brown 
336f349c1fdSJed Brown   PetscFunctionBegin;
337f349c1fdSJed Brown   PetscValidPointer(list,3);
338f349c1fdSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,4);
339f349c1fdSJed Brown   PetscValidIntPointer(found,5);
340f349c1fdSJed Brown   *found = -1;
341f349c1fdSJed Brown   for (i=0; i<n; i++) {
342207556f9SJed Brown     if (!list[i]) continue;
343f349c1fdSJed Brown     ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr);
344f349c1fdSJed Brown     if (flg) {
345f349c1fdSJed Brown       *found = i;
346f349c1fdSJed Brown       PetscFunctionReturn(0);
347f349c1fdSJed Brown     }
348f349c1fdSJed Brown   }
349f349c1fdSJed Brown   SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Could not find index set");
350f349c1fdSJed Brown   PetscFunctionReturn(0);
351f349c1fdSJed Brown }
352f349c1fdSJed Brown 
353f349c1fdSJed Brown #undef __FUNCT__
354*8188e55aSJed Brown #define __FUNCT__ "MatNestGetRow"
355*8188e55aSJed Brown /* Get a block row as a new MatNest */
356*8188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B)
357*8188e55aSJed Brown {
358*8188e55aSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
359*8188e55aSJed Brown   Mat            C;
360*8188e55aSJed Brown   char           keyname[256];
361*8188e55aSJed Brown   PetscErrorCode ierr;
362*8188e55aSJed Brown 
363*8188e55aSJed Brown   PetscFunctionBegin;
364*8188e55aSJed Brown   *B = PETSC_NULL;
365*8188e55aSJed Brown   ierr = PetscSNPrintf(keyname,sizeof keyname,"NestRow_%D",row);CHKERRQ(ierr);
366*8188e55aSJed Brown   ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr);
367*8188e55aSJed Brown   if (*B) PetscFunctionReturn(0);
368*8188e55aSJed Brown 
369*8188e55aSJed Brown   ierr = MatCreateNest(((PetscObject)A)->comm,1,PETSC_NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr);
370*8188e55aSJed Brown   (*B)->assembled = A->assembled;
371*8188e55aSJed Brown   ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr);
372*8188e55aSJed Brown   ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */
373*8188e55aSJed Brown   PetscFunctionReturn(0);
374*8188e55aSJed Brown }
375*8188e55aSJed Brown 
376*8188e55aSJed Brown #undef __FUNCT__
377f349c1fdSJed Brown #define __FUNCT__ "MatNestFindSubMat"
378f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
379f349c1fdSJed Brown {
380f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
381*8188e55aSJed Brown   PetscErrorCode ierr;
382*8188e55aSJed Brown   PetscInt       row,col,i;
383*8188e55aSJed Brown   IS             *iscopy;
384*8188e55aSJed Brown   Mat            *matcopy;
385*8188e55aSJed Brown   PetscBool      same,isFullCol;
386f349c1fdSJed Brown 
387f349c1fdSJed Brown   PetscFunctionBegin;
388*8188e55aSJed Brown   /* Check if full column space. This is a hack */
389*8188e55aSJed Brown   isFullCol = PETSC_FALSE;
390*8188e55aSJed Brown   ierr = PetscTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr);
391*8188e55aSJed Brown   if (same) {
392*8188e55aSJed Brown     PetscInt n,first,step;
393*8188e55aSJed Brown     ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
394*8188e55aSJed Brown     ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
395*8188e55aSJed Brown     if (A->cmap->n == n && A->cmap->rstart == first && step == 1) isFullCol = PETSC_TRUE;
396*8188e55aSJed Brown   }
397*8188e55aSJed Brown 
398*8188e55aSJed Brown   if (isFullCol) {
399*8188e55aSJed Brown     PetscInt row;
400*8188e55aSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
401*8188e55aSJed Brown     ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr);
402*8188e55aSJed Brown   } else {
403f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
404f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr);
405f349c1fdSJed Brown     *B = vs->m[row][col];
406*8188e55aSJed Brown   }
407f349c1fdSJed Brown   PetscFunctionReturn(0);
408f349c1fdSJed Brown }
409f349c1fdSJed Brown 
410f349c1fdSJed Brown #undef __FUNCT__
411f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest"
412f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
413f349c1fdSJed Brown {
414f349c1fdSJed Brown   PetscErrorCode ierr;
415f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
416f349c1fdSJed Brown   Mat            sub;
417f349c1fdSJed Brown 
418f349c1fdSJed Brown   PetscFunctionBegin;
419f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr);
420f349c1fdSJed Brown   switch (reuse) {
421f349c1fdSJed Brown   case MAT_INITIAL_MATRIX:
4227874fa86SDave May     if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); }
423f349c1fdSJed Brown     *B = sub;
424f349c1fdSJed Brown     break;
425f349c1fdSJed Brown   case MAT_REUSE_MATRIX:
426f349c1fdSJed Brown     if (sub != *B) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
427f349c1fdSJed Brown     break;
428f349c1fdSJed Brown   case MAT_IGNORE_MATRIX:       /* Nothing to do */
429f349c1fdSJed Brown     break;
430f349c1fdSJed Brown   }
431f349c1fdSJed Brown   PetscFunctionReturn(0);
432f349c1fdSJed Brown }
433f349c1fdSJed Brown 
434f349c1fdSJed Brown #undef __FUNCT__
435f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest"
436f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
437f349c1fdSJed Brown {
438f349c1fdSJed Brown   PetscErrorCode ierr;
439f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
440f349c1fdSJed Brown   Mat            sub;
441f349c1fdSJed Brown 
442f349c1fdSJed Brown   PetscFunctionBegin;
443f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
444f349c1fdSJed Brown   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
445f349c1fdSJed Brown   if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);}
446f349c1fdSJed Brown   *B = sub;
447d8588912SDave May   PetscFunctionReturn(0);
448d8588912SDave May }
449d8588912SDave May 
450d8588912SDave May #undef __FUNCT__
451d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest"
452207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
453d8588912SDave May {
454d8588912SDave May   PetscErrorCode ierr;
455f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
456f349c1fdSJed Brown   Mat            sub;
457d8588912SDave May 
458d8588912SDave May   PetscFunctionBegin;
459f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
460f349c1fdSJed Brown   if (*B != sub) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
461f349c1fdSJed Brown   if (sub) {
462f349c1fdSJed Brown     if (((PetscObject)sub)->refct <= 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
4636bf464f9SBarry Smith     ierr = MatDestroy(B);CHKERRQ(ierr);
464d8588912SDave May   }
465d8588912SDave May   PetscFunctionReturn(0);
466d8588912SDave May }
467d8588912SDave May 
468d8588912SDave May #undef __FUNCT__
4697874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest"
4707874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
4717874fa86SDave May {
4727874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
4737874fa86SDave May   Vec            *bdiag;
4747874fa86SDave May   PetscInt       i;
4757874fa86SDave May   PetscErrorCode ierr;
4767874fa86SDave May 
4777874fa86SDave May   PetscFunctionBegin;
478a0769a71SSatish Balay   /*  ierr = MatGetVecs_Nest(A,&diag,PETSC_NULL);CHKERRQ(ierr); */
479a0769a71SSatish Balay   /*  ierr = VecNestGetSubVecs(diag,PETSC_NULL,&bdiag);CHKERRQ(ierr); */
4807874fa86SDave May   ierr = VecNestGetSubVecs(v,PETSC_NULL,&bdiag);CHKERRQ(ierr);
4817874fa86SDave May   for (i=0; i<bA->nr; i++) {
4827874fa86SDave May     if (bA->m[i][i]) {
4837874fa86SDave May       ierr = MatGetDiagonal(bA->m[i][i],bdiag[i]);CHKERRQ(ierr);
4847874fa86SDave May     } else {
4857874fa86SDave May       ierr = VecSet(bdiag[i],1.0);CHKERRQ(ierr);
4867874fa86SDave May     }
4877874fa86SDave May   }
4887874fa86SDave May   PetscFunctionReturn(0);
4897874fa86SDave May }
4907874fa86SDave May 
4917874fa86SDave May #undef __FUNCT__
4927874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest2"
4937874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest2(Mat A,Vec v)
4947874fa86SDave May {
4957874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
4967874fa86SDave May   Vec            diag,*bdiag;
4977874fa86SDave May   VecScatter     *vscat;
4987874fa86SDave May   PetscInt       i;
4997874fa86SDave May   PetscErrorCode ierr;
5007874fa86SDave May 
5017874fa86SDave May   PetscFunctionBegin;
5027874fa86SDave May   ierr = MatGetVecs_Nest(A,&diag,PETSC_NULL);CHKERRQ(ierr);
5037874fa86SDave May   ierr = MatGetDiagonal_Nest(A,diag);CHKERRQ(ierr);
5047874fa86SDave May 
5057874fa86SDave May   /* scatter diag into v */
5067874fa86SDave May   ierr = VecNestGetSubVecs(diag,PETSC_NULL,&bdiag);CHKERRQ(ierr);
5077874fa86SDave May   ierr = PetscMalloc( sizeof(VecScatter) * bA->nr, &vscat );CHKERRQ(ierr);
5087874fa86SDave May   for (i=0; i<bA->nr; i++) {
5097874fa86SDave May     ierr = VecScatterCreate(v,bA->isglobal.row[i], bdiag[i],PETSC_NULL,&vscat[i]);CHKERRQ(ierr);
5107874fa86SDave May     ierr = VecScatterBegin(vscat[i],bdiag[i],v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5117874fa86SDave May   }
5127874fa86SDave May   for (i=0; i<bA->nr; i++) {
5137874fa86SDave May     ierr = VecScatterEnd(vscat[i],bdiag[i],v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5147874fa86SDave May   }
5157874fa86SDave May   for (i=0; i<bA->nr; i++) {
5166bf464f9SBarry Smith     ierr = VecScatterDestroy(&vscat[i]);CHKERRQ(ierr);
5177874fa86SDave May   }
5187874fa86SDave May   ierr = PetscFree(vscat);CHKERRQ(ierr);
5196bf464f9SBarry Smith   ierr = VecDestroy(&diag);CHKERRQ(ierr);
5207874fa86SDave May   PetscFunctionReturn(0);
5217874fa86SDave May }
5227874fa86SDave May 
5237874fa86SDave May #undef __FUNCT__
5247874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest"
5257874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
5267874fa86SDave May {
5277874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
5287874fa86SDave May   Vec            *bl,*br;
5297874fa86SDave May   PetscInt       i,j;
5307874fa86SDave May   PetscErrorCode ierr;
5317874fa86SDave May 
5327874fa86SDave May   PetscFunctionBegin;
5337874fa86SDave May   ierr = VecNestGetSubVecs(l,PETSC_NULL,&bl);CHKERRQ(ierr);
5347874fa86SDave May   ierr = VecNestGetSubVecs(r,PETSC_NULL,&br);CHKERRQ(ierr);
5357874fa86SDave May   for (i=0; i<bA->nr; i++) {
5367874fa86SDave May     for (j=0; j<bA->nc; j++) {
5377874fa86SDave May       if (bA->m[i][j]) {
5387874fa86SDave May         ierr = MatDiagonalScale(bA->m[i][j],bl[i],br[j]);CHKERRQ(ierr);
5397874fa86SDave May       }
5407874fa86SDave May     }
5417874fa86SDave May   }
5427874fa86SDave May   PetscFunctionReturn(0);
5437874fa86SDave May }
5447874fa86SDave May 
5457874fa86SDave May #undef __FUNCT__
5467874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest2"
5477874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest2(Mat A,Vec l,Vec r)
5487874fa86SDave May {
5497874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
5507874fa86SDave May   Vec            bl,br,*ble,*bre;
5517874fa86SDave May   VecScatter     *vscatl,*vscatr;
5527874fa86SDave May   PetscInt       i;
5537874fa86SDave May   PetscErrorCode ierr;
5547874fa86SDave May 
5557874fa86SDave May   PetscFunctionBegin;
5567874fa86SDave May   /* scatter l,r into bl,br */
5577874fa86SDave May   ierr = MatGetVecs_Nest(A,&bl,&br);CHKERRQ(ierr);
5587874fa86SDave May   ierr = VecNestGetSubVecs(bl,PETSC_NULL,&ble);CHKERRQ(ierr);
5597874fa86SDave May   ierr = VecNestGetSubVecs(br,PETSC_NULL,&bre);CHKERRQ(ierr);
5607874fa86SDave May 
5617874fa86SDave May   /* row */
5627874fa86SDave May   ierr = PetscMalloc( sizeof(VecScatter) * bA->nr, &vscatl );CHKERRQ(ierr);
5637874fa86SDave May   for (i=0; i<bA->nr; i++) {
5647874fa86SDave May     ierr = VecScatterCreate(l,bA->isglobal.row[i], ble[i],PETSC_NULL,&vscatl[i]);CHKERRQ(ierr);
5657874fa86SDave May     ierr = VecScatterBegin(vscatl[i],l,ble[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5667874fa86SDave May   }
5677874fa86SDave May   for (i=0; i<bA->nr; i++) {
5687874fa86SDave May     ierr = VecScatterEnd(vscatl[i],l,ble[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5697874fa86SDave May   }
5707874fa86SDave May   /* col */
5717874fa86SDave May   ierr = PetscMalloc( sizeof(VecScatter) * bA->nc, &vscatr );CHKERRQ(ierr);
5727874fa86SDave May   for (i=0; i<bA->nc; i++) {
5737874fa86SDave May     ierr = VecScatterCreate(r,bA->isglobal.col[i], bre[i],PETSC_NULL,&vscatr[i]);CHKERRQ(ierr);
5747874fa86SDave May     ierr = VecScatterBegin(vscatr[i],l,bre[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5757874fa86SDave May   }
5767874fa86SDave May   for (i=0; i<bA->nc; i++) {
5777874fa86SDave May     ierr = VecScatterEnd(vscatr[i],r,bre[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5787874fa86SDave May   }
5797874fa86SDave May 
5807874fa86SDave May   ierr = MatDiagonalScale_Nest(A,bl,br);CHKERRQ(ierr);
5817874fa86SDave May 
5827874fa86SDave May   for (i=0; i<bA->nr; i++) {
5836bf464f9SBarry Smith     ierr = VecScatterDestroy(&vscatl[i]);CHKERRQ(ierr);
5847874fa86SDave May   }
5857874fa86SDave May   for (i=0; i<bA->nc; i++) {
5866bf464f9SBarry Smith     ierr = VecScatterDestroy(&vscatr[i]);CHKERRQ(ierr);
5877874fa86SDave May   }
5887874fa86SDave May   ierr = PetscFree(vscatl);CHKERRQ(ierr);
5897874fa86SDave May   ierr = PetscFree(vscatr);CHKERRQ(ierr);
5906bf464f9SBarry Smith   ierr = VecDestroy(&bl);CHKERRQ(ierr);
5916bf464f9SBarry Smith   ierr = VecDestroy(&br);CHKERRQ(ierr);
5927874fa86SDave May   PetscFunctionReturn(0);
5937874fa86SDave May }
5947874fa86SDave May 
5957874fa86SDave May #undef __FUNCT__
596d8588912SDave May #define __FUNCT__ "MatGetVecs_Nest"
597207556f9SJed Brown static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left)
598d8588912SDave May {
599d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
600d8588912SDave May   Vec            *L,*R;
601d8588912SDave May   MPI_Comm       comm;
602d8588912SDave May   PetscInt       i,j;
603d8588912SDave May   PetscErrorCode ierr;
604d8588912SDave May 
605d8588912SDave May   PetscFunctionBegin;
606d8588912SDave May   comm = ((PetscObject)A)->comm;
607d8588912SDave May   if (right) {
608d8588912SDave May     /* allocate R */
609d8588912SDave May     ierr = PetscMalloc( sizeof(Vec) * bA->nc, &R );CHKERRQ(ierr);
610d8588912SDave May     /* Create the right vectors */
611d8588912SDave May     for (j=0; j<bA->nc; j++) {
612d8588912SDave May       for (i=0; i<bA->nr; i++) {
613d8588912SDave May         if (bA->m[i][j]) {
614d8588912SDave May           ierr = MatGetVecs(bA->m[i][j],&R[j],PETSC_NULL);CHKERRQ(ierr);
615d8588912SDave May           break;
616d8588912SDave May         }
617d8588912SDave May       }
618d8588912SDave May       if (i==bA->nr) {
619d8588912SDave May         /* have an empty column */
620d8588912SDave May         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
621d8588912SDave May       }
622d8588912SDave May     }
623f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
624d8588912SDave May     /* hand back control to the nest vector */
625d8588912SDave May     for (j=0; j<bA->nc; j++) {
6266bf464f9SBarry Smith       ierr = VecDestroy(&R[j]);CHKERRQ(ierr);
627d8588912SDave May     }
628d8588912SDave May     ierr = PetscFree(R);CHKERRQ(ierr);
629d8588912SDave May   }
630d8588912SDave May 
631d8588912SDave May   if (left) {
632d8588912SDave May     /* allocate L */
633d8588912SDave May     ierr = PetscMalloc( sizeof(Vec) * bA->nr, &L );CHKERRQ(ierr);
634d8588912SDave May     /* Create the left vectors */
635d8588912SDave May     for (i=0; i<bA->nr; i++) {
636d8588912SDave May       for (j=0; j<bA->nc; j++) {
637d8588912SDave May         if (bA->m[i][j]) {
638d8588912SDave May           ierr = MatGetVecs(bA->m[i][j],PETSC_NULL,&L[i]);CHKERRQ(ierr);
639d8588912SDave May           break;
640d8588912SDave May         }
641d8588912SDave May       }
642d8588912SDave May       if (j==bA->nc) {
643d8588912SDave May         /* have an empty row */
644d8588912SDave May         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
645d8588912SDave May       }
646d8588912SDave May     }
647d8588912SDave May 
648f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
649d8588912SDave May     for (i=0; i<bA->nr; i++) {
6506bf464f9SBarry Smith       ierr = VecDestroy(&L[i]);CHKERRQ(ierr);
651d8588912SDave May     }
652d8588912SDave May 
653d8588912SDave May     ierr = PetscFree(L);CHKERRQ(ierr);
654d8588912SDave May   }
655d8588912SDave May   PetscFunctionReturn(0);
656d8588912SDave May }
657d8588912SDave May 
658d8588912SDave May #undef __FUNCT__
659d8588912SDave May #define __FUNCT__ "MatView_Nest"
660207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
661d8588912SDave May {
662d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
663d8588912SDave May   PetscBool      isascii;
664d8588912SDave May   PetscInt       i,j;
665d8588912SDave May   PetscErrorCode ierr;
666d8588912SDave May 
667d8588912SDave May   PetscFunctionBegin;
668d8588912SDave May   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
669d8588912SDave May   if (isascii) {
670d8588912SDave May 
671d8588912SDave May     PetscViewerASCIIPrintf(viewer,"Matrix object: \n" );
672d8588912SDave May     PetscViewerASCIIPushTab( viewer );    /* push0 */
673d8588912SDave May     PetscViewerASCIIPrintf( viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc);
674d8588912SDave May 
675d8588912SDave May     PetscViewerASCIIPrintf(viewer,"MatNest structure: \n" );
676d8588912SDave May     for (i=0; i<bA->nr; i++) {
677d8588912SDave May       for (j=0; j<bA->nc; j++) {
678d8588912SDave May         const MatType type;
679d8588912SDave May         const char *name;
680d8588912SDave May         PetscInt NR,NC;
681d8588912SDave May         PetscBool isNest = PETSC_FALSE;
682d8588912SDave May 
683d8588912SDave May         if (!bA->m[i][j]) {
684d8588912SDave May           PetscViewerASCIIPrintf( viewer, "(%D,%D) : PETSC_NULL \n",i,j);
685d8588912SDave May           continue;
686d8588912SDave May         }
687d8588912SDave May         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
688d8588912SDave May         ierr = MatGetType( bA->m[i][j], &type );CHKERRQ(ierr);
689d8588912SDave May         name = ((PetscObject)bA->m[i][j])->prefix;
690d8588912SDave May         ierr = PetscTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
691d8588912SDave May 
692d8588912SDave May         PetscViewerASCIIPrintf( viewer, "(%D,%D) : name=\"%s\", type=%s, rows=%D, cols=%D \n",i,j,name,type,NR,NC);
693d8588912SDave May 
694d8588912SDave May         if (isNest) {
695d8588912SDave May           PetscViewerASCIIPushTab( viewer );  /* push1 */
696d8588912SDave May           ierr = MatView( bA->m[i][j], viewer );CHKERRQ(ierr);
697d8588912SDave May           PetscViewerASCIIPopTab(viewer);    /* pop1 */
698d8588912SDave May         }
699d8588912SDave May       }
700d8588912SDave May     }
701d8588912SDave May     PetscViewerASCIIPopTab(viewer);    /* pop0 */
702d8588912SDave May   }
703d8588912SDave May   PetscFunctionReturn(0);
704d8588912SDave May }
705d8588912SDave May 
706d8588912SDave May #undef __FUNCT__
707d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest"
708207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A)
709d8588912SDave May {
710d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
711d8588912SDave May   PetscInt       i,j;
712d8588912SDave May   PetscErrorCode ierr;
713d8588912SDave May 
714d8588912SDave May   PetscFunctionBegin;
715d8588912SDave May   for (i=0; i<bA->nr; i++) {
716d8588912SDave May     for (j=0; j<bA->nc; j++) {
717d8588912SDave May       if (!bA->m[i][j]) continue;
718d8588912SDave May       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
719d8588912SDave May     }
720d8588912SDave May   }
721d8588912SDave May   PetscFunctionReturn(0);
722d8588912SDave May }
723d8588912SDave May 
724d8588912SDave May #undef __FUNCT__
725d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest"
726207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
727d8588912SDave May {
728d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
729841e96a3SJed Brown   Mat            *b;
730841e96a3SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
731d8588912SDave May   PetscErrorCode ierr;
732d8588912SDave May 
733d8588912SDave May   PetscFunctionBegin;
734841e96a3SJed Brown   ierr = PetscMalloc(nr*nc*sizeof(Mat),&b);CHKERRQ(ierr);
735841e96a3SJed Brown   for (i=0; i<nr; i++) {
736841e96a3SJed Brown     for (j=0; j<nc; j++) {
737841e96a3SJed Brown       if (bA->m[i][j]) {
738841e96a3SJed Brown         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
739841e96a3SJed Brown       } else {
740841e96a3SJed Brown         b[i*nc+j] = PETSC_NULL;
741d8588912SDave May       }
742d8588912SDave May     }
743d8588912SDave May   }
744f349c1fdSJed Brown   ierr = MatCreateNest(((PetscObject)A)->comm,nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
745841e96a3SJed Brown   /* Give the new MatNest exclusive ownership */
746841e96a3SJed Brown   for (i=0; i<nr*nc; i++) {
7476bf464f9SBarry Smith     ierr = MatDestroy(&b[i]);CHKERRQ(ierr);
748d8588912SDave May   }
749d8588912SDave May   ierr = PetscFree(b);CHKERRQ(ierr);
750d8588912SDave May 
751841e96a3SJed Brown   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
752841e96a3SJed Brown   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
753d8588912SDave May   PetscFunctionReturn(0);
754d8588912SDave May }
755d8588912SDave May 
756d8588912SDave May /* nest api */
757d8588912SDave May EXTERN_C_BEGIN
758d8588912SDave May #undef __FUNCT__
759d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest"
760d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
761d8588912SDave May {
762d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
763d8588912SDave May   PetscFunctionBegin;
764d8588912SDave May   if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
765d8588912SDave May   if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
766d8588912SDave May   *mat = bA->m[idxm][jdxm];
767d8588912SDave May   PetscFunctionReturn(0);
768d8588912SDave May }
769d8588912SDave May EXTERN_C_END
770d8588912SDave May 
771d8588912SDave May #undef __FUNCT__
772d8588912SDave May #define __FUNCT__ "MatNestGetSubMat"
773d8588912SDave May /*@C
774d8588912SDave May  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
775d8588912SDave May 
776d8588912SDave May  Not collective
777d8588912SDave May 
778d8588912SDave May  Input Parameters:
779d8588912SDave May  .  A  - nest matrix
780d8588912SDave May  .  idxm - index of the matrix within the nest matrix
781d8588912SDave May  .  jdxm - index of the matrix within the nest matrix
782d8588912SDave May 
783d8588912SDave May  Output Parameter:
784d8588912SDave May  .  sub - matrix at index idxm,jdxm within the nest matrix
785d8588912SDave May 
786d8588912SDave May  Notes:
787d8588912SDave May 
788d8588912SDave May  Level: developer
789d8588912SDave May 
790d8588912SDave May  .seealso: MatNestGetSize(), MatNestGetSubMats()
791d8588912SDave May @*/
7927087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
793d8588912SDave May {
794699a902aSJed Brown   PetscErrorCode ierr;
795d8588912SDave May 
796d8588912SDave May   PetscFunctionBegin;
797699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
798d8588912SDave May   PetscFunctionReturn(0);
799d8588912SDave May }
800d8588912SDave May 
801d8588912SDave May EXTERN_C_BEGIN
802d8588912SDave May #undef __FUNCT__
803d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest"
804d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
805d8588912SDave May {
806d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
807d8588912SDave May   PetscFunctionBegin;
808d8588912SDave May   if (M)   { *M   = bA->nr; }
809d8588912SDave May   if (N)   { *N   = bA->nc; }
810d8588912SDave May   if (mat) { *mat = bA->m;  }
811d8588912SDave May   PetscFunctionReturn(0);
812d8588912SDave May }
813d8588912SDave May EXTERN_C_END
814d8588912SDave May 
815d8588912SDave May #undef __FUNCT__
816d8588912SDave May #define __FUNCT__ "MatNestGetSubMats"
817d8588912SDave May /*@C
818d8588912SDave May  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
819d8588912SDave May 
820d8588912SDave May  Not collective
821d8588912SDave May 
822d8588912SDave May  Input Parameters:
823d8588912SDave May  .  A  - nest matri
824d8588912SDave May 
825d8588912SDave May  Output Parameter:
826d8588912SDave May  .  M - number of rows in the nest matrix
827d8588912SDave May  .  N - number of cols in the nest matrix
828d8588912SDave May  .  mat - 2d array of matrices
829d8588912SDave May 
830d8588912SDave May  Notes:
831d8588912SDave May 
832d8588912SDave May  The user should not free the array mat.
833d8588912SDave May 
834d8588912SDave May  Level: developer
835d8588912SDave May 
836d8588912SDave May  .seealso: MatNestGetSize(), MatNestGetSubMat()
837d8588912SDave May @*/
8387087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
839d8588912SDave May {
840699a902aSJed Brown   PetscErrorCode ierr;
841d8588912SDave May 
842d8588912SDave May   PetscFunctionBegin;
843699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
844d8588912SDave May   PetscFunctionReturn(0);
845d8588912SDave May }
846d8588912SDave May 
847d8588912SDave May EXTERN_C_BEGIN
848d8588912SDave May #undef __FUNCT__
849d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest"
8507087cfbeSBarry Smith PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
851d8588912SDave May {
852d8588912SDave May   Mat_Nest  *bA = (Mat_Nest*)A->data;
853d8588912SDave May 
854d8588912SDave May   PetscFunctionBegin;
855d8588912SDave May   if (M) { *M  = bA->nr; }
856d8588912SDave May   if (N) { *N  = bA->nc; }
857d8588912SDave May   PetscFunctionReturn(0);
858d8588912SDave May }
859d8588912SDave May EXTERN_C_END
860d8588912SDave May 
861d8588912SDave May #undef __FUNCT__
862d8588912SDave May #define __FUNCT__ "MatNestGetSize"
863d8588912SDave May /*@C
864d8588912SDave May  MatNestGetSize - Returns the size of the nest matrix.
865d8588912SDave May 
866d8588912SDave May  Not collective
867d8588912SDave May 
868d8588912SDave May  Input Parameters:
869d8588912SDave May  .  A  - nest matrix
870d8588912SDave May 
871d8588912SDave May  Output Parameter:
872d8588912SDave May  .  M - number of rows in the nested mat
873d8588912SDave May  .  N - number of cols in the nested mat
874d8588912SDave May 
875d8588912SDave May  Notes:
876d8588912SDave May 
877d8588912SDave May  Level: developer
878d8588912SDave May 
879d8588912SDave May  .seealso: MatNestGetSubMat(), MatNestGetSubMats()
880d8588912SDave May @*/
8817087cfbeSBarry Smith PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
882d8588912SDave May {
883699a902aSJed Brown   PetscErrorCode ierr;
884d8588912SDave May 
885d8588912SDave May   PetscFunctionBegin;
886699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
887d8588912SDave May   PetscFunctionReturn(0);
888d8588912SDave May }
889d8588912SDave May 
890207556f9SJed Brown EXTERN_C_BEGIN
891207556f9SJed Brown #undef __FUNCT__
892207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest"
8937087cfbeSBarry Smith PetscErrorCode  MatNestSetVecType_Nest(Mat A,const VecType vtype)
894207556f9SJed Brown {
895207556f9SJed Brown   PetscErrorCode ierr;
896207556f9SJed Brown   PetscBool      flg;
897207556f9SJed Brown 
898207556f9SJed Brown   PetscFunctionBegin;
899207556f9SJed Brown   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
900207556f9SJed Brown   /* In reality, this only distinguishes VECNEST and "other" */
901207556f9SJed Brown   A->ops->getvecs       = flg ? MatGetVecs_Nest       : 0;
9027874fa86SDave May   A->ops->getdiagonal   = flg ? MatGetDiagonal_Nest   : 0;
9037874fa86SDave May   A->ops->diagonalscale = flg ? MatDiagonalScale_Nest : 0;
9047874fa86SDave May 
905207556f9SJed Brown  PetscFunctionReturn(0);
906207556f9SJed Brown }
907207556f9SJed Brown EXTERN_C_END
908207556f9SJed Brown 
909207556f9SJed Brown #undef __FUNCT__
910207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType"
911207556f9SJed Brown /*@C
912207556f9SJed Brown  MatNestSetVecType - Sets the type of Vec returned by MatGetVecs()
913207556f9SJed Brown 
914207556f9SJed Brown  Not collective
915207556f9SJed Brown 
916207556f9SJed Brown  Input Parameters:
917207556f9SJed Brown +  A  - nest matrix
918207556f9SJed Brown -  vtype - type to use for creating vectors
919207556f9SJed Brown 
920207556f9SJed Brown  Notes:
921207556f9SJed Brown 
922207556f9SJed Brown  Level: developer
923207556f9SJed Brown 
924207556f9SJed Brown  .seealso: MatGetVecs()
925207556f9SJed Brown @*/
9267087cfbeSBarry Smith PetscErrorCode  MatNestSetVecType(Mat A,const VecType vtype)
927207556f9SJed Brown {
928207556f9SJed Brown   PetscErrorCode ierr;
929207556f9SJed Brown 
930207556f9SJed Brown   PetscFunctionBegin;
931207556f9SJed Brown   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,const VecType),(A,vtype));CHKERRQ(ierr);
932207556f9SJed Brown   PetscFunctionReturn(0);
933207556f9SJed Brown }
934207556f9SJed Brown 
935c8883902SJed Brown EXTERN_C_BEGIN
936d8588912SDave May #undef __FUNCT__
937c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats_Nest"
938c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
939d8588912SDave May {
940c8883902SJed Brown   Mat_Nest       *s = (Mat_Nest*)A->data;
941c8883902SJed Brown   PetscInt       i,j,m,n,M,N;
942d8588912SDave May   PetscErrorCode ierr;
943d8588912SDave May 
944d8588912SDave May   PetscFunctionBegin;
945c8883902SJed Brown   s->nr = nr;
946c8883902SJed Brown   s->nc = nc;
947d8588912SDave May 
948c8883902SJed Brown   /* Create space for submatrices */
949c8883902SJed Brown   ierr = PetscMalloc(sizeof(Mat*)*nr,&s->m);CHKERRQ(ierr);
950c8883902SJed Brown   for (i=0; i<nr; i++) {
951c8883902SJed Brown     ierr = PetscMalloc(sizeof(Mat)*nc,&s->m[i]);CHKERRQ(ierr);
952d8588912SDave May   }
953c8883902SJed Brown   for (i=0; i<nr; i++) {
954c8883902SJed Brown     for (j=0; j<nc; j++) {
955c8883902SJed Brown       s->m[i][j] = a[i*nc+j];
956c8883902SJed Brown       if (a[i*nc+j]) {
957c8883902SJed Brown         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
958d8588912SDave May       }
959d8588912SDave May     }
960d8588912SDave May   }
961d8588912SDave May 
962*8188e55aSJed Brown   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
963d8588912SDave May 
964c8883902SJed Brown   ierr = PetscMalloc(sizeof(PetscInt)*nr,&s->row_len);CHKERRQ(ierr);
965c8883902SJed Brown   ierr = PetscMalloc(sizeof(PetscInt)*nc,&s->col_len);CHKERRQ(ierr);
966c8883902SJed Brown   for (i=0; i<nr; i++) s->row_len[i]=-1;
967c8883902SJed Brown   for (j=0; j<nc; j++) s->col_len[j]=-1;
968d8588912SDave May 
969*8188e55aSJed Brown   ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr);
970d8588912SDave May 
971c8883902SJed Brown   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
972c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
973c8883902SJed Brown   ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr);
974c8883902SJed Brown   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
975c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
976c8883902SJed Brown   ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr);
977c8883902SJed Brown 
978c8883902SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
979c8883902SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
980c8883902SJed Brown 
981c8883902SJed Brown   ierr = PetscMalloc2(nr,Vec,&s->left,nc,Vec,&s->right);CHKERRQ(ierr);
982c8883902SJed Brown   ierr = PetscMemzero(s->left,nr*sizeof(Vec));CHKERRQ(ierr);
983c8883902SJed Brown   ierr = PetscMemzero(s->right,nc*sizeof(Vec));CHKERRQ(ierr);
984d8588912SDave May   PetscFunctionReturn(0);
985d8588912SDave May }
986c8883902SJed Brown EXTERN_C_END
987d8588912SDave May 
988c8883902SJed Brown #undef __FUNCT__
989c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats"
990c8883902SJed Brown /*@
991c8883902SJed Brown    MatNestSetSubMats - Sets the nested submatrices
992c8883902SJed Brown 
993c8883902SJed Brown    Collective on Mat
994c8883902SJed Brown 
995c8883902SJed Brown    Input Parameter:
996c8883902SJed Brown +  N - nested matrix
997c8883902SJed Brown .  nr - number of nested row blocks
998c8883902SJed Brown .  is_row - index sets for each nested row block, or PETSC_NULL to make contiguous
999c8883902SJed Brown .  nc - number of nested column blocks
1000c8883902SJed Brown .  is_col - index sets for each nested column block, or PETSC_NULL to make contiguous
1001c8883902SJed Brown -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL
1002c8883902SJed Brown 
1003c8883902SJed Brown    Level: advanced
1004c8883902SJed Brown 
1005c8883902SJed Brown .seealso: MatCreateNest(), MATNEST
1006c8883902SJed Brown @*/
1007c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1008c8883902SJed Brown {
1009c8883902SJed Brown   PetscErrorCode ierr;
1010c8883902SJed Brown   PetscInt i;
1011c8883902SJed Brown 
1012c8883902SJed Brown   PetscFunctionBegin;
1013c8883902SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1014c8883902SJed Brown   if (nr < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
1015c8883902SJed Brown   if (nr && is_row) {
1016c8883902SJed Brown     PetscValidPointer(is_row,3);
1017c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
1018c8883902SJed Brown   }
1019c8883902SJed Brown   if (nc < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
1020c8883902SJed Brown   if (nc && is_row) {
1021c8883902SJed Brown     PetscValidPointer(is_col,5);
1022c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
1023c8883902SJed Brown   }
1024c8883902SJed Brown   if (nr*nc) PetscValidPointer(a,6);
1025c8883902SJed 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);
1026c8883902SJed Brown   PetscFunctionReturn(0);
1027c8883902SJed Brown }
1028d8588912SDave May 
1029d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1030d8588912SDave May /*
1031d8588912SDave May   nprocessors = NP
1032d8588912SDave May   Nest x^T = ( (g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1) )
1033d8588912SDave May        proc 0: => (g_0,h_0,)
1034d8588912SDave May        proc 1: => (g_1,h_1,)
1035d8588912SDave May        ...
1036d8588912SDave May        proc nprocs-1: => (g_NP-1,h_NP-1,)
1037d8588912SDave May 
1038d8588912SDave May             proc 0:                      proc 1:                    proc nprocs-1:
1039d8588912SDave May     is[0] = ( 0,1,2,...,nlocal(g_0)-1 )  ( 0,1,...,nlocal(g_1)-1 )  ( 0,1,...,nlocal(g_NP-1) )
1040d8588912SDave May 
1041d8588912SDave May             proc 0:
1042d8588912SDave May     is[1] = ( nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1 )
1043d8588912SDave May             proc 1:
1044d8588912SDave May     is[1] = ( nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1 )
1045d8588912SDave May 
1046d8588912SDave May             proc NP-1:
1047d8588912SDave May     is[1] = ( nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1 )
1048d8588912SDave May */
1049d8588912SDave May #undef __FUNCT__
1050d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private"
1051841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1052d8588912SDave May {
1053e2d7f03fSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
1054*8188e55aSJed Brown   PetscInt       i,j,offset,n,nsum,bs;
1055d8588912SDave May   PetscErrorCode ierr;
10562ae74bdbSJed Brown   Mat            sub;
1057d8588912SDave May 
1058d8588912SDave May   PetscFunctionBegin;
1059*8188e55aSJed Brown   ierr = PetscMalloc(sizeof(IS)*nr,&vs->isglobal.row);CHKERRQ(ierr);
1060*8188e55aSJed Brown   ierr = PetscMalloc(sizeof(IS)*nc,&vs->isglobal.col);CHKERRQ(ierr);
1061d8588912SDave May   if (is_row) { /* valid IS is passed in */
1062d8588912SDave May     /* refs on is[] are incremeneted */
1063e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1064d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
1065e2d7f03fSJed Brown       vs->isglobal.row[i] = is_row[i];
1066d8588912SDave May     }
10672ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
1068*8188e55aSJed Brown     nsum = 0;
1069*8188e55aSJed Brown     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
1070*8188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1071*8188e55aSJed Brown       if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
1072*8188e55aSJed Brown       ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr);
1073*8188e55aSJed Brown       if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
1074*8188e55aSJed Brown       nsum += n;
1075*8188e55aSJed Brown     }
1076*8188e55aSJed Brown     offset = 0;
1077*8188e55aSJed Brown     ierr = MPI_Exscan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
1078e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1079f349c1fdSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
10802ae74bdbSJed Brown       ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr);
10812ae74bdbSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1082e2d7f03fSJed Brown       ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
1083e2d7f03fSJed Brown       ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
10842ae74bdbSJed Brown       offset += n;
1085d8588912SDave May     }
1086d8588912SDave May   }
1087d8588912SDave May 
1088d8588912SDave May   if (is_col) { /* valid IS is passed in */
1089d8588912SDave May     /* refs on is[] are incremeneted */
1090e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1091d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
1092e2d7f03fSJed Brown       vs->isglobal.col[j] = is_col[j];
1093d8588912SDave May     }
10942ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
10952ae74bdbSJed Brown     offset = A->cmap->rstart;
1096*8188e55aSJed Brown     nsum = 0;
1097*8188e55aSJed Brown     for (j=0; j<vs->nc; j++) {
1098*8188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
1099*8188e55aSJed Brown       if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
1100*8188e55aSJed Brown       ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr);
1101*8188e55aSJed Brown       if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
1102*8188e55aSJed Brown       nsum += n;
1103*8188e55aSJed Brown     }
1104*8188e55aSJed Brown     offset = 0;
1105*8188e55aSJed Brown     ierr = MPI_Exscan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
1106e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1107f349c1fdSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
11082ae74bdbSJed Brown       ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr);
11092ae74bdbSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1110e2d7f03fSJed Brown       ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
1111e2d7f03fSJed Brown       ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
11122ae74bdbSJed Brown       offset += n;
1113d8588912SDave May     }
1114d8588912SDave May   }
1115e2d7f03fSJed Brown 
1116e2d7f03fSJed Brown   /* Set up the local ISs */
1117e2d7f03fSJed Brown   ierr = PetscMalloc(vs->nr*sizeof(IS),&vs->islocal.row);CHKERRQ(ierr);
1118e2d7f03fSJed Brown   ierr = PetscMalloc(vs->nc*sizeof(IS),&vs->islocal.col);CHKERRQ(ierr);
1119e2d7f03fSJed Brown   for (i=0,offset=0; i<vs->nr; i++) {
1120e2d7f03fSJed Brown     IS                     isloc;
1121*8188e55aSJed Brown     ISLocalToGlobalMapping rmap = PETSC_NULL;
1122e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1123e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1124*8188e55aSJed Brown     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,PETSC_NULL);CHKERRQ(ierr);}
1125207556f9SJed Brown     if (rmap) {
1126e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1127e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
1128e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1129e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1130207556f9SJed Brown     } else {
1131207556f9SJed Brown       nlocal = 0;
1132207556f9SJed Brown       isloc  = PETSC_NULL;
1133207556f9SJed Brown     }
1134e2d7f03fSJed Brown     vs->islocal.row[i] = isloc;
1135e2d7f03fSJed Brown     offset += nlocal;
1136e2d7f03fSJed Brown   }
1137*8188e55aSJed Brown   for (i=0,offset=0; i<vs->nc; i++) {
1138e2d7f03fSJed Brown     IS                     isloc;
1139*8188e55aSJed Brown     ISLocalToGlobalMapping cmap = PETSC_NULL;
1140e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1141e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
1142*8188e55aSJed Brown     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,PETSC_NULL,&cmap);CHKERRQ(ierr);}
1143207556f9SJed Brown     if (cmap) {
1144e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1145e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr);
1146e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1147e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1148207556f9SJed Brown     } else {
1149207556f9SJed Brown       nlocal = 0;
1150207556f9SJed Brown       isloc  = PETSC_NULL;
1151207556f9SJed Brown     }
1152e2d7f03fSJed Brown     vs->islocal.col[i] = isloc;
1153e2d7f03fSJed Brown     offset += nlocal;
1154e2d7f03fSJed Brown   }
1155d8588912SDave May   PetscFunctionReturn(0);
1156d8588912SDave May }
1157d8588912SDave May 
1158d8588912SDave May #undef __FUNCT__
1159d8588912SDave May #define __FUNCT__ "MatCreateNest"
1160659c6bb0SJed Brown /*@
1161659c6bb0SJed Brown    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
1162659c6bb0SJed Brown 
1163659c6bb0SJed Brown    Collective on Mat
1164659c6bb0SJed Brown 
1165659c6bb0SJed Brown    Input Parameter:
1166659c6bb0SJed Brown +  comm - Communicator for the new Mat
1167659c6bb0SJed Brown .  nr - number of nested row blocks
1168659c6bb0SJed Brown .  is_row - index sets for each nested row block, or PETSC_NULL to make contiguous
1169659c6bb0SJed Brown .  nc - number of nested column blocks
1170659c6bb0SJed Brown .  is_col - index sets for each nested column block, or PETSC_NULL to make contiguous
1171659c6bb0SJed Brown -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL
1172659c6bb0SJed Brown 
1173659c6bb0SJed Brown    Output Parameter:
1174659c6bb0SJed Brown .  B - new matrix
1175659c6bb0SJed Brown 
1176659c6bb0SJed Brown    Level: advanced
1177659c6bb0SJed Brown 
1178659c6bb0SJed Brown .seealso: MatCreate(), VecCreateNest(), DMGetMatrix(), MATNEST
1179659c6bb0SJed Brown @*/
11807087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1181d8588912SDave May {
1182d8588912SDave May   Mat            A;
1183d8588912SDave May   PetscErrorCode ierr;
1184d8588912SDave May 
1185d8588912SDave May   PetscFunctionBegin;
1186c8883902SJed Brown   *B = 0;
1187d8588912SDave May   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1188c8883902SJed Brown   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
1189c8883902SJed Brown   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1190d8588912SDave May   *B = A;
1191d8588912SDave May   PetscFunctionReturn(0);
1192d8588912SDave May }
1193659c6bb0SJed Brown 
1194659c6bb0SJed Brown /*MC
1195659c6bb0SJed Brown   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
1196659c6bb0SJed Brown 
1197659c6bb0SJed Brown   Level: intermediate
1198659c6bb0SJed Brown 
1199659c6bb0SJed Brown   Notes:
1200659c6bb0SJed Brown   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1201659c6bb0SJed Brown   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1202659c6bb0SJed Brown   It is usually used with DMComposite and DMGetMatrix()
1203659c6bb0SJed Brown 
1204659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest()
1205659c6bb0SJed Brown M*/
1206c8883902SJed Brown EXTERN_C_BEGIN
1207c8883902SJed Brown #undef __FUNCT__
1208c8883902SJed Brown #define __FUNCT__ "MatCreate_Nest"
1209c8883902SJed Brown PetscErrorCode MatCreate_Nest(Mat A)
1210c8883902SJed Brown {
1211c8883902SJed Brown   Mat_Nest       *s;
1212c8883902SJed Brown   PetscErrorCode ierr;
1213c8883902SJed Brown 
1214c8883902SJed Brown   PetscFunctionBegin;
1215c8883902SJed Brown   ierr = PetscNewLog(A,Mat_Nest,&s);CHKERRQ(ierr);
1216c8883902SJed Brown   A->data         = (void*)s;
1217c8883902SJed Brown   s->nr = s->nc   = -1;
1218c8883902SJed Brown   s->m            = PETSC_NULL;
1219c8883902SJed Brown 
1220c8883902SJed Brown   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
1221c8883902SJed Brown   A->ops->mult                  = MatMult_Nest;
1222c8883902SJed Brown   A->ops->multadd               = 0;
1223c8883902SJed Brown   A->ops->multtranspose         = MatMultTranspose_Nest;
1224c8883902SJed Brown   A->ops->multtransposeadd      = 0;
1225c8883902SJed Brown   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1226c8883902SJed Brown   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1227c8883902SJed Brown   A->ops->zeroentries           = MatZeroEntries_Nest;
1228c8883902SJed Brown   A->ops->duplicate             = MatDuplicate_Nest;
1229c8883902SJed Brown   A->ops->getsubmatrix          = MatGetSubMatrix_Nest;
1230c8883902SJed Brown   A->ops->destroy               = MatDestroy_Nest;
1231c8883902SJed Brown   A->ops->view                  = MatView_Nest;
1232c8883902SJed Brown   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1233c8883902SJed Brown   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1234c8883902SJed Brown   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
12357874fa86SDave May   A->ops->getdiagonal           = MatGetDiagonal_Nest2; /* VECNEST version activated by calling MatNestSetVecType(A,VECNEST) */
12367874fa86SDave May   A->ops->diagonalscale         = MatDiagonalScale_Nest2; /* VECNEST version activated by calling MatNestSetVecType(A,VECNEST) */
1237c8883902SJed Brown 
1238c8883902SJed Brown   A->spptr        = 0;
1239c8883902SJed Brown   A->same_nonzero = PETSC_FALSE;
1240c8883902SJed Brown   A->assembled    = PETSC_FALSE;
1241c8883902SJed Brown 
1242c8883902SJed Brown   /* expose Nest api's */
1243c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C",   "MatNestGetSubMat_Nest",   MatNestGetSubMat_Nest);CHKERRQ(ierr);
1244c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C",  "MatNestGetSubMats_Nest",  MatNestGetSubMats_Nest);CHKERRQ(ierr);
1245c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C",     "MatNestGetSize_Nest",     MatNestGetSize_Nest);CHKERRQ(ierr);
1246c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C",  "MatNestSetVecType_Nest",  MatNestSetVecType_Nest);CHKERRQ(ierr);
1247c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C",  "MatNestSetSubMats_Nest",  MatNestSetSubMats_Nest);CHKERRQ(ierr);
1248c8883902SJed Brown 
1249c8883902SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
1250c8883902SJed Brown   PetscFunctionReturn(0);
1251c8883902SJed Brown }
125286e80854SHong Zhang EXTERN_C_END
1253