xref: /petsc/src/mat/impls/nest/matnest.c (revision 13135bc6d86bc2722b0d134ace4ca52edf6ae00b)
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__
148e2d7f03fSJed Brown #define __FUNCT__ "MatNestDestroyISList"
149e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list)
150e2d7f03fSJed Brown {
151e2d7f03fSJed Brown   PetscErrorCode ierr;
152e2d7f03fSJed Brown   IS             *lst = *list;
153e2d7f03fSJed Brown   PetscInt       i;
154e2d7f03fSJed Brown 
155e2d7f03fSJed Brown   PetscFunctionBegin;
156e2d7f03fSJed Brown   if (!lst) PetscFunctionReturn(0);
1576bf464f9SBarry Smith   for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);}
158e2d7f03fSJed Brown   ierr  = PetscFree(lst);CHKERRQ(ierr);
1590298fd71SBarry Smith   *list = NULL;
160e2d7f03fSJed Brown   PetscFunctionReturn(0);
161e2d7f03fSJed Brown }
162e2d7f03fSJed Brown 
163e2d7f03fSJed Brown #undef __FUNCT__
164d8588912SDave May #define __FUNCT__ "MatDestroy_Nest"
165207556f9SJed Brown static PetscErrorCode MatDestroy_Nest(Mat A)
166d8588912SDave May {
167d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
168d8588912SDave May   PetscInt       i,j;
169d8588912SDave May   PetscErrorCode ierr;
170d8588912SDave May 
171d8588912SDave May   PetscFunctionBegin;
172d8588912SDave May   /* release the matrices and the place holders */
173e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr);
174e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr);
175e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
176e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
177d8588912SDave May 
178d8588912SDave May   ierr = PetscFree(vs->row_len);CHKERRQ(ierr);
179d8588912SDave May   ierr = PetscFree(vs->col_len);CHKERRQ(ierr);
180d8588912SDave May 
181207556f9SJed Brown   ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr);
182207556f9SJed Brown 
183d8588912SDave May   /* release the matrices and the place holders */
184d8588912SDave May   if (vs->m) {
185d8588912SDave May     for (i=0; i<vs->nr; i++) {
186d8588912SDave May       for (j=0; j<vs->nc; j++) {
1876bf464f9SBarry Smith         ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr);
188d8588912SDave May       }
189d8588912SDave May       ierr = PetscFree(vs->m[i]);CHKERRQ(ierr);
190d8588912SDave May     }
191d8588912SDave May     ierr = PetscFree(vs->m);CHKERRQ(ierr);
192d8588912SDave May   }
193bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
194d8588912SDave May 
195bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",0);CHKERRQ(ierr);
196bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",0);CHKERRQ(ierr);
197bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",0);CHKERRQ(ierr);
198bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",0);CHKERRQ(ierr);
199bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",0);CHKERRQ(ierr);
200bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",0);CHKERRQ(ierr);
201bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",0);CHKERRQ(ierr);
202bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",0);CHKERRQ(ierr);
203d8588912SDave May   PetscFunctionReturn(0);
204d8588912SDave May }
205d8588912SDave May 
206d8588912SDave May #undef __FUNCT__
207d8588912SDave May #define __FUNCT__ "MatAssemblyBegin_Nest"
208207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
209d8588912SDave May {
210d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
211d8588912SDave May   PetscInt       i,j;
212d8588912SDave May   PetscErrorCode ierr;
213d8588912SDave May 
214d8588912SDave May   PetscFunctionBegin;
215d8588912SDave May   for (i=0; i<vs->nr; i++) {
216d8588912SDave May     for (j=0; j<vs->nc; j++) {
217e7c19651SJed Brown       if (vs->m[i][j]) {
218e7c19651SJed Brown         ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr);
219e7c19651SJed Brown         if (!vs->splitassembly) {
220e7c19651SJed Brown           /* Note: split assembly will fail if the same block appears more than once (even indirectly through a nested
221e7c19651SJed 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
222e7c19651SJed Brown            * already performing an assembly, but the result would by more complicated and appears to offer less
223e7c19651SJed Brown            * potential for diagnostics and correctness checking. Split assembly should be fixed once there is an
224e7c19651SJed Brown            * interface for libraries to make asynchronous progress in "user-defined non-blocking collectives".
225e7c19651SJed Brown            */
226e7c19651SJed Brown           ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr);
227e7c19651SJed Brown         }
228e7c19651SJed Brown       }
229d8588912SDave May     }
230d8588912SDave May   }
231d8588912SDave May   PetscFunctionReturn(0);
232d8588912SDave May }
233d8588912SDave May 
234d8588912SDave May #undef __FUNCT__
235d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest"
236207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
237d8588912SDave May {
238d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
239d8588912SDave May   PetscInt       i,j;
240d8588912SDave May   PetscErrorCode ierr;
241d8588912SDave May 
242d8588912SDave May   PetscFunctionBegin;
243d8588912SDave May   for (i=0; i<vs->nr; i++) {
244d8588912SDave May     for (j=0; j<vs->nc; j++) {
245e7c19651SJed Brown       if (vs->m[i][j]) {
246e7c19651SJed Brown         if (vs->splitassembly) {
247e7c19651SJed Brown           ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr);
248e7c19651SJed Brown         }
249e7c19651SJed Brown       }
250d8588912SDave May     }
251d8588912SDave May   }
252d8588912SDave May   PetscFunctionReturn(0);
253d8588912SDave May }
254d8588912SDave May 
255d8588912SDave May #undef __FUNCT__
256f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatRow"
257f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
258d8588912SDave May {
259207556f9SJed Brown   PetscErrorCode ierr;
260f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
261f349c1fdSJed Brown   PetscInt       j;
262f349c1fdSJed Brown   Mat            sub;
263d8588912SDave May 
264d8588912SDave May   PetscFunctionBegin;
2650298fd71SBarry Smith   sub = (row < vs->nc) ? vs->m[row][row] : (Mat)NULL; /* Prefer to find on the diagonal */
266f349c1fdSJed Brown   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
2674994cf47SJed Brown   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
268f349c1fdSJed Brown   *B = sub;
269f349c1fdSJed Brown   PetscFunctionReturn(0);
270d8588912SDave May }
271d8588912SDave May 
272f349c1fdSJed Brown #undef __FUNCT__
273f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatCol"
274f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
275f349c1fdSJed Brown {
276207556f9SJed Brown   PetscErrorCode ierr;
277f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
278f349c1fdSJed Brown   PetscInt       i;
279f349c1fdSJed Brown   Mat            sub;
280f349c1fdSJed Brown 
281f349c1fdSJed Brown   PetscFunctionBegin;
2820298fd71SBarry Smith   sub = (col < vs->nr) ? vs->m[col][col] : (Mat)NULL; /* Prefer to find on the diagonal */
283f349c1fdSJed Brown   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
2844994cf47SJed Brown   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
285f349c1fdSJed Brown   *B = sub;
286f349c1fdSJed Brown   PetscFunctionReturn(0);
287d8588912SDave May }
288d8588912SDave May 
289f349c1fdSJed Brown #undef __FUNCT__
290f349c1fdSJed Brown #define __FUNCT__ "MatNestFindIS"
291f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
292f349c1fdSJed Brown {
293f349c1fdSJed Brown   PetscErrorCode ierr;
294f349c1fdSJed Brown   PetscInt       i;
295f349c1fdSJed Brown   PetscBool      flg;
296f349c1fdSJed Brown 
297f349c1fdSJed Brown   PetscFunctionBegin;
298f349c1fdSJed Brown   PetscValidPointer(list,3);
299f349c1fdSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,4);
300f349c1fdSJed Brown   PetscValidIntPointer(found,5);
301f349c1fdSJed Brown   *found = -1;
302f349c1fdSJed Brown   for (i=0; i<n; i++) {
303207556f9SJed Brown     if (!list[i]) continue;
304f349c1fdSJed Brown     ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr);
305f349c1fdSJed Brown     if (flg) {
306f349c1fdSJed Brown       *found = i;
307f349c1fdSJed Brown       PetscFunctionReturn(0);
308f349c1fdSJed Brown     }
309f349c1fdSJed Brown   }
310ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Could not find index set");
311f349c1fdSJed Brown   PetscFunctionReturn(0);
312f349c1fdSJed Brown }
313f349c1fdSJed Brown 
314f349c1fdSJed Brown #undef __FUNCT__
3158188e55aSJed Brown #define __FUNCT__ "MatNestGetRow"
3168188e55aSJed Brown /* Get a block row as a new MatNest */
3178188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B)
3188188e55aSJed Brown {
3198188e55aSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
3208188e55aSJed Brown   char           keyname[256];
3218188e55aSJed Brown   PetscErrorCode ierr;
3228188e55aSJed Brown 
3238188e55aSJed Brown   PetscFunctionBegin;
3240298fd71SBarry Smith   *B   = NULL;
3258caf3d72SBarry Smith   ierr = PetscSNPrintf(keyname,sizeof(keyname),"NestRow_%D",row);CHKERRQ(ierr);
3268188e55aSJed Brown   ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr);
3278188e55aSJed Brown   if (*B) PetscFunctionReturn(0);
3288188e55aSJed Brown 
329ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),1,NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr);
33026fbe8dcSKarl Rupp 
3318188e55aSJed Brown   (*B)->assembled = A->assembled;
33226fbe8dcSKarl Rupp 
3338188e55aSJed Brown   ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr);
3348188e55aSJed Brown   ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */
3358188e55aSJed Brown   PetscFunctionReturn(0);
3368188e55aSJed Brown }
3378188e55aSJed Brown 
3388188e55aSJed Brown #undef __FUNCT__
339f349c1fdSJed Brown #define __FUNCT__ "MatNestFindSubMat"
340f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
341f349c1fdSJed Brown {
342f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
3438188e55aSJed Brown   PetscErrorCode ierr;
3446b3a5b13SJed Brown   PetscInt       row,col;
345e072481dSJed Brown   PetscBool      same,isFullCol,isFullColGlobal;
346f349c1fdSJed Brown 
347f349c1fdSJed Brown   PetscFunctionBegin;
3488188e55aSJed Brown   /* Check if full column space. This is a hack */
3498188e55aSJed Brown   isFullCol = PETSC_FALSE;
350251f4c67SDmitry Karpeev   ierr      = PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr);
3518188e55aSJed Brown   if (same) {
35277019fcaSJed Brown     PetscInt n,first,step,i,an,am,afirst,astep;
3538188e55aSJed Brown     ierr      = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
3548188e55aSJed Brown     ierr      = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
35577019fcaSJed Brown     isFullCol = PETSC_TRUE;
35605ce4453SJed Brown     for (i=0,an=A->cmap->rstart; i<vs->nc; i++) {
35777019fcaSJed Brown       ierr = ISStrideGetInfo(is->col[i],&afirst,&astep);CHKERRQ(ierr);
35877019fcaSJed Brown       ierr = ISGetLocalSize(is->col[i],&am);CHKERRQ(ierr);
35977019fcaSJed Brown       if (afirst != an || astep != step) isFullCol = PETSC_FALSE;
36077019fcaSJed Brown       an += am;
36177019fcaSJed Brown     }
36205ce4453SJed Brown     if (an != A->cmap->rstart+n) isFullCol = PETSC_FALSE;
3638188e55aSJed Brown   }
364b2566f29SBarry Smith   ierr = MPIU_Allreduce(&isFullCol,&isFullColGlobal,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)iscol));CHKERRQ(ierr);
3658188e55aSJed Brown 
366e072481dSJed Brown   if (isFullColGlobal) {
3678188e55aSJed Brown     PetscInt row;
3688188e55aSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
3698188e55aSJed Brown     ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr);
3708188e55aSJed Brown   } else {
371f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
372f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr);
373b6480e04SStefano Zampini     if (!vs->m[row][col]) {
374b6480e04SStefano Zampini       PetscInt lr,lc;
375b6480e04SStefano Zampini 
376b6480e04SStefano Zampini       ierr = MatCreate(PetscObjectComm((PetscObject)A),&vs->m[row][col]);CHKERRQ(ierr);
377b6480e04SStefano Zampini       ierr = ISGetLocalSize(vs->isglobal.row[row],&lr);CHKERRQ(ierr);
378b6480e04SStefano Zampini       ierr = ISGetLocalSize(vs->isglobal.col[col],&lc);CHKERRQ(ierr);
379b6480e04SStefano Zampini       ierr = MatSetSizes(vs->m[row][col],lr,lc,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
380b6480e04SStefano Zampini       ierr = MatSetUp(vs->m[row][col]);CHKERRQ(ierr);
381b6480e04SStefano Zampini       ierr = MatAssemblyBegin(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
382b6480e04SStefano Zampini       ierr = MatAssemblyEnd(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
383b6480e04SStefano Zampini     }
384f349c1fdSJed Brown     *B = vs->m[row][col];
3858188e55aSJed Brown   }
386f349c1fdSJed Brown   PetscFunctionReturn(0);
387f349c1fdSJed Brown }
388f349c1fdSJed Brown 
389f349c1fdSJed Brown #undef __FUNCT__
390f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest"
391f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
392f349c1fdSJed Brown {
393f349c1fdSJed Brown   PetscErrorCode ierr;
394f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
395f349c1fdSJed Brown   Mat            sub;
396f349c1fdSJed Brown 
397f349c1fdSJed Brown   PetscFunctionBegin;
398f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr);
399f349c1fdSJed Brown   switch (reuse) {
400f349c1fdSJed Brown   case MAT_INITIAL_MATRIX:
4017874fa86SDave May     if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); }
402f349c1fdSJed Brown     *B = sub;
403f349c1fdSJed Brown     break;
404f349c1fdSJed Brown   case MAT_REUSE_MATRIX:
405ce94432eSBarry Smith     if (sub != *B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
406f349c1fdSJed Brown     break;
407f349c1fdSJed Brown   case MAT_IGNORE_MATRIX:       /* Nothing to do */
408f349c1fdSJed Brown     break;
409511c6705SHong Zhang   case MAT_INPLACE_MATRIX:       /* Nothing to do */
410511c6705SHong Zhang     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX is not supported yet");
411511c6705SHong Zhang     break;
412f349c1fdSJed Brown   }
413f349c1fdSJed Brown   PetscFunctionReturn(0);
414f349c1fdSJed Brown }
415f349c1fdSJed Brown 
416f349c1fdSJed Brown #undef __FUNCT__
417f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest"
418f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
419f349c1fdSJed Brown {
420f349c1fdSJed Brown   PetscErrorCode ierr;
421f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
422f349c1fdSJed Brown   Mat            sub;
423f349c1fdSJed Brown 
424f349c1fdSJed Brown   PetscFunctionBegin;
425f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
426f349c1fdSJed Brown   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
427f349c1fdSJed Brown   if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);}
428f349c1fdSJed Brown   *B = sub;
429d8588912SDave May   PetscFunctionReturn(0);
430d8588912SDave May }
431d8588912SDave May 
432d8588912SDave May #undef __FUNCT__
433d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest"
434207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
435d8588912SDave May {
436d8588912SDave May   PetscErrorCode ierr;
437f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
438f349c1fdSJed Brown   Mat            sub;
439d8588912SDave May 
440d8588912SDave May   PetscFunctionBegin;
441f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
442ce94432eSBarry Smith   if (*B != sub) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
443f349c1fdSJed Brown   if (sub) {
444ce94432eSBarry Smith     if (((PetscObject)sub)->refct <= 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
4456bf464f9SBarry Smith     ierr = MatDestroy(B);CHKERRQ(ierr);
446d8588912SDave May   }
447d8588912SDave May   PetscFunctionReturn(0);
448d8588912SDave May }
449d8588912SDave May 
450d8588912SDave May #undef __FUNCT__
4517874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest"
4527874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
4537874fa86SDave May {
4547874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
4557874fa86SDave May   PetscInt       i;
4567874fa86SDave May   PetscErrorCode ierr;
4577874fa86SDave May 
4587874fa86SDave May   PetscFunctionBegin;
4597874fa86SDave May   for (i=0; i<bA->nr; i++) {
460429bac76SJed Brown     Vec bv;
461429bac76SJed Brown     ierr = VecGetSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
4627874fa86SDave May     if (bA->m[i][i]) {
463429bac76SJed Brown       ierr = MatGetDiagonal(bA->m[i][i],bv);CHKERRQ(ierr);
4647874fa86SDave May     } else {
4655159a857SMatthew G. Knepley       ierr = VecSet(bv,0.0);CHKERRQ(ierr);
4667874fa86SDave May     }
467429bac76SJed Brown     ierr = VecRestoreSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
4687874fa86SDave May   }
4697874fa86SDave May   PetscFunctionReturn(0);
4707874fa86SDave May }
4717874fa86SDave May 
4727874fa86SDave May #undef __FUNCT__
4737874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest"
4747874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
4757874fa86SDave May {
4767874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
477429bac76SJed Brown   Vec            bl,*br;
4787874fa86SDave May   PetscInt       i,j;
4797874fa86SDave May   PetscErrorCode ierr;
4807874fa86SDave May 
4817874fa86SDave May   PetscFunctionBegin;
4823f800ebeSJed Brown   ierr = PetscCalloc1(bA->nc,&br);CHKERRQ(ierr);
4832e6472ebSElliott Sales de Andrade   if (r) {
484429bac76SJed Brown     for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
4852e6472ebSElliott Sales de Andrade   }
4862e6472ebSElliott Sales de Andrade   bl = NULL;
4877874fa86SDave May   for (i=0; i<bA->nr; i++) {
4882e6472ebSElliott Sales de Andrade     if (l) {
489429bac76SJed Brown       ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
4902e6472ebSElliott Sales de Andrade     }
4917874fa86SDave May     for (j=0; j<bA->nc; j++) {
4927874fa86SDave May       if (bA->m[i][j]) {
493429bac76SJed Brown         ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr);
4947874fa86SDave May       }
4957874fa86SDave May     }
4962e6472ebSElliott Sales de Andrade     if (l) {
497a061e289SJed Brown       ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
4987874fa86SDave May     }
4992e6472ebSElliott Sales de Andrade   }
5002e6472ebSElliott Sales de Andrade   if (r) {
501429bac76SJed Brown     for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
5022e6472ebSElliott Sales de Andrade   }
503429bac76SJed Brown   ierr = PetscFree(br);CHKERRQ(ierr);
5047874fa86SDave May   PetscFunctionReturn(0);
5057874fa86SDave May }
5067874fa86SDave May 
5077874fa86SDave May #undef __FUNCT__
508a061e289SJed Brown #define __FUNCT__ "MatScale_Nest"
509a061e289SJed Brown static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a)
510a061e289SJed Brown {
511a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
512a061e289SJed Brown   PetscInt       i,j;
513a061e289SJed Brown   PetscErrorCode ierr;
514a061e289SJed Brown 
515a061e289SJed Brown   PetscFunctionBegin;
516a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
517a061e289SJed Brown     for (j=0; j<bA->nc; j++) {
518a061e289SJed Brown       if (bA->m[i][j]) {
519a061e289SJed Brown         ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr);
520a061e289SJed Brown       }
521a061e289SJed Brown     }
522a061e289SJed Brown   }
523a061e289SJed Brown   PetscFunctionReturn(0);
524a061e289SJed Brown }
525a061e289SJed Brown 
526a061e289SJed Brown #undef __FUNCT__
527a061e289SJed Brown #define __FUNCT__ "MatShift_Nest"
528a061e289SJed Brown static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a)
529a061e289SJed Brown {
530a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
531a061e289SJed Brown   PetscInt       i;
532a061e289SJed Brown   PetscErrorCode ierr;
533a061e289SJed Brown 
534a061e289SJed Brown   PetscFunctionBegin;
535a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
536ce94432eSBarry 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);
537a061e289SJed Brown     ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr);
538a061e289SJed Brown   }
539a061e289SJed Brown   PetscFunctionReturn(0);
540a061e289SJed Brown }
541a061e289SJed Brown 
542a061e289SJed Brown #undef __FUNCT__
543*13135bc6SAlex Fikl #define __FUNCT__ "MatDiagonalSet_Nest"
544*13135bc6SAlex Fikl static PetscErrorCode MatDiagonalSet_Nest(Mat A,Vec D,InsertMode is)
545*13135bc6SAlex Fikl {
546*13135bc6SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data;
547*13135bc6SAlex Fikl   PetscInt       i;
548*13135bc6SAlex Fikl   PetscErrorCode ierr;
549*13135bc6SAlex Fikl 
550*13135bc6SAlex Fikl   PetscFunctionBegin;
551*13135bc6SAlex Fikl   for (i=0; i<bA->nr; i++) {
552*13135bc6SAlex Fikl     Vec bv;
553*13135bc6SAlex Fikl     ierr = VecGetSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
554*13135bc6SAlex Fikl     if (bA->m[i][i]) {
555*13135bc6SAlex Fikl       ierr = MatDiagonalSet(bA->m[i][i],bv,is);CHKERRQ(ierr);
556*13135bc6SAlex Fikl     }
557*13135bc6SAlex Fikl     ierr = VecRestoreSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
558*13135bc6SAlex Fikl   }
559*13135bc6SAlex Fikl   PetscFunctionReturn(0);
560*13135bc6SAlex Fikl }
561*13135bc6SAlex Fikl 
562*13135bc6SAlex Fikl #undef __FUNCT__
5632a7a6963SBarry Smith #define __FUNCT__ "MatCreateVecs_Nest"
5642a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left)
565d8588912SDave May {
566d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
567d8588912SDave May   Vec            *L,*R;
568d8588912SDave May   MPI_Comm       comm;
569d8588912SDave May   PetscInt       i,j;
570d8588912SDave May   PetscErrorCode ierr;
571d8588912SDave May 
572d8588912SDave May   PetscFunctionBegin;
573ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
574d8588912SDave May   if (right) {
575d8588912SDave May     /* allocate R */
576854ce69bSBarry Smith     ierr = PetscMalloc1(bA->nc, &R);CHKERRQ(ierr);
577d8588912SDave May     /* Create the right vectors */
578d8588912SDave May     for (j=0; j<bA->nc; j++) {
579d8588912SDave May       for (i=0; i<bA->nr; i++) {
580d8588912SDave May         if (bA->m[i][j]) {
5812a7a6963SBarry Smith           ierr = MatCreateVecs(bA->m[i][j],&R[j],NULL);CHKERRQ(ierr);
582d8588912SDave May           break;
583d8588912SDave May         }
584d8588912SDave May       }
5856c4ed002SBarry Smith       if (i==bA->nr) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
586d8588912SDave May     }
587f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
588d8588912SDave May     /* hand back control to the nest vector */
589d8588912SDave May     for (j=0; j<bA->nc; j++) {
5906bf464f9SBarry Smith       ierr = VecDestroy(&R[j]);CHKERRQ(ierr);
591d8588912SDave May     }
592d8588912SDave May     ierr = PetscFree(R);CHKERRQ(ierr);
593d8588912SDave May   }
594d8588912SDave May 
595d8588912SDave May   if (left) {
596d8588912SDave May     /* allocate L */
597854ce69bSBarry Smith     ierr = PetscMalloc1(bA->nr, &L);CHKERRQ(ierr);
598d8588912SDave May     /* Create the left vectors */
599d8588912SDave May     for (i=0; i<bA->nr; i++) {
600d8588912SDave May       for (j=0; j<bA->nc; j++) {
601d8588912SDave May         if (bA->m[i][j]) {
6022a7a6963SBarry Smith           ierr = MatCreateVecs(bA->m[i][j],NULL,&L[i]);CHKERRQ(ierr);
603d8588912SDave May           break;
604d8588912SDave May         }
605d8588912SDave May       }
6066c4ed002SBarry Smith       if (j==bA->nc) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
607d8588912SDave May     }
608d8588912SDave May 
609f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
610d8588912SDave May     for (i=0; i<bA->nr; i++) {
6116bf464f9SBarry Smith       ierr = VecDestroy(&L[i]);CHKERRQ(ierr);
612d8588912SDave May     }
613d8588912SDave May 
614d8588912SDave May     ierr = PetscFree(L);CHKERRQ(ierr);
615d8588912SDave May   }
616d8588912SDave May   PetscFunctionReturn(0);
617d8588912SDave May }
618d8588912SDave May 
619d8588912SDave May #undef __FUNCT__
620d8588912SDave May #define __FUNCT__ "MatView_Nest"
621207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
622d8588912SDave May {
623d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
624d8588912SDave May   PetscBool      isascii;
625d8588912SDave May   PetscInt       i,j;
626d8588912SDave May   PetscErrorCode ierr;
627d8588912SDave May 
628d8588912SDave May   PetscFunctionBegin;
629251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
630d8588912SDave May   if (isascii) {
631d8588912SDave May 
632d8588912SDave May     PetscViewerASCIIPrintf(viewer,"Matrix object: \n");
633d8588912SDave May     PetscViewerASCIIPushTab(viewer);    /* push0 */
634d8588912SDave May     PetscViewerASCIIPrintf(viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc);
635d8588912SDave May 
636d8588912SDave May     PetscViewerASCIIPrintf(viewer,"MatNest structure: \n");
637d8588912SDave May     for (i=0; i<bA->nr; i++) {
638d8588912SDave May       for (j=0; j<bA->nc; j++) {
63919fd82e9SBarry Smith         MatType   type;
640270f95d7SJed Brown         char      name[256] = "",prefix[256] = "";
641d8588912SDave May         PetscInt  NR,NC;
642d8588912SDave May         PetscBool isNest = PETSC_FALSE;
643d8588912SDave May 
644d8588912SDave May         if (!bA->m[i][j]) {
6450298fd71SBarry Smith           PetscViewerASCIIPrintf(viewer, "(%D,%D) : NULL \n",i,j);
646d8588912SDave May           continue;
647d8588912SDave May         }
648d8588912SDave May         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
649d8588912SDave May         ierr = MatGetType(bA->m[i][j], &type);CHKERRQ(ierr);
6508caf3d72SBarry Smith         if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);}
6518caf3d72SBarry Smith         if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);}
652251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
653d8588912SDave May 
654270f95d7SJed Brown         ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr);
655d8588912SDave May 
656d8588912SDave May         if (isNest) {
657270f95d7SJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);  /* push1 */
658d8588912SDave May           ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr);
659270f95d7SJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop1 */
660d8588912SDave May         }
661d8588912SDave May       }
662d8588912SDave May     }
663d8588912SDave May     PetscViewerASCIIPopTab(viewer);    /* pop0 */
664d8588912SDave May   }
665d8588912SDave May   PetscFunctionReturn(0);
666d8588912SDave May }
667d8588912SDave May 
668d8588912SDave May #undef __FUNCT__
669d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest"
670207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A)
671d8588912SDave May {
672d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
673d8588912SDave May   PetscInt       i,j;
674d8588912SDave May   PetscErrorCode ierr;
675d8588912SDave May 
676d8588912SDave May   PetscFunctionBegin;
677d8588912SDave May   for (i=0; i<bA->nr; i++) {
678d8588912SDave May     for (j=0; j<bA->nc; j++) {
679d8588912SDave May       if (!bA->m[i][j]) continue;
680d8588912SDave May       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
681d8588912SDave May     }
682d8588912SDave May   }
683d8588912SDave May   PetscFunctionReturn(0);
684d8588912SDave May }
685d8588912SDave May 
686d8588912SDave May #undef __FUNCT__
687c222c20dSDavid Ham #define __FUNCT__ "MatCopy_Nest"
688c222c20dSDavid Ham static PetscErrorCode MatCopy_Nest(Mat A,Mat B,MatStructure str)
689c222c20dSDavid Ham {
690c222c20dSDavid Ham   Mat_Nest       *bA = (Mat_Nest*)A->data,*bB = (Mat_Nest*)B->data;
691c222c20dSDavid Ham   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
692c222c20dSDavid Ham   PetscErrorCode ierr;
693c222c20dSDavid Ham 
694c222c20dSDavid Ham   PetscFunctionBegin;
695c222c20dSDavid 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);
696c222c20dSDavid Ham   for (i=0; i<nr; i++) {
697c222c20dSDavid Ham     for (j=0; j<nc; j++) {
69846a2b97cSJed Brown       if (bA->m[i][j] && bB->m[i][j]) {
699c222c20dSDavid Ham         ierr = MatCopy(bA->m[i][j],bB->m[i][j],str);CHKERRQ(ierr);
70046a2b97cSJed 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);
701c222c20dSDavid Ham     }
702c222c20dSDavid Ham   }
703c222c20dSDavid Ham   PetscFunctionReturn(0);
704c222c20dSDavid Ham }
705c222c20dSDavid Ham 
706c222c20dSDavid Ham #undef __FUNCT__
707d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest"
708207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
709d8588912SDave May {
710d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
711841e96a3SJed Brown   Mat            *b;
712841e96a3SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
713d8588912SDave May   PetscErrorCode ierr;
714d8588912SDave May 
715d8588912SDave May   PetscFunctionBegin;
716785e854fSJed Brown   ierr = PetscMalloc1(nr*nc,&b);CHKERRQ(ierr);
717841e96a3SJed Brown   for (i=0; i<nr; i++) {
718841e96a3SJed Brown     for (j=0; j<nc; j++) {
719841e96a3SJed Brown       if (bA->m[i][j]) {
720841e96a3SJed Brown         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
721841e96a3SJed Brown       } else {
7220298fd71SBarry Smith         b[i*nc+j] = NULL;
723d8588912SDave May       }
724d8588912SDave May     }
725d8588912SDave May   }
726ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
727841e96a3SJed Brown   /* Give the new MatNest exclusive ownership */
728841e96a3SJed Brown   for (i=0; i<nr*nc; i++) {
7296bf464f9SBarry Smith     ierr = MatDestroy(&b[i]);CHKERRQ(ierr);
730d8588912SDave May   }
731d8588912SDave May   ierr = PetscFree(b);CHKERRQ(ierr);
732d8588912SDave May 
733841e96a3SJed Brown   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
734841e96a3SJed Brown   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
735d8588912SDave May   PetscFunctionReturn(0);
736d8588912SDave May }
737d8588912SDave May 
738d8588912SDave May /* nest api */
739d8588912SDave May #undef __FUNCT__
740d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest"
741d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
742d8588912SDave May {
743d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
7445fd66863SKarl Rupp 
745d8588912SDave May   PetscFunctionBegin;
746ce94432eSBarry Smith   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
747ce94432eSBarry Smith   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
748d8588912SDave May   *mat = bA->m[idxm][jdxm];
749d8588912SDave May   PetscFunctionReturn(0);
750d8588912SDave May }
751d8588912SDave May 
752d8588912SDave May #undef __FUNCT__
753d8588912SDave May #define __FUNCT__ "MatNestGetSubMat"
7549ba0d327SJed Brown /*@
755d8588912SDave May  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
756d8588912SDave May 
757d8588912SDave May  Not collective
758d8588912SDave May 
759d8588912SDave May  Input Parameters:
760629881c0SJed Brown +   A  - nest matrix
761d8588912SDave May .   idxm - index of the matrix within the nest matrix
762629881c0SJed Brown -   jdxm - index of the matrix within the nest matrix
763d8588912SDave May 
764d8588912SDave May  Output Parameter:
765d8588912SDave May .   sub - matrix at index idxm,jdxm within the nest matrix
766d8588912SDave May 
767d8588912SDave May  Level: developer
768d8588912SDave May 
769d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMats()
770d8588912SDave May @*/
7717087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
772d8588912SDave May {
773699a902aSJed Brown   PetscErrorCode ierr;
774d8588912SDave May 
775d8588912SDave May   PetscFunctionBegin;
776699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
777d8588912SDave May   PetscFunctionReturn(0);
778d8588912SDave May }
779d8588912SDave May 
780d8588912SDave May #undef __FUNCT__
7810782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat_Nest"
7820782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat)
7830782ca92SJed Brown {
7840782ca92SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
7850782ca92SJed Brown   PetscInt       m,n,M,N,mi,ni,Mi,Ni;
7860782ca92SJed Brown   PetscErrorCode ierr;
7870782ca92SJed Brown 
7880782ca92SJed Brown   PetscFunctionBegin;
789ce94432eSBarry Smith   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
790ce94432eSBarry Smith   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
7910782ca92SJed Brown   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
7920782ca92SJed Brown   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
7930782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr);
7940782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr);
7950782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr);
7960782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr);
797ce94432eSBarry 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);
798ce94432eSBarry 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);
79926fbe8dcSKarl Rupp 
8000782ca92SJed Brown   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
8010782ca92SJed Brown   ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr);
8020782ca92SJed Brown   bA->m[idxm][jdxm] = mat;
8030782ca92SJed Brown   PetscFunctionReturn(0);
8040782ca92SJed Brown }
8050782ca92SJed Brown 
8060782ca92SJed Brown #undef __FUNCT__
8070782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat"
8089ba0d327SJed Brown /*@
8090782ca92SJed Brown  MatNestSetSubMat - Set a single submatrix in the nest matrix.
8100782ca92SJed Brown 
8110782ca92SJed Brown  Logically collective on the submatrix communicator
8120782ca92SJed Brown 
8130782ca92SJed Brown  Input Parameters:
8140782ca92SJed Brown +   A  - nest matrix
8150782ca92SJed Brown .   idxm - index of the matrix within the nest matrix
8160782ca92SJed Brown .   jdxm - index of the matrix within the nest matrix
8170782ca92SJed Brown -   sub - matrix at index idxm,jdxm within the nest matrix
8180782ca92SJed Brown 
8190782ca92SJed Brown  Notes:
8200782ca92SJed Brown  The new submatrix must have the same size and communicator as that block of the nest.
8210782ca92SJed Brown 
8220782ca92SJed Brown  This increments the reference count of the submatrix.
8230782ca92SJed Brown 
8240782ca92SJed Brown  Level: developer
8250782ca92SJed Brown 
8260782ca92SJed Brown .seealso: MatNestSetSubMats(), MatNestGetSubMat()
8270782ca92SJed Brown @*/
8280782ca92SJed Brown PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
8290782ca92SJed Brown {
8300782ca92SJed Brown   PetscErrorCode ierr;
8310782ca92SJed Brown 
8320782ca92SJed Brown   PetscFunctionBegin;
8330782ca92SJed Brown   ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr);
8340782ca92SJed Brown   PetscFunctionReturn(0);
8350782ca92SJed Brown }
8360782ca92SJed Brown 
8370782ca92SJed Brown #undef __FUNCT__
838d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest"
839d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
840d8588912SDave May {
841d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
8425fd66863SKarl Rupp 
843d8588912SDave May   PetscFunctionBegin;
84426fbe8dcSKarl Rupp   if (M)   *M   = bA->nr;
84526fbe8dcSKarl Rupp   if (N)   *N   = bA->nc;
84626fbe8dcSKarl Rupp   if (mat) *mat = bA->m;
847d8588912SDave May   PetscFunctionReturn(0);
848d8588912SDave May }
849d8588912SDave May 
850d8588912SDave May #undef __FUNCT__
851d8588912SDave May #define __FUNCT__ "MatNestGetSubMats"
852d8588912SDave May /*@C
853d8588912SDave May  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
854d8588912SDave May 
855d8588912SDave May  Not collective
856d8588912SDave May 
857d8588912SDave May  Input Parameters:
858629881c0SJed Brown .   A  - nest matrix
859d8588912SDave May 
860d8588912SDave May  Output Parameter:
861629881c0SJed Brown +   M - number of rows in the nest matrix
862d8588912SDave May .   N - number of cols in the nest matrix
863629881c0SJed Brown -   mat - 2d array of matrices
864d8588912SDave May 
865d8588912SDave May  Notes:
866d8588912SDave May 
867d8588912SDave May  The user should not free the array mat.
868d8588912SDave May 
869351962e3SVincent Le Chenadec  In Fortran, this routine has a calling sequence
870351962e3SVincent Le Chenadec $   call MatNestGetSubMats(A, M, N, mat, ierr)
871351962e3SVincent Le Chenadec  where the space allocated for the optional argument mat is assumed large enough (if provided).
872351962e3SVincent Le Chenadec 
873d8588912SDave May  Level: developer
874d8588912SDave May 
875d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMat()
876d8588912SDave May @*/
8777087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
878d8588912SDave May {
879699a902aSJed Brown   PetscErrorCode ierr;
880d8588912SDave May 
881d8588912SDave May   PetscFunctionBegin;
882699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
883d8588912SDave May   PetscFunctionReturn(0);
884d8588912SDave May }
885d8588912SDave May 
886d8588912SDave May #undef __FUNCT__
887d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest"
8887087cfbeSBarry Smith PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
889d8588912SDave May {
890d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
891d8588912SDave May 
892d8588912SDave May   PetscFunctionBegin;
89326fbe8dcSKarl Rupp   if (M) *M = bA->nr;
89426fbe8dcSKarl Rupp   if (N) *N = bA->nc;
895d8588912SDave May   PetscFunctionReturn(0);
896d8588912SDave May }
897d8588912SDave May 
898d8588912SDave May #undef __FUNCT__
899d8588912SDave May #define __FUNCT__ "MatNestGetSize"
9009ba0d327SJed Brown /*@
901d8588912SDave May  MatNestGetSize - Returns the size of the nest matrix.
902d8588912SDave May 
903d8588912SDave May  Not collective
904d8588912SDave May 
905d8588912SDave May  Input Parameters:
906d8588912SDave May .   A  - nest matrix
907d8588912SDave May 
908d8588912SDave May  Output Parameter:
909629881c0SJed Brown +   M - number of rows in the nested mat
910629881c0SJed Brown -   N - number of cols in the nested mat
911d8588912SDave May 
912d8588912SDave May  Notes:
913d8588912SDave May 
914d8588912SDave May  Level: developer
915d8588912SDave May 
916d8588912SDave May .seealso: MatNestGetSubMat(), MatNestGetSubMats()
917d8588912SDave May @*/
9187087cfbeSBarry Smith PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
919d8588912SDave May {
920699a902aSJed Brown   PetscErrorCode ierr;
921d8588912SDave May 
922d8588912SDave May   PetscFunctionBegin;
923699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
924d8588912SDave May   PetscFunctionReturn(0);
925d8588912SDave May }
926d8588912SDave May 
927900e7ff2SJed Brown #undef __FUNCT__
928900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs_Nest"
929f7a08781SBarry Smith static PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[])
930900e7ff2SJed Brown {
931900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest*)A->data;
932900e7ff2SJed Brown   PetscInt i;
933900e7ff2SJed Brown 
934900e7ff2SJed Brown   PetscFunctionBegin;
935900e7ff2SJed Brown   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i];
936900e7ff2SJed Brown   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i];
937900e7ff2SJed Brown   PetscFunctionReturn(0);
938900e7ff2SJed Brown }
939900e7ff2SJed Brown 
940900e7ff2SJed Brown #undef __FUNCT__
941900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs"
9423a4d7b9aSSatish Balay /*@C
943900e7ff2SJed Brown  MatNestGetISs - Returns the index sets partitioning the row and column spaces
944900e7ff2SJed Brown 
945900e7ff2SJed Brown  Not collective
946900e7ff2SJed Brown 
947900e7ff2SJed Brown  Input Parameters:
948900e7ff2SJed Brown .   A  - nest matrix
949900e7ff2SJed Brown 
950900e7ff2SJed Brown  Output Parameter:
951900e7ff2SJed Brown +   rows - array of row index sets
952900e7ff2SJed Brown -   cols - array of column index sets
953900e7ff2SJed Brown 
954900e7ff2SJed Brown  Level: advanced
955900e7ff2SJed Brown 
956900e7ff2SJed Brown  Notes:
957900e7ff2SJed Brown  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
958900e7ff2SJed Brown 
959900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs()
960900e7ff2SJed Brown @*/
961900e7ff2SJed Brown PetscErrorCode  MatNestGetISs(Mat A,IS rows[],IS cols[])
962900e7ff2SJed Brown {
963900e7ff2SJed Brown   PetscErrorCode ierr;
964900e7ff2SJed Brown 
965900e7ff2SJed Brown   PetscFunctionBegin;
966900e7ff2SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
967900e7ff2SJed Brown   ierr = PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
968900e7ff2SJed Brown   PetscFunctionReturn(0);
969900e7ff2SJed Brown }
970900e7ff2SJed Brown 
971900e7ff2SJed Brown #undef __FUNCT__
972900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs_Nest"
973f7a08781SBarry Smith static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[])
974900e7ff2SJed Brown {
975900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest*)A->data;
976900e7ff2SJed Brown   PetscInt i;
977900e7ff2SJed Brown 
978900e7ff2SJed Brown   PetscFunctionBegin;
979900e7ff2SJed Brown   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i];
980900e7ff2SJed Brown   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i];
981900e7ff2SJed Brown   PetscFunctionReturn(0);
982900e7ff2SJed Brown }
983900e7ff2SJed Brown 
984900e7ff2SJed Brown #undef __FUNCT__
985900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs"
986900e7ff2SJed Brown /*@C
987900e7ff2SJed Brown  MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces
988900e7ff2SJed Brown 
989900e7ff2SJed Brown  Not collective
990900e7ff2SJed Brown 
991900e7ff2SJed Brown  Input Parameters:
992900e7ff2SJed Brown .   A  - nest matrix
993900e7ff2SJed Brown 
994900e7ff2SJed Brown  Output Parameter:
9950298fd71SBarry Smith +   rows - array of row index sets (or NULL to ignore)
9960298fd71SBarry Smith -   cols - array of column index sets (or NULL to ignore)
997900e7ff2SJed Brown 
998900e7ff2SJed Brown  Level: advanced
999900e7ff2SJed Brown 
1000900e7ff2SJed Brown  Notes:
1001900e7ff2SJed Brown  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
1002900e7ff2SJed Brown 
1003900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs()
1004900e7ff2SJed Brown @*/
1005900e7ff2SJed Brown PetscErrorCode  MatNestGetLocalISs(Mat A,IS rows[],IS cols[])
1006900e7ff2SJed Brown {
1007900e7ff2SJed Brown   PetscErrorCode ierr;
1008900e7ff2SJed Brown 
1009900e7ff2SJed Brown   PetscFunctionBegin;
1010900e7ff2SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1011900e7ff2SJed Brown   ierr = PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
1012900e7ff2SJed Brown   PetscFunctionReturn(0);
1013900e7ff2SJed Brown }
1014900e7ff2SJed Brown 
1015207556f9SJed Brown #undef __FUNCT__
1016207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest"
101719fd82e9SBarry Smith PetscErrorCode  MatNestSetVecType_Nest(Mat A,VecType vtype)
1018207556f9SJed Brown {
1019207556f9SJed Brown   PetscErrorCode ierr;
1020207556f9SJed Brown   PetscBool      flg;
1021207556f9SJed Brown 
1022207556f9SJed Brown   PetscFunctionBegin;
1023207556f9SJed Brown   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
1024207556f9SJed Brown   /* In reality, this only distinguishes VECNEST and "other" */
10252a7a6963SBarry Smith   if (flg) A->ops->getvecs = MatCreateVecs_Nest;
102612b53f24SSatish Balay   else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0;
1027207556f9SJed Brown   PetscFunctionReturn(0);
1028207556f9SJed Brown }
1029207556f9SJed Brown 
1030207556f9SJed Brown #undef __FUNCT__
1031207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType"
1032207556f9SJed Brown /*@C
10332a7a6963SBarry Smith  MatNestSetVecType - Sets the type of Vec returned by MatCreateVecs()
1034207556f9SJed Brown 
1035207556f9SJed Brown  Not collective
1036207556f9SJed Brown 
1037207556f9SJed Brown  Input Parameters:
1038207556f9SJed Brown +  A  - nest matrix
1039207556f9SJed Brown -  vtype - type to use for creating vectors
1040207556f9SJed Brown 
1041207556f9SJed Brown  Notes:
1042207556f9SJed Brown 
1043207556f9SJed Brown  Level: developer
1044207556f9SJed Brown 
10452a7a6963SBarry Smith .seealso: MatCreateVecs()
1046207556f9SJed Brown @*/
104719fd82e9SBarry Smith PetscErrorCode  MatNestSetVecType(Mat A,VecType vtype)
1048207556f9SJed Brown {
1049207556f9SJed Brown   PetscErrorCode ierr;
1050207556f9SJed Brown 
1051207556f9SJed Brown   PetscFunctionBegin;
105219fd82e9SBarry Smith   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));CHKERRQ(ierr);
1053207556f9SJed Brown   PetscFunctionReturn(0);
1054207556f9SJed Brown }
1055207556f9SJed Brown 
1056d8588912SDave May #undef __FUNCT__
1057c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats_Nest"
1058c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1059d8588912SDave May {
1060c8883902SJed Brown   Mat_Nest       *s = (Mat_Nest*)A->data;
1061c8883902SJed Brown   PetscInt       i,j,m,n,M,N;
1062d8588912SDave May   PetscErrorCode ierr;
1063d8588912SDave May 
1064d8588912SDave May   PetscFunctionBegin;
1065c8883902SJed Brown   s->nr = nr;
1066c8883902SJed Brown   s->nc = nc;
1067d8588912SDave May 
1068c8883902SJed Brown   /* Create space for submatrices */
1069854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&s->m);CHKERRQ(ierr);
1070c8883902SJed Brown   for (i=0; i<nr; i++) {
1071854ce69bSBarry Smith     ierr = PetscMalloc1(nc,&s->m[i]);CHKERRQ(ierr);
1072d8588912SDave May   }
1073c8883902SJed Brown   for (i=0; i<nr; i++) {
1074c8883902SJed Brown     for (j=0; j<nc; j++) {
1075c8883902SJed Brown       s->m[i][j] = a[i*nc+j];
1076c8883902SJed Brown       if (a[i*nc+j]) {
1077c8883902SJed Brown         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
1078d8588912SDave May       }
1079d8588912SDave May     }
1080d8588912SDave May   }
1081d8588912SDave May 
10828188e55aSJed Brown   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
1083d8588912SDave May 
1084854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&s->row_len);CHKERRQ(ierr);
1085854ce69bSBarry Smith   ierr = PetscMalloc1(nc,&s->col_len);CHKERRQ(ierr);
1086c8883902SJed Brown   for (i=0; i<nr; i++) s->row_len[i]=-1;
1087c8883902SJed Brown   for (j=0; j<nc; j++) s->col_len[j]=-1;
1088d8588912SDave May 
10898188e55aSJed Brown   ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr);
1090d8588912SDave May 
1091c8883902SJed Brown   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
1092c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
1093c8883902SJed Brown   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
1094c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
1095c8883902SJed Brown 
1096c8883902SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1097c8883902SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1098c8883902SJed Brown 
10991795a4d1SJed Brown   ierr = PetscCalloc2(nr,&s->left,nc,&s->right);CHKERRQ(ierr);
1100d8588912SDave May   PetscFunctionReturn(0);
1101d8588912SDave May }
1102d8588912SDave May 
1103c8883902SJed Brown #undef __FUNCT__
1104c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats"
1105c8883902SJed Brown /*@
1106c8883902SJed Brown    MatNestSetSubMats - Sets the nested submatrices
1107c8883902SJed Brown 
1108c8883902SJed Brown    Collective on Mat
1109c8883902SJed Brown 
1110c8883902SJed Brown    Input Parameter:
1111c8883902SJed Brown +  N - nested matrix
1112c8883902SJed Brown .  nr - number of nested row blocks
11130298fd71SBarry Smith .  is_row - index sets for each nested row block, or NULL to make contiguous
1114c8883902SJed Brown .  nc - number of nested column blocks
11150298fd71SBarry Smith .  is_col - index sets for each nested column block, or NULL to make contiguous
11160298fd71SBarry Smith -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1117c8883902SJed Brown 
1118c8883902SJed Brown    Level: advanced
1119c8883902SJed Brown 
1120c8883902SJed Brown .seealso: MatCreateNest(), MATNEST
1121c8883902SJed Brown @*/
1122c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1123c8883902SJed Brown {
1124c8883902SJed Brown   PetscErrorCode ierr;
1125c8883902SJed Brown   PetscInt       i;
1126c8883902SJed Brown 
1127c8883902SJed Brown   PetscFunctionBegin;
1128c8883902SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1129ce94432eSBarry Smith   if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
1130c8883902SJed Brown   if (nr && is_row) {
1131c8883902SJed Brown     PetscValidPointer(is_row,3);
1132c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
1133c8883902SJed Brown   }
1134ce94432eSBarry Smith   if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
11351664e352SJed Brown   if (nc && is_col) {
1136c8883902SJed Brown     PetscValidPointer(is_col,5);
11379b30a8f6SBarry Smith     for (i=0; i<nc; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
1138c8883902SJed Brown   }
1139c8883902SJed Brown   if (nr*nc) PetscValidPointer(a,6);
1140c8883902SJed 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);
1141c8883902SJed Brown   PetscFunctionReturn(0);
1142c8883902SJed Brown }
1143d8588912SDave May 
114477019fcaSJed Brown #undef __FUNCT__
114577019fcaSJed Brown #define __FUNCT__ "MatNestCreateAggregateL2G_Private"
114645b6f7e9SBarry Smith static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog)
114777019fcaSJed Brown {
114877019fcaSJed Brown   PetscErrorCode ierr;
114977019fcaSJed Brown   PetscBool      flg;
115077019fcaSJed Brown   PetscInt       i,j,m,mi,*ix;
115177019fcaSJed Brown 
115277019fcaSJed Brown   PetscFunctionBegin;
115377019fcaSJed Brown   for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) {
115477019fcaSJed Brown     if (islocal[i]) {
115577019fcaSJed Brown       ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
115677019fcaSJed Brown       flg  = PETSC_TRUE;      /* We found a non-trivial entry */
115777019fcaSJed Brown     } else {
115877019fcaSJed Brown       ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
115977019fcaSJed Brown     }
116077019fcaSJed Brown     m += mi;
116177019fcaSJed Brown   }
116277019fcaSJed Brown   if (flg) {
1163785e854fSJed Brown     ierr = PetscMalloc1(m,&ix);CHKERRQ(ierr);
116477019fcaSJed Brown     for (i=0,n=0; i<n; i++) {
11650298fd71SBarry Smith       ISLocalToGlobalMapping smap = NULL;
116677019fcaSJed Brown       VecScatter             scat;
116777019fcaSJed Brown       IS                     isreq;
116877019fcaSJed Brown       Vec                    lvec,gvec;
11693361c9a7SJed Brown       union {char padding[sizeof(PetscScalar)]; PetscInt integer;} *x;
117077019fcaSJed Brown       Mat sub;
117177019fcaSJed Brown 
1172ce94432eSBarry Smith       if (sizeof(*x) != sizeof(PetscScalar)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support when scalars smaller than integers");
117377019fcaSJed Brown       if (colflg) {
117477019fcaSJed Brown         ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
117577019fcaSJed Brown       } else {
117677019fcaSJed Brown         ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
117777019fcaSJed Brown       }
11780298fd71SBarry Smith       if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&smap,NULL);CHKERRQ(ierr);}
117977019fcaSJed Brown       if (islocal[i]) {
118077019fcaSJed Brown         ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
118177019fcaSJed Brown       } else {
118277019fcaSJed Brown         ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
118377019fcaSJed Brown       }
118477019fcaSJed Brown       for (j=0; j<mi; j++) ix[m+j] = j;
118577019fcaSJed Brown       if (smap) {ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr);}
118677019fcaSJed Brown       /*
118777019fcaSJed Brown         Now we need to extract the monolithic global indices that correspond to the given split global indices.
118877019fcaSJed 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.
118977019fcaSJed Brown         The approach here is ugly because it uses VecScatter to move indices.
119077019fcaSJed Brown        */
119177019fcaSJed Brown       ierr = VecCreateSeq(PETSC_COMM_SELF,mi,&lvec);CHKERRQ(ierr);
119277019fcaSJed Brown       ierr = VecCreateMPI(((PetscObject)isglobal[i])->comm,mi,PETSC_DECIDE,&gvec);CHKERRQ(ierr);
119377019fcaSJed Brown       ierr = ISCreateGeneral(((PetscObject)isglobal[i])->comm,mi,ix+m,PETSC_COPY_VALUES,&isreq);CHKERRQ(ierr);
11940298fd71SBarry Smith       ierr = VecScatterCreate(gvec,isreq,lvec,NULL,&scat);CHKERRQ(ierr);
119577019fcaSJed Brown       ierr = VecGetArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
119677019fcaSJed Brown       for (j=0; j<mi; j++) x[j].integer = ix[m+j];
119777019fcaSJed Brown       ierr = VecRestoreArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
119877019fcaSJed Brown       ierr = VecScatterBegin(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
119977019fcaSJed Brown       ierr = VecScatterEnd(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
120077019fcaSJed Brown       ierr = VecGetArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
120177019fcaSJed Brown       for (j=0; j<mi; j++) ix[m+j] = x[j].integer;
120277019fcaSJed Brown       ierr = VecRestoreArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
120377019fcaSJed Brown       ierr = VecDestroy(&lvec);CHKERRQ(ierr);
120477019fcaSJed Brown       ierr = VecDestroy(&gvec);CHKERRQ(ierr);
120577019fcaSJed Brown       ierr = ISDestroy(&isreq);CHKERRQ(ierr);
120677019fcaSJed Brown       ierr = VecScatterDestroy(&scat);CHKERRQ(ierr);
120777019fcaSJed Brown       m   += mi;
120877019fcaSJed Brown     }
1209f0413b6fSBarry Smith     ierr   = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),1,m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
121077019fcaSJed Brown   } else {
12110298fd71SBarry Smith     *ltog  = NULL;
121277019fcaSJed Brown   }
121377019fcaSJed Brown   PetscFunctionReturn(0);
121477019fcaSJed Brown }
121577019fcaSJed Brown 
121677019fcaSJed Brown 
1217d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1218d8588912SDave May /*
1219d8588912SDave May   nprocessors = NP
1220d8588912SDave May   Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1))
1221d8588912SDave May        proc 0: => (g_0,h_0,)
1222d8588912SDave May        proc 1: => (g_1,h_1,)
1223d8588912SDave May        ...
1224d8588912SDave May        proc nprocs-1: => (g_NP-1,h_NP-1,)
1225d8588912SDave May 
1226d8588912SDave May             proc 0:                      proc 1:                    proc nprocs-1:
1227d8588912SDave May     is[0] = (0,1,2,...,nlocal(g_0)-1)  (0,1,...,nlocal(g_1)-1)  (0,1,...,nlocal(g_NP-1))
1228d8588912SDave May 
1229d8588912SDave May             proc 0:
1230d8588912SDave May     is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1)
1231d8588912SDave May             proc 1:
1232d8588912SDave May     is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1)
1233d8588912SDave May 
1234d8588912SDave May             proc NP-1:
1235d8588912SDave May     is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1)
1236d8588912SDave May */
1237d8588912SDave May #undef __FUNCT__
1238d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private"
1239841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1240d8588912SDave May {
1241e2d7f03fSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
12428188e55aSJed Brown   PetscInt       i,j,offset,n,nsum,bs;
1243d8588912SDave May   PetscErrorCode ierr;
12440298fd71SBarry Smith   Mat            sub = NULL;
1245d8588912SDave May 
1246d8588912SDave May   PetscFunctionBegin;
1247854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&vs->isglobal.row);CHKERRQ(ierr);
1248854ce69bSBarry Smith   ierr = PetscMalloc1(nc,&vs->isglobal.col);CHKERRQ(ierr);
1249d8588912SDave May   if (is_row) { /* valid IS is passed in */
1250d8588912SDave May     /* refs on is[] are incremeneted */
1251e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1252d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
125326fbe8dcSKarl Rupp 
1254e2d7f03fSJed Brown       vs->isglobal.row[i] = is_row[i];
1255d8588912SDave May     }
12562ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
12578188e55aSJed Brown     nsum = 0;
12588188e55aSJed Brown     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
12598188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1260ce94432eSBarry Smith       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
12610298fd71SBarry Smith       ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
1262ce94432eSBarry Smith       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
12638188e55aSJed Brown       nsum += n;
12648188e55aSJed Brown     }
1265ce94432eSBarry Smith     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
126630bc264bSJed Brown     offset -= nsum;
1267e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1268f349c1fdSJed Brown       ierr    = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
12690298fd71SBarry Smith       ierr    = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
12702ae74bdbSJed Brown       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1271ce94432eSBarry Smith       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
1272e2d7f03fSJed Brown       ierr    = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
12732ae74bdbSJed Brown       offset += n;
1274d8588912SDave May     }
1275d8588912SDave May   }
1276d8588912SDave May 
1277d8588912SDave May   if (is_col) { /* valid IS is passed in */
1278d8588912SDave May     /* refs on is[] are incremeneted */
1279e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1280d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
128126fbe8dcSKarl Rupp 
1282e2d7f03fSJed Brown       vs->isglobal.col[j] = is_col[j];
1283d8588912SDave May     }
12842ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
12852ae74bdbSJed Brown     offset = A->cmap->rstart;
12868188e55aSJed Brown     nsum   = 0;
12878188e55aSJed Brown     for (j=0; j<vs->nc; j++) {
12888188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
1289ce94432eSBarry Smith       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
12900298fd71SBarry Smith       ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
1291ce94432eSBarry Smith       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
12928188e55aSJed Brown       nsum += n;
12938188e55aSJed Brown     }
1294ce94432eSBarry Smith     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
129530bc264bSJed Brown     offset -= nsum;
1296e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1297f349c1fdSJed Brown       ierr    = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
12980298fd71SBarry Smith       ierr    = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
12992ae74bdbSJed Brown       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1300ce94432eSBarry Smith       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
1301e2d7f03fSJed Brown       ierr    = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
13022ae74bdbSJed Brown       offset += n;
1303d8588912SDave May     }
1304d8588912SDave May   }
1305e2d7f03fSJed Brown 
1306e2d7f03fSJed Brown   /* Set up the local ISs */
1307785e854fSJed Brown   ierr = PetscMalloc1(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
1308785e854fSJed Brown   ierr = PetscMalloc1(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
1309e2d7f03fSJed Brown   for (i=0,offset=0; i<vs->nr; i++) {
1310e2d7f03fSJed Brown     IS                     isloc;
13110298fd71SBarry Smith     ISLocalToGlobalMapping rmap = NULL;
1312e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1313e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
13140298fd71SBarry Smith     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,NULL);CHKERRQ(ierr);}
1315207556f9SJed Brown     if (rmap) {
1316e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1317e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
1318e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1319e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1320207556f9SJed Brown     } else {
1321207556f9SJed Brown       nlocal = 0;
13220298fd71SBarry Smith       isloc  = NULL;
1323207556f9SJed Brown     }
1324e2d7f03fSJed Brown     vs->islocal.row[i] = isloc;
1325e2d7f03fSJed Brown     offset            += nlocal;
1326e2d7f03fSJed Brown   }
13278188e55aSJed Brown   for (i=0,offset=0; i<vs->nc; i++) {
1328e2d7f03fSJed Brown     IS                     isloc;
13290298fd71SBarry Smith     ISLocalToGlobalMapping cmap = NULL;
1330e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1331e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
13320298fd71SBarry Smith     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,NULL,&cmap);CHKERRQ(ierr);}
1333207556f9SJed Brown     if (cmap) {
1334e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1335e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr);
1336e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1337e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1338207556f9SJed Brown     } else {
1339207556f9SJed Brown       nlocal = 0;
13400298fd71SBarry Smith       isloc  = NULL;
1341207556f9SJed Brown     }
1342e2d7f03fSJed Brown     vs->islocal.col[i] = isloc;
1343e2d7f03fSJed Brown     offset            += nlocal;
1344e2d7f03fSJed Brown   }
13450189643fSJed Brown 
134677019fcaSJed Brown   /* Set up the aggregate ISLocalToGlobalMapping */
134777019fcaSJed Brown   {
134845b6f7e9SBarry Smith     ISLocalToGlobalMapping rmap,cmap;
134945b6f7e9SBarry Smith     ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap);CHKERRQ(ierr);
135045b6f7e9SBarry Smith     ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap);CHKERRQ(ierr);
135177019fcaSJed Brown     if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);}
135277019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
135377019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
135477019fcaSJed Brown   }
135577019fcaSJed Brown 
13560189643fSJed Brown #if defined(PETSC_USE_DEBUG)
13570189643fSJed Brown   for (i=0; i<vs->nr; i++) {
13580189643fSJed Brown     for (j=0; j<vs->nc; j++) {
13590189643fSJed Brown       PetscInt m,n,M,N,mi,ni,Mi,Ni;
13600189643fSJed Brown       Mat      B = vs->m[i][j];
13610189643fSJed Brown       if (!B) continue;
13620189643fSJed Brown       ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr);
13630189643fSJed Brown       ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr);
13640189643fSJed Brown       ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr);
13650189643fSJed Brown       ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr);
13660189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr);
13670189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr);
1368ce94432eSBarry 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);
1369ce94432eSBarry 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);
13700189643fSJed Brown     }
13710189643fSJed Brown   }
13720189643fSJed Brown #endif
1373a061e289SJed Brown 
1374a061e289SJed Brown   /* Set A->assembled if all non-null blocks are currently assembled */
1375a061e289SJed Brown   for (i=0; i<vs->nr; i++) {
1376a061e289SJed Brown     for (j=0; j<vs->nc; j++) {
1377a061e289SJed Brown       if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0);
1378a061e289SJed Brown     }
1379a061e289SJed Brown   }
1380a061e289SJed Brown   A->assembled = PETSC_TRUE;
1381d8588912SDave May   PetscFunctionReturn(0);
1382d8588912SDave May }
1383d8588912SDave May 
1384d8588912SDave May #undef __FUNCT__
1385d8588912SDave May #define __FUNCT__ "MatCreateNest"
138645c38901SJed Brown /*@C
1387659c6bb0SJed Brown    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
1388659c6bb0SJed Brown 
1389659c6bb0SJed Brown    Collective on Mat
1390659c6bb0SJed Brown 
1391659c6bb0SJed Brown    Input Parameter:
1392659c6bb0SJed Brown +  comm - Communicator for the new Mat
1393659c6bb0SJed Brown .  nr - number of nested row blocks
13940298fd71SBarry Smith .  is_row - index sets for each nested row block, or NULL to make contiguous
1395659c6bb0SJed Brown .  nc - number of nested column blocks
13960298fd71SBarry Smith .  is_col - index sets for each nested column block, or NULL to make contiguous
13970298fd71SBarry Smith -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1398659c6bb0SJed Brown 
1399659c6bb0SJed Brown    Output Parameter:
1400659c6bb0SJed Brown .  B - new matrix
1401659c6bb0SJed Brown 
1402659c6bb0SJed Brown    Level: advanced
1403659c6bb0SJed Brown 
1404950540a4SJed Brown .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST
1405659c6bb0SJed Brown @*/
14067087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1407d8588912SDave May {
1408d8588912SDave May   Mat            A;
1409d8588912SDave May   PetscErrorCode ierr;
1410d8588912SDave May 
1411d8588912SDave May   PetscFunctionBegin;
1412c8883902SJed Brown   *B   = 0;
1413d8588912SDave May   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1414c8883902SJed Brown   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
141591a28eb3SBarry Smith   A->preallocated = PETSC_TRUE;
1416c8883902SJed Brown   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1417d8588912SDave May   *B   = A;
1418d8588912SDave May   PetscFunctionReturn(0);
1419d8588912SDave May }
1420659c6bb0SJed Brown 
1421629c3df2SDmitry Karpeev #undef __FUNCT__
1422629c3df2SDmitry Karpeev #define __FUNCT__ "MatConvert_Nest_AIJ"
1423cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1424629c3df2SDmitry Karpeev {
1425629c3df2SDmitry Karpeev   PetscErrorCode ierr;
1426629c3df2SDmitry Karpeev   Mat_Nest       *nest = (Mat_Nest*)A->data;
142783b1a929SMark Adams   PetscInt       m,n,M,N,i,j,k,*dnnz,*onnz,rstart;
1428649b366bSFande Kong   PetscInt       cstart,cend;
1429629c3df2SDmitry Karpeev   Mat            C;
1430629c3df2SDmitry Karpeev 
1431629c3df2SDmitry Karpeev   PetscFunctionBegin;
1432629c3df2SDmitry Karpeev   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
1433629c3df2SDmitry Karpeev   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
1434649b366bSFande Kong   ierr = MatGetOwnershipRangeColumn(A,&cstart,&cend);CHKERRQ(ierr);
1435629c3df2SDmitry Karpeev   switch (reuse) {
1436629c3df2SDmitry Karpeev   case MAT_INITIAL_MATRIX:
1437ce94432eSBarry Smith     ierr    = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
1438629c3df2SDmitry Karpeev     ierr    = MatSetType(C,newtype);CHKERRQ(ierr);
1439629c3df2SDmitry Karpeev     ierr    = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr);
1440629c3df2SDmitry Karpeev     *newmat = C;
1441629c3df2SDmitry Karpeev     break;
1442629c3df2SDmitry Karpeev   case MAT_REUSE_MATRIX:
1443629c3df2SDmitry Karpeev     C = *newmat;
1444629c3df2SDmitry Karpeev     break;
1445ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse");
1446629c3df2SDmitry Karpeev   }
1447785e854fSJed Brown   ierr = PetscMalloc1(2*m,&dnnz);CHKERRQ(ierr);
1448629c3df2SDmitry Karpeev   onnz = dnnz + m;
1449629c3df2SDmitry Karpeev   for (k=0; k<m; k++) {
1450629c3df2SDmitry Karpeev     dnnz[k] = 0;
1451629c3df2SDmitry Karpeev     onnz[k] = 0;
1452629c3df2SDmitry Karpeev   }
1453629c3df2SDmitry Karpeev   for (j=0; j<nest->nc; ++j) {
1454629c3df2SDmitry Karpeev     IS             bNis;
1455629c3df2SDmitry Karpeev     PetscInt       bN;
1456629c3df2SDmitry Karpeev     const PetscInt *bNindices;
1457629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
1458629c3df2SDmitry Karpeev     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1459629c3df2SDmitry Karpeev     ierr = ISGetSize(bNis, &bN);CHKERRQ(ierr);
1460629c3df2SDmitry Karpeev     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1461629c3df2SDmitry Karpeev     for (i=0; i<nest->nr; ++i) {
1462629c3df2SDmitry Karpeev       PetscSF        bmsf;
1463649b366bSFande Kong       PetscSFNode    *iremote;
1464629c3df2SDmitry Karpeev       Mat            B;
1465649b366bSFande Kong       PetscInt       bm, *sub_dnnz,*sub_onnz, br;
1466629c3df2SDmitry Karpeev       const PetscInt *bmindices;
1467629c3df2SDmitry Karpeev       B = nest->m[i][j];
1468629c3df2SDmitry Karpeev       if (!B) continue;
1469629c3df2SDmitry Karpeev       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1470629c3df2SDmitry Karpeev       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1471ce94432eSBarry Smith       ierr = PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);CHKERRQ(ierr);
1472649b366bSFande Kong       ierr = PetscMalloc1(bm,&iremote);CHKERRQ(ierr);
1473649b366bSFande Kong       ierr = PetscMalloc1(bm,&sub_dnnz);CHKERRQ(ierr);
1474649b366bSFande Kong       ierr = PetscMalloc1(bm,&sub_onnz);CHKERRQ(ierr);
1475649b366bSFande Kong       for (k = 0; k < bm; ++k){
1476649b366bSFande Kong     	sub_dnnz[k] = 0;
1477649b366bSFande Kong     	sub_onnz[k] = 0;
1478649b366bSFande Kong       }
1479629c3df2SDmitry Karpeev       /*
1480629c3df2SDmitry Karpeev        Locate the owners for all of the locally-owned global row indices for this row block.
1481629c3df2SDmitry Karpeev        These determine the roots of PetscSF used to communicate preallocation data to row owners.
1482629c3df2SDmitry Karpeev        The roots correspond to the dnnz and onnz entries; thus, there are two roots per row.
1483629c3df2SDmitry Karpeev        */
148483b1a929SMark Adams       ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr);
1485629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
1486649b366bSFande Kong         PetscInt       row = bmindices[br], rowowner = 0, brncols, col;
1487629c3df2SDmitry Karpeev         const PetscInt *brcols;
1488a4b3d3acSMatthew G Knepley         PetscInt       rowrel = 0; /* row's relative index on its owner rank */
1489629c3df2SDmitry Karpeev         ierr      = PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);CHKERRQ(ierr);
1490649b366bSFande Kong         /* how many roots  */
1491649b366bSFande Kong         iremote[br].rank = rowowner; iremote[br].index = rowrel;           /* edge from bmdnnz to dnnz */
1492649b366bSFande Kong         /* get nonzero pattern */
149383b1a929SMark Adams         ierr = MatGetRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr);
1494629c3df2SDmitry Karpeev         for (k=0; k<brncols; k++) {
1495629c3df2SDmitry Karpeev           col  = bNindices[brcols[k]];
1496649b366bSFande Kong           if(col>=A->cmap->range[rowowner] && col<A->cmap->range[rowowner+1]){
1497649b366bSFande Kong         	sub_dnnz[br]++;
1498649b366bSFande Kong           }else{
1499649b366bSFande Kong         	sub_onnz[br]++;
1500649b366bSFande Kong           }
1501629c3df2SDmitry Karpeev         }
150283b1a929SMark Adams         ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr);
1503629c3df2SDmitry Karpeev       }
1504629c3df2SDmitry Karpeev       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1505629c3df2SDmitry Karpeev       /* bsf will have to take care of disposing of bedges. */
1506649b366bSFande Kong       ierr = PetscSFSetGraph(bmsf,m,bm,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1507649b366bSFande Kong       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr);
1508649b366bSFande Kong       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr);
1509649b366bSFande Kong       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr);
1510649b366bSFande Kong       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr);
1511649b366bSFande Kong       ierr = PetscFree(sub_dnnz);CHKERRQ(ierr);
1512649b366bSFande Kong       ierr = PetscFree(sub_onnz);CHKERRQ(ierr);
1513629c3df2SDmitry Karpeev       ierr = PetscSFDestroy(&bmsf);CHKERRQ(ierr);
1514629c3df2SDmitry Karpeev     }
151522d28d08SBarry Smith     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
1516629c3df2SDmitry Karpeev     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
1517629c3df2SDmitry Karpeev   }
1518629c3df2SDmitry Karpeev   ierr = MatSeqAIJSetPreallocation(C,0,dnnz);CHKERRQ(ierr);
1519629c3df2SDmitry Karpeev   ierr = MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);CHKERRQ(ierr);
1520629c3df2SDmitry Karpeev   ierr = PetscFree(dnnz);CHKERRQ(ierr);
1521629c3df2SDmitry Karpeev 
1522629c3df2SDmitry Karpeev   /* Fill by row */
1523629c3df2SDmitry Karpeev   for (j=0; j<nest->nc; ++j) {
1524629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
1525629c3df2SDmitry Karpeev     IS             bNis;
1526629c3df2SDmitry Karpeev     PetscInt       bN;
1527629c3df2SDmitry Karpeev     const PetscInt *bNindices;
1528629c3df2SDmitry Karpeev     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1529629c3df2SDmitry Karpeev     ierr = ISGetSize(bNis,&bN);CHKERRQ(ierr);
1530629c3df2SDmitry Karpeev     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1531629c3df2SDmitry Karpeev     for (i=0; i<nest->nr; ++i) {
1532629c3df2SDmitry Karpeev       Mat            B;
1533629c3df2SDmitry Karpeev       PetscInt       bm, br;
1534629c3df2SDmitry Karpeev       const PetscInt *bmindices;
1535629c3df2SDmitry Karpeev       B = nest->m[i][j];
1536629c3df2SDmitry Karpeev       if (!B) continue;
1537629c3df2SDmitry Karpeev       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1538629c3df2SDmitry Karpeev       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
153983b1a929SMark Adams       ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr);
1540629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
1541629c3df2SDmitry Karpeev         PetscInt          row = bmindices[br], brncols,  *cols;
1542629c3df2SDmitry Karpeev         const PetscInt    *brcols;
1543629c3df2SDmitry Karpeev         const PetscScalar *brcoldata;
154483b1a929SMark Adams         ierr = MatGetRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
1545785e854fSJed Brown         ierr = PetscMalloc1(brncols,&cols);CHKERRQ(ierr);
154626fbe8dcSKarl Rupp         for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]];
1547629c3df2SDmitry Karpeev         /*
1548629c3df2SDmitry Karpeev           Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match.
1549629c3df2SDmitry Karpeev           Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES.
1550629c3df2SDmitry Karpeev          */
1551a2ea699eSBarry Smith         ierr = MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);CHKERRQ(ierr);
155283b1a929SMark Adams         ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
1553629c3df2SDmitry Karpeev         ierr = PetscFree(cols);CHKERRQ(ierr);
1554629c3df2SDmitry Karpeev       }
1555629c3df2SDmitry Karpeev       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1556629c3df2SDmitry Karpeev     }
1557a2ea699eSBarry Smith     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
1558629c3df2SDmitry Karpeev     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
1559629c3df2SDmitry Karpeev   }
1560629c3df2SDmitry Karpeev   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1561629c3df2SDmitry Karpeev   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1562629c3df2SDmitry Karpeev   PetscFunctionReturn(0);
1563629c3df2SDmitry Karpeev }
1564629c3df2SDmitry Karpeev 
1565659c6bb0SJed Brown /*MC
1566659c6bb0SJed Brown   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
1567659c6bb0SJed Brown 
1568659c6bb0SJed Brown   Level: intermediate
1569659c6bb0SJed Brown 
1570659c6bb0SJed Brown   Notes:
1571659c6bb0SJed Brown   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1572659c6bb0SJed Brown   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1573950540a4SJed Brown   It is usually used with DMComposite and DMCreateMatrix()
1574659c6bb0SJed Brown 
1575659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest()
1576659c6bb0SJed Brown M*/
1577c8883902SJed Brown #undef __FUNCT__
1578c8883902SJed Brown #define __FUNCT__ "MatCreate_Nest"
15798cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A)
1580c8883902SJed Brown {
1581c8883902SJed Brown   Mat_Nest       *s;
1582c8883902SJed Brown   PetscErrorCode ierr;
1583c8883902SJed Brown 
1584c8883902SJed Brown   PetscFunctionBegin;
1585b00a9115SJed Brown   ierr    = PetscNewLog(A,&s);CHKERRQ(ierr);
1586c8883902SJed Brown   A->data = (void*)s;
1587e7c19651SJed Brown 
1588e7c19651SJed Brown   s->nr            = -1;
1589e7c19651SJed Brown   s->nc            = -1;
15900298fd71SBarry Smith   s->m             = NULL;
1591e7c19651SJed Brown   s->splitassembly = PETSC_FALSE;
1592c8883902SJed Brown 
1593c8883902SJed Brown   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
159426fbe8dcSKarl Rupp 
1595c8883902SJed Brown   A->ops->mult                  = MatMult_Nest;
15969194d70fSJed Brown   A->ops->multadd               = MatMultAdd_Nest;
1597c8883902SJed Brown   A->ops->multtranspose         = MatMultTranspose_Nest;
15989194d70fSJed Brown   A->ops->multtransposeadd      = MatMultTransposeAdd_Nest;
1599c8883902SJed Brown   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1600c8883902SJed Brown   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1601c8883902SJed Brown   A->ops->zeroentries           = MatZeroEntries_Nest;
1602c222c20dSDavid Ham   A->ops->copy                  = MatCopy_Nest;
1603c8883902SJed Brown   A->ops->duplicate             = MatDuplicate_Nest;
1604c8883902SJed Brown   A->ops->getsubmatrix          = MatGetSubMatrix_Nest;
1605c8883902SJed Brown   A->ops->destroy               = MatDestroy_Nest;
1606c8883902SJed Brown   A->ops->view                  = MatView_Nest;
1607c8883902SJed Brown   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1608c8883902SJed Brown   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1609c8883902SJed Brown   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
1610429bac76SJed Brown   A->ops->getdiagonal           = MatGetDiagonal_Nest;
1611429bac76SJed Brown   A->ops->diagonalscale         = MatDiagonalScale_Nest;
1612a061e289SJed Brown   A->ops->scale                 = MatScale_Nest;
1613a061e289SJed Brown   A->ops->shift                 = MatShift_Nest;
1614*13135bc6SAlex Fikl   A->ops->diagonalset           = MatDiagonalSet_Nest;
1615c8883902SJed Brown 
1616c8883902SJed Brown   A->spptr        = 0;
1617c8883902SJed Brown   A->assembled    = PETSC_FALSE;
1618c8883902SJed Brown 
1619c8883902SJed Brown   /* expose Nest api's */
1620bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",   MatNestGetSubMat_Nest);CHKERRQ(ierr);
1621bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",   MatNestSetSubMat_Nest);CHKERRQ(ierr);
1622bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",  MatNestGetSubMats_Nest);CHKERRQ(ierr);
1623bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",     MatNestGetSize_Nest);CHKERRQ(ierr);
1624bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",      MatNestGetISs_Nest);CHKERRQ(ierr);
1625bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C", MatNestGetLocalISs_Nest);CHKERRQ(ierr);
1626bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",  MatNestSetVecType_Nest);CHKERRQ(ierr);
1627bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",  MatNestSetSubMats_Nest);CHKERRQ(ierr);
162883b1a929SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",MatConvert_Nest_AIJ);CHKERRQ(ierr);
1629c8883902SJed Brown 
1630c8883902SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
1631c8883902SJed Brown   PetscFunctionReturn(0);
1632c8883902SJed Brown }
1633