xref: /petsc/src/mat/impls/nest/matnest.c (revision f8170845fd294a162af8260e56f6608dda617098)
1d8588912SDave May 
2aaa7dc30SBarry Smith #include <../src/mat/impls/nest/matnestimpl.h> /*I   "petscmat.h"   I*/
30c312b8eSJed Brown #include <petscsf.h>
4d8588912SDave May 
5c8883902SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]);
62a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left);
7c8883902SJed Brown 
8d8588912SDave May /* private functions */
9d8588912SDave May #undef __FUNCT__
108188e55aSJed Brown #define __FUNCT__ "MatNestGetSizes_Private"
118188e55aSJed Brown static PetscErrorCode MatNestGetSizes_Private(Mat A,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N)
12d8588912SDave May {
13d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
148188e55aSJed Brown   PetscInt       i,j;
15d8588912SDave May   PetscErrorCode ierr;
16d8588912SDave May 
17d8588912SDave May   PetscFunctionBegin;
188188e55aSJed Brown   *m = *n = *M = *N = 0;
198188e55aSJed Brown   for (i=0; i<bA->nr; i++) {  /* rows */
208188e55aSJed Brown     PetscInt sm,sM;
218188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.row[i],&sm);CHKERRQ(ierr);
228188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.row[i],&sM);CHKERRQ(ierr);
238188e55aSJed Brown     *m  += sm;
248188e55aSJed Brown     *M  += sM;
25d8588912SDave May   }
268188e55aSJed Brown   for (j=0; j<bA->nc; j++) {  /* cols */
278188e55aSJed Brown     PetscInt sn,sN;
288188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.col[j],&sn);CHKERRQ(ierr);
298188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.col[j],&sN);CHKERRQ(ierr);
308188e55aSJed Brown     *n  += sn;
318188e55aSJed Brown     *N  += sN;
32d8588912SDave May   }
33d8588912SDave May   PetscFunctionReturn(0);
34d8588912SDave May }
35d8588912SDave May 
36d8588912SDave May /* operations */
37d8588912SDave May #undef __FUNCT__
38d8588912SDave May #define __FUNCT__ "MatMult_Nest"
39207556f9SJed Brown static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y)
40d8588912SDave May {
41d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
42207556f9SJed Brown   Vec            *bx = bA->right,*by = bA->left;
43207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
44d8588912SDave May   PetscErrorCode ierr;
45d8588912SDave May 
46d8588912SDave May   PetscFunctionBegin;
47207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
48207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
49207556f9SJed Brown   for (i=0; i<nr; i++) {
50d8588912SDave May     ierr = VecZeroEntries(by[i]);CHKERRQ(ierr);
51207556f9SJed Brown     for (j=0; j<nc; j++) {
52207556f9SJed Brown       if (!bA->m[i][j]) continue;
53d8588912SDave May       /* y[i] <- y[i] + A[i][j] * x[j] */
54d8588912SDave May       ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr);
55d8588912SDave May     }
56d8588912SDave May   }
57207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
58207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
59d8588912SDave May   PetscFunctionReturn(0);
60d8588912SDave May }
61d8588912SDave May 
62d8588912SDave May #undef __FUNCT__
639194d70fSJed Brown #define __FUNCT__ "MatMultAdd_Nest"
649194d70fSJed Brown static PetscErrorCode MatMultAdd_Nest(Mat A,Vec x,Vec y,Vec z)
659194d70fSJed Brown {
669194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
679194d70fSJed Brown   Vec            *bx = bA->right,*bz = bA->left;
689194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
699194d70fSJed Brown   PetscErrorCode ierr;
709194d70fSJed Brown 
719194d70fSJed Brown   PetscFunctionBegin;
729194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
739194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
749194d70fSJed Brown   for (i=0; i<nr; i++) {
759194d70fSJed Brown     if (y != z) {
769194d70fSJed Brown       Vec by;
779194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
789194d70fSJed Brown       ierr = VecCopy(by,bz[i]);CHKERRQ(ierr);
79336d21e7SJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
809194d70fSJed Brown     }
819194d70fSJed Brown     for (j=0; j<nc; j++) {
829194d70fSJed Brown       if (!bA->m[i][j]) continue;
839194d70fSJed Brown       /* y[i] <- y[i] + A[i][j] * x[j] */
849194d70fSJed Brown       ierr = MatMultAdd(bA->m[i][j],bx[j],bz[i],bz[i]);CHKERRQ(ierr);
859194d70fSJed Brown     }
869194d70fSJed Brown   }
879194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
889194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
899194d70fSJed Brown   PetscFunctionReturn(0);
909194d70fSJed Brown }
919194d70fSJed Brown 
929194d70fSJed Brown #undef __FUNCT__
93d8588912SDave May #define __FUNCT__ "MatMultTranspose_Nest"
94207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y)
95d8588912SDave May {
96d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
97207556f9SJed Brown   Vec            *bx = bA->left,*by = bA->right;
98207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
99d8588912SDave May   PetscErrorCode ierr;
100d8588912SDave May 
101d8588912SDave May   PetscFunctionBegin;
102609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
103609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
104207556f9SJed Brown   for (j=0; j<nc; j++) {
105609e31cbSJed Brown     ierr = VecZeroEntries(by[j]);CHKERRQ(ierr);
106609e31cbSJed Brown     for (i=0; i<nr; i++) {
1076c75ac25SJed Brown       if (!bA->m[i][j]) continue;
108609e31cbSJed Brown       /* y[j] <- y[j] + (A[i][j])^T * x[i] */
109609e31cbSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],by[j],by[j]);CHKERRQ(ierr);
110d8588912SDave May     }
111d8588912SDave May   }
112609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
113609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
114d8588912SDave May   PetscFunctionReturn(0);
115d8588912SDave May }
116d8588912SDave May 
117d8588912SDave May #undef __FUNCT__
1189194d70fSJed Brown #define __FUNCT__ "MatMultTransposeAdd_Nest"
1199194d70fSJed Brown static PetscErrorCode MatMultTransposeAdd_Nest(Mat A,Vec x,Vec y,Vec z)
1209194d70fSJed Brown {
1219194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
1229194d70fSJed Brown   Vec            *bx = bA->left,*bz = bA->right;
1239194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
1249194d70fSJed Brown   PetscErrorCode ierr;
1259194d70fSJed Brown 
1269194d70fSJed Brown   PetscFunctionBegin;
1279194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
1289194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
1299194d70fSJed Brown   for (j=0; j<nc; j++) {
1309194d70fSJed Brown     if (y != z) {
1319194d70fSJed Brown       Vec by;
1329194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
1339194d70fSJed Brown       ierr = VecCopy(by,bz[j]);CHKERRQ(ierr);
1349194d70fSJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
1359194d70fSJed Brown     }
1369194d70fSJed Brown     for (i=0; i<nr; i++) {
1376c75ac25SJed Brown       if (!bA->m[i][j]) continue;
1389194d70fSJed Brown       /* z[j] <- y[j] + (A[i][j])^T * x[i] */
1399194d70fSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],bz[j],bz[j]);CHKERRQ(ierr);
1409194d70fSJed Brown     }
1419194d70fSJed Brown   }
1429194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
1439194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
1449194d70fSJed Brown   PetscFunctionReturn(0);
1459194d70fSJed Brown }
1469194d70fSJed Brown 
1479194d70fSJed Brown #undef __FUNCT__
148*f8170845SAlex Fikl #define __FUNCT__ "MatTranspose_Nest"
149*f8170845SAlex Fikl static PetscErrorCode MatTranspose_Nest(Mat A,MatReuse reuse,Mat *B)
150*f8170845SAlex Fikl {
151*f8170845SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data, *bC;
152*f8170845SAlex Fikl   Mat            C;
153*f8170845SAlex Fikl   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
154*f8170845SAlex Fikl   PetscErrorCode ierr;
155*f8170845SAlex Fikl 
156*f8170845SAlex Fikl   PetscFunctionBegin;
157*f8170845SAlex Fikl   if (reuse == MAT_REUSE_MATRIX && A == *B && nr != nc) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Square nested matrix only for in-place");
158*f8170845SAlex Fikl 
159*f8170845SAlex Fikl   if (reuse == MAT_INITIAL_MATRIX || A == *B) {
160*f8170845SAlex Fikl     Mat *subs;
161*f8170845SAlex Fikl     IS  *is_row,*is_col;
162*f8170845SAlex Fikl 
163*f8170845SAlex Fikl     ierr = PetscCalloc1(nr * nc,&subs);CHKERRQ(ierr);
164*f8170845SAlex Fikl     ierr = PetscMalloc2(nr,&is_row,nc,&is_col);CHKERRQ(ierr);
165*f8170845SAlex Fikl     ierr = MatNestGetISs(A,is_row,is_col);CHKERRQ(ierr);
166*f8170845SAlex Fikl     ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nc,is_col,nr,is_row,subs,&C);CHKERRQ(ierr);
167*f8170845SAlex Fikl     ierr = PetscFree(subs);CHKERRQ(ierr);
168*f8170845SAlex Fikl     ierr = PetscFree(is_row);CHKERRQ(ierr);
169*f8170845SAlex Fikl     ierr = PetscFree(is_col);CHKERRQ(ierr);
170*f8170845SAlex Fikl   } else {
171*f8170845SAlex Fikl     C = *B;
172*f8170845SAlex Fikl   }
173*f8170845SAlex Fikl 
174*f8170845SAlex Fikl   bC = (Mat_Nest*)C->data;
175*f8170845SAlex Fikl   for (i=0; i<nr; i++) {
176*f8170845SAlex Fikl     for (j=0; j<nc; j++) {
177*f8170845SAlex Fikl       if (bA->m[i][j]) {
178*f8170845SAlex Fikl         ierr = MatTranspose(bA->m[i][j], reuse, &(bC->m[j][i]));CHKERRQ(ierr);
179*f8170845SAlex Fikl       } else {
180*f8170845SAlex Fikl         bC->m[j][i] = NULL;
181*f8170845SAlex Fikl       }
182*f8170845SAlex Fikl     }
183*f8170845SAlex Fikl   }
184*f8170845SAlex Fikl 
185*f8170845SAlex Fikl   if (reuse == MAT_INITIAL_MATRIX || A != *B) {
186*f8170845SAlex Fikl     *B = C;
187*f8170845SAlex Fikl   } else {
188*f8170845SAlex Fikl     ierr = MatHeaderMerge(A, &C);CHKERRQ(ierr);
189*f8170845SAlex Fikl   }
190*f8170845SAlex Fikl   PetscFunctionReturn(0);
191*f8170845SAlex Fikl }
192*f8170845SAlex Fikl 
193*f8170845SAlex Fikl #undef __FUNCT__
194e2d7f03fSJed Brown #define __FUNCT__ "MatNestDestroyISList"
195e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list)
196e2d7f03fSJed Brown {
197e2d7f03fSJed Brown   PetscErrorCode ierr;
198e2d7f03fSJed Brown   IS             *lst = *list;
199e2d7f03fSJed Brown   PetscInt       i;
200e2d7f03fSJed Brown 
201e2d7f03fSJed Brown   PetscFunctionBegin;
202e2d7f03fSJed Brown   if (!lst) PetscFunctionReturn(0);
2036bf464f9SBarry Smith   for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);}
204e2d7f03fSJed Brown   ierr  = PetscFree(lst);CHKERRQ(ierr);
2050298fd71SBarry Smith   *list = NULL;
206e2d7f03fSJed Brown   PetscFunctionReturn(0);
207e2d7f03fSJed Brown }
208e2d7f03fSJed Brown 
209e2d7f03fSJed Brown #undef __FUNCT__
210d8588912SDave May #define __FUNCT__ "MatDestroy_Nest"
211207556f9SJed Brown static PetscErrorCode MatDestroy_Nest(Mat A)
212d8588912SDave May {
213d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
214d8588912SDave May   PetscInt       i,j;
215d8588912SDave May   PetscErrorCode ierr;
216d8588912SDave May 
217d8588912SDave May   PetscFunctionBegin;
218d8588912SDave May   /* release the matrices and the place holders */
219e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr);
220e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr);
221e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
222e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
223d8588912SDave May 
224d8588912SDave May   ierr = PetscFree(vs->row_len);CHKERRQ(ierr);
225d8588912SDave May   ierr = PetscFree(vs->col_len);CHKERRQ(ierr);
226d8588912SDave May 
227207556f9SJed Brown   ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr);
228207556f9SJed Brown 
229d8588912SDave May   /* release the matrices and the place holders */
230d8588912SDave May   if (vs->m) {
231d8588912SDave May     for (i=0; i<vs->nr; i++) {
232d8588912SDave May       for (j=0; j<vs->nc; j++) {
2336bf464f9SBarry Smith         ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr);
234d8588912SDave May       }
235d8588912SDave May       ierr = PetscFree(vs->m[i]);CHKERRQ(ierr);
236d8588912SDave May     }
237d8588912SDave May     ierr = PetscFree(vs->m);CHKERRQ(ierr);
238d8588912SDave May   }
239bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
240d8588912SDave May 
241bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",0);CHKERRQ(ierr);
242bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",0);CHKERRQ(ierr);
243bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",0);CHKERRQ(ierr);
244bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",0);CHKERRQ(ierr);
245bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",0);CHKERRQ(ierr);
246bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",0);CHKERRQ(ierr);
247bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",0);CHKERRQ(ierr);
248bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",0);CHKERRQ(ierr);
249d8588912SDave May   PetscFunctionReturn(0);
250d8588912SDave May }
251d8588912SDave May 
252d8588912SDave May #undef __FUNCT__
253d8588912SDave May #define __FUNCT__ "MatAssemblyBegin_Nest"
254207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
255d8588912SDave May {
256d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
257d8588912SDave May   PetscInt       i,j;
258d8588912SDave May   PetscErrorCode ierr;
259d8588912SDave May 
260d8588912SDave May   PetscFunctionBegin;
261d8588912SDave May   for (i=0; i<vs->nr; i++) {
262d8588912SDave May     for (j=0; j<vs->nc; j++) {
263e7c19651SJed Brown       if (vs->m[i][j]) {
264e7c19651SJed Brown         ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr);
265e7c19651SJed Brown         if (!vs->splitassembly) {
266e7c19651SJed Brown           /* Note: split assembly will fail if the same block appears more than once (even indirectly through a nested
267e7c19651SJed Brown            * sub-block). This could be fixed by adding a flag to Mat so that there was a way to check if a Mat was
268e7c19651SJed Brown            * already performing an assembly, but the result would by more complicated and appears to offer less
269e7c19651SJed Brown            * potential for diagnostics and correctness checking. Split assembly should be fixed once there is an
270e7c19651SJed Brown            * interface for libraries to make asynchronous progress in "user-defined non-blocking collectives".
271e7c19651SJed Brown            */
272e7c19651SJed Brown           ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr);
273e7c19651SJed Brown         }
274e7c19651SJed Brown       }
275d8588912SDave May     }
276d8588912SDave May   }
277d8588912SDave May   PetscFunctionReturn(0);
278d8588912SDave May }
279d8588912SDave May 
280d8588912SDave May #undef __FUNCT__
281d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest"
282207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
283d8588912SDave May {
284d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
285d8588912SDave May   PetscInt       i,j;
286d8588912SDave May   PetscErrorCode ierr;
287d8588912SDave May 
288d8588912SDave May   PetscFunctionBegin;
289d8588912SDave May   for (i=0; i<vs->nr; i++) {
290d8588912SDave May     for (j=0; j<vs->nc; j++) {
291e7c19651SJed Brown       if (vs->m[i][j]) {
292e7c19651SJed Brown         if (vs->splitassembly) {
293e7c19651SJed Brown           ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr);
294e7c19651SJed Brown         }
295e7c19651SJed Brown       }
296d8588912SDave May     }
297d8588912SDave May   }
298d8588912SDave May   PetscFunctionReturn(0);
299d8588912SDave May }
300d8588912SDave May 
301d8588912SDave May #undef __FUNCT__
302f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatRow"
303f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
304d8588912SDave May {
305207556f9SJed Brown   PetscErrorCode ierr;
306f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
307f349c1fdSJed Brown   PetscInt       j;
308f349c1fdSJed Brown   Mat            sub;
309d8588912SDave May 
310d8588912SDave May   PetscFunctionBegin;
3110298fd71SBarry Smith   sub = (row < vs->nc) ? vs->m[row][row] : (Mat)NULL; /* Prefer to find on the diagonal */
312f349c1fdSJed Brown   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
3134994cf47SJed Brown   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
314f349c1fdSJed Brown   *B = sub;
315f349c1fdSJed Brown   PetscFunctionReturn(0);
316d8588912SDave May }
317d8588912SDave May 
318f349c1fdSJed Brown #undef __FUNCT__
319f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatCol"
320f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
321f349c1fdSJed Brown {
322207556f9SJed Brown   PetscErrorCode ierr;
323f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
324f349c1fdSJed Brown   PetscInt       i;
325f349c1fdSJed Brown   Mat            sub;
326f349c1fdSJed Brown 
327f349c1fdSJed Brown   PetscFunctionBegin;
3280298fd71SBarry Smith   sub = (col < vs->nr) ? vs->m[col][col] : (Mat)NULL; /* Prefer to find on the diagonal */
329f349c1fdSJed Brown   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
3304994cf47SJed Brown   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
331f349c1fdSJed Brown   *B = sub;
332f349c1fdSJed Brown   PetscFunctionReturn(0);
333d8588912SDave May }
334d8588912SDave May 
335f349c1fdSJed Brown #undef __FUNCT__
336f349c1fdSJed Brown #define __FUNCT__ "MatNestFindIS"
337f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
338f349c1fdSJed Brown {
339f349c1fdSJed Brown   PetscErrorCode ierr;
340f349c1fdSJed Brown   PetscInt       i;
341f349c1fdSJed Brown   PetscBool      flg;
342f349c1fdSJed Brown 
343f349c1fdSJed Brown   PetscFunctionBegin;
344f349c1fdSJed Brown   PetscValidPointer(list,3);
345f349c1fdSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,4);
346f349c1fdSJed Brown   PetscValidIntPointer(found,5);
347f349c1fdSJed Brown   *found = -1;
348f349c1fdSJed Brown   for (i=0; i<n; i++) {
349207556f9SJed Brown     if (!list[i]) continue;
350f349c1fdSJed Brown     ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr);
351f349c1fdSJed Brown     if (flg) {
352f349c1fdSJed Brown       *found = i;
353f349c1fdSJed Brown       PetscFunctionReturn(0);
354f349c1fdSJed Brown     }
355f349c1fdSJed Brown   }
356ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Could not find index set");
357f349c1fdSJed Brown   PetscFunctionReturn(0);
358f349c1fdSJed Brown }
359f349c1fdSJed Brown 
360f349c1fdSJed Brown #undef __FUNCT__
3618188e55aSJed Brown #define __FUNCT__ "MatNestGetRow"
3628188e55aSJed Brown /* Get a block row as a new MatNest */
3638188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B)
3648188e55aSJed Brown {
3658188e55aSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
3668188e55aSJed Brown   char           keyname[256];
3678188e55aSJed Brown   PetscErrorCode ierr;
3688188e55aSJed Brown 
3698188e55aSJed Brown   PetscFunctionBegin;
3700298fd71SBarry Smith   *B   = NULL;
3718caf3d72SBarry Smith   ierr = PetscSNPrintf(keyname,sizeof(keyname),"NestRow_%D",row);CHKERRQ(ierr);
3728188e55aSJed Brown   ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr);
3738188e55aSJed Brown   if (*B) PetscFunctionReturn(0);
3748188e55aSJed Brown 
375ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),1,NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr);
37626fbe8dcSKarl Rupp 
3778188e55aSJed Brown   (*B)->assembled = A->assembled;
37826fbe8dcSKarl Rupp 
3798188e55aSJed Brown   ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr);
3808188e55aSJed Brown   ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */
3818188e55aSJed Brown   PetscFunctionReturn(0);
3828188e55aSJed Brown }
3838188e55aSJed Brown 
3848188e55aSJed Brown #undef __FUNCT__
385f349c1fdSJed Brown #define __FUNCT__ "MatNestFindSubMat"
386f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
387f349c1fdSJed Brown {
388f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
3898188e55aSJed Brown   PetscErrorCode ierr;
3906b3a5b13SJed Brown   PetscInt       row,col;
391e072481dSJed Brown   PetscBool      same,isFullCol,isFullColGlobal;
392f349c1fdSJed Brown 
393f349c1fdSJed Brown   PetscFunctionBegin;
3948188e55aSJed Brown   /* Check if full column space. This is a hack */
3958188e55aSJed Brown   isFullCol = PETSC_FALSE;
396251f4c67SDmitry Karpeev   ierr      = PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr);
3978188e55aSJed Brown   if (same) {
39877019fcaSJed Brown     PetscInt n,first,step,i,an,am,afirst,astep;
3998188e55aSJed Brown     ierr      = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
4008188e55aSJed Brown     ierr      = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
40177019fcaSJed Brown     isFullCol = PETSC_TRUE;
40205ce4453SJed Brown     for (i=0,an=A->cmap->rstart; i<vs->nc; i++) {
40377019fcaSJed Brown       ierr = ISStrideGetInfo(is->col[i],&afirst,&astep);CHKERRQ(ierr);
40477019fcaSJed Brown       ierr = ISGetLocalSize(is->col[i],&am);CHKERRQ(ierr);
40577019fcaSJed Brown       if (afirst != an || astep != step) isFullCol = PETSC_FALSE;
40677019fcaSJed Brown       an += am;
40777019fcaSJed Brown     }
40805ce4453SJed Brown     if (an != A->cmap->rstart+n) isFullCol = PETSC_FALSE;
4098188e55aSJed Brown   }
410b2566f29SBarry Smith   ierr = MPIU_Allreduce(&isFullCol,&isFullColGlobal,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)iscol));CHKERRQ(ierr);
4118188e55aSJed Brown 
412e072481dSJed Brown   if (isFullColGlobal) {
4138188e55aSJed Brown     PetscInt row;
4148188e55aSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
4158188e55aSJed Brown     ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr);
4168188e55aSJed Brown   } else {
417f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
418f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr);
419b6480e04SStefano Zampini     if (!vs->m[row][col]) {
420b6480e04SStefano Zampini       PetscInt lr,lc;
421b6480e04SStefano Zampini 
422b6480e04SStefano Zampini       ierr = MatCreate(PetscObjectComm((PetscObject)A),&vs->m[row][col]);CHKERRQ(ierr);
423b6480e04SStefano Zampini       ierr = ISGetLocalSize(vs->isglobal.row[row],&lr);CHKERRQ(ierr);
424b6480e04SStefano Zampini       ierr = ISGetLocalSize(vs->isglobal.col[col],&lc);CHKERRQ(ierr);
425b6480e04SStefano Zampini       ierr = MatSetSizes(vs->m[row][col],lr,lc,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
426b6480e04SStefano Zampini       ierr = MatSetUp(vs->m[row][col]);CHKERRQ(ierr);
427b6480e04SStefano Zampini       ierr = MatAssemblyBegin(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
428b6480e04SStefano Zampini       ierr = MatAssemblyEnd(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
429b6480e04SStefano Zampini     }
430f349c1fdSJed Brown     *B = vs->m[row][col];
4318188e55aSJed Brown   }
432f349c1fdSJed Brown   PetscFunctionReturn(0);
433f349c1fdSJed Brown }
434f349c1fdSJed Brown 
435f349c1fdSJed Brown #undef __FUNCT__
436f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest"
437f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
438f349c1fdSJed Brown {
439f349c1fdSJed Brown   PetscErrorCode ierr;
440f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
441f349c1fdSJed Brown   Mat            sub;
442f349c1fdSJed Brown 
443f349c1fdSJed Brown   PetscFunctionBegin;
444f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr);
445f349c1fdSJed Brown   switch (reuse) {
446f349c1fdSJed Brown   case MAT_INITIAL_MATRIX:
4477874fa86SDave May     if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); }
448f349c1fdSJed Brown     *B = sub;
449f349c1fdSJed Brown     break;
450f349c1fdSJed Brown   case MAT_REUSE_MATRIX:
451ce94432eSBarry Smith     if (sub != *B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
452f349c1fdSJed Brown     break;
453f349c1fdSJed Brown   case MAT_IGNORE_MATRIX:       /* Nothing to do */
454f349c1fdSJed Brown     break;
455511c6705SHong Zhang   case MAT_INPLACE_MATRIX:       /* Nothing to do */
456511c6705SHong Zhang     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX is not supported yet");
457511c6705SHong Zhang     break;
458f349c1fdSJed Brown   }
459f349c1fdSJed Brown   PetscFunctionReturn(0);
460f349c1fdSJed Brown }
461f349c1fdSJed Brown 
462f349c1fdSJed Brown #undef __FUNCT__
463f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest"
464f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
465f349c1fdSJed Brown {
466f349c1fdSJed Brown   PetscErrorCode ierr;
467f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
468f349c1fdSJed Brown   Mat            sub;
469f349c1fdSJed Brown 
470f349c1fdSJed Brown   PetscFunctionBegin;
471f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
472f349c1fdSJed Brown   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
473f349c1fdSJed Brown   if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);}
474f349c1fdSJed Brown   *B = sub;
475d8588912SDave May   PetscFunctionReturn(0);
476d8588912SDave May }
477d8588912SDave May 
478d8588912SDave May #undef __FUNCT__
479d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest"
480207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
481d8588912SDave May {
482d8588912SDave May   PetscErrorCode ierr;
483f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
484f349c1fdSJed Brown   Mat            sub;
485d8588912SDave May 
486d8588912SDave May   PetscFunctionBegin;
487f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
488ce94432eSBarry Smith   if (*B != sub) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
489f349c1fdSJed Brown   if (sub) {
490ce94432eSBarry Smith     if (((PetscObject)sub)->refct <= 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
4916bf464f9SBarry Smith     ierr = MatDestroy(B);CHKERRQ(ierr);
492d8588912SDave May   }
493d8588912SDave May   PetscFunctionReturn(0);
494d8588912SDave May }
495d8588912SDave May 
496d8588912SDave May #undef __FUNCT__
4977874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest"
4987874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
4997874fa86SDave May {
5007874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
5017874fa86SDave May   PetscInt       i;
5027874fa86SDave May   PetscErrorCode ierr;
5037874fa86SDave May 
5047874fa86SDave May   PetscFunctionBegin;
5057874fa86SDave May   for (i=0; i<bA->nr; i++) {
506429bac76SJed Brown     Vec bv;
507429bac76SJed Brown     ierr = VecGetSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
5087874fa86SDave May     if (bA->m[i][i]) {
509429bac76SJed Brown       ierr = MatGetDiagonal(bA->m[i][i],bv);CHKERRQ(ierr);
5107874fa86SDave May     } else {
5115159a857SMatthew G. Knepley       ierr = VecSet(bv,0.0);CHKERRQ(ierr);
5127874fa86SDave May     }
513429bac76SJed Brown     ierr = VecRestoreSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
5147874fa86SDave May   }
5157874fa86SDave May   PetscFunctionReturn(0);
5167874fa86SDave May }
5177874fa86SDave May 
5187874fa86SDave May #undef __FUNCT__
5197874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest"
5207874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
5217874fa86SDave May {
5227874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
523429bac76SJed Brown   Vec            bl,*br;
5247874fa86SDave May   PetscInt       i,j;
5257874fa86SDave May   PetscErrorCode ierr;
5267874fa86SDave May 
5277874fa86SDave May   PetscFunctionBegin;
5283f800ebeSJed Brown   ierr = PetscCalloc1(bA->nc,&br);CHKERRQ(ierr);
5292e6472ebSElliott Sales de Andrade   if (r) {
530429bac76SJed Brown     for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
5312e6472ebSElliott Sales de Andrade   }
5322e6472ebSElliott Sales de Andrade   bl = NULL;
5337874fa86SDave May   for (i=0; i<bA->nr; i++) {
5342e6472ebSElliott Sales de Andrade     if (l) {
535429bac76SJed Brown       ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
5362e6472ebSElliott Sales de Andrade     }
5377874fa86SDave May     for (j=0; j<bA->nc; j++) {
5387874fa86SDave May       if (bA->m[i][j]) {
539429bac76SJed Brown         ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr);
5407874fa86SDave May       }
5417874fa86SDave May     }
5422e6472ebSElliott Sales de Andrade     if (l) {
543a061e289SJed Brown       ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
5447874fa86SDave May     }
5452e6472ebSElliott Sales de Andrade   }
5462e6472ebSElliott Sales de Andrade   if (r) {
547429bac76SJed Brown     for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
5482e6472ebSElliott Sales de Andrade   }
549429bac76SJed Brown   ierr = PetscFree(br);CHKERRQ(ierr);
5507874fa86SDave May   PetscFunctionReturn(0);
5517874fa86SDave May }
5527874fa86SDave May 
5537874fa86SDave May #undef __FUNCT__
554a061e289SJed Brown #define __FUNCT__ "MatScale_Nest"
555a061e289SJed Brown static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a)
556a061e289SJed Brown {
557a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
558a061e289SJed Brown   PetscInt       i,j;
559a061e289SJed Brown   PetscErrorCode ierr;
560a061e289SJed Brown 
561a061e289SJed Brown   PetscFunctionBegin;
562a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
563a061e289SJed Brown     for (j=0; j<bA->nc; j++) {
564a061e289SJed Brown       if (bA->m[i][j]) {
565a061e289SJed Brown         ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr);
566a061e289SJed Brown       }
567a061e289SJed Brown     }
568a061e289SJed Brown   }
569a061e289SJed Brown   PetscFunctionReturn(0);
570a061e289SJed Brown }
571a061e289SJed Brown 
572a061e289SJed Brown #undef __FUNCT__
573a061e289SJed Brown #define __FUNCT__ "MatShift_Nest"
574a061e289SJed Brown static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a)
575a061e289SJed Brown {
576a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
577a061e289SJed Brown   PetscInt       i;
578a061e289SJed Brown   PetscErrorCode ierr;
579a061e289SJed Brown 
580a061e289SJed Brown   PetscFunctionBegin;
581a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
582ce94432eSBarry Smith     if (!bA->m[i][i]) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for shifting an empty diagonal block, insert a matrix in block (%D,%D)",i,i);
583a061e289SJed Brown     ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr);
584a061e289SJed Brown   }
585a061e289SJed Brown   PetscFunctionReturn(0);
586a061e289SJed Brown }
587a061e289SJed Brown 
588a061e289SJed Brown #undef __FUNCT__
58913135bc6SAlex Fikl #define __FUNCT__ "MatDiagonalSet_Nest"
59013135bc6SAlex Fikl static PetscErrorCode MatDiagonalSet_Nest(Mat A,Vec D,InsertMode is)
59113135bc6SAlex Fikl {
59213135bc6SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data;
59313135bc6SAlex Fikl   PetscInt       i;
59413135bc6SAlex Fikl   PetscErrorCode ierr;
59513135bc6SAlex Fikl 
59613135bc6SAlex Fikl   PetscFunctionBegin;
59713135bc6SAlex Fikl   for (i=0; i<bA->nr; i++) {
59813135bc6SAlex Fikl     Vec bv;
59913135bc6SAlex Fikl     ierr = VecGetSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
60013135bc6SAlex Fikl     if (bA->m[i][i]) {
60113135bc6SAlex Fikl       ierr = MatDiagonalSet(bA->m[i][i],bv,is);CHKERRQ(ierr);
60213135bc6SAlex Fikl     }
60313135bc6SAlex Fikl     ierr = VecRestoreSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
60413135bc6SAlex Fikl   }
60513135bc6SAlex Fikl   PetscFunctionReturn(0);
60613135bc6SAlex Fikl }
60713135bc6SAlex Fikl 
60813135bc6SAlex Fikl #undef __FUNCT__
609*f8170845SAlex Fikl #define __FUNCT__ "MatSetRandom_Nest"
610*f8170845SAlex Fikl static PetscErrorCode MatSetRandom_Nest(Mat A,PetscRandom rctx)
611*f8170845SAlex Fikl {
612*f8170845SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data;
613*f8170845SAlex Fikl   PetscInt       i,j;
614*f8170845SAlex Fikl   PetscErrorCode ierr;
615*f8170845SAlex Fikl 
616*f8170845SAlex Fikl   PetscFunctionBegin;
617*f8170845SAlex Fikl   for (i=0; i<bA->nr; i++) {
618*f8170845SAlex Fikl     for (j=0; j<bA->nc; j++) {
619*f8170845SAlex Fikl       if (bA->m[i][j]) {
620*f8170845SAlex Fikl         ierr = MatSetRandom(bA->m[i][j],rctx);CHKERRQ(ierr);
621*f8170845SAlex Fikl       }
622*f8170845SAlex Fikl     }
623*f8170845SAlex Fikl   }
624*f8170845SAlex Fikl   PetscFunctionReturn(0);
625*f8170845SAlex Fikl }
626*f8170845SAlex Fikl 
627*f8170845SAlex Fikl #undef __FUNCT__
6282a7a6963SBarry Smith #define __FUNCT__ "MatCreateVecs_Nest"
6292a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left)
630d8588912SDave May {
631d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
632d8588912SDave May   Vec            *L,*R;
633d8588912SDave May   MPI_Comm       comm;
634d8588912SDave May   PetscInt       i,j;
635d8588912SDave May   PetscErrorCode ierr;
636d8588912SDave May 
637d8588912SDave May   PetscFunctionBegin;
638ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
639d8588912SDave May   if (right) {
640d8588912SDave May     /* allocate R */
641854ce69bSBarry Smith     ierr = PetscMalloc1(bA->nc, &R);CHKERRQ(ierr);
642d8588912SDave May     /* Create the right vectors */
643d8588912SDave May     for (j=0; j<bA->nc; j++) {
644d8588912SDave May       for (i=0; i<bA->nr; i++) {
645d8588912SDave May         if (bA->m[i][j]) {
6462a7a6963SBarry Smith           ierr = MatCreateVecs(bA->m[i][j],&R[j],NULL);CHKERRQ(ierr);
647d8588912SDave May           break;
648d8588912SDave May         }
649d8588912SDave May       }
6506c4ed002SBarry Smith       if (i==bA->nr) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
651d8588912SDave May     }
652f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
653d8588912SDave May     /* hand back control to the nest vector */
654d8588912SDave May     for (j=0; j<bA->nc; j++) {
6556bf464f9SBarry Smith       ierr = VecDestroy(&R[j]);CHKERRQ(ierr);
656d8588912SDave May     }
657d8588912SDave May     ierr = PetscFree(R);CHKERRQ(ierr);
658d8588912SDave May   }
659d8588912SDave May 
660d8588912SDave May   if (left) {
661d8588912SDave May     /* allocate L */
662854ce69bSBarry Smith     ierr = PetscMalloc1(bA->nr, &L);CHKERRQ(ierr);
663d8588912SDave May     /* Create the left vectors */
664d8588912SDave May     for (i=0; i<bA->nr; i++) {
665d8588912SDave May       for (j=0; j<bA->nc; j++) {
666d8588912SDave May         if (bA->m[i][j]) {
6672a7a6963SBarry Smith           ierr = MatCreateVecs(bA->m[i][j],NULL,&L[i]);CHKERRQ(ierr);
668d8588912SDave May           break;
669d8588912SDave May         }
670d8588912SDave May       }
6716c4ed002SBarry Smith       if (j==bA->nc) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
672d8588912SDave May     }
673d8588912SDave May 
674f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
675d8588912SDave May     for (i=0; i<bA->nr; i++) {
6766bf464f9SBarry Smith       ierr = VecDestroy(&L[i]);CHKERRQ(ierr);
677d8588912SDave May     }
678d8588912SDave May 
679d8588912SDave May     ierr = PetscFree(L);CHKERRQ(ierr);
680d8588912SDave May   }
681d8588912SDave May   PetscFunctionReturn(0);
682d8588912SDave May }
683d8588912SDave May 
684d8588912SDave May #undef __FUNCT__
685d8588912SDave May #define __FUNCT__ "MatView_Nest"
686207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
687d8588912SDave May {
688d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
689d8588912SDave May   PetscBool      isascii;
690d8588912SDave May   PetscInt       i,j;
691d8588912SDave May   PetscErrorCode ierr;
692d8588912SDave May 
693d8588912SDave May   PetscFunctionBegin;
694251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
695d8588912SDave May   if (isascii) {
696d8588912SDave May 
697d8588912SDave May     PetscViewerASCIIPrintf(viewer,"Matrix object: \n");
698d8588912SDave May     PetscViewerASCIIPushTab(viewer);    /* push0 */
699d8588912SDave May     PetscViewerASCIIPrintf(viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc);
700d8588912SDave May 
701d8588912SDave May     PetscViewerASCIIPrintf(viewer,"MatNest structure: \n");
702d8588912SDave May     for (i=0; i<bA->nr; i++) {
703d8588912SDave May       for (j=0; j<bA->nc; j++) {
70419fd82e9SBarry Smith         MatType   type;
705270f95d7SJed Brown         char      name[256] = "",prefix[256] = "";
706d8588912SDave May         PetscInt  NR,NC;
707d8588912SDave May         PetscBool isNest = PETSC_FALSE;
708d8588912SDave May 
709d8588912SDave May         if (!bA->m[i][j]) {
7100298fd71SBarry Smith           PetscViewerASCIIPrintf(viewer, "(%D,%D) : NULL \n",i,j);
711d8588912SDave May           continue;
712d8588912SDave May         }
713d8588912SDave May         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
714d8588912SDave May         ierr = MatGetType(bA->m[i][j], &type);CHKERRQ(ierr);
7158caf3d72SBarry Smith         if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);}
7168caf3d72SBarry Smith         if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);}
717251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
718d8588912SDave May 
719270f95d7SJed Brown         ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr);
720d8588912SDave May 
721d8588912SDave May         if (isNest) {
722270f95d7SJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);  /* push1 */
723d8588912SDave May           ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr);
724270f95d7SJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop1 */
725d8588912SDave May         }
726d8588912SDave May       }
727d8588912SDave May     }
728d8588912SDave May     PetscViewerASCIIPopTab(viewer);    /* pop0 */
729d8588912SDave May   }
730d8588912SDave May   PetscFunctionReturn(0);
731d8588912SDave May }
732d8588912SDave May 
733d8588912SDave May #undef __FUNCT__
734d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest"
735207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A)
736d8588912SDave May {
737d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
738d8588912SDave May   PetscInt       i,j;
739d8588912SDave May   PetscErrorCode ierr;
740d8588912SDave May 
741d8588912SDave May   PetscFunctionBegin;
742d8588912SDave May   for (i=0; i<bA->nr; i++) {
743d8588912SDave May     for (j=0; j<bA->nc; j++) {
744d8588912SDave May       if (!bA->m[i][j]) continue;
745d8588912SDave May       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
746d8588912SDave May     }
747d8588912SDave May   }
748d8588912SDave May   PetscFunctionReturn(0);
749d8588912SDave May }
750d8588912SDave May 
751d8588912SDave May #undef __FUNCT__
752c222c20dSDavid Ham #define __FUNCT__ "MatCopy_Nest"
753c222c20dSDavid Ham static PetscErrorCode MatCopy_Nest(Mat A,Mat B,MatStructure str)
754c222c20dSDavid Ham {
755c222c20dSDavid Ham   Mat_Nest       *bA = (Mat_Nest*)A->data,*bB = (Mat_Nest*)B->data;
756c222c20dSDavid Ham   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
757c222c20dSDavid Ham   PetscErrorCode ierr;
758c222c20dSDavid Ham 
759c222c20dSDavid Ham   PetscFunctionBegin;
760c222c20dSDavid Ham   if (nr != bB->nr || nc != bB->nc) SETERRQ4(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Cannot copy a Mat_Nest of block size (%D,%D) to a Mat_Nest of block size (%D,%D)",bB->nr,bB->nc,nr,nc);
761c222c20dSDavid Ham   for (i=0; i<nr; i++) {
762c222c20dSDavid Ham     for (j=0; j<nc; j++) {
76346a2b97cSJed Brown       if (bA->m[i][j] && bB->m[i][j]) {
764c222c20dSDavid Ham         ierr = MatCopy(bA->m[i][j],bB->m[i][j],str);CHKERRQ(ierr);
76546a2b97cSJed Brown       } else if (bA->m[i][j] || bB->m[i][j]) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Matrix block does not exist at %D,%D",i,j);
766c222c20dSDavid Ham     }
767c222c20dSDavid Ham   }
768c222c20dSDavid Ham   PetscFunctionReturn(0);
769c222c20dSDavid Ham }
770c222c20dSDavid Ham 
771c222c20dSDavid Ham #undef __FUNCT__
772d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest"
773207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
774d8588912SDave May {
775d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
776841e96a3SJed Brown   Mat            *b;
777841e96a3SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
778d8588912SDave May   PetscErrorCode ierr;
779d8588912SDave May 
780d8588912SDave May   PetscFunctionBegin;
781785e854fSJed Brown   ierr = PetscMalloc1(nr*nc,&b);CHKERRQ(ierr);
782841e96a3SJed Brown   for (i=0; i<nr; i++) {
783841e96a3SJed Brown     for (j=0; j<nc; j++) {
784841e96a3SJed Brown       if (bA->m[i][j]) {
785841e96a3SJed Brown         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
786841e96a3SJed Brown       } else {
7870298fd71SBarry Smith         b[i*nc+j] = NULL;
788d8588912SDave May       }
789d8588912SDave May     }
790d8588912SDave May   }
791ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
792841e96a3SJed Brown   /* Give the new MatNest exclusive ownership */
793841e96a3SJed Brown   for (i=0; i<nr*nc; i++) {
7946bf464f9SBarry Smith     ierr = MatDestroy(&b[i]);CHKERRQ(ierr);
795d8588912SDave May   }
796d8588912SDave May   ierr = PetscFree(b);CHKERRQ(ierr);
797d8588912SDave May 
798841e96a3SJed Brown   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
799841e96a3SJed Brown   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
800d8588912SDave May   PetscFunctionReturn(0);
801d8588912SDave May }
802d8588912SDave May 
803d8588912SDave May /* nest api */
804d8588912SDave May #undef __FUNCT__
805d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest"
806d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
807d8588912SDave May {
808d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
8095fd66863SKarl Rupp 
810d8588912SDave May   PetscFunctionBegin;
811ce94432eSBarry Smith   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
812ce94432eSBarry Smith   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
813d8588912SDave May   *mat = bA->m[idxm][jdxm];
814d8588912SDave May   PetscFunctionReturn(0);
815d8588912SDave May }
816d8588912SDave May 
817d8588912SDave May #undef __FUNCT__
818d8588912SDave May #define __FUNCT__ "MatNestGetSubMat"
8199ba0d327SJed Brown /*@
820d8588912SDave May  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
821d8588912SDave May 
822d8588912SDave May  Not collective
823d8588912SDave May 
824d8588912SDave May  Input Parameters:
825629881c0SJed Brown +   A  - nest matrix
826d8588912SDave May .   idxm - index of the matrix within the nest matrix
827629881c0SJed Brown -   jdxm - index of the matrix within the nest matrix
828d8588912SDave May 
829d8588912SDave May  Output Parameter:
830d8588912SDave May .   sub - matrix at index idxm,jdxm within the nest matrix
831d8588912SDave May 
832d8588912SDave May  Level: developer
833d8588912SDave May 
834d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMats()
835d8588912SDave May @*/
8367087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
837d8588912SDave May {
838699a902aSJed Brown   PetscErrorCode ierr;
839d8588912SDave May 
840d8588912SDave May   PetscFunctionBegin;
841699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
842d8588912SDave May   PetscFunctionReturn(0);
843d8588912SDave May }
844d8588912SDave May 
845d8588912SDave May #undef __FUNCT__
8460782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat_Nest"
8470782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat)
8480782ca92SJed Brown {
8490782ca92SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
8500782ca92SJed Brown   PetscInt       m,n,M,N,mi,ni,Mi,Ni;
8510782ca92SJed Brown   PetscErrorCode ierr;
8520782ca92SJed Brown 
8530782ca92SJed Brown   PetscFunctionBegin;
854ce94432eSBarry Smith   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
855ce94432eSBarry Smith   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
8560782ca92SJed Brown   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
8570782ca92SJed Brown   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
8580782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr);
8590782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr);
8600782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr);
8610782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr);
862ce94432eSBarry Smith   if (M != Mi || N != Ni) SETERRQ4(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"Submatrix dimension (%D,%D) incompatible with nest block (%D,%D)",M,N,Mi,Ni);
863ce94432eSBarry Smith   if (m != mi || n != ni) SETERRQ4(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"Submatrix local dimension (%D,%D) incompatible with nest block (%D,%D)",m,n,mi,ni);
86426fbe8dcSKarl Rupp 
8650782ca92SJed Brown   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
8660782ca92SJed Brown   ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr);
8670782ca92SJed Brown   bA->m[idxm][jdxm] = mat;
8680782ca92SJed Brown   PetscFunctionReturn(0);
8690782ca92SJed Brown }
8700782ca92SJed Brown 
8710782ca92SJed Brown #undef __FUNCT__
8720782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat"
8739ba0d327SJed Brown /*@
8740782ca92SJed Brown  MatNestSetSubMat - Set a single submatrix in the nest matrix.
8750782ca92SJed Brown 
8760782ca92SJed Brown  Logically collective on the submatrix communicator
8770782ca92SJed Brown 
8780782ca92SJed Brown  Input Parameters:
8790782ca92SJed Brown +   A  - nest matrix
8800782ca92SJed Brown .   idxm - index of the matrix within the nest matrix
8810782ca92SJed Brown .   jdxm - index of the matrix within the nest matrix
8820782ca92SJed Brown -   sub - matrix at index idxm,jdxm within the nest matrix
8830782ca92SJed Brown 
8840782ca92SJed Brown  Notes:
8850782ca92SJed Brown  The new submatrix must have the same size and communicator as that block of the nest.
8860782ca92SJed Brown 
8870782ca92SJed Brown  This increments the reference count of the submatrix.
8880782ca92SJed Brown 
8890782ca92SJed Brown  Level: developer
8900782ca92SJed Brown 
8910782ca92SJed Brown .seealso: MatNestSetSubMats(), MatNestGetSubMat()
8920782ca92SJed Brown @*/
8930782ca92SJed Brown PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
8940782ca92SJed Brown {
8950782ca92SJed Brown   PetscErrorCode ierr;
8960782ca92SJed Brown 
8970782ca92SJed Brown   PetscFunctionBegin;
8980782ca92SJed Brown   ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr);
8990782ca92SJed Brown   PetscFunctionReturn(0);
9000782ca92SJed Brown }
9010782ca92SJed Brown 
9020782ca92SJed Brown #undef __FUNCT__
903d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest"
904d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
905d8588912SDave May {
906d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
9075fd66863SKarl Rupp 
908d8588912SDave May   PetscFunctionBegin;
90926fbe8dcSKarl Rupp   if (M)   *M   = bA->nr;
91026fbe8dcSKarl Rupp   if (N)   *N   = bA->nc;
91126fbe8dcSKarl Rupp   if (mat) *mat = bA->m;
912d8588912SDave May   PetscFunctionReturn(0);
913d8588912SDave May }
914d8588912SDave May 
915d8588912SDave May #undef __FUNCT__
916d8588912SDave May #define __FUNCT__ "MatNestGetSubMats"
917d8588912SDave May /*@C
918d8588912SDave May  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
919d8588912SDave May 
920d8588912SDave May  Not collective
921d8588912SDave May 
922d8588912SDave May  Input Parameters:
923629881c0SJed Brown .   A  - nest matrix
924d8588912SDave May 
925d8588912SDave May  Output Parameter:
926629881c0SJed Brown +   M - number of rows in the nest matrix
927d8588912SDave May .   N - number of cols in the nest matrix
928629881c0SJed Brown -   mat - 2d array of matrices
929d8588912SDave May 
930d8588912SDave May  Notes:
931d8588912SDave May 
932d8588912SDave May  The user should not free the array mat.
933d8588912SDave May 
934351962e3SVincent Le Chenadec  In Fortran, this routine has a calling sequence
935351962e3SVincent Le Chenadec $   call MatNestGetSubMats(A, M, N, mat, ierr)
936351962e3SVincent Le Chenadec  where the space allocated for the optional argument mat is assumed large enough (if provided).
937351962e3SVincent Le Chenadec 
938d8588912SDave May  Level: developer
939d8588912SDave May 
940d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMat()
941d8588912SDave May @*/
9427087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
943d8588912SDave May {
944699a902aSJed Brown   PetscErrorCode ierr;
945d8588912SDave May 
946d8588912SDave May   PetscFunctionBegin;
947699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
948d8588912SDave May   PetscFunctionReturn(0);
949d8588912SDave May }
950d8588912SDave May 
951d8588912SDave May #undef __FUNCT__
952d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest"
9537087cfbeSBarry Smith PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
954d8588912SDave May {
955d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
956d8588912SDave May 
957d8588912SDave May   PetscFunctionBegin;
95826fbe8dcSKarl Rupp   if (M) *M = bA->nr;
95926fbe8dcSKarl Rupp   if (N) *N = bA->nc;
960d8588912SDave May   PetscFunctionReturn(0);
961d8588912SDave May }
962d8588912SDave May 
963d8588912SDave May #undef __FUNCT__
964d8588912SDave May #define __FUNCT__ "MatNestGetSize"
9659ba0d327SJed Brown /*@
966d8588912SDave May  MatNestGetSize - Returns the size of the nest matrix.
967d8588912SDave May 
968d8588912SDave May  Not collective
969d8588912SDave May 
970d8588912SDave May  Input Parameters:
971d8588912SDave May .   A  - nest matrix
972d8588912SDave May 
973d8588912SDave May  Output Parameter:
974629881c0SJed Brown +   M - number of rows in the nested mat
975629881c0SJed Brown -   N - number of cols in the nested mat
976d8588912SDave May 
977d8588912SDave May  Notes:
978d8588912SDave May 
979d8588912SDave May  Level: developer
980d8588912SDave May 
981d8588912SDave May .seealso: MatNestGetSubMat(), MatNestGetSubMats()
982d8588912SDave May @*/
9837087cfbeSBarry Smith PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
984d8588912SDave May {
985699a902aSJed Brown   PetscErrorCode ierr;
986d8588912SDave May 
987d8588912SDave May   PetscFunctionBegin;
988699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
989d8588912SDave May   PetscFunctionReturn(0);
990d8588912SDave May }
991d8588912SDave May 
992900e7ff2SJed Brown #undef __FUNCT__
993900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs_Nest"
994f7a08781SBarry Smith static PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[])
995900e7ff2SJed Brown {
996900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest*)A->data;
997900e7ff2SJed Brown   PetscInt i;
998900e7ff2SJed Brown 
999900e7ff2SJed Brown   PetscFunctionBegin;
1000900e7ff2SJed Brown   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i];
1001900e7ff2SJed Brown   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i];
1002900e7ff2SJed Brown   PetscFunctionReturn(0);
1003900e7ff2SJed Brown }
1004900e7ff2SJed Brown 
1005900e7ff2SJed Brown #undef __FUNCT__
1006900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs"
10073a4d7b9aSSatish Balay /*@C
1008900e7ff2SJed Brown  MatNestGetISs - Returns the index sets partitioning the row and column spaces
1009900e7ff2SJed Brown 
1010900e7ff2SJed Brown  Not collective
1011900e7ff2SJed Brown 
1012900e7ff2SJed Brown  Input Parameters:
1013900e7ff2SJed Brown .   A  - nest matrix
1014900e7ff2SJed Brown 
1015900e7ff2SJed Brown  Output Parameter:
1016900e7ff2SJed Brown +   rows - array of row index sets
1017900e7ff2SJed Brown -   cols - array of column index sets
1018900e7ff2SJed Brown 
1019900e7ff2SJed Brown  Level: advanced
1020900e7ff2SJed Brown 
1021900e7ff2SJed Brown  Notes:
1022900e7ff2SJed Brown  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
1023900e7ff2SJed Brown 
1024900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs()
1025900e7ff2SJed Brown @*/
1026900e7ff2SJed Brown PetscErrorCode  MatNestGetISs(Mat A,IS rows[],IS cols[])
1027900e7ff2SJed Brown {
1028900e7ff2SJed Brown   PetscErrorCode ierr;
1029900e7ff2SJed Brown 
1030900e7ff2SJed Brown   PetscFunctionBegin;
1031900e7ff2SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1032900e7ff2SJed Brown   ierr = PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
1033900e7ff2SJed Brown   PetscFunctionReturn(0);
1034900e7ff2SJed Brown }
1035900e7ff2SJed Brown 
1036900e7ff2SJed Brown #undef __FUNCT__
1037900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs_Nest"
1038f7a08781SBarry Smith static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[])
1039900e7ff2SJed Brown {
1040900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest*)A->data;
1041900e7ff2SJed Brown   PetscInt i;
1042900e7ff2SJed Brown 
1043900e7ff2SJed Brown   PetscFunctionBegin;
1044900e7ff2SJed Brown   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i];
1045900e7ff2SJed Brown   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i];
1046900e7ff2SJed Brown   PetscFunctionReturn(0);
1047900e7ff2SJed Brown }
1048900e7ff2SJed Brown 
1049900e7ff2SJed Brown #undef __FUNCT__
1050900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs"
1051900e7ff2SJed Brown /*@C
1052900e7ff2SJed Brown  MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces
1053900e7ff2SJed Brown 
1054900e7ff2SJed Brown  Not collective
1055900e7ff2SJed Brown 
1056900e7ff2SJed Brown  Input Parameters:
1057900e7ff2SJed Brown .   A  - nest matrix
1058900e7ff2SJed Brown 
1059900e7ff2SJed Brown  Output Parameter:
10600298fd71SBarry Smith +   rows - array of row index sets (or NULL to ignore)
10610298fd71SBarry Smith -   cols - array of column index sets (or NULL to ignore)
1062900e7ff2SJed Brown 
1063900e7ff2SJed Brown  Level: advanced
1064900e7ff2SJed Brown 
1065900e7ff2SJed Brown  Notes:
1066900e7ff2SJed Brown  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
1067900e7ff2SJed Brown 
1068900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs()
1069900e7ff2SJed Brown @*/
1070900e7ff2SJed Brown PetscErrorCode  MatNestGetLocalISs(Mat A,IS rows[],IS cols[])
1071900e7ff2SJed Brown {
1072900e7ff2SJed Brown   PetscErrorCode ierr;
1073900e7ff2SJed Brown 
1074900e7ff2SJed Brown   PetscFunctionBegin;
1075900e7ff2SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1076900e7ff2SJed Brown   ierr = PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
1077900e7ff2SJed Brown   PetscFunctionReturn(0);
1078900e7ff2SJed Brown }
1079900e7ff2SJed Brown 
1080207556f9SJed Brown #undef __FUNCT__
1081207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest"
108219fd82e9SBarry Smith PetscErrorCode  MatNestSetVecType_Nest(Mat A,VecType vtype)
1083207556f9SJed Brown {
1084207556f9SJed Brown   PetscErrorCode ierr;
1085207556f9SJed Brown   PetscBool      flg;
1086207556f9SJed Brown 
1087207556f9SJed Brown   PetscFunctionBegin;
1088207556f9SJed Brown   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
1089207556f9SJed Brown   /* In reality, this only distinguishes VECNEST and "other" */
10902a7a6963SBarry Smith   if (flg) A->ops->getvecs = MatCreateVecs_Nest;
109112b53f24SSatish Balay   else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0;
1092207556f9SJed Brown   PetscFunctionReturn(0);
1093207556f9SJed Brown }
1094207556f9SJed Brown 
1095207556f9SJed Brown #undef __FUNCT__
1096207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType"
1097207556f9SJed Brown /*@C
10982a7a6963SBarry Smith  MatNestSetVecType - Sets the type of Vec returned by MatCreateVecs()
1099207556f9SJed Brown 
1100207556f9SJed Brown  Not collective
1101207556f9SJed Brown 
1102207556f9SJed Brown  Input Parameters:
1103207556f9SJed Brown +  A  - nest matrix
1104207556f9SJed Brown -  vtype - type to use for creating vectors
1105207556f9SJed Brown 
1106207556f9SJed Brown  Notes:
1107207556f9SJed Brown 
1108207556f9SJed Brown  Level: developer
1109207556f9SJed Brown 
11102a7a6963SBarry Smith .seealso: MatCreateVecs()
1111207556f9SJed Brown @*/
111219fd82e9SBarry Smith PetscErrorCode  MatNestSetVecType(Mat A,VecType vtype)
1113207556f9SJed Brown {
1114207556f9SJed Brown   PetscErrorCode ierr;
1115207556f9SJed Brown 
1116207556f9SJed Brown   PetscFunctionBegin;
111719fd82e9SBarry Smith   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));CHKERRQ(ierr);
1118207556f9SJed Brown   PetscFunctionReturn(0);
1119207556f9SJed Brown }
1120207556f9SJed Brown 
1121d8588912SDave May #undef __FUNCT__
1122c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats_Nest"
1123c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1124d8588912SDave May {
1125c8883902SJed Brown   Mat_Nest       *s = (Mat_Nest*)A->data;
1126c8883902SJed Brown   PetscInt       i,j,m,n,M,N;
1127d8588912SDave May   PetscErrorCode ierr;
1128d8588912SDave May 
1129d8588912SDave May   PetscFunctionBegin;
1130c8883902SJed Brown   s->nr = nr;
1131c8883902SJed Brown   s->nc = nc;
1132d8588912SDave May 
1133c8883902SJed Brown   /* Create space for submatrices */
1134854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&s->m);CHKERRQ(ierr);
1135c8883902SJed Brown   for (i=0; i<nr; i++) {
1136854ce69bSBarry Smith     ierr = PetscMalloc1(nc,&s->m[i]);CHKERRQ(ierr);
1137d8588912SDave May   }
1138c8883902SJed Brown   for (i=0; i<nr; i++) {
1139c8883902SJed Brown     for (j=0; j<nc; j++) {
1140c8883902SJed Brown       s->m[i][j] = a[i*nc+j];
1141c8883902SJed Brown       if (a[i*nc+j]) {
1142c8883902SJed Brown         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
1143d8588912SDave May       }
1144d8588912SDave May     }
1145d8588912SDave May   }
1146d8588912SDave May 
11478188e55aSJed Brown   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
1148d8588912SDave May 
1149854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&s->row_len);CHKERRQ(ierr);
1150854ce69bSBarry Smith   ierr = PetscMalloc1(nc,&s->col_len);CHKERRQ(ierr);
1151c8883902SJed Brown   for (i=0; i<nr; i++) s->row_len[i]=-1;
1152c8883902SJed Brown   for (j=0; j<nc; j++) s->col_len[j]=-1;
1153d8588912SDave May 
11548188e55aSJed Brown   ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr);
1155d8588912SDave May 
1156c8883902SJed Brown   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
1157c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
1158c8883902SJed Brown   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
1159c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
1160c8883902SJed Brown 
1161c8883902SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1162c8883902SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1163c8883902SJed Brown 
11641795a4d1SJed Brown   ierr = PetscCalloc2(nr,&s->left,nc,&s->right);CHKERRQ(ierr);
1165d8588912SDave May   PetscFunctionReturn(0);
1166d8588912SDave May }
1167d8588912SDave May 
1168c8883902SJed Brown #undef __FUNCT__
1169c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats"
1170c8883902SJed Brown /*@
1171c8883902SJed Brown    MatNestSetSubMats - Sets the nested submatrices
1172c8883902SJed Brown 
1173c8883902SJed Brown    Collective on Mat
1174c8883902SJed Brown 
1175c8883902SJed Brown    Input Parameter:
1176c8883902SJed Brown +  N - nested matrix
1177c8883902SJed Brown .  nr - number of nested row blocks
11780298fd71SBarry Smith .  is_row - index sets for each nested row block, or NULL to make contiguous
1179c8883902SJed Brown .  nc - number of nested column blocks
11800298fd71SBarry Smith .  is_col - index sets for each nested column block, or NULL to make contiguous
11810298fd71SBarry Smith -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1182c8883902SJed Brown 
1183c8883902SJed Brown    Level: advanced
1184c8883902SJed Brown 
1185c8883902SJed Brown .seealso: MatCreateNest(), MATNEST
1186c8883902SJed Brown @*/
1187c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1188c8883902SJed Brown {
1189c8883902SJed Brown   PetscErrorCode ierr;
1190c8883902SJed Brown   PetscInt       i;
1191c8883902SJed Brown 
1192c8883902SJed Brown   PetscFunctionBegin;
1193c8883902SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1194ce94432eSBarry Smith   if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
1195c8883902SJed Brown   if (nr && is_row) {
1196c8883902SJed Brown     PetscValidPointer(is_row,3);
1197c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
1198c8883902SJed Brown   }
1199ce94432eSBarry Smith   if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
12001664e352SJed Brown   if (nc && is_col) {
1201c8883902SJed Brown     PetscValidPointer(is_col,5);
12029b30a8f6SBarry Smith     for (i=0; i<nc; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
1203c8883902SJed Brown   }
1204c8883902SJed Brown   if (nr*nc) PetscValidPointer(a,6);
1205c8883902SJed 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);
1206c8883902SJed Brown   PetscFunctionReturn(0);
1207c8883902SJed Brown }
1208d8588912SDave May 
120977019fcaSJed Brown #undef __FUNCT__
121077019fcaSJed Brown #define __FUNCT__ "MatNestCreateAggregateL2G_Private"
121145b6f7e9SBarry Smith static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog)
121277019fcaSJed Brown {
121377019fcaSJed Brown   PetscErrorCode ierr;
121477019fcaSJed Brown   PetscBool      flg;
121577019fcaSJed Brown   PetscInt       i,j,m,mi,*ix;
121677019fcaSJed Brown 
121777019fcaSJed Brown   PetscFunctionBegin;
121877019fcaSJed Brown   for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) {
121977019fcaSJed Brown     if (islocal[i]) {
122077019fcaSJed Brown       ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
122177019fcaSJed Brown       flg  = PETSC_TRUE;      /* We found a non-trivial entry */
122277019fcaSJed Brown     } else {
122377019fcaSJed Brown       ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
122477019fcaSJed Brown     }
122577019fcaSJed Brown     m += mi;
122677019fcaSJed Brown   }
122777019fcaSJed Brown   if (flg) {
1228785e854fSJed Brown     ierr = PetscMalloc1(m,&ix);CHKERRQ(ierr);
122977019fcaSJed Brown     for (i=0,n=0; i<n; i++) {
12300298fd71SBarry Smith       ISLocalToGlobalMapping smap = NULL;
123177019fcaSJed Brown       VecScatter             scat;
123277019fcaSJed Brown       IS                     isreq;
123377019fcaSJed Brown       Vec                    lvec,gvec;
12343361c9a7SJed Brown       union {char padding[sizeof(PetscScalar)]; PetscInt integer;} *x;
123577019fcaSJed Brown       Mat sub;
123677019fcaSJed Brown 
1237ce94432eSBarry Smith       if (sizeof(*x) != sizeof(PetscScalar)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support when scalars smaller than integers");
123877019fcaSJed Brown       if (colflg) {
123977019fcaSJed Brown         ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
124077019fcaSJed Brown       } else {
124177019fcaSJed Brown         ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
124277019fcaSJed Brown       }
12430298fd71SBarry Smith       if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&smap,NULL);CHKERRQ(ierr);}
124477019fcaSJed Brown       if (islocal[i]) {
124577019fcaSJed Brown         ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
124677019fcaSJed Brown       } else {
124777019fcaSJed Brown         ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
124877019fcaSJed Brown       }
124977019fcaSJed Brown       for (j=0; j<mi; j++) ix[m+j] = j;
125077019fcaSJed Brown       if (smap) {ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr);}
125177019fcaSJed Brown       /*
125277019fcaSJed Brown         Now we need to extract the monolithic global indices that correspond to the given split global indices.
125377019fcaSJed Brown         In many/most cases, we only want MatGetLocalSubMatrix() to work, in which case we only need to know the size of the local spaces.
125477019fcaSJed Brown         The approach here is ugly because it uses VecScatter to move indices.
125577019fcaSJed Brown        */
125677019fcaSJed Brown       ierr = VecCreateSeq(PETSC_COMM_SELF,mi,&lvec);CHKERRQ(ierr);
125777019fcaSJed Brown       ierr = VecCreateMPI(((PetscObject)isglobal[i])->comm,mi,PETSC_DECIDE,&gvec);CHKERRQ(ierr);
125877019fcaSJed Brown       ierr = ISCreateGeneral(((PetscObject)isglobal[i])->comm,mi,ix+m,PETSC_COPY_VALUES,&isreq);CHKERRQ(ierr);
12590298fd71SBarry Smith       ierr = VecScatterCreate(gvec,isreq,lvec,NULL,&scat);CHKERRQ(ierr);
126077019fcaSJed Brown       ierr = VecGetArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
126177019fcaSJed Brown       for (j=0; j<mi; j++) x[j].integer = ix[m+j];
126277019fcaSJed Brown       ierr = VecRestoreArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
126377019fcaSJed Brown       ierr = VecScatterBegin(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
126477019fcaSJed Brown       ierr = VecScatterEnd(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
126577019fcaSJed Brown       ierr = VecGetArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
126677019fcaSJed Brown       for (j=0; j<mi; j++) ix[m+j] = x[j].integer;
126777019fcaSJed Brown       ierr = VecRestoreArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
126877019fcaSJed Brown       ierr = VecDestroy(&lvec);CHKERRQ(ierr);
126977019fcaSJed Brown       ierr = VecDestroy(&gvec);CHKERRQ(ierr);
127077019fcaSJed Brown       ierr = ISDestroy(&isreq);CHKERRQ(ierr);
127177019fcaSJed Brown       ierr = VecScatterDestroy(&scat);CHKERRQ(ierr);
127277019fcaSJed Brown       m   += mi;
127377019fcaSJed Brown     }
1274f0413b6fSBarry Smith     ierr   = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),1,m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
127577019fcaSJed Brown   } else {
12760298fd71SBarry Smith     *ltog  = NULL;
127777019fcaSJed Brown   }
127877019fcaSJed Brown   PetscFunctionReturn(0);
127977019fcaSJed Brown }
128077019fcaSJed Brown 
128177019fcaSJed Brown 
1282d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1283d8588912SDave May /*
1284d8588912SDave May   nprocessors = NP
1285d8588912SDave May   Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1))
1286d8588912SDave May        proc 0: => (g_0,h_0,)
1287d8588912SDave May        proc 1: => (g_1,h_1,)
1288d8588912SDave May        ...
1289d8588912SDave May        proc nprocs-1: => (g_NP-1,h_NP-1,)
1290d8588912SDave May 
1291d8588912SDave May             proc 0:                      proc 1:                    proc nprocs-1:
1292d8588912SDave May     is[0] = (0,1,2,...,nlocal(g_0)-1)  (0,1,...,nlocal(g_1)-1)  (0,1,...,nlocal(g_NP-1))
1293d8588912SDave May 
1294d8588912SDave May             proc 0:
1295d8588912SDave May     is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1)
1296d8588912SDave May             proc 1:
1297d8588912SDave May     is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1)
1298d8588912SDave May 
1299d8588912SDave May             proc NP-1:
1300d8588912SDave May     is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1)
1301d8588912SDave May */
1302d8588912SDave May #undef __FUNCT__
1303d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private"
1304841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1305d8588912SDave May {
1306e2d7f03fSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
13078188e55aSJed Brown   PetscInt       i,j,offset,n,nsum,bs;
1308d8588912SDave May   PetscErrorCode ierr;
13090298fd71SBarry Smith   Mat            sub = NULL;
1310d8588912SDave May 
1311d8588912SDave May   PetscFunctionBegin;
1312854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&vs->isglobal.row);CHKERRQ(ierr);
1313854ce69bSBarry Smith   ierr = PetscMalloc1(nc,&vs->isglobal.col);CHKERRQ(ierr);
1314d8588912SDave May   if (is_row) { /* valid IS is passed in */
1315d8588912SDave May     /* refs on is[] are incremeneted */
1316e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1317d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
131826fbe8dcSKarl Rupp 
1319e2d7f03fSJed Brown       vs->isglobal.row[i] = is_row[i];
1320d8588912SDave May     }
13212ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
13228188e55aSJed Brown     nsum = 0;
13238188e55aSJed Brown     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
13248188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1325ce94432eSBarry Smith       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
13260298fd71SBarry Smith       ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
1327ce94432eSBarry Smith       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
13288188e55aSJed Brown       nsum += n;
13298188e55aSJed Brown     }
1330ce94432eSBarry Smith     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
133130bc264bSJed Brown     offset -= nsum;
1332e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1333f349c1fdSJed Brown       ierr    = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
13340298fd71SBarry Smith       ierr    = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
13352ae74bdbSJed Brown       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1336ce94432eSBarry Smith       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
1337e2d7f03fSJed Brown       ierr    = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
13382ae74bdbSJed Brown       offset += n;
1339d8588912SDave May     }
1340d8588912SDave May   }
1341d8588912SDave May 
1342d8588912SDave May   if (is_col) { /* valid IS is passed in */
1343d8588912SDave May     /* refs on is[] are incremeneted */
1344e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1345d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
134626fbe8dcSKarl Rupp 
1347e2d7f03fSJed Brown       vs->isglobal.col[j] = is_col[j];
1348d8588912SDave May     }
13492ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
13502ae74bdbSJed Brown     offset = A->cmap->rstart;
13518188e55aSJed Brown     nsum   = 0;
13528188e55aSJed Brown     for (j=0; j<vs->nc; j++) {
13538188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
1354ce94432eSBarry Smith       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
13550298fd71SBarry Smith       ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
1356ce94432eSBarry Smith       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
13578188e55aSJed Brown       nsum += n;
13588188e55aSJed Brown     }
1359ce94432eSBarry Smith     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
136030bc264bSJed Brown     offset -= nsum;
1361e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1362f349c1fdSJed Brown       ierr    = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
13630298fd71SBarry Smith       ierr    = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
13642ae74bdbSJed Brown       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1365ce94432eSBarry Smith       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
1366e2d7f03fSJed Brown       ierr    = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
13672ae74bdbSJed Brown       offset += n;
1368d8588912SDave May     }
1369d8588912SDave May   }
1370e2d7f03fSJed Brown 
1371e2d7f03fSJed Brown   /* Set up the local ISs */
1372785e854fSJed Brown   ierr = PetscMalloc1(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
1373785e854fSJed Brown   ierr = PetscMalloc1(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
1374e2d7f03fSJed Brown   for (i=0,offset=0; i<vs->nr; i++) {
1375e2d7f03fSJed Brown     IS                     isloc;
13760298fd71SBarry Smith     ISLocalToGlobalMapping rmap = NULL;
1377e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1378e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
13790298fd71SBarry Smith     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,NULL);CHKERRQ(ierr);}
1380207556f9SJed Brown     if (rmap) {
1381e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1382e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
1383e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1384e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1385207556f9SJed Brown     } else {
1386207556f9SJed Brown       nlocal = 0;
13870298fd71SBarry Smith       isloc  = NULL;
1388207556f9SJed Brown     }
1389e2d7f03fSJed Brown     vs->islocal.row[i] = isloc;
1390e2d7f03fSJed Brown     offset            += nlocal;
1391e2d7f03fSJed Brown   }
13928188e55aSJed Brown   for (i=0,offset=0; i<vs->nc; i++) {
1393e2d7f03fSJed Brown     IS                     isloc;
13940298fd71SBarry Smith     ISLocalToGlobalMapping cmap = NULL;
1395e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1396e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
13970298fd71SBarry Smith     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,NULL,&cmap);CHKERRQ(ierr);}
1398207556f9SJed Brown     if (cmap) {
1399e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1400e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr);
1401e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1402e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1403207556f9SJed Brown     } else {
1404207556f9SJed Brown       nlocal = 0;
14050298fd71SBarry Smith       isloc  = NULL;
1406207556f9SJed Brown     }
1407e2d7f03fSJed Brown     vs->islocal.col[i] = isloc;
1408e2d7f03fSJed Brown     offset            += nlocal;
1409e2d7f03fSJed Brown   }
14100189643fSJed Brown 
141177019fcaSJed Brown   /* Set up the aggregate ISLocalToGlobalMapping */
141277019fcaSJed Brown   {
141345b6f7e9SBarry Smith     ISLocalToGlobalMapping rmap,cmap;
141445b6f7e9SBarry Smith     ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap);CHKERRQ(ierr);
141545b6f7e9SBarry Smith     ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap);CHKERRQ(ierr);
141677019fcaSJed Brown     if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);}
141777019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
141877019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
141977019fcaSJed Brown   }
142077019fcaSJed Brown 
14210189643fSJed Brown #if defined(PETSC_USE_DEBUG)
14220189643fSJed Brown   for (i=0; i<vs->nr; i++) {
14230189643fSJed Brown     for (j=0; j<vs->nc; j++) {
14240189643fSJed Brown       PetscInt m,n,M,N,mi,ni,Mi,Ni;
14250189643fSJed Brown       Mat      B = vs->m[i][j];
14260189643fSJed Brown       if (!B) continue;
14270189643fSJed Brown       ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr);
14280189643fSJed Brown       ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr);
14290189643fSJed Brown       ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr);
14300189643fSJed Brown       ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr);
14310189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr);
14320189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr);
1433ce94432eSBarry Smith       if (M != Mi || N != Ni) SETERRQ6(PetscObjectComm((PetscObject)sub),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);
1434ce94432eSBarry Smith       if (m != mi || n != ni) SETERRQ6(PetscObjectComm((PetscObject)sub),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);
14350189643fSJed Brown     }
14360189643fSJed Brown   }
14370189643fSJed Brown #endif
1438a061e289SJed Brown 
1439a061e289SJed Brown   /* Set A->assembled if all non-null blocks are currently assembled */
1440a061e289SJed Brown   for (i=0; i<vs->nr; i++) {
1441a061e289SJed Brown     for (j=0; j<vs->nc; j++) {
1442a061e289SJed Brown       if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0);
1443a061e289SJed Brown     }
1444a061e289SJed Brown   }
1445a061e289SJed Brown   A->assembled = PETSC_TRUE;
1446d8588912SDave May   PetscFunctionReturn(0);
1447d8588912SDave May }
1448d8588912SDave May 
1449d8588912SDave May #undef __FUNCT__
1450d8588912SDave May #define __FUNCT__ "MatCreateNest"
145145c38901SJed Brown /*@C
1452659c6bb0SJed Brown    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
1453659c6bb0SJed Brown 
1454659c6bb0SJed Brown    Collective on Mat
1455659c6bb0SJed Brown 
1456659c6bb0SJed Brown    Input Parameter:
1457659c6bb0SJed Brown +  comm - Communicator for the new Mat
1458659c6bb0SJed Brown .  nr - number of nested row blocks
14590298fd71SBarry Smith .  is_row - index sets for each nested row block, or NULL to make contiguous
1460659c6bb0SJed Brown .  nc - number of nested column blocks
14610298fd71SBarry Smith .  is_col - index sets for each nested column block, or NULL to make contiguous
14620298fd71SBarry Smith -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1463659c6bb0SJed Brown 
1464659c6bb0SJed Brown    Output Parameter:
1465659c6bb0SJed Brown .  B - new matrix
1466659c6bb0SJed Brown 
1467659c6bb0SJed Brown    Level: advanced
1468659c6bb0SJed Brown 
1469950540a4SJed Brown .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST
1470659c6bb0SJed Brown @*/
14717087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1472d8588912SDave May {
1473d8588912SDave May   Mat            A;
1474d8588912SDave May   PetscErrorCode ierr;
1475d8588912SDave May 
1476d8588912SDave May   PetscFunctionBegin;
1477c8883902SJed Brown   *B   = 0;
1478d8588912SDave May   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1479c8883902SJed Brown   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
148091a28eb3SBarry Smith   A->preallocated = PETSC_TRUE;
1481c8883902SJed Brown   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1482d8588912SDave May   *B   = A;
1483d8588912SDave May   PetscFunctionReturn(0);
1484d8588912SDave May }
1485659c6bb0SJed Brown 
1486629c3df2SDmitry Karpeev #undef __FUNCT__
1487629c3df2SDmitry Karpeev #define __FUNCT__ "MatConvert_Nest_AIJ"
1488cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1489629c3df2SDmitry Karpeev {
1490629c3df2SDmitry Karpeev   PetscErrorCode ierr;
1491629c3df2SDmitry Karpeev   Mat_Nest       *nest = (Mat_Nest*)A->data;
149283b1a929SMark Adams   PetscInt       m,n,M,N,i,j,k,*dnnz,*onnz,rstart;
1493649b366bSFande Kong   PetscInt       cstart,cend;
1494629c3df2SDmitry Karpeev   Mat            C;
1495629c3df2SDmitry Karpeev 
1496629c3df2SDmitry Karpeev   PetscFunctionBegin;
1497629c3df2SDmitry Karpeev   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
1498629c3df2SDmitry Karpeev   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
1499649b366bSFande Kong   ierr = MatGetOwnershipRangeColumn(A,&cstart,&cend);CHKERRQ(ierr);
1500629c3df2SDmitry Karpeev   switch (reuse) {
1501629c3df2SDmitry Karpeev   case MAT_INITIAL_MATRIX:
1502ce94432eSBarry Smith     ierr    = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
1503629c3df2SDmitry Karpeev     ierr    = MatSetType(C,newtype);CHKERRQ(ierr);
1504629c3df2SDmitry Karpeev     ierr    = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr);
1505629c3df2SDmitry Karpeev     *newmat = C;
1506629c3df2SDmitry Karpeev     break;
1507629c3df2SDmitry Karpeev   case MAT_REUSE_MATRIX:
1508629c3df2SDmitry Karpeev     C = *newmat;
1509629c3df2SDmitry Karpeev     break;
1510ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse");
1511629c3df2SDmitry Karpeev   }
1512785e854fSJed Brown   ierr = PetscMalloc1(2*m,&dnnz);CHKERRQ(ierr);
1513629c3df2SDmitry Karpeev   onnz = dnnz + m;
1514629c3df2SDmitry Karpeev   for (k=0; k<m; k++) {
1515629c3df2SDmitry Karpeev     dnnz[k] = 0;
1516629c3df2SDmitry Karpeev     onnz[k] = 0;
1517629c3df2SDmitry Karpeev   }
1518629c3df2SDmitry Karpeev   for (j=0; j<nest->nc; ++j) {
1519629c3df2SDmitry Karpeev     IS             bNis;
1520629c3df2SDmitry Karpeev     PetscInt       bN;
1521629c3df2SDmitry Karpeev     const PetscInt *bNindices;
1522629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
1523629c3df2SDmitry Karpeev     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1524629c3df2SDmitry Karpeev     ierr = ISGetSize(bNis, &bN);CHKERRQ(ierr);
1525629c3df2SDmitry Karpeev     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1526629c3df2SDmitry Karpeev     for (i=0; i<nest->nr; ++i) {
1527629c3df2SDmitry Karpeev       PetscSF        bmsf;
1528649b366bSFande Kong       PetscSFNode    *iremote;
1529629c3df2SDmitry Karpeev       Mat            B;
1530649b366bSFande Kong       PetscInt       bm, *sub_dnnz,*sub_onnz, br;
1531629c3df2SDmitry Karpeev       const PetscInt *bmindices;
1532629c3df2SDmitry Karpeev       B = nest->m[i][j];
1533629c3df2SDmitry Karpeev       if (!B) continue;
1534629c3df2SDmitry Karpeev       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1535629c3df2SDmitry Karpeev       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1536ce94432eSBarry Smith       ierr = PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);CHKERRQ(ierr);
1537649b366bSFande Kong       ierr = PetscMalloc1(bm,&iremote);CHKERRQ(ierr);
1538649b366bSFande Kong       ierr = PetscMalloc1(bm,&sub_dnnz);CHKERRQ(ierr);
1539649b366bSFande Kong       ierr = PetscMalloc1(bm,&sub_onnz);CHKERRQ(ierr);
1540649b366bSFande Kong       for (k = 0; k < bm; ++k){
1541649b366bSFande Kong     	sub_dnnz[k] = 0;
1542649b366bSFande Kong     	sub_onnz[k] = 0;
1543649b366bSFande Kong       }
1544629c3df2SDmitry Karpeev       /*
1545629c3df2SDmitry Karpeev        Locate the owners for all of the locally-owned global row indices for this row block.
1546629c3df2SDmitry Karpeev        These determine the roots of PetscSF used to communicate preallocation data to row owners.
1547629c3df2SDmitry Karpeev        The roots correspond to the dnnz and onnz entries; thus, there are two roots per row.
1548629c3df2SDmitry Karpeev        */
154983b1a929SMark Adams       ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr);
1550629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
1551649b366bSFande Kong         PetscInt       row = bmindices[br], rowowner = 0, brncols, col;
1552629c3df2SDmitry Karpeev         const PetscInt *brcols;
1553a4b3d3acSMatthew G Knepley         PetscInt       rowrel = 0; /* row's relative index on its owner rank */
1554629c3df2SDmitry Karpeev         ierr      = PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);CHKERRQ(ierr);
1555649b366bSFande Kong         /* how many roots  */
1556649b366bSFande Kong         iremote[br].rank = rowowner; iremote[br].index = rowrel;           /* edge from bmdnnz to dnnz */
1557649b366bSFande Kong         /* get nonzero pattern */
155883b1a929SMark Adams         ierr = MatGetRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr);
1559629c3df2SDmitry Karpeev         for (k=0; k<brncols; k++) {
1560629c3df2SDmitry Karpeev           col  = bNindices[brcols[k]];
1561649b366bSFande Kong           if(col>=A->cmap->range[rowowner] && col<A->cmap->range[rowowner+1]){
1562649b366bSFande Kong         	sub_dnnz[br]++;
1563649b366bSFande Kong           }else{
1564649b366bSFande Kong         	sub_onnz[br]++;
1565649b366bSFande Kong           }
1566629c3df2SDmitry Karpeev         }
156783b1a929SMark Adams         ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr);
1568629c3df2SDmitry Karpeev       }
1569629c3df2SDmitry Karpeev       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1570629c3df2SDmitry Karpeev       /* bsf will have to take care of disposing of bedges. */
1571649b366bSFande Kong       ierr = PetscSFSetGraph(bmsf,m,bm,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1572649b366bSFande Kong       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr);
1573649b366bSFande Kong       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr);
1574649b366bSFande Kong       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr);
1575649b366bSFande Kong       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr);
1576649b366bSFande Kong       ierr = PetscFree(sub_dnnz);CHKERRQ(ierr);
1577649b366bSFande Kong       ierr = PetscFree(sub_onnz);CHKERRQ(ierr);
1578629c3df2SDmitry Karpeev       ierr = PetscSFDestroy(&bmsf);CHKERRQ(ierr);
1579629c3df2SDmitry Karpeev     }
158022d28d08SBarry Smith     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
1581629c3df2SDmitry Karpeev     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
1582629c3df2SDmitry Karpeev   }
1583629c3df2SDmitry Karpeev   ierr = MatSeqAIJSetPreallocation(C,0,dnnz);CHKERRQ(ierr);
1584629c3df2SDmitry Karpeev   ierr = MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);CHKERRQ(ierr);
1585629c3df2SDmitry Karpeev   ierr = PetscFree(dnnz);CHKERRQ(ierr);
1586629c3df2SDmitry Karpeev 
1587629c3df2SDmitry Karpeev   /* Fill by row */
1588629c3df2SDmitry Karpeev   for (j=0; j<nest->nc; ++j) {
1589629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
1590629c3df2SDmitry Karpeev     IS             bNis;
1591629c3df2SDmitry Karpeev     PetscInt       bN;
1592629c3df2SDmitry Karpeev     const PetscInt *bNindices;
1593629c3df2SDmitry Karpeev     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1594629c3df2SDmitry Karpeev     ierr = ISGetSize(bNis,&bN);CHKERRQ(ierr);
1595629c3df2SDmitry Karpeev     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1596629c3df2SDmitry Karpeev     for (i=0; i<nest->nr; ++i) {
1597629c3df2SDmitry Karpeev       Mat            B;
1598629c3df2SDmitry Karpeev       PetscInt       bm, br;
1599629c3df2SDmitry Karpeev       const PetscInt *bmindices;
1600629c3df2SDmitry Karpeev       B = nest->m[i][j];
1601629c3df2SDmitry Karpeev       if (!B) continue;
1602629c3df2SDmitry Karpeev       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1603629c3df2SDmitry Karpeev       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
160483b1a929SMark Adams       ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr);
1605629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
1606629c3df2SDmitry Karpeev         PetscInt          row = bmindices[br], brncols,  *cols;
1607629c3df2SDmitry Karpeev         const PetscInt    *brcols;
1608629c3df2SDmitry Karpeev         const PetscScalar *brcoldata;
160983b1a929SMark Adams         ierr = MatGetRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
1610785e854fSJed Brown         ierr = PetscMalloc1(brncols,&cols);CHKERRQ(ierr);
161126fbe8dcSKarl Rupp         for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]];
1612629c3df2SDmitry Karpeev         /*
1613629c3df2SDmitry Karpeev           Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match.
1614629c3df2SDmitry Karpeev           Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES.
1615629c3df2SDmitry Karpeev          */
1616a2ea699eSBarry Smith         ierr = MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);CHKERRQ(ierr);
161783b1a929SMark Adams         ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
1618629c3df2SDmitry Karpeev         ierr = PetscFree(cols);CHKERRQ(ierr);
1619629c3df2SDmitry Karpeev       }
1620629c3df2SDmitry Karpeev       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1621629c3df2SDmitry Karpeev     }
1622a2ea699eSBarry Smith     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
1623629c3df2SDmitry Karpeev     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
1624629c3df2SDmitry Karpeev   }
1625629c3df2SDmitry Karpeev   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1626629c3df2SDmitry Karpeev   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1627629c3df2SDmitry Karpeev   PetscFunctionReturn(0);
1628629c3df2SDmitry Karpeev }
1629629c3df2SDmitry Karpeev 
1630659c6bb0SJed Brown /*MC
1631659c6bb0SJed Brown   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
1632659c6bb0SJed Brown 
1633659c6bb0SJed Brown   Level: intermediate
1634659c6bb0SJed Brown 
1635659c6bb0SJed Brown   Notes:
1636659c6bb0SJed Brown   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1637659c6bb0SJed Brown   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1638950540a4SJed Brown   It is usually used with DMComposite and DMCreateMatrix()
1639659c6bb0SJed Brown 
1640659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest()
1641659c6bb0SJed Brown M*/
1642c8883902SJed Brown #undef __FUNCT__
1643c8883902SJed Brown #define __FUNCT__ "MatCreate_Nest"
16448cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A)
1645c8883902SJed Brown {
1646c8883902SJed Brown   Mat_Nest       *s;
1647c8883902SJed Brown   PetscErrorCode ierr;
1648c8883902SJed Brown 
1649c8883902SJed Brown   PetscFunctionBegin;
1650b00a9115SJed Brown   ierr    = PetscNewLog(A,&s);CHKERRQ(ierr);
1651c8883902SJed Brown   A->data = (void*)s;
1652e7c19651SJed Brown 
1653e7c19651SJed Brown   s->nr            = -1;
1654e7c19651SJed Brown   s->nc            = -1;
16550298fd71SBarry Smith   s->m             = NULL;
1656e7c19651SJed Brown   s->splitassembly = PETSC_FALSE;
1657c8883902SJed Brown 
1658c8883902SJed Brown   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
165926fbe8dcSKarl Rupp 
1660c8883902SJed Brown   A->ops->mult                  = MatMult_Nest;
16619194d70fSJed Brown   A->ops->multadd               = MatMultAdd_Nest;
1662c8883902SJed Brown   A->ops->multtranspose         = MatMultTranspose_Nest;
16639194d70fSJed Brown   A->ops->multtransposeadd      = MatMultTransposeAdd_Nest;
1664*f8170845SAlex Fikl   A->ops->transpose             = MatTranspose_Nest;
1665c8883902SJed Brown   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1666c8883902SJed Brown   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1667c8883902SJed Brown   A->ops->zeroentries           = MatZeroEntries_Nest;
1668c222c20dSDavid Ham   A->ops->copy                  = MatCopy_Nest;
1669c8883902SJed Brown   A->ops->duplicate             = MatDuplicate_Nest;
1670c8883902SJed Brown   A->ops->getsubmatrix          = MatGetSubMatrix_Nest;
1671c8883902SJed Brown   A->ops->destroy               = MatDestroy_Nest;
1672c8883902SJed Brown   A->ops->view                  = MatView_Nest;
1673c8883902SJed Brown   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1674c8883902SJed Brown   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1675c8883902SJed Brown   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
1676429bac76SJed Brown   A->ops->getdiagonal           = MatGetDiagonal_Nest;
1677429bac76SJed Brown   A->ops->diagonalscale         = MatDiagonalScale_Nest;
1678a061e289SJed Brown   A->ops->scale                 = MatScale_Nest;
1679a061e289SJed Brown   A->ops->shift                 = MatShift_Nest;
168013135bc6SAlex Fikl   A->ops->diagonalset           = MatDiagonalSet_Nest;
1681*f8170845SAlex Fikl   A->ops->setrandom             = MatSetRandom_Nest;
1682c8883902SJed Brown 
1683c8883902SJed Brown   A->spptr        = 0;
1684c8883902SJed Brown   A->assembled    = PETSC_FALSE;
1685c8883902SJed Brown 
1686c8883902SJed Brown   /* expose Nest api's */
1687bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",   MatNestGetSubMat_Nest);CHKERRQ(ierr);
1688bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",   MatNestSetSubMat_Nest);CHKERRQ(ierr);
1689bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",  MatNestGetSubMats_Nest);CHKERRQ(ierr);
1690bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",     MatNestGetSize_Nest);CHKERRQ(ierr);
1691bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",      MatNestGetISs_Nest);CHKERRQ(ierr);
1692bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C", MatNestGetLocalISs_Nest);CHKERRQ(ierr);
1693bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",  MatNestSetVecType_Nest);CHKERRQ(ierr);
1694bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",  MatNestSetSubMats_Nest);CHKERRQ(ierr);
169583b1a929SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",MatConvert_Nest_AIJ);CHKERRQ(ierr);
1696c8883902SJed Brown 
1697c8883902SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
1698c8883902SJed Brown   PetscFunctionReturn(0);
1699c8883902SJed Brown }
1700