xref: /petsc/src/mat/impls/nest/matnest.c (revision 9194d70fb74746bb29895a4f46eca3d615b126ad)
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__
98188e55aSJed Brown #define __FUNCT__ "MatNestGetSizes_Private"
108188e55aSJed 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;
138188e55aSJed Brown   PetscInt       i,j;
14d8588912SDave May   PetscErrorCode ierr;
15d8588912SDave May 
16d8588912SDave May   PetscFunctionBegin;
178188e55aSJed Brown   *m = *n = *M = *N = 0;
188188e55aSJed Brown   for (i=0; i<bA->nr; i++) {  /* rows */
198188e55aSJed Brown     PetscInt sm,sM;
208188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.row[i],&sm);CHKERRQ(ierr);
218188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.row[i],&sM);CHKERRQ(ierr);
228188e55aSJed Brown     *m += sm;
238188e55aSJed Brown     *M += sM;
24d8588912SDave May   }
258188e55aSJed Brown   for (j=0; j<bA->nc; j++) {  /* cols */
268188e55aSJed Brown     PetscInt sn,sN;
278188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.col[j],&sn);CHKERRQ(ierr);
288188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.col[j],&sN);CHKERRQ(ierr);
298188e55aSJed Brown     *n += sn;
308188e55aSJed 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__
176*9194d70fSJed Brown #define __FUNCT__ "MatMultAdd_Nest"
177*9194d70fSJed Brown static PetscErrorCode MatMultAdd_Nest(Mat A,Vec x,Vec y,Vec z)
178*9194d70fSJed Brown {
179*9194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
180*9194d70fSJed Brown   Vec            *bx = bA->right,*bz = bA->left;
181*9194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
182*9194d70fSJed Brown   PetscErrorCode ierr;
183*9194d70fSJed Brown 
184*9194d70fSJed Brown   PetscFunctionBegin;
185*9194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
186*9194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
187*9194d70fSJed Brown   for (i=0; i<nr; i++) {
188*9194d70fSJed Brown     if (y != z) {
189*9194d70fSJed Brown       Vec by;
190*9194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
191*9194d70fSJed Brown       ierr = VecCopy(by,bz[i]);CHKERRQ(ierr);
192*9194d70fSJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
193*9194d70fSJed Brown     }
194*9194d70fSJed Brown     for (j=0; j<nc; j++) {
195*9194d70fSJed Brown       if (!bA->m[i][j]) continue;
196*9194d70fSJed Brown       /* y[i] <- y[i] + A[i][j] * x[j] */
197*9194d70fSJed Brown       ierr = MatMultAdd(bA->m[i][j],bx[j],bz[i],bz[i]);CHKERRQ(ierr);
198*9194d70fSJed Brown     }
199*9194d70fSJed Brown   }
200*9194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
201*9194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
202*9194d70fSJed Brown   PetscFunctionReturn(0);
203*9194d70fSJed Brown }
204*9194d70fSJed Brown 
205*9194d70fSJed Brown #undef __FUNCT__
206d8588912SDave May #define __FUNCT__ "MatMultTranspose_Nest"
207207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y)
208d8588912SDave May {
209d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
210207556f9SJed Brown   Vec            *bx = bA->left,*by = bA->right;
211207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
212d8588912SDave May   PetscErrorCode ierr;
213d8588912SDave May 
214d8588912SDave May   PetscFunctionBegin;
215609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
216609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
217207556f9SJed Brown   for (j=0; j<nc; j++) {
218609e31cbSJed Brown     ierr = VecZeroEntries(by[j]);CHKERRQ(ierr);
219609e31cbSJed Brown     for (i=0; i<nr; i++) {
220207556f9SJed Brown       if (!bA->m[j][i]) continue;
221609e31cbSJed Brown       /* y[j] <- y[j] + (A[i][j])^T * x[i] */
222609e31cbSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],by[j],by[j]);CHKERRQ(ierr);
223d8588912SDave May     }
224d8588912SDave May   }
225609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
226609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
227d8588912SDave May   PetscFunctionReturn(0);
228d8588912SDave May }
229d8588912SDave May 
230d8588912SDave May #undef __FUNCT__
231*9194d70fSJed Brown #define __FUNCT__ "MatMultTransposeAdd_Nest"
232*9194d70fSJed Brown static PetscErrorCode MatMultTransposeAdd_Nest(Mat A,Vec x,Vec y,Vec z)
233*9194d70fSJed Brown {
234*9194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
235*9194d70fSJed Brown   Vec            *bx = bA->left,*bz = bA->right;
236*9194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
237*9194d70fSJed Brown   PetscErrorCode ierr;
238*9194d70fSJed Brown 
239*9194d70fSJed Brown   PetscFunctionBegin;
240*9194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
241*9194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
242*9194d70fSJed Brown   for (j=0; j<nc; j++) {
243*9194d70fSJed Brown     if (y != z) {
244*9194d70fSJed Brown       Vec by;
245*9194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
246*9194d70fSJed Brown       ierr = VecCopy(by,bz[j]);CHKERRQ(ierr);
247*9194d70fSJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
248*9194d70fSJed Brown     }
249*9194d70fSJed Brown     for (i=0; i<nr; i++) {
250*9194d70fSJed Brown       if (!bA->m[j][i]) continue;
251*9194d70fSJed Brown       /* z[j] <- y[j] + (A[i][j])^T * x[i] */
252*9194d70fSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],bz[j],bz[j]);CHKERRQ(ierr);
253*9194d70fSJed Brown     }
254*9194d70fSJed Brown   }
255*9194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
256*9194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
257*9194d70fSJed Brown   PetscFunctionReturn(0);
258*9194d70fSJed Brown }
259*9194d70fSJed Brown 
260*9194d70fSJed Brown #undef __FUNCT__
261e2d7f03fSJed Brown #define __FUNCT__ "MatNestDestroyISList"
262e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list)
263e2d7f03fSJed Brown {
264e2d7f03fSJed Brown   PetscErrorCode ierr;
265e2d7f03fSJed Brown   IS             *lst = *list;
266e2d7f03fSJed Brown   PetscInt       i;
267e2d7f03fSJed Brown 
268e2d7f03fSJed Brown   PetscFunctionBegin;
269e2d7f03fSJed Brown   if (!lst) PetscFunctionReturn(0);
2706bf464f9SBarry Smith   for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);}
271e2d7f03fSJed Brown   ierr = PetscFree(lst);CHKERRQ(ierr);
272e2d7f03fSJed Brown   *list = PETSC_NULL;
273e2d7f03fSJed Brown   PetscFunctionReturn(0);
274e2d7f03fSJed Brown }
275e2d7f03fSJed Brown 
276e2d7f03fSJed Brown #undef __FUNCT__
277d8588912SDave May #define __FUNCT__ "MatDestroy_Nest"
278207556f9SJed Brown static PetscErrorCode MatDestroy_Nest(Mat A)
279d8588912SDave May {
280d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
281d8588912SDave May   PetscInt       i,j;
282d8588912SDave May   PetscErrorCode ierr;
283d8588912SDave May 
284d8588912SDave May   PetscFunctionBegin;
285d8588912SDave May   /* release the matrices and the place holders */
286e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr);
287e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr);
288e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
289e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
290d8588912SDave May 
291d8588912SDave May   ierr = PetscFree(vs->row_len);CHKERRQ(ierr);
292d8588912SDave May   ierr = PetscFree(vs->col_len);CHKERRQ(ierr);
293d8588912SDave May 
294207556f9SJed Brown   ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr);
295207556f9SJed Brown 
296d8588912SDave May   /* release the matrices and the place holders */
297d8588912SDave May   if (vs->m) {
298d8588912SDave May     for (i=0; i<vs->nr; i++) {
299d8588912SDave May       for (j=0; j<vs->nc; j++) {
3006bf464f9SBarry Smith         ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr);
301d8588912SDave May       }
302d8588912SDave May       ierr = PetscFree( vs->m[i] );CHKERRQ(ierr);
303d8588912SDave May     }
304d8588912SDave May     ierr = PetscFree(vs->m);CHKERRQ(ierr);
305d8588912SDave May   }
306bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
307d8588912SDave May 
308c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C",   "",0);CHKERRQ(ierr);
3090782ca92SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMat_C",   "",0);CHKERRQ(ierr);
310c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C",  "",0);CHKERRQ(ierr);
311c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C",     "",0);CHKERRQ(ierr);
312c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C",  "",0);CHKERRQ(ierr);
313c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C",   "",0);CHKERRQ(ierr);
314d8588912SDave May   PetscFunctionReturn(0);
315d8588912SDave May }
316d8588912SDave May 
317d8588912SDave May #undef __FUNCT__
318d8588912SDave May #define __FUNCT__ "MatAssemblyBegin_Nest"
319207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
320d8588912SDave May {
321d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
322d8588912SDave May   PetscInt       i,j;
323d8588912SDave May   PetscErrorCode ierr;
324d8588912SDave May 
325d8588912SDave May   PetscFunctionBegin;
326d8588912SDave May   for (i=0; i<vs->nr; i++) {
327d8588912SDave May     for (j=0; j<vs->nc; j++) {
328d8588912SDave May       if (vs->m[i][j]) { ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); }
329d8588912SDave May     }
330d8588912SDave May   }
331d8588912SDave May   PetscFunctionReturn(0);
332d8588912SDave May }
333d8588912SDave May 
334d8588912SDave May #undef __FUNCT__
335d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest"
336207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
337d8588912SDave May {
338d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
339d8588912SDave May   PetscInt       i,j;
340d8588912SDave May   PetscErrorCode ierr;
341d8588912SDave May 
342d8588912SDave May   PetscFunctionBegin;
343d8588912SDave May   for (i=0; i<vs->nr; i++) {
344d8588912SDave May     for (j=0; j<vs->nc; j++) {
345d8588912SDave May       if (vs->m[i][j]) { ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); }
346d8588912SDave May     }
347d8588912SDave May   }
348d8588912SDave May   PetscFunctionReturn(0);
349d8588912SDave May }
350d8588912SDave May 
351d8588912SDave May #undef __FUNCT__
352f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatRow"
353f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
354d8588912SDave May {
355207556f9SJed Brown   PetscErrorCode ierr;
356f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
357f349c1fdSJed Brown   PetscInt       j;
358f349c1fdSJed Brown   Mat            sub;
359d8588912SDave May 
360d8588912SDave May   PetscFunctionBegin;
3618188e55aSJed Brown   sub = (row < vs->nc) ? vs->m[row][row] : PETSC_NULL; /* Prefer to find on the diagonal */
362f349c1fdSJed Brown   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
3638188e55aSJed Brown   if (sub) {ierr = MatPreallocated(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */
364f349c1fdSJed Brown   *B = sub;
365f349c1fdSJed Brown   PetscFunctionReturn(0);
366d8588912SDave May }
367d8588912SDave May 
368f349c1fdSJed Brown #undef __FUNCT__
369f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatCol"
370f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
371f349c1fdSJed Brown {
372207556f9SJed Brown   PetscErrorCode ierr;
373f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
374f349c1fdSJed Brown   PetscInt       i;
375f349c1fdSJed Brown   Mat            sub;
376f349c1fdSJed Brown 
377f349c1fdSJed Brown   PetscFunctionBegin;
3788188e55aSJed Brown   sub = (col < vs->nr) ? vs->m[col][col] : PETSC_NULL; /* Prefer to find on the diagonal */
379f349c1fdSJed Brown   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
3808188e55aSJed Brown   if (sub) {ierr = MatPreallocated(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */
381f349c1fdSJed Brown   *B = sub;
382f349c1fdSJed Brown   PetscFunctionReturn(0);
383d8588912SDave May }
384d8588912SDave May 
385f349c1fdSJed Brown #undef __FUNCT__
386f349c1fdSJed Brown #define __FUNCT__ "MatNestFindIS"
387f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
388f349c1fdSJed Brown {
389f349c1fdSJed Brown   PetscErrorCode ierr;
390f349c1fdSJed Brown   PetscInt       i;
391f349c1fdSJed Brown   PetscBool      flg;
392f349c1fdSJed Brown 
393f349c1fdSJed Brown   PetscFunctionBegin;
394f349c1fdSJed Brown   PetscValidPointer(list,3);
395f349c1fdSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,4);
396f349c1fdSJed Brown   PetscValidIntPointer(found,5);
397f349c1fdSJed Brown   *found = -1;
398f349c1fdSJed Brown   for (i=0; i<n; i++) {
399207556f9SJed Brown     if (!list[i]) continue;
400f349c1fdSJed Brown     ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr);
401f349c1fdSJed Brown     if (flg) {
402f349c1fdSJed Brown       *found = i;
403f349c1fdSJed Brown       PetscFunctionReturn(0);
404f349c1fdSJed Brown     }
405f349c1fdSJed Brown   }
406f349c1fdSJed Brown   SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Could not find index set");
407f349c1fdSJed Brown   PetscFunctionReturn(0);
408f349c1fdSJed Brown }
409f349c1fdSJed Brown 
410f349c1fdSJed Brown #undef __FUNCT__
4118188e55aSJed Brown #define __FUNCT__ "MatNestGetRow"
4128188e55aSJed Brown /* Get a block row as a new MatNest */
4138188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B)
4148188e55aSJed Brown {
4158188e55aSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
4168188e55aSJed Brown   Mat            C;
4178188e55aSJed Brown   char           keyname[256];
4188188e55aSJed Brown   PetscErrorCode ierr;
4198188e55aSJed Brown 
4208188e55aSJed Brown   PetscFunctionBegin;
4218188e55aSJed Brown   *B = PETSC_NULL;
4228188e55aSJed Brown   ierr = PetscSNPrintf(keyname,sizeof keyname,"NestRow_%D",row);CHKERRQ(ierr);
4238188e55aSJed Brown   ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr);
4248188e55aSJed Brown   if (*B) PetscFunctionReturn(0);
4258188e55aSJed Brown 
4268188e55aSJed Brown   ierr = MatCreateNest(((PetscObject)A)->comm,1,PETSC_NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr);
4278188e55aSJed Brown   (*B)->assembled = A->assembled;
4288188e55aSJed Brown   ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr);
4298188e55aSJed Brown   ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */
4308188e55aSJed Brown   PetscFunctionReturn(0);
4318188e55aSJed Brown }
4328188e55aSJed Brown 
4338188e55aSJed Brown #undef __FUNCT__
434f349c1fdSJed Brown #define __FUNCT__ "MatNestFindSubMat"
435f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
436f349c1fdSJed Brown {
437f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
4388188e55aSJed Brown   PetscErrorCode ierr;
4398188e55aSJed Brown   PetscInt       row,col,i;
4408188e55aSJed Brown   IS             *iscopy;
4418188e55aSJed Brown   Mat            *matcopy;
4428188e55aSJed Brown   PetscBool      same,isFullCol;
443f349c1fdSJed Brown 
444f349c1fdSJed Brown   PetscFunctionBegin;
4458188e55aSJed Brown   /* Check if full column space. This is a hack */
4468188e55aSJed Brown   isFullCol = PETSC_FALSE;
4478188e55aSJed Brown   ierr = PetscTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr);
4488188e55aSJed Brown   if (same) {
4498188e55aSJed Brown     PetscInt n,first,step;
4508188e55aSJed Brown     ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
4518188e55aSJed Brown     ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
4528188e55aSJed Brown     if (A->cmap->n == n && A->cmap->rstart == first && step == 1) isFullCol = PETSC_TRUE;
4538188e55aSJed Brown   }
4548188e55aSJed Brown 
4558188e55aSJed Brown   if (isFullCol) {
4568188e55aSJed Brown     PetscInt row;
4578188e55aSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
4588188e55aSJed Brown     ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr);
4598188e55aSJed Brown   } else {
460f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
461f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr);
462f349c1fdSJed Brown     *B = vs->m[row][col];
4638188e55aSJed Brown   }
464f349c1fdSJed Brown   PetscFunctionReturn(0);
465f349c1fdSJed Brown }
466f349c1fdSJed Brown 
467f349c1fdSJed Brown #undef __FUNCT__
468f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest"
469f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
470f349c1fdSJed Brown {
471f349c1fdSJed Brown   PetscErrorCode ierr;
472f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
473f349c1fdSJed Brown   Mat            sub;
474f349c1fdSJed Brown 
475f349c1fdSJed Brown   PetscFunctionBegin;
476f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr);
477f349c1fdSJed Brown   switch (reuse) {
478f349c1fdSJed Brown   case MAT_INITIAL_MATRIX:
4797874fa86SDave May     if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); }
480f349c1fdSJed Brown     *B = sub;
481f349c1fdSJed Brown     break;
482f349c1fdSJed Brown   case MAT_REUSE_MATRIX:
483f349c1fdSJed Brown     if (sub != *B) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
484f349c1fdSJed Brown     break;
485f349c1fdSJed Brown   case MAT_IGNORE_MATRIX:       /* Nothing to do */
486f349c1fdSJed Brown     break;
487f349c1fdSJed Brown   }
488f349c1fdSJed Brown   PetscFunctionReturn(0);
489f349c1fdSJed Brown }
490f349c1fdSJed Brown 
491f349c1fdSJed Brown #undef __FUNCT__
492f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest"
493f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
494f349c1fdSJed Brown {
495f349c1fdSJed Brown   PetscErrorCode ierr;
496f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
497f349c1fdSJed Brown   Mat            sub;
498f349c1fdSJed Brown 
499f349c1fdSJed Brown   PetscFunctionBegin;
500f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
501f349c1fdSJed Brown   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
502f349c1fdSJed Brown   if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);}
503f349c1fdSJed Brown   *B = sub;
504d8588912SDave May   PetscFunctionReturn(0);
505d8588912SDave May }
506d8588912SDave May 
507d8588912SDave May #undef __FUNCT__
508d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest"
509207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
510d8588912SDave May {
511d8588912SDave May   PetscErrorCode ierr;
512f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
513f349c1fdSJed Brown   Mat            sub;
514d8588912SDave May 
515d8588912SDave May   PetscFunctionBegin;
516f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
517f349c1fdSJed Brown   if (*B != sub) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
518f349c1fdSJed Brown   if (sub) {
519f349c1fdSJed Brown     if (((PetscObject)sub)->refct <= 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
5206bf464f9SBarry Smith     ierr = MatDestroy(B);CHKERRQ(ierr);
521d8588912SDave May   }
522d8588912SDave May   PetscFunctionReturn(0);
523d8588912SDave May }
524d8588912SDave May 
525d8588912SDave May #undef __FUNCT__
5267874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest"
5277874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
5287874fa86SDave May {
5297874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
5307874fa86SDave May   Vec            *bdiag;
5317874fa86SDave May   PetscInt       i;
5327874fa86SDave May   PetscErrorCode ierr;
5337874fa86SDave May 
5347874fa86SDave May   PetscFunctionBegin;
535a0769a71SSatish Balay   /*  ierr = MatGetVecs_Nest(A,&diag,PETSC_NULL);CHKERRQ(ierr); */
536a0769a71SSatish Balay   /*  ierr = VecNestGetSubVecs(diag,PETSC_NULL,&bdiag);CHKERRQ(ierr); */
5377874fa86SDave May   ierr = VecNestGetSubVecs(v,PETSC_NULL,&bdiag);CHKERRQ(ierr);
5387874fa86SDave May   for (i=0; i<bA->nr; i++) {
5397874fa86SDave May     if (bA->m[i][i]) {
5407874fa86SDave May       ierr = MatGetDiagonal(bA->m[i][i],bdiag[i]);CHKERRQ(ierr);
5417874fa86SDave May     } else {
5427874fa86SDave May       ierr = VecSet(bdiag[i],1.0);CHKERRQ(ierr);
5437874fa86SDave May     }
5447874fa86SDave May   }
5457874fa86SDave May   PetscFunctionReturn(0);
5467874fa86SDave May }
5477874fa86SDave May 
5487874fa86SDave May #undef __FUNCT__
5497874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest2"
5507874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest2(Mat A,Vec v)
5517874fa86SDave May {
5527874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
5537874fa86SDave May   Vec            diag,*bdiag;
5547874fa86SDave May   VecScatter     *vscat;
5557874fa86SDave May   PetscInt       i;
5567874fa86SDave May   PetscErrorCode ierr;
5577874fa86SDave May 
5587874fa86SDave May   PetscFunctionBegin;
5597874fa86SDave May   ierr = MatGetVecs_Nest(A,&diag,PETSC_NULL);CHKERRQ(ierr);
5607874fa86SDave May   ierr = MatGetDiagonal_Nest(A,diag);CHKERRQ(ierr);
5617874fa86SDave May 
5627874fa86SDave May   /* scatter diag into v */
5637874fa86SDave May   ierr = VecNestGetSubVecs(diag,PETSC_NULL,&bdiag);CHKERRQ(ierr);
5647874fa86SDave May   ierr = PetscMalloc( sizeof(VecScatter) * bA->nr, &vscat );CHKERRQ(ierr);
5657874fa86SDave May   for (i=0; i<bA->nr; i++) {
5667874fa86SDave May     ierr = VecScatterCreate(v,bA->isglobal.row[i], bdiag[i],PETSC_NULL,&vscat[i]);CHKERRQ(ierr);
5677874fa86SDave May     ierr = VecScatterBegin(vscat[i],bdiag[i],v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5687874fa86SDave May   }
5697874fa86SDave May   for (i=0; i<bA->nr; i++) {
5707874fa86SDave May     ierr = VecScatterEnd(vscat[i],bdiag[i],v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
5717874fa86SDave May   }
5727874fa86SDave May   for (i=0; i<bA->nr; i++) {
5736bf464f9SBarry Smith     ierr = VecScatterDestroy(&vscat[i]);CHKERRQ(ierr);
5747874fa86SDave May   }
5757874fa86SDave May   ierr = PetscFree(vscat);CHKERRQ(ierr);
5766bf464f9SBarry Smith   ierr = VecDestroy(&diag);CHKERRQ(ierr);
5777874fa86SDave May   PetscFunctionReturn(0);
5787874fa86SDave May }
5797874fa86SDave May 
5807874fa86SDave May #undef __FUNCT__
5817874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest"
5827874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
5837874fa86SDave May {
5847874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
5857874fa86SDave May   Vec            *bl,*br;
5867874fa86SDave May   PetscInt       i,j;
5877874fa86SDave May   PetscErrorCode ierr;
5887874fa86SDave May 
5897874fa86SDave May   PetscFunctionBegin;
5907874fa86SDave May   ierr = VecNestGetSubVecs(l,PETSC_NULL,&bl);CHKERRQ(ierr);
5917874fa86SDave May   ierr = VecNestGetSubVecs(r,PETSC_NULL,&br);CHKERRQ(ierr);
5927874fa86SDave May   for (i=0; i<bA->nr; i++) {
5937874fa86SDave May     for (j=0; j<bA->nc; j++) {
5947874fa86SDave May       if (bA->m[i][j]) {
5957874fa86SDave May         ierr = MatDiagonalScale(bA->m[i][j],bl[i],br[j]);CHKERRQ(ierr);
5967874fa86SDave May       }
5977874fa86SDave May     }
5987874fa86SDave May   }
5997874fa86SDave May   PetscFunctionReturn(0);
6007874fa86SDave May }
6017874fa86SDave May 
6027874fa86SDave May #undef __FUNCT__
6037874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest2"
6047874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest2(Mat A,Vec l,Vec r)
6057874fa86SDave May {
6067874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
6077874fa86SDave May   Vec            bl,br,*ble,*bre;
6087874fa86SDave May   VecScatter     *vscatl,*vscatr;
6097874fa86SDave May   PetscInt       i;
6107874fa86SDave May   PetscErrorCode ierr;
6117874fa86SDave May 
6127874fa86SDave May   PetscFunctionBegin;
6137874fa86SDave May   /* scatter l,r into bl,br */
6147874fa86SDave May   ierr = MatGetVecs_Nest(A,&bl,&br);CHKERRQ(ierr);
6157874fa86SDave May   ierr = VecNestGetSubVecs(bl,PETSC_NULL,&ble);CHKERRQ(ierr);
6167874fa86SDave May   ierr = VecNestGetSubVecs(br,PETSC_NULL,&bre);CHKERRQ(ierr);
6177874fa86SDave May 
6187874fa86SDave May   /* row */
6197874fa86SDave May   ierr = PetscMalloc( sizeof(VecScatter) * bA->nr, &vscatl );CHKERRQ(ierr);
6207874fa86SDave May   for (i=0; i<bA->nr; i++) {
6217874fa86SDave May     ierr = VecScatterCreate(l,bA->isglobal.row[i], ble[i],PETSC_NULL,&vscatl[i]);CHKERRQ(ierr);
6227874fa86SDave May     ierr = VecScatterBegin(vscatl[i],l,ble[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6237874fa86SDave May   }
6247874fa86SDave May   for (i=0; i<bA->nr; i++) {
6257874fa86SDave May     ierr = VecScatterEnd(vscatl[i],l,ble[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6267874fa86SDave May   }
6277874fa86SDave May   /* col */
6287874fa86SDave May   ierr = PetscMalloc( sizeof(VecScatter) * bA->nc, &vscatr );CHKERRQ(ierr);
6297874fa86SDave May   for (i=0; i<bA->nc; i++) {
6307874fa86SDave May     ierr = VecScatterCreate(r,bA->isglobal.col[i], bre[i],PETSC_NULL,&vscatr[i]);CHKERRQ(ierr);
6317874fa86SDave May     ierr = VecScatterBegin(vscatr[i],l,bre[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6327874fa86SDave May   }
6337874fa86SDave May   for (i=0; i<bA->nc; i++) {
6347874fa86SDave May     ierr = VecScatterEnd(vscatr[i],r,bre[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
6357874fa86SDave May   }
6367874fa86SDave May 
6377874fa86SDave May   ierr = MatDiagonalScale_Nest(A,bl,br);CHKERRQ(ierr);
6387874fa86SDave May 
6397874fa86SDave May   for (i=0; i<bA->nr; i++) {
6406bf464f9SBarry Smith     ierr = VecScatterDestroy(&vscatl[i]);CHKERRQ(ierr);
6417874fa86SDave May   }
6427874fa86SDave May   for (i=0; i<bA->nc; i++) {
6436bf464f9SBarry Smith     ierr = VecScatterDestroy(&vscatr[i]);CHKERRQ(ierr);
6447874fa86SDave May   }
6457874fa86SDave May   ierr = PetscFree(vscatl);CHKERRQ(ierr);
6467874fa86SDave May   ierr = PetscFree(vscatr);CHKERRQ(ierr);
6476bf464f9SBarry Smith   ierr = VecDestroy(&bl);CHKERRQ(ierr);
6486bf464f9SBarry Smith   ierr = VecDestroy(&br);CHKERRQ(ierr);
6497874fa86SDave May   PetscFunctionReturn(0);
6507874fa86SDave May }
6517874fa86SDave May 
6527874fa86SDave May #undef __FUNCT__
653d8588912SDave May #define __FUNCT__ "MatGetVecs_Nest"
654207556f9SJed Brown static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left)
655d8588912SDave May {
656d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
657d8588912SDave May   Vec            *L,*R;
658d8588912SDave May   MPI_Comm       comm;
659d8588912SDave May   PetscInt       i,j;
660d8588912SDave May   PetscErrorCode ierr;
661d8588912SDave May 
662d8588912SDave May   PetscFunctionBegin;
663d8588912SDave May   comm = ((PetscObject)A)->comm;
664d8588912SDave May   if (right) {
665d8588912SDave May     /* allocate R */
666d8588912SDave May     ierr = PetscMalloc( sizeof(Vec) * bA->nc, &R );CHKERRQ(ierr);
667d8588912SDave May     /* Create the right vectors */
668d8588912SDave May     for (j=0; j<bA->nc; j++) {
669d8588912SDave May       for (i=0; i<bA->nr; i++) {
670d8588912SDave May         if (bA->m[i][j]) {
671d8588912SDave May           ierr = MatGetVecs(bA->m[i][j],&R[j],PETSC_NULL);CHKERRQ(ierr);
672d8588912SDave May           break;
673d8588912SDave May         }
674d8588912SDave May       }
675d8588912SDave May       if (i==bA->nr) {
676d8588912SDave May         /* have an empty column */
677d8588912SDave May         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
678d8588912SDave May       }
679d8588912SDave May     }
680f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
681d8588912SDave May     /* hand back control to the nest vector */
682d8588912SDave May     for (j=0; j<bA->nc; j++) {
6836bf464f9SBarry Smith       ierr = VecDestroy(&R[j]);CHKERRQ(ierr);
684d8588912SDave May     }
685d8588912SDave May     ierr = PetscFree(R);CHKERRQ(ierr);
686d8588912SDave May   }
687d8588912SDave May 
688d8588912SDave May   if (left) {
689d8588912SDave May     /* allocate L */
690d8588912SDave May     ierr = PetscMalloc( sizeof(Vec) * bA->nr, &L );CHKERRQ(ierr);
691d8588912SDave May     /* Create the left vectors */
692d8588912SDave May     for (i=0; i<bA->nr; i++) {
693d8588912SDave May       for (j=0; j<bA->nc; j++) {
694d8588912SDave May         if (bA->m[i][j]) {
695d8588912SDave May           ierr = MatGetVecs(bA->m[i][j],PETSC_NULL,&L[i]);CHKERRQ(ierr);
696d8588912SDave May           break;
697d8588912SDave May         }
698d8588912SDave May       }
699d8588912SDave May       if (j==bA->nc) {
700d8588912SDave May         /* have an empty row */
701d8588912SDave May         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
702d8588912SDave May       }
703d8588912SDave May     }
704d8588912SDave May 
705f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
706d8588912SDave May     for (i=0; i<bA->nr; i++) {
7076bf464f9SBarry Smith       ierr = VecDestroy(&L[i]);CHKERRQ(ierr);
708d8588912SDave May     }
709d8588912SDave May 
710d8588912SDave May     ierr = PetscFree(L);CHKERRQ(ierr);
711d8588912SDave May   }
712d8588912SDave May   PetscFunctionReturn(0);
713d8588912SDave May }
714d8588912SDave May 
715d8588912SDave May #undef __FUNCT__
716d8588912SDave May #define __FUNCT__ "MatView_Nest"
717207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
718d8588912SDave May {
719d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
720d8588912SDave May   PetscBool      isascii;
721d8588912SDave May   PetscInt       i,j;
722d8588912SDave May   PetscErrorCode ierr;
723d8588912SDave May 
724d8588912SDave May   PetscFunctionBegin;
725d8588912SDave May   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
726d8588912SDave May   if (isascii) {
727d8588912SDave May 
728d8588912SDave May     PetscViewerASCIIPrintf(viewer,"Matrix object: \n" );
729d8588912SDave May     PetscViewerASCIIPushTab( viewer );    /* push0 */
730d8588912SDave May     PetscViewerASCIIPrintf( viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc);
731d8588912SDave May 
732d8588912SDave May     PetscViewerASCIIPrintf(viewer,"MatNest structure: \n" );
733d8588912SDave May     for (i=0; i<bA->nr; i++) {
734d8588912SDave May       for (j=0; j<bA->nc; j++) {
735d8588912SDave May         const MatType type;
736d8588912SDave May         const char *name;
737d8588912SDave May         PetscInt NR,NC;
738d8588912SDave May         PetscBool isNest = PETSC_FALSE;
739d8588912SDave May 
740d8588912SDave May         if (!bA->m[i][j]) {
741d8588912SDave May           PetscViewerASCIIPrintf( viewer, "(%D,%D) : PETSC_NULL \n",i,j);
742d8588912SDave May           continue;
743d8588912SDave May         }
744d8588912SDave May         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
745d8588912SDave May         ierr = MatGetType( bA->m[i][j], &type );CHKERRQ(ierr);
746d8588912SDave May         name = ((PetscObject)bA->m[i][j])->prefix;
747d8588912SDave May         ierr = PetscTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
748d8588912SDave May 
749d8588912SDave May         PetscViewerASCIIPrintf( viewer, "(%D,%D) : name=\"%s\", type=%s, rows=%D, cols=%D \n",i,j,name,type,NR,NC);
750d8588912SDave May 
751d8588912SDave May         if (isNest) {
752d8588912SDave May           PetscViewerASCIIPushTab( viewer );  /* push1 */
753d8588912SDave May           ierr = MatView( bA->m[i][j], viewer );CHKERRQ(ierr);
754d8588912SDave May           PetscViewerASCIIPopTab(viewer);    /* pop1 */
755d8588912SDave May         }
756d8588912SDave May       }
757d8588912SDave May     }
758d8588912SDave May     PetscViewerASCIIPopTab(viewer);    /* pop0 */
759d8588912SDave May   }
760d8588912SDave May   PetscFunctionReturn(0);
761d8588912SDave May }
762d8588912SDave May 
763d8588912SDave May #undef __FUNCT__
764d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest"
765207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A)
766d8588912SDave May {
767d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
768d8588912SDave May   PetscInt       i,j;
769d8588912SDave May   PetscErrorCode ierr;
770d8588912SDave May 
771d8588912SDave May   PetscFunctionBegin;
772d8588912SDave May   for (i=0; i<bA->nr; i++) {
773d8588912SDave May     for (j=0; j<bA->nc; j++) {
774d8588912SDave May       if (!bA->m[i][j]) continue;
775d8588912SDave May       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
776d8588912SDave May     }
777d8588912SDave May   }
778d8588912SDave May   PetscFunctionReturn(0);
779d8588912SDave May }
780d8588912SDave May 
781d8588912SDave May #undef __FUNCT__
782d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest"
783207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
784d8588912SDave May {
785d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
786841e96a3SJed Brown   Mat            *b;
787841e96a3SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
788d8588912SDave May   PetscErrorCode ierr;
789d8588912SDave May 
790d8588912SDave May   PetscFunctionBegin;
791841e96a3SJed Brown   ierr = PetscMalloc(nr*nc*sizeof(Mat),&b);CHKERRQ(ierr);
792841e96a3SJed Brown   for (i=0; i<nr; i++) {
793841e96a3SJed Brown     for (j=0; j<nc; j++) {
794841e96a3SJed Brown       if (bA->m[i][j]) {
795841e96a3SJed Brown         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
796841e96a3SJed Brown       } else {
797841e96a3SJed Brown         b[i*nc+j] = PETSC_NULL;
798d8588912SDave May       }
799d8588912SDave May     }
800d8588912SDave May   }
801f349c1fdSJed Brown   ierr = MatCreateNest(((PetscObject)A)->comm,nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
802841e96a3SJed Brown   /* Give the new MatNest exclusive ownership */
803841e96a3SJed Brown   for (i=0; i<nr*nc; i++) {
8046bf464f9SBarry Smith     ierr = MatDestroy(&b[i]);CHKERRQ(ierr);
805d8588912SDave May   }
806d8588912SDave May   ierr = PetscFree(b);CHKERRQ(ierr);
807d8588912SDave May 
808841e96a3SJed Brown   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
809841e96a3SJed Brown   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
810d8588912SDave May   PetscFunctionReturn(0);
811d8588912SDave May }
812d8588912SDave May 
813d8588912SDave May /* nest api */
814d8588912SDave May EXTERN_C_BEGIN
815d8588912SDave May #undef __FUNCT__
816d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest"
817d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
818d8588912SDave May {
819d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
820d8588912SDave May   PetscFunctionBegin;
821d8588912SDave May   if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
822d8588912SDave May   if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
823d8588912SDave May   *mat = bA->m[idxm][jdxm];
824d8588912SDave May   PetscFunctionReturn(0);
825d8588912SDave May }
826d8588912SDave May EXTERN_C_END
827d8588912SDave May 
828d8588912SDave May #undef __FUNCT__
829d8588912SDave May #define __FUNCT__ "MatNestGetSubMat"
830d8588912SDave May /*@C
831d8588912SDave May  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
832d8588912SDave May 
833d8588912SDave May  Not collective
834d8588912SDave May 
835d8588912SDave May  Input Parameters:
836629881c0SJed Brown +   A  - nest matrix
837d8588912SDave May .   idxm - index of the matrix within the nest matrix
838629881c0SJed Brown -   jdxm - index of the matrix within the nest matrix
839d8588912SDave May 
840d8588912SDave May  Output Parameter:
841d8588912SDave May .   sub - matrix at index idxm,jdxm within the nest matrix
842d8588912SDave May 
843d8588912SDave May  Level: developer
844d8588912SDave May 
845d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMats()
846d8588912SDave May @*/
8477087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
848d8588912SDave May {
849699a902aSJed Brown   PetscErrorCode ierr;
850d8588912SDave May 
851d8588912SDave May   PetscFunctionBegin;
852699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
853d8588912SDave May   PetscFunctionReturn(0);
854d8588912SDave May }
855d8588912SDave May 
856d8588912SDave May EXTERN_C_BEGIN
857d8588912SDave May #undef __FUNCT__
8580782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat_Nest"
8590782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat)
8600782ca92SJed Brown {
8610782ca92SJed Brown   Mat_Nest *bA = (Mat_Nest*)A->data;
8620782ca92SJed Brown   PetscInt m,n,M,N,mi,ni,Mi,Ni;
8630782ca92SJed Brown   PetscErrorCode ierr;
8640782ca92SJed Brown 
8650782ca92SJed Brown   PetscFunctionBegin;
8660782ca92SJed Brown   if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
8670782ca92SJed Brown   if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
8680782ca92SJed Brown   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
8690782ca92SJed Brown   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
8700782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr);
8710782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr);
8720782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr);
8730782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr);
8740782ca92SJed Brown   if (M != Mi || N != Ni) SETERRQ4(((PetscObject)mat)->comm,PETSC_ERR_ARG_INCOMP,"Submatrix dimension (%D,%D) incompatible with nest block (%D,%D)",M,N,Mi,Ni);
8750782ca92SJed Brown   if (m != mi || n != ni) SETERRQ4(((PetscObject)mat)->comm,PETSC_ERR_ARG_INCOMP,"Submatrix local dimension (%D,%D) incompatible with nest block (%D,%D)",m,n,mi,ni);
8760782ca92SJed Brown   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
8770782ca92SJed Brown   ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr);
8780782ca92SJed Brown   bA->m[idxm][jdxm] = mat;
8790782ca92SJed Brown   PetscFunctionReturn(0);
8800782ca92SJed Brown }
8810782ca92SJed Brown EXTERN_C_END
8820782ca92SJed Brown 
8830782ca92SJed Brown #undef __FUNCT__
8840782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat"
8850782ca92SJed Brown /*@C
8860782ca92SJed Brown  MatNestSetSubMat - Set a single submatrix in the nest matrix.
8870782ca92SJed Brown 
8880782ca92SJed Brown  Logically collective on the submatrix communicator
8890782ca92SJed Brown 
8900782ca92SJed Brown  Input Parameters:
8910782ca92SJed Brown +   A  - nest matrix
8920782ca92SJed Brown .   idxm - index of the matrix within the nest matrix
8930782ca92SJed Brown .   jdxm - index of the matrix within the nest matrix
8940782ca92SJed Brown -   sub - matrix at index idxm,jdxm within the nest matrix
8950782ca92SJed Brown 
8960782ca92SJed Brown  Notes:
8970782ca92SJed Brown  The new submatrix must have the same size and communicator as that block of the nest.
8980782ca92SJed Brown 
8990782ca92SJed Brown  This increments the reference count of the submatrix.
9000782ca92SJed Brown 
9010782ca92SJed Brown  Level: developer
9020782ca92SJed Brown 
9030782ca92SJed Brown .seealso: MatNestSetSubMats(), MatNestGetSubMat()
9040782ca92SJed Brown @*/
9050782ca92SJed Brown PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
9060782ca92SJed Brown {
9070782ca92SJed Brown   PetscErrorCode ierr;
9080782ca92SJed Brown 
9090782ca92SJed Brown   PetscFunctionBegin;
9100782ca92SJed Brown   ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr);
9110782ca92SJed Brown   PetscFunctionReturn(0);
9120782ca92SJed Brown }
9130782ca92SJed Brown 
9140782ca92SJed Brown EXTERN_C_BEGIN
9150782ca92SJed Brown #undef __FUNCT__
916d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest"
917d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
918d8588912SDave May {
919d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
920d8588912SDave May   PetscFunctionBegin;
921d8588912SDave May   if (M)   { *M   = bA->nr; }
922d8588912SDave May   if (N)   { *N   = bA->nc; }
923d8588912SDave May   if (mat) { *mat = bA->m;  }
924d8588912SDave May   PetscFunctionReturn(0);
925d8588912SDave May }
926d8588912SDave May EXTERN_C_END
927d8588912SDave May 
928d8588912SDave May #undef __FUNCT__
929d8588912SDave May #define __FUNCT__ "MatNestGetSubMats"
930d8588912SDave May /*@C
931d8588912SDave May  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
932d8588912SDave May 
933d8588912SDave May  Not collective
934d8588912SDave May 
935d8588912SDave May  Input Parameters:
936629881c0SJed Brown .   A  - nest matrix
937d8588912SDave May 
938d8588912SDave May  Output Parameter:
939629881c0SJed Brown +   M - number of rows in the nest matrix
940d8588912SDave May .   N - number of cols in the nest matrix
941629881c0SJed Brown -   mat - 2d array of matrices
942d8588912SDave May 
943d8588912SDave May  Notes:
944d8588912SDave May 
945d8588912SDave May  The user should not free the array mat.
946d8588912SDave May 
947d8588912SDave May  Level: developer
948d8588912SDave May 
949d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMat()
950d8588912SDave May @*/
9517087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
952d8588912SDave May {
953699a902aSJed Brown   PetscErrorCode ierr;
954d8588912SDave May 
955d8588912SDave May   PetscFunctionBegin;
956699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
957d8588912SDave May   PetscFunctionReturn(0);
958d8588912SDave May }
959d8588912SDave May 
960d8588912SDave May EXTERN_C_BEGIN
961d8588912SDave May #undef __FUNCT__
962d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest"
9637087cfbeSBarry Smith PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
964d8588912SDave May {
965d8588912SDave May   Mat_Nest  *bA = (Mat_Nest*)A->data;
966d8588912SDave May 
967d8588912SDave May   PetscFunctionBegin;
968d8588912SDave May   if (M) { *M  = bA->nr; }
969d8588912SDave May   if (N) { *N  = bA->nc; }
970d8588912SDave May   PetscFunctionReturn(0);
971d8588912SDave May }
972d8588912SDave May EXTERN_C_END
973d8588912SDave May 
974d8588912SDave May #undef __FUNCT__
975d8588912SDave May #define __FUNCT__ "MatNestGetSize"
976d8588912SDave May /*@C
977d8588912SDave May  MatNestGetSize - Returns the size of the nest matrix.
978d8588912SDave May 
979d8588912SDave May  Not collective
980d8588912SDave May 
981d8588912SDave May  Input Parameters:
982d8588912SDave May .   A  - nest matrix
983d8588912SDave May 
984d8588912SDave May  Output Parameter:
985629881c0SJed Brown +   M - number of rows in the nested mat
986629881c0SJed Brown -   N - number of cols in the nested mat
987d8588912SDave May 
988d8588912SDave May  Notes:
989d8588912SDave May 
990d8588912SDave May  Level: developer
991d8588912SDave May 
992d8588912SDave May .seealso: MatNestGetSubMat(), MatNestGetSubMats()
993d8588912SDave May @*/
9947087cfbeSBarry Smith PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
995d8588912SDave May {
996699a902aSJed Brown   PetscErrorCode ierr;
997d8588912SDave May 
998d8588912SDave May   PetscFunctionBegin;
999699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
1000d8588912SDave May   PetscFunctionReturn(0);
1001d8588912SDave May }
1002d8588912SDave May 
1003207556f9SJed Brown EXTERN_C_BEGIN
1004207556f9SJed Brown #undef __FUNCT__
1005207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest"
10067087cfbeSBarry Smith PetscErrorCode  MatNestSetVecType_Nest(Mat A,const VecType vtype)
1007207556f9SJed Brown {
1008207556f9SJed Brown   PetscErrorCode ierr;
1009207556f9SJed Brown   PetscBool      flg;
1010207556f9SJed Brown 
1011207556f9SJed Brown   PetscFunctionBegin;
1012207556f9SJed Brown   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
1013207556f9SJed Brown   /* In reality, this only distinguishes VECNEST and "other" */
1014207556f9SJed Brown   A->ops->getvecs       = flg ? MatGetVecs_Nest       : 0;
10157874fa86SDave May   A->ops->getdiagonal   = flg ? MatGetDiagonal_Nest   : 0;
10167874fa86SDave May   A->ops->diagonalscale = flg ? MatDiagonalScale_Nest : 0;
10177874fa86SDave May 
1018207556f9SJed Brown  PetscFunctionReturn(0);
1019207556f9SJed Brown }
1020207556f9SJed Brown EXTERN_C_END
1021207556f9SJed Brown 
1022207556f9SJed Brown #undef __FUNCT__
1023207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType"
1024207556f9SJed Brown /*@C
1025207556f9SJed Brown  MatNestSetVecType - Sets the type of Vec returned by MatGetVecs()
1026207556f9SJed Brown 
1027207556f9SJed Brown  Not collective
1028207556f9SJed Brown 
1029207556f9SJed Brown  Input Parameters:
1030207556f9SJed Brown +  A  - nest matrix
1031207556f9SJed Brown -  vtype - type to use for creating vectors
1032207556f9SJed Brown 
1033207556f9SJed Brown  Notes:
1034207556f9SJed Brown 
1035207556f9SJed Brown  Level: developer
1036207556f9SJed Brown 
1037207556f9SJed Brown .seealso: MatGetVecs()
1038207556f9SJed Brown @*/
10397087cfbeSBarry Smith PetscErrorCode  MatNestSetVecType(Mat A,const VecType vtype)
1040207556f9SJed Brown {
1041207556f9SJed Brown   PetscErrorCode ierr;
1042207556f9SJed Brown 
1043207556f9SJed Brown   PetscFunctionBegin;
1044207556f9SJed Brown   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,const VecType),(A,vtype));CHKERRQ(ierr);
1045207556f9SJed Brown   PetscFunctionReturn(0);
1046207556f9SJed Brown }
1047207556f9SJed Brown 
1048c8883902SJed Brown EXTERN_C_BEGIN
1049d8588912SDave May #undef __FUNCT__
1050c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats_Nest"
1051c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1052d8588912SDave May {
1053c8883902SJed Brown   Mat_Nest       *s = (Mat_Nest*)A->data;
1054c8883902SJed Brown   PetscInt       i,j,m,n,M,N;
1055d8588912SDave May   PetscErrorCode ierr;
1056d8588912SDave May 
1057d8588912SDave May   PetscFunctionBegin;
1058c8883902SJed Brown   s->nr = nr;
1059c8883902SJed Brown   s->nc = nc;
1060d8588912SDave May 
1061c8883902SJed Brown   /* Create space for submatrices */
1062c8883902SJed Brown   ierr = PetscMalloc(sizeof(Mat*)*nr,&s->m);CHKERRQ(ierr);
1063c8883902SJed Brown   for (i=0; i<nr; i++) {
1064c8883902SJed Brown     ierr = PetscMalloc(sizeof(Mat)*nc,&s->m[i]);CHKERRQ(ierr);
1065d8588912SDave May   }
1066c8883902SJed Brown   for (i=0; i<nr; i++) {
1067c8883902SJed Brown     for (j=0; j<nc; j++) {
1068c8883902SJed Brown       s->m[i][j] = a[i*nc+j];
1069c8883902SJed Brown       if (a[i*nc+j]) {
1070c8883902SJed Brown         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
1071d8588912SDave May       }
1072d8588912SDave May     }
1073d8588912SDave May   }
1074d8588912SDave May 
10758188e55aSJed Brown   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
1076d8588912SDave May 
1077c8883902SJed Brown   ierr = PetscMalloc(sizeof(PetscInt)*nr,&s->row_len);CHKERRQ(ierr);
1078c8883902SJed Brown   ierr = PetscMalloc(sizeof(PetscInt)*nc,&s->col_len);CHKERRQ(ierr);
1079c8883902SJed Brown   for (i=0; i<nr; i++) s->row_len[i]=-1;
1080c8883902SJed Brown   for (j=0; j<nc; j++) s->col_len[j]=-1;
1081d8588912SDave May 
10828188e55aSJed Brown   ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr);
1083d8588912SDave May 
1084c8883902SJed Brown   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
1085c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
1086c8883902SJed Brown   ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr);
1087c8883902SJed Brown   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
1088c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
1089c8883902SJed Brown   ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr);
1090c8883902SJed Brown 
1091c8883902SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1092c8883902SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1093c8883902SJed Brown 
1094c8883902SJed Brown   ierr = PetscMalloc2(nr,Vec,&s->left,nc,Vec,&s->right);CHKERRQ(ierr);
1095c8883902SJed Brown   ierr = PetscMemzero(s->left,nr*sizeof(Vec));CHKERRQ(ierr);
1096c8883902SJed Brown   ierr = PetscMemzero(s->right,nc*sizeof(Vec));CHKERRQ(ierr);
1097d8588912SDave May   PetscFunctionReturn(0);
1098d8588912SDave May }
1099c8883902SJed Brown EXTERN_C_END
1100d8588912SDave May 
1101c8883902SJed Brown #undef __FUNCT__
1102c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats"
1103c8883902SJed Brown /*@
1104c8883902SJed Brown    MatNestSetSubMats - Sets the nested submatrices
1105c8883902SJed Brown 
1106c8883902SJed Brown    Collective on Mat
1107c8883902SJed Brown 
1108c8883902SJed Brown    Input Parameter:
1109c8883902SJed Brown +  N - nested matrix
1110c8883902SJed Brown .  nr - number of nested row blocks
1111c8883902SJed Brown .  is_row - index sets for each nested row block, or PETSC_NULL to make contiguous
1112c8883902SJed Brown .  nc - number of nested column blocks
1113c8883902SJed Brown .  is_col - index sets for each nested column block, or PETSC_NULL to make contiguous
1114c8883902SJed Brown -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL
1115c8883902SJed Brown 
1116c8883902SJed Brown    Level: advanced
1117c8883902SJed Brown 
1118c8883902SJed Brown .seealso: MatCreateNest(), MATNEST
1119c8883902SJed Brown @*/
1120c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1121c8883902SJed Brown {
1122c8883902SJed Brown   PetscErrorCode ierr;
1123c8883902SJed Brown   PetscInt i;
1124c8883902SJed Brown 
1125c8883902SJed Brown   PetscFunctionBegin;
1126c8883902SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1127c8883902SJed Brown   if (nr < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
1128c8883902SJed Brown   if (nr && is_row) {
1129c8883902SJed Brown     PetscValidPointer(is_row,3);
1130c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
1131c8883902SJed Brown   }
1132c8883902SJed Brown   if (nc < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
11331664e352SJed Brown   if (nc && is_col) {
1134c8883902SJed Brown     PetscValidPointer(is_col,5);
1135c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
1136c8883902SJed Brown   }
1137c8883902SJed Brown   if (nr*nc) PetscValidPointer(a,6);
1138c8883902SJed 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);
1139c8883902SJed Brown   PetscFunctionReturn(0);
1140c8883902SJed Brown }
1141d8588912SDave May 
1142d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1143d8588912SDave May /*
1144d8588912SDave May   nprocessors = NP
1145d8588912SDave May   Nest x^T = ( (g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1) )
1146d8588912SDave May        proc 0: => (g_0,h_0,)
1147d8588912SDave May        proc 1: => (g_1,h_1,)
1148d8588912SDave May        ...
1149d8588912SDave May        proc nprocs-1: => (g_NP-1,h_NP-1,)
1150d8588912SDave May 
1151d8588912SDave May             proc 0:                      proc 1:                    proc nprocs-1:
1152d8588912SDave May     is[0] = ( 0,1,2,...,nlocal(g_0)-1 )  ( 0,1,...,nlocal(g_1)-1 )  ( 0,1,...,nlocal(g_NP-1) )
1153d8588912SDave May 
1154d8588912SDave May             proc 0:
1155d8588912SDave May     is[1] = ( nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1 )
1156d8588912SDave May             proc 1:
1157d8588912SDave May     is[1] = ( nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1 )
1158d8588912SDave May 
1159d8588912SDave May             proc NP-1:
1160d8588912SDave May     is[1] = ( nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1 )
1161d8588912SDave May */
1162d8588912SDave May #undef __FUNCT__
1163d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private"
1164841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1165d8588912SDave May {
1166e2d7f03fSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
11678188e55aSJed Brown   PetscInt       i,j,offset,n,nsum,bs;
1168d8588912SDave May   PetscErrorCode ierr;
11692ae74bdbSJed Brown   Mat            sub;
1170d8588912SDave May 
1171d8588912SDave May   PetscFunctionBegin;
11728188e55aSJed Brown   ierr = PetscMalloc(sizeof(IS)*nr,&vs->isglobal.row);CHKERRQ(ierr);
11738188e55aSJed Brown   ierr = PetscMalloc(sizeof(IS)*nc,&vs->isglobal.col);CHKERRQ(ierr);
1174d8588912SDave May   if (is_row) { /* valid IS is passed in */
1175d8588912SDave May     /* refs on is[] are incremeneted */
1176e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1177d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
1178e2d7f03fSJed Brown       vs->isglobal.row[i] = is_row[i];
1179d8588912SDave May     }
11802ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
11818188e55aSJed Brown     nsum = 0;
11828188e55aSJed Brown     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
11838188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
11848188e55aSJed Brown       if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
11858188e55aSJed Brown       ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr);
11868188e55aSJed Brown       if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
11878188e55aSJed Brown       nsum += n;
11888188e55aSJed Brown     }
11898188e55aSJed Brown     offset = 0;
11908188e55aSJed Brown     ierr = MPI_Exscan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
1191e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1192f349c1fdSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
11932ae74bdbSJed Brown       ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr);
11942ae74bdbSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1195e2d7f03fSJed Brown       ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
1196e2d7f03fSJed Brown       ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
11972ae74bdbSJed Brown       offset += n;
1198d8588912SDave May     }
1199d8588912SDave May   }
1200d8588912SDave May 
1201d8588912SDave May   if (is_col) { /* valid IS is passed in */
1202d8588912SDave May     /* refs on is[] are incremeneted */
1203e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1204d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
1205e2d7f03fSJed Brown       vs->isglobal.col[j] = is_col[j];
1206d8588912SDave May     }
12072ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
12082ae74bdbSJed Brown     offset = A->cmap->rstart;
12098188e55aSJed Brown     nsum = 0;
12108188e55aSJed Brown     for (j=0; j<vs->nc; j++) {
12118188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
12128188e55aSJed Brown       if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
12138188e55aSJed Brown       ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr);
12148188e55aSJed Brown       if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
12158188e55aSJed Brown       nsum += n;
12168188e55aSJed Brown     }
12178188e55aSJed Brown     offset = 0;
12188188e55aSJed Brown     ierr = MPI_Exscan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
1219e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1220f349c1fdSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
12212ae74bdbSJed Brown       ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr);
12222ae74bdbSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1223e2d7f03fSJed Brown       ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
1224e2d7f03fSJed Brown       ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
12252ae74bdbSJed Brown       offset += n;
1226d8588912SDave May     }
1227d8588912SDave May   }
1228e2d7f03fSJed Brown 
1229e2d7f03fSJed Brown   /* Set up the local ISs */
1230e2d7f03fSJed Brown   ierr = PetscMalloc(vs->nr*sizeof(IS),&vs->islocal.row);CHKERRQ(ierr);
1231e2d7f03fSJed Brown   ierr = PetscMalloc(vs->nc*sizeof(IS),&vs->islocal.col);CHKERRQ(ierr);
1232e2d7f03fSJed Brown   for (i=0,offset=0; i<vs->nr; i++) {
1233e2d7f03fSJed Brown     IS                     isloc;
12348188e55aSJed Brown     ISLocalToGlobalMapping rmap = PETSC_NULL;
1235e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1236e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
12378188e55aSJed Brown     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,PETSC_NULL);CHKERRQ(ierr);}
1238207556f9SJed Brown     if (rmap) {
1239e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1240e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
1241e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1242e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1243207556f9SJed Brown     } else {
1244207556f9SJed Brown       nlocal = 0;
1245207556f9SJed Brown       isloc  = PETSC_NULL;
1246207556f9SJed Brown     }
1247e2d7f03fSJed Brown     vs->islocal.row[i] = isloc;
1248e2d7f03fSJed Brown     offset += nlocal;
1249e2d7f03fSJed Brown   }
12508188e55aSJed Brown   for (i=0,offset=0; i<vs->nc; i++) {
1251e2d7f03fSJed Brown     IS                     isloc;
12528188e55aSJed Brown     ISLocalToGlobalMapping cmap = PETSC_NULL;
1253e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1254e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
12558188e55aSJed Brown     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,PETSC_NULL,&cmap);CHKERRQ(ierr);}
1256207556f9SJed Brown     if (cmap) {
1257e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1258e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr);
1259e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1260e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1261207556f9SJed Brown     } else {
1262207556f9SJed Brown       nlocal = 0;
1263207556f9SJed Brown       isloc  = PETSC_NULL;
1264207556f9SJed Brown     }
1265e2d7f03fSJed Brown     vs->islocal.col[i] = isloc;
1266e2d7f03fSJed Brown     offset += nlocal;
1267e2d7f03fSJed Brown   }
12680189643fSJed Brown 
12690189643fSJed Brown #if defined(PETSC_USE_DEBUG)
12700189643fSJed Brown   for (i=0; i<vs->nr; i++) {
12710189643fSJed Brown     for (j=0; j<vs->nc; j++) {
12720189643fSJed Brown       PetscInt m,n,M,N,mi,ni,Mi,Ni;
12730189643fSJed Brown       Mat B = vs->m[i][j];
12740189643fSJed Brown       if (!B) continue;
12750189643fSJed Brown       ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr);
12760189643fSJed Brown       ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr);
12770189643fSJed Brown       ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr);
12780189643fSJed Brown       ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr);
12790189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr);
12800189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr);
12810189643fSJed Brown       if (M != Mi || N != Ni) SETERRQ6(((PetscObject)sub)->comm,PETSC_ERR_ARG_INCOMP,"Global sizes (%D,%D) of nested submatrix (%D,%D) do not agree with space defined by index sets (%D,%D)",M,N,i,j,Mi,Ni);
12820189643fSJed Brown       if (m != mi || n != ni) SETERRQ6(((PetscObject)sub)->comm,PETSC_ERR_ARG_INCOMP,"Local sizes (%D,%D) of nested submatrix (%D,%D) do not agree with space defined by index sets (%D,%D)",m,n,i,j,mi,ni);
12830189643fSJed Brown     }
12840189643fSJed Brown   }
12850189643fSJed Brown #endif
1286d8588912SDave May   PetscFunctionReturn(0);
1287d8588912SDave May }
1288d8588912SDave May 
1289d8588912SDave May #undef __FUNCT__
1290d8588912SDave May #define __FUNCT__ "MatCreateNest"
1291659c6bb0SJed Brown /*@
1292659c6bb0SJed Brown    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
1293659c6bb0SJed Brown 
1294659c6bb0SJed Brown    Collective on Mat
1295659c6bb0SJed Brown 
1296659c6bb0SJed Brown    Input Parameter:
1297659c6bb0SJed Brown +  comm - Communicator for the new Mat
1298659c6bb0SJed Brown .  nr - number of nested row blocks
1299659c6bb0SJed Brown .  is_row - index sets for each nested row block, or PETSC_NULL to make contiguous
1300659c6bb0SJed Brown .  nc - number of nested column blocks
1301659c6bb0SJed Brown .  is_col - index sets for each nested column block, or PETSC_NULL to make contiguous
1302659c6bb0SJed Brown -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL
1303659c6bb0SJed Brown 
1304659c6bb0SJed Brown    Output Parameter:
1305659c6bb0SJed Brown .  B - new matrix
1306659c6bb0SJed Brown 
1307659c6bb0SJed Brown    Level: advanced
1308659c6bb0SJed Brown 
1309659c6bb0SJed Brown .seealso: MatCreate(), VecCreateNest(), DMGetMatrix(), MATNEST
1310659c6bb0SJed Brown @*/
13117087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1312d8588912SDave May {
1313d8588912SDave May   Mat            A;
1314d8588912SDave May   PetscErrorCode ierr;
1315d8588912SDave May 
1316d8588912SDave May   PetscFunctionBegin;
1317c8883902SJed Brown   *B = 0;
1318d8588912SDave May   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1319c8883902SJed Brown   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
1320c8883902SJed Brown   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1321d8588912SDave May   *B = A;
1322d8588912SDave May   PetscFunctionReturn(0);
1323d8588912SDave May }
1324659c6bb0SJed Brown 
1325659c6bb0SJed Brown /*MC
1326659c6bb0SJed Brown   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
1327659c6bb0SJed Brown 
1328659c6bb0SJed Brown   Level: intermediate
1329659c6bb0SJed Brown 
1330659c6bb0SJed Brown   Notes:
1331659c6bb0SJed Brown   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1332659c6bb0SJed Brown   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1333659c6bb0SJed Brown   It is usually used with DMComposite and DMGetMatrix()
1334659c6bb0SJed Brown 
1335659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest()
1336659c6bb0SJed Brown M*/
1337c8883902SJed Brown EXTERN_C_BEGIN
1338c8883902SJed Brown #undef __FUNCT__
1339c8883902SJed Brown #define __FUNCT__ "MatCreate_Nest"
1340c8883902SJed Brown PetscErrorCode MatCreate_Nest(Mat A)
1341c8883902SJed Brown {
1342c8883902SJed Brown   Mat_Nest       *s;
1343c8883902SJed Brown   PetscErrorCode ierr;
1344c8883902SJed Brown 
1345c8883902SJed Brown   PetscFunctionBegin;
1346c8883902SJed Brown   ierr = PetscNewLog(A,Mat_Nest,&s);CHKERRQ(ierr);
1347c8883902SJed Brown   A->data         = (void*)s;
1348c8883902SJed Brown   s->nr = s->nc   = -1;
1349c8883902SJed Brown   s->m            = PETSC_NULL;
1350c8883902SJed Brown 
1351c8883902SJed Brown   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
1352c8883902SJed Brown   A->ops->mult                  = MatMult_Nest;
1353*9194d70fSJed Brown   A->ops->multadd               = MatMultAdd_Nest;
1354c8883902SJed Brown   A->ops->multtranspose         = MatMultTranspose_Nest;
1355*9194d70fSJed Brown   A->ops->multtransposeadd      = MatMultTransposeAdd_Nest;
1356c8883902SJed Brown   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1357c8883902SJed Brown   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1358c8883902SJed Brown   A->ops->zeroentries           = MatZeroEntries_Nest;
1359c8883902SJed Brown   A->ops->duplicate             = MatDuplicate_Nest;
1360c8883902SJed Brown   A->ops->getsubmatrix          = MatGetSubMatrix_Nest;
1361c8883902SJed Brown   A->ops->destroy               = MatDestroy_Nest;
1362c8883902SJed Brown   A->ops->view                  = MatView_Nest;
1363c8883902SJed Brown   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1364c8883902SJed Brown   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1365c8883902SJed Brown   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
13667874fa86SDave May   A->ops->getdiagonal           = MatGetDiagonal_Nest2; /* VECNEST version activated by calling MatNestSetVecType(A,VECNEST) */
13677874fa86SDave May   A->ops->diagonalscale         = MatDiagonalScale_Nest2; /* VECNEST version activated by calling MatNestSetVecType(A,VECNEST) */
1368c8883902SJed Brown 
1369c8883902SJed Brown   A->spptr        = 0;
1370c8883902SJed Brown   A->same_nonzero = PETSC_FALSE;
1371c8883902SJed Brown   A->assembled    = PETSC_FALSE;
1372c8883902SJed Brown 
1373c8883902SJed Brown   /* expose Nest api's */
1374c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C",   "MatNestGetSubMat_Nest",   MatNestGetSubMat_Nest);CHKERRQ(ierr);
13750782ca92SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMat_C",   "MatNestSetSubMat_Nest",   MatNestSetSubMat_Nest);CHKERRQ(ierr);
1376c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C",  "MatNestGetSubMats_Nest",  MatNestGetSubMats_Nest);CHKERRQ(ierr);
1377c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C",     "MatNestGetSize_Nest",     MatNestGetSize_Nest);CHKERRQ(ierr);
1378c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C",  "MatNestSetVecType_Nest",  MatNestSetVecType_Nest);CHKERRQ(ierr);
1379c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C",  "MatNestSetSubMats_Nest",  MatNestSetSubMats_Nest);CHKERRQ(ierr);
1380c8883902SJed Brown 
1381c8883902SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
1382c8883902SJed Brown   PetscFunctionReturn(0);
1383c8883902SJed Brown }
138486e80854SHong Zhang EXTERN_C_END
1385