xref: /petsc/src/mat/impls/nest/matnest.c (revision d5dfb6943b42ca1cf8a1a2498d9ef78f9cad9cdc)
1d8588912SDave May 
2aaa7dc30SBarry Smith #include <../src/mat/impls/nest/matnestimpl.h> /*I   "petscmat.h"   I*/
3b68353e5Sstefano_zampini #include <../src/mat/impls/aij/seq/aij.h>
40c312b8eSJed Brown #include <petscsf.h>
5d8588912SDave May 
6c8883902SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]);
72a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left);
85e3038f0Sstefano_zampini PETSC_INTERN PetscErrorCode MatConvert_Nest_IS(Mat,MatType,MatReuse,Mat*);
9c8883902SJed Brown 
10d8588912SDave May /* private functions */
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 */
37207556f9SJed Brown static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y)
38d8588912SDave May {
39d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
40207556f9SJed Brown   Vec            *bx = bA->right,*by = bA->left;
41207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
42d8588912SDave May   PetscErrorCode ierr;
43d8588912SDave May 
44d8588912SDave May   PetscFunctionBegin;
45207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
46207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
47207556f9SJed Brown   for (i=0; i<nr; i++) {
48d8588912SDave May     ierr = VecZeroEntries(by[i]);CHKERRQ(ierr);
49207556f9SJed Brown     for (j=0; j<nc; j++) {
50207556f9SJed Brown       if (!bA->m[i][j]) continue;
51d8588912SDave May       /* y[i] <- y[i] + A[i][j] * x[j] */
52d8588912SDave May       ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr);
53d8588912SDave May     }
54d8588912SDave May   }
55207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
56207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
57d8588912SDave May   PetscFunctionReturn(0);
58d8588912SDave May }
59d8588912SDave May 
609194d70fSJed Brown static PetscErrorCode MatMultAdd_Nest(Mat A,Vec x,Vec y,Vec z)
619194d70fSJed Brown {
629194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
639194d70fSJed Brown   Vec            *bx = bA->right,*bz = bA->left;
649194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
659194d70fSJed Brown   PetscErrorCode ierr;
669194d70fSJed Brown 
679194d70fSJed Brown   PetscFunctionBegin;
689194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
699194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
709194d70fSJed Brown   for (i=0; i<nr; i++) {
719194d70fSJed Brown     if (y != z) {
729194d70fSJed Brown       Vec by;
739194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
749194d70fSJed Brown       ierr = VecCopy(by,bz[i]);CHKERRQ(ierr);
75336d21e7SJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
769194d70fSJed Brown     }
779194d70fSJed Brown     for (j=0; j<nc; j++) {
789194d70fSJed Brown       if (!bA->m[i][j]) continue;
799194d70fSJed Brown       /* y[i] <- y[i] + A[i][j] * x[j] */
809194d70fSJed Brown       ierr = MatMultAdd(bA->m[i][j],bx[j],bz[i],bz[i]);CHKERRQ(ierr);
819194d70fSJed Brown     }
829194d70fSJed Brown   }
839194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
849194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
859194d70fSJed Brown   PetscFunctionReturn(0);
869194d70fSJed Brown }
879194d70fSJed Brown 
88207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y)
89d8588912SDave May {
90d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
91207556f9SJed Brown   Vec            *bx = bA->left,*by = bA->right;
92207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
93d8588912SDave May   PetscErrorCode ierr;
94d8588912SDave May 
95d8588912SDave May   PetscFunctionBegin;
96609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
97609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
98207556f9SJed Brown   for (j=0; j<nc; j++) {
99609e31cbSJed Brown     ierr = VecZeroEntries(by[j]);CHKERRQ(ierr);
100609e31cbSJed Brown     for (i=0; i<nr; i++) {
1016c75ac25SJed Brown       if (!bA->m[i][j]) continue;
102609e31cbSJed Brown       /* y[j] <- y[j] + (A[i][j])^T * x[i] */
103609e31cbSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],by[j],by[j]);CHKERRQ(ierr);
104d8588912SDave May     }
105d8588912SDave May   }
106609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
107609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
108d8588912SDave May   PetscFunctionReturn(0);
109d8588912SDave May }
110d8588912SDave May 
1119194d70fSJed Brown static PetscErrorCode MatMultTransposeAdd_Nest(Mat A,Vec x,Vec y,Vec z)
1129194d70fSJed Brown {
1139194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
1149194d70fSJed Brown   Vec            *bx = bA->left,*bz = bA->right;
1159194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
1169194d70fSJed Brown   PetscErrorCode ierr;
1179194d70fSJed Brown 
1189194d70fSJed Brown   PetscFunctionBegin;
1199194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
1209194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
1219194d70fSJed Brown   for (j=0; j<nc; j++) {
1229194d70fSJed Brown     if (y != z) {
1239194d70fSJed Brown       Vec by;
1249194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
1259194d70fSJed Brown       ierr = VecCopy(by,bz[j]);CHKERRQ(ierr);
1269194d70fSJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
1279194d70fSJed Brown     }
1289194d70fSJed Brown     for (i=0; i<nr; i++) {
1296c75ac25SJed Brown       if (!bA->m[i][j]) continue;
1309194d70fSJed Brown       /* z[j] <- y[j] + (A[i][j])^T * x[i] */
1319194d70fSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],bz[j],bz[j]);CHKERRQ(ierr);
1329194d70fSJed Brown     }
1339194d70fSJed Brown   }
1349194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
1359194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
1369194d70fSJed Brown   PetscFunctionReturn(0);
1379194d70fSJed Brown }
1389194d70fSJed Brown 
139f8170845SAlex Fikl static PetscErrorCode MatTranspose_Nest(Mat A,MatReuse reuse,Mat *B)
140f8170845SAlex Fikl {
141f8170845SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data, *bC;
142f8170845SAlex Fikl   Mat            C;
143f8170845SAlex Fikl   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
144f8170845SAlex Fikl   PetscErrorCode ierr;
145f8170845SAlex Fikl 
146f8170845SAlex Fikl   PetscFunctionBegin;
147f8170845SAlex Fikl   if (reuse == MAT_REUSE_MATRIX && A == *B && nr != nc) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Square nested matrix only for in-place");
148f8170845SAlex Fikl 
149f8170845SAlex Fikl   if (reuse == MAT_INITIAL_MATRIX || A == *B) {
150f8170845SAlex Fikl     Mat *subs;
151f8170845SAlex Fikl     IS  *is_row,*is_col;
152f8170845SAlex Fikl 
153f8170845SAlex Fikl     ierr = PetscCalloc1(nr * nc,&subs);CHKERRQ(ierr);
154f8170845SAlex Fikl     ierr = PetscMalloc2(nr,&is_row,nc,&is_col);CHKERRQ(ierr);
155f8170845SAlex Fikl     ierr = MatNestGetISs(A,is_row,is_col);CHKERRQ(ierr);
156ddeb9bd8SAlex Fikl     if (reuse == MAT_REUSE_MATRIX) {
157ddeb9bd8SAlex Fikl       for (i=0; i<nr; i++) {
158ddeb9bd8SAlex Fikl         for (j=0; j<nc; j++) {
159ddeb9bd8SAlex Fikl           subs[i + nr * j] = bA->m[i][j];
160ddeb9bd8SAlex Fikl         }
161ddeb9bd8SAlex Fikl       }
162ddeb9bd8SAlex Fikl     }
163ddeb9bd8SAlex Fikl 
164f8170845SAlex Fikl     ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nc,is_col,nr,is_row,subs,&C);CHKERRQ(ierr);
165f8170845SAlex Fikl     ierr = PetscFree(subs);CHKERRQ(ierr);
1663d994f23SBarry Smith     ierr = PetscFree2(is_row,is_col);CHKERRQ(ierr);
167f8170845SAlex Fikl   } else {
168f8170845SAlex Fikl     C = *B;
169f8170845SAlex Fikl   }
170f8170845SAlex Fikl 
171f8170845SAlex Fikl   bC = (Mat_Nest*)C->data;
172f8170845SAlex Fikl   for (i=0; i<nr; i++) {
173f8170845SAlex Fikl     for (j=0; j<nc; j++) {
174f8170845SAlex Fikl       if (bA->m[i][j]) {
175f8170845SAlex Fikl         ierr = MatTranspose(bA->m[i][j], reuse, &(bC->m[j][i]));CHKERRQ(ierr);
176f8170845SAlex Fikl       } else {
177f8170845SAlex Fikl         bC->m[j][i] = NULL;
178f8170845SAlex Fikl       }
179f8170845SAlex Fikl     }
180f8170845SAlex Fikl   }
181f8170845SAlex Fikl 
182f8170845SAlex Fikl   if (reuse == MAT_INITIAL_MATRIX || A != *B) {
183f8170845SAlex Fikl     *B = C;
184f8170845SAlex Fikl   } else {
185f8170845SAlex Fikl     ierr = MatHeaderMerge(A, &C);CHKERRQ(ierr);
186f8170845SAlex Fikl   }
187f8170845SAlex Fikl   PetscFunctionReturn(0);
188f8170845SAlex Fikl }
189f8170845SAlex Fikl 
190e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list)
191e2d7f03fSJed Brown {
192e2d7f03fSJed Brown   PetscErrorCode ierr;
193e2d7f03fSJed Brown   IS             *lst = *list;
194e2d7f03fSJed Brown   PetscInt       i;
195e2d7f03fSJed Brown 
196e2d7f03fSJed Brown   PetscFunctionBegin;
197e2d7f03fSJed Brown   if (!lst) PetscFunctionReturn(0);
1986bf464f9SBarry Smith   for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);}
199e2d7f03fSJed Brown   ierr  = PetscFree(lst);CHKERRQ(ierr);
2000298fd71SBarry Smith   *list = NULL;
201e2d7f03fSJed Brown   PetscFunctionReturn(0);
202e2d7f03fSJed Brown }
203e2d7f03fSJed Brown 
204207556f9SJed Brown static PetscErrorCode MatDestroy_Nest(Mat A)
205d8588912SDave May {
206d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
207d8588912SDave May   PetscInt       i,j;
208d8588912SDave May   PetscErrorCode ierr;
209d8588912SDave May 
210d8588912SDave May   PetscFunctionBegin;
211d8588912SDave May   /* release the matrices and the place holders */
212e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr);
213e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr);
214e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
215e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
216d8588912SDave May 
217d8588912SDave May   ierr = PetscFree(vs->row_len);CHKERRQ(ierr);
218d8588912SDave May   ierr = PetscFree(vs->col_len);CHKERRQ(ierr);
219d8588912SDave May 
220207556f9SJed Brown   ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr);
221207556f9SJed Brown 
222d8588912SDave May   /* release the matrices and the place holders */
223d8588912SDave May   if (vs->m) {
224d8588912SDave May     for (i=0; i<vs->nr; i++) {
225d8588912SDave May       for (j=0; j<vs->nc; j++) {
2266bf464f9SBarry Smith         ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr);
227d8588912SDave May       }
228d8588912SDave May       ierr = PetscFree(vs->m[i]);CHKERRQ(ierr);
229d8588912SDave May     }
230d8588912SDave May     ierr = PetscFree(vs->m);CHKERRQ(ierr);
231d8588912SDave May   }
232bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
233d8588912SDave May 
234bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",0);CHKERRQ(ierr);
235bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",0);CHKERRQ(ierr);
236bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",0);CHKERRQ(ierr);
237bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",0);CHKERRQ(ierr);
238bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",0);CHKERRQ(ierr);
239bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",0);CHKERRQ(ierr);
240bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",0);CHKERRQ(ierr);
241bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",0);CHKERRQ(ierr);
2425e3038f0Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",0);CHKERRQ(ierr);
2435e3038f0Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_is_C",0);CHKERRQ(ierr);
244d8588912SDave May   PetscFunctionReturn(0);
245d8588912SDave May }
246d8588912SDave May 
247207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
248d8588912SDave May {
249d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
250d8588912SDave May   PetscInt       i,j;
251d8588912SDave May   PetscErrorCode ierr;
252d8588912SDave May 
253d8588912SDave May   PetscFunctionBegin;
254d8588912SDave May   for (i=0; i<vs->nr; i++) {
255d8588912SDave May     for (j=0; j<vs->nc; j++) {
256e7c19651SJed Brown       if (vs->m[i][j]) {
257e7c19651SJed Brown         ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr);
258e7c19651SJed Brown         if (!vs->splitassembly) {
259e7c19651SJed Brown           /* Note: split assembly will fail if the same block appears more than once (even indirectly through a nested
260e7c19651SJed 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
261e7c19651SJed Brown            * already performing an assembly, but the result would by more complicated and appears to offer less
262e7c19651SJed Brown            * potential for diagnostics and correctness checking. Split assembly should be fixed once there is an
263e7c19651SJed Brown            * interface for libraries to make asynchronous progress in "user-defined non-blocking collectives".
264e7c19651SJed Brown            */
265e7c19651SJed Brown           ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr);
266e7c19651SJed Brown         }
267e7c19651SJed Brown       }
268d8588912SDave May     }
269d8588912SDave May   }
270d8588912SDave May   PetscFunctionReturn(0);
271d8588912SDave May }
272d8588912SDave May 
273207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
274d8588912SDave May {
275d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
276d8588912SDave May   PetscInt       i,j;
277d8588912SDave May   PetscErrorCode ierr;
278d8588912SDave May 
279d8588912SDave May   PetscFunctionBegin;
280d8588912SDave May   for (i=0; i<vs->nr; i++) {
281d8588912SDave May     for (j=0; j<vs->nc; j++) {
282e7c19651SJed Brown       if (vs->m[i][j]) {
283e7c19651SJed Brown         if (vs->splitassembly) {
284e7c19651SJed Brown           ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr);
285e7c19651SJed Brown         }
286e7c19651SJed Brown       }
287d8588912SDave May     }
288d8588912SDave May   }
289d8588912SDave May   PetscFunctionReturn(0);
290d8588912SDave May }
291d8588912SDave May 
292f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
293d8588912SDave May {
294207556f9SJed Brown   PetscErrorCode ierr;
295f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
296f349c1fdSJed Brown   PetscInt       j;
297f349c1fdSJed Brown   Mat            sub;
298d8588912SDave May 
299d8588912SDave May   PetscFunctionBegin;
3000298fd71SBarry Smith   sub = (row < vs->nc) ? vs->m[row][row] : (Mat)NULL; /* Prefer to find on the diagonal */
301f349c1fdSJed Brown   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
3024994cf47SJed Brown   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
303f349c1fdSJed Brown   *B = sub;
304f349c1fdSJed Brown   PetscFunctionReturn(0);
305d8588912SDave May }
306d8588912SDave May 
307f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
308f349c1fdSJed Brown {
309207556f9SJed Brown   PetscErrorCode ierr;
310f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
311f349c1fdSJed Brown   PetscInt       i;
312f349c1fdSJed Brown   Mat            sub;
313f349c1fdSJed Brown 
314f349c1fdSJed Brown   PetscFunctionBegin;
3150298fd71SBarry Smith   sub = (col < vs->nr) ? vs->m[col][col] : (Mat)NULL; /* Prefer to find on the diagonal */
316f349c1fdSJed Brown   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
3174994cf47SJed Brown   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
318f349c1fdSJed Brown   *B = sub;
319f349c1fdSJed Brown   PetscFunctionReturn(0);
320d8588912SDave May }
321d8588912SDave May 
322f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
323f349c1fdSJed Brown {
324f349c1fdSJed Brown   PetscErrorCode ierr;
325f349c1fdSJed Brown   PetscInt       i;
326f349c1fdSJed Brown   PetscBool      flg;
327f349c1fdSJed Brown 
328f349c1fdSJed Brown   PetscFunctionBegin;
329f349c1fdSJed Brown   PetscValidPointer(list,3);
330f349c1fdSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,4);
331f349c1fdSJed Brown   PetscValidIntPointer(found,5);
332f349c1fdSJed Brown   *found = -1;
333f349c1fdSJed Brown   for (i=0; i<n; i++) {
334207556f9SJed Brown     if (!list[i]) continue;
335f349c1fdSJed Brown     ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr);
336f349c1fdSJed Brown     if (flg) {
337f349c1fdSJed Brown       *found = i;
338f349c1fdSJed Brown       PetscFunctionReturn(0);
339f349c1fdSJed Brown     }
340f349c1fdSJed Brown   }
341ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Could not find index set");
342f349c1fdSJed Brown   PetscFunctionReturn(0);
343f349c1fdSJed Brown }
344f349c1fdSJed Brown 
3458188e55aSJed Brown /* Get a block row as a new MatNest */
3468188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B)
3478188e55aSJed Brown {
3488188e55aSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
3498188e55aSJed Brown   char           keyname[256];
3508188e55aSJed Brown   PetscErrorCode ierr;
3518188e55aSJed Brown 
3528188e55aSJed Brown   PetscFunctionBegin;
3530298fd71SBarry Smith   *B   = NULL;
3548caf3d72SBarry Smith   ierr = PetscSNPrintf(keyname,sizeof(keyname),"NestRow_%D",row);CHKERRQ(ierr);
3558188e55aSJed Brown   ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr);
3568188e55aSJed Brown   if (*B) PetscFunctionReturn(0);
3578188e55aSJed Brown 
358ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),1,NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr);
35926fbe8dcSKarl Rupp 
3608188e55aSJed Brown   (*B)->assembled = A->assembled;
36126fbe8dcSKarl Rupp 
3628188e55aSJed Brown   ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr);
3638188e55aSJed Brown   ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */
3648188e55aSJed Brown   PetscFunctionReturn(0);
3658188e55aSJed Brown }
3668188e55aSJed Brown 
367f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
368f349c1fdSJed Brown {
369f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
3708188e55aSJed Brown   PetscErrorCode ierr;
3716b3a5b13SJed Brown   PetscInt       row,col;
372e072481dSJed Brown   PetscBool      same,isFullCol,isFullColGlobal;
373f349c1fdSJed Brown 
374f349c1fdSJed Brown   PetscFunctionBegin;
3758188e55aSJed Brown   /* Check if full column space. This is a hack */
3768188e55aSJed Brown   isFullCol = PETSC_FALSE;
377251f4c67SDmitry Karpeev   ierr      = PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr);
3788188e55aSJed Brown   if (same) {
37977019fcaSJed Brown     PetscInt n,first,step,i,an,am,afirst,astep;
3808188e55aSJed Brown     ierr      = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
3818188e55aSJed Brown     ierr      = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
38277019fcaSJed Brown     isFullCol = PETSC_TRUE;
38305ce4453SJed Brown     for (i=0,an=A->cmap->rstart; i<vs->nc; i++) {
38477019fcaSJed Brown       ierr = ISStrideGetInfo(is->col[i],&afirst,&astep);CHKERRQ(ierr);
38577019fcaSJed Brown       ierr = ISGetLocalSize(is->col[i],&am);CHKERRQ(ierr);
38677019fcaSJed Brown       if (afirst != an || astep != step) isFullCol = PETSC_FALSE;
38777019fcaSJed Brown       an += am;
38877019fcaSJed Brown     }
38905ce4453SJed Brown     if (an != A->cmap->rstart+n) isFullCol = PETSC_FALSE;
3908188e55aSJed Brown   }
391b2566f29SBarry Smith   ierr = MPIU_Allreduce(&isFullCol,&isFullColGlobal,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)iscol));CHKERRQ(ierr);
3928188e55aSJed Brown 
393e072481dSJed Brown   if (isFullColGlobal) {
3948188e55aSJed Brown     PetscInt row;
3958188e55aSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
3968188e55aSJed Brown     ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr);
3978188e55aSJed Brown   } else {
398f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
399f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr);
400b6480e04SStefano Zampini     if (!vs->m[row][col]) {
401b6480e04SStefano Zampini       PetscInt lr,lc;
402b6480e04SStefano Zampini 
403b6480e04SStefano Zampini       ierr = MatCreate(PetscObjectComm((PetscObject)A),&vs->m[row][col]);CHKERRQ(ierr);
404b6480e04SStefano Zampini       ierr = ISGetLocalSize(vs->isglobal.row[row],&lr);CHKERRQ(ierr);
405b6480e04SStefano Zampini       ierr = ISGetLocalSize(vs->isglobal.col[col],&lc);CHKERRQ(ierr);
406b6480e04SStefano Zampini       ierr = MatSetSizes(vs->m[row][col],lr,lc,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
407b6480e04SStefano Zampini       ierr = MatSetUp(vs->m[row][col]);CHKERRQ(ierr);
408b6480e04SStefano Zampini       ierr = MatAssemblyBegin(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
409b6480e04SStefano Zampini       ierr = MatAssemblyEnd(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
410b6480e04SStefano Zampini     }
411f349c1fdSJed Brown     *B = vs->m[row][col];
4128188e55aSJed Brown   }
413f349c1fdSJed Brown   PetscFunctionReturn(0);
414f349c1fdSJed Brown }
415f349c1fdSJed Brown 
416f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
417f349c1fdSJed Brown {
418f349c1fdSJed Brown   PetscErrorCode ierr;
419f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
420f349c1fdSJed Brown   Mat            sub;
421f349c1fdSJed Brown 
422f349c1fdSJed Brown   PetscFunctionBegin;
423f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr);
424f349c1fdSJed Brown   switch (reuse) {
425f349c1fdSJed Brown   case MAT_INITIAL_MATRIX:
4267874fa86SDave May     if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); }
427f349c1fdSJed Brown     *B = sub;
428f349c1fdSJed Brown     break;
429f349c1fdSJed Brown   case MAT_REUSE_MATRIX:
430ce94432eSBarry Smith     if (sub != *B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
431f349c1fdSJed Brown     break;
432f349c1fdSJed Brown   case MAT_IGNORE_MATRIX:       /* Nothing to do */
433f349c1fdSJed Brown     break;
434511c6705SHong Zhang   case MAT_INPLACE_MATRIX:       /* Nothing to do */
435511c6705SHong Zhang     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX is not supported yet");
436511c6705SHong Zhang     break;
437f349c1fdSJed Brown   }
438f349c1fdSJed Brown   PetscFunctionReturn(0);
439f349c1fdSJed Brown }
440f349c1fdSJed Brown 
441f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
442f349c1fdSJed Brown {
443f349c1fdSJed Brown   PetscErrorCode ierr;
444f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
445f349c1fdSJed Brown   Mat            sub;
446f349c1fdSJed Brown 
447f349c1fdSJed Brown   PetscFunctionBegin;
448f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
449f349c1fdSJed Brown   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
450f349c1fdSJed Brown   if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);}
451f349c1fdSJed Brown   *B = sub;
452d8588912SDave May   PetscFunctionReturn(0);
453d8588912SDave May }
454d8588912SDave May 
455207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
456d8588912SDave May {
457d8588912SDave May   PetscErrorCode ierr;
458f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
459f349c1fdSJed Brown   Mat            sub;
460d8588912SDave May 
461d8588912SDave May   PetscFunctionBegin;
462f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
463ce94432eSBarry Smith   if (*B != sub) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
464f349c1fdSJed Brown   if (sub) {
465ce94432eSBarry Smith     if (((PetscObject)sub)->refct <= 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
4666bf464f9SBarry Smith     ierr = MatDestroy(B);CHKERRQ(ierr);
467d8588912SDave May   }
468d8588912SDave May   PetscFunctionReturn(0);
469d8588912SDave May }
470d8588912SDave May 
4717874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
4727874fa86SDave May {
4737874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
4747874fa86SDave May   PetscInt       i;
4757874fa86SDave May   PetscErrorCode ierr;
4767874fa86SDave May 
4777874fa86SDave May   PetscFunctionBegin;
4787874fa86SDave May   for (i=0; i<bA->nr; i++) {
479429bac76SJed Brown     Vec bv;
480429bac76SJed Brown     ierr = VecGetSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
4817874fa86SDave May     if (bA->m[i][i]) {
482429bac76SJed Brown       ierr = MatGetDiagonal(bA->m[i][i],bv);CHKERRQ(ierr);
4837874fa86SDave May     } else {
4845159a857SMatthew G. Knepley       ierr = VecSet(bv,0.0);CHKERRQ(ierr);
4857874fa86SDave May     }
486429bac76SJed Brown     ierr = VecRestoreSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
4877874fa86SDave May   }
4887874fa86SDave May   PetscFunctionReturn(0);
4897874fa86SDave May }
4907874fa86SDave May 
4917874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
4927874fa86SDave May {
4937874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
494429bac76SJed Brown   Vec            bl,*br;
4957874fa86SDave May   PetscInt       i,j;
4967874fa86SDave May   PetscErrorCode ierr;
4977874fa86SDave May 
4987874fa86SDave May   PetscFunctionBegin;
4993f800ebeSJed Brown   ierr = PetscCalloc1(bA->nc,&br);CHKERRQ(ierr);
5002e6472ebSElliott Sales de Andrade   if (r) {
501429bac76SJed Brown     for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
5022e6472ebSElliott Sales de Andrade   }
5032e6472ebSElliott Sales de Andrade   bl = NULL;
5047874fa86SDave May   for (i=0; i<bA->nr; i++) {
5052e6472ebSElliott Sales de Andrade     if (l) {
506429bac76SJed Brown       ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
5072e6472ebSElliott Sales de Andrade     }
5087874fa86SDave May     for (j=0; j<bA->nc; j++) {
5097874fa86SDave May       if (bA->m[i][j]) {
510429bac76SJed Brown         ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr);
5117874fa86SDave May       }
5127874fa86SDave May     }
5132e6472ebSElliott Sales de Andrade     if (l) {
514a061e289SJed Brown       ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
5157874fa86SDave May     }
5162e6472ebSElliott Sales de Andrade   }
5172e6472ebSElliott Sales de Andrade   if (r) {
518429bac76SJed Brown     for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
5192e6472ebSElliott Sales de Andrade   }
520429bac76SJed Brown   ierr = PetscFree(br);CHKERRQ(ierr);
5217874fa86SDave May   PetscFunctionReturn(0);
5227874fa86SDave May }
5237874fa86SDave May 
524a061e289SJed Brown static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a)
525a061e289SJed Brown {
526a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
527a061e289SJed Brown   PetscInt       i,j;
528a061e289SJed Brown   PetscErrorCode ierr;
529a061e289SJed Brown 
530a061e289SJed Brown   PetscFunctionBegin;
531a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
532a061e289SJed Brown     for (j=0; j<bA->nc; j++) {
533a061e289SJed Brown       if (bA->m[i][j]) {
534a061e289SJed Brown         ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr);
535a061e289SJed Brown       }
536a061e289SJed Brown     }
537a061e289SJed Brown   }
538a061e289SJed Brown   PetscFunctionReturn(0);
539a061e289SJed Brown }
540a061e289SJed Brown 
541a061e289SJed Brown static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a)
542a061e289SJed Brown {
543a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
544a061e289SJed Brown   PetscInt       i;
545a061e289SJed Brown   PetscErrorCode ierr;
546a061e289SJed Brown 
547a061e289SJed Brown   PetscFunctionBegin;
548a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
549ce94432eSBarry 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);
550a061e289SJed Brown     ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr);
551a061e289SJed Brown   }
552a061e289SJed Brown   PetscFunctionReturn(0);
553a061e289SJed Brown }
554a061e289SJed Brown 
55513135bc6SAlex Fikl static PetscErrorCode MatDiagonalSet_Nest(Mat A,Vec D,InsertMode is)
55613135bc6SAlex Fikl {
55713135bc6SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data;
55813135bc6SAlex Fikl   PetscInt       i;
55913135bc6SAlex Fikl   PetscErrorCode ierr;
56013135bc6SAlex Fikl 
56113135bc6SAlex Fikl   PetscFunctionBegin;
56213135bc6SAlex Fikl   for (i=0; i<bA->nr; i++) {
56313135bc6SAlex Fikl     Vec bv;
56413135bc6SAlex Fikl     ierr = VecGetSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
56513135bc6SAlex Fikl     if (bA->m[i][i]) {
56613135bc6SAlex Fikl       ierr = MatDiagonalSet(bA->m[i][i],bv,is);CHKERRQ(ierr);
56713135bc6SAlex Fikl     }
56813135bc6SAlex Fikl     ierr = VecRestoreSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
56913135bc6SAlex Fikl   }
57013135bc6SAlex Fikl   PetscFunctionReturn(0);
57113135bc6SAlex Fikl }
57213135bc6SAlex Fikl 
573f8170845SAlex Fikl static PetscErrorCode MatSetRandom_Nest(Mat A,PetscRandom rctx)
574f8170845SAlex Fikl {
575f8170845SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data;
576f8170845SAlex Fikl   PetscInt       i,j;
577f8170845SAlex Fikl   PetscErrorCode ierr;
578f8170845SAlex Fikl 
579f8170845SAlex Fikl   PetscFunctionBegin;
580f8170845SAlex Fikl   for (i=0; i<bA->nr; i++) {
581f8170845SAlex Fikl     for (j=0; j<bA->nc; j++) {
582f8170845SAlex Fikl       if (bA->m[i][j]) {
583f8170845SAlex Fikl         ierr = MatSetRandom(bA->m[i][j],rctx);CHKERRQ(ierr);
584f8170845SAlex Fikl       }
585f8170845SAlex Fikl     }
586f8170845SAlex Fikl   }
587f8170845SAlex Fikl   PetscFunctionReturn(0);
588f8170845SAlex Fikl }
589f8170845SAlex Fikl 
5902a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left)
591d8588912SDave May {
592d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
593d8588912SDave May   Vec            *L,*R;
594d8588912SDave May   MPI_Comm       comm;
595d8588912SDave May   PetscInt       i,j;
596d8588912SDave May   PetscErrorCode ierr;
597d8588912SDave May 
598d8588912SDave May   PetscFunctionBegin;
599ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
600d8588912SDave May   if (right) {
601d8588912SDave May     /* allocate R */
602854ce69bSBarry Smith     ierr = PetscMalloc1(bA->nc, &R);CHKERRQ(ierr);
603d8588912SDave May     /* Create the right vectors */
604d8588912SDave May     for (j=0; j<bA->nc; j++) {
605d8588912SDave May       for (i=0; i<bA->nr; i++) {
606d8588912SDave May         if (bA->m[i][j]) {
6072a7a6963SBarry Smith           ierr = MatCreateVecs(bA->m[i][j],&R[j],NULL);CHKERRQ(ierr);
608d8588912SDave May           break;
609d8588912SDave May         }
610d8588912SDave May       }
6116c4ed002SBarry Smith       if (i==bA->nr) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
612d8588912SDave May     }
613f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
614d8588912SDave May     /* hand back control to the nest vector */
615d8588912SDave May     for (j=0; j<bA->nc; j++) {
6166bf464f9SBarry Smith       ierr = VecDestroy(&R[j]);CHKERRQ(ierr);
617d8588912SDave May     }
618d8588912SDave May     ierr = PetscFree(R);CHKERRQ(ierr);
619d8588912SDave May   }
620d8588912SDave May 
621d8588912SDave May   if (left) {
622d8588912SDave May     /* allocate L */
623854ce69bSBarry Smith     ierr = PetscMalloc1(bA->nr, &L);CHKERRQ(ierr);
624d8588912SDave May     /* Create the left vectors */
625d8588912SDave May     for (i=0; i<bA->nr; i++) {
626d8588912SDave May       for (j=0; j<bA->nc; j++) {
627d8588912SDave May         if (bA->m[i][j]) {
6282a7a6963SBarry Smith           ierr = MatCreateVecs(bA->m[i][j],NULL,&L[i]);CHKERRQ(ierr);
629d8588912SDave May           break;
630d8588912SDave May         }
631d8588912SDave May       }
6326c4ed002SBarry Smith       if (j==bA->nc) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
633d8588912SDave May     }
634d8588912SDave May 
635f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
636d8588912SDave May     for (i=0; i<bA->nr; i++) {
6376bf464f9SBarry Smith       ierr = VecDestroy(&L[i]);CHKERRQ(ierr);
638d8588912SDave May     }
639d8588912SDave May 
640d8588912SDave May     ierr = PetscFree(L);CHKERRQ(ierr);
641d8588912SDave May   }
642d8588912SDave May   PetscFunctionReturn(0);
643d8588912SDave May }
644d8588912SDave May 
645207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
646d8588912SDave May {
647d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
648d8588912SDave May   PetscBool      isascii;
649d8588912SDave May   PetscInt       i,j;
650d8588912SDave May   PetscErrorCode ierr;
651d8588912SDave May 
652d8588912SDave May   PetscFunctionBegin;
653251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
654d8588912SDave May   if (isascii) {
655d8588912SDave May 
656d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Matrix object: \n");CHKERRQ(ierr);
657d86155a6SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
658d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer, "type=nest, rows=%D, cols=%D \n",bA->nr,bA->nc);CHKERRQ(ierr);
659d8588912SDave May 
660d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"MatNest structure: \n");CHKERRQ(ierr);
661d8588912SDave May     for (i=0; i<bA->nr; i++) {
662d8588912SDave May       for (j=0; j<bA->nc; j++) {
66319fd82e9SBarry Smith         MatType   type;
664270f95d7SJed Brown         char      name[256] = "",prefix[256] = "";
665d8588912SDave May         PetscInt  NR,NC;
666d8588912SDave May         PetscBool isNest = PETSC_FALSE;
667d8588912SDave May 
668d8588912SDave May         if (!bA->m[i][j]) {
669d86155a6SBarry Smith           CHKERRQ(ierr);PetscViewerASCIIPrintf(viewer, "(%D,%D) : NULL \n",i,j);CHKERRQ(ierr);
670d8588912SDave May           continue;
671d8588912SDave May         }
672d8588912SDave May         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
673d8588912SDave May         ierr = MatGetType(bA->m[i][j], &type);CHKERRQ(ierr);
6748caf3d72SBarry Smith         if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);}
6758caf3d72SBarry Smith         if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);}
676251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
677d8588912SDave May 
678270f95d7SJed Brown         ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr);
679d8588912SDave May 
680d8588912SDave May         if (isNest) {
681270f95d7SJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);  /* push1 */
682d8588912SDave May           ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr);
683270f95d7SJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop1 */
684d8588912SDave May         }
685d8588912SDave May       }
686d8588912SDave May     }
687d86155a6SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop0 */
688d8588912SDave May   }
689d8588912SDave May   PetscFunctionReturn(0);
690d8588912SDave May }
691d8588912SDave May 
692207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A)
693d8588912SDave May {
694d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
695d8588912SDave May   PetscInt       i,j;
696d8588912SDave May   PetscErrorCode ierr;
697d8588912SDave May 
698d8588912SDave May   PetscFunctionBegin;
699d8588912SDave May   for (i=0; i<bA->nr; i++) {
700d8588912SDave May     for (j=0; j<bA->nc; j++) {
701d8588912SDave May       if (!bA->m[i][j]) continue;
702d8588912SDave May       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
703d8588912SDave May     }
704d8588912SDave May   }
705d8588912SDave May   PetscFunctionReturn(0);
706d8588912SDave May }
707d8588912SDave May 
708c222c20dSDavid Ham static PetscErrorCode MatCopy_Nest(Mat A,Mat B,MatStructure str)
709c222c20dSDavid Ham {
710c222c20dSDavid Ham   Mat_Nest       *bA = (Mat_Nest*)A->data,*bB = (Mat_Nest*)B->data;
711c222c20dSDavid Ham   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
712c222c20dSDavid Ham   PetscErrorCode ierr;
713c222c20dSDavid Ham 
714c222c20dSDavid Ham   PetscFunctionBegin;
715c222c20dSDavid 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);
716c222c20dSDavid Ham   for (i=0; i<nr; i++) {
717c222c20dSDavid Ham     for (j=0; j<nc; j++) {
71846a2b97cSJed Brown       if (bA->m[i][j] && bB->m[i][j]) {
719c222c20dSDavid Ham         ierr = MatCopy(bA->m[i][j],bB->m[i][j],str);CHKERRQ(ierr);
72046a2b97cSJed 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);
721c222c20dSDavid Ham     }
722c222c20dSDavid Ham   }
723c222c20dSDavid Ham   PetscFunctionReturn(0);
724c222c20dSDavid Ham }
725c222c20dSDavid Ham 
726207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
727d8588912SDave May {
728d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
729841e96a3SJed Brown   Mat            *b;
730841e96a3SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
731d8588912SDave May   PetscErrorCode ierr;
732d8588912SDave May 
733d8588912SDave May   PetscFunctionBegin;
734785e854fSJed Brown   ierr = PetscMalloc1(nr*nc,&b);CHKERRQ(ierr);
735841e96a3SJed Brown   for (i=0; i<nr; i++) {
736841e96a3SJed Brown     for (j=0; j<nc; j++) {
737841e96a3SJed Brown       if (bA->m[i][j]) {
738841e96a3SJed Brown         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
739841e96a3SJed Brown       } else {
7400298fd71SBarry Smith         b[i*nc+j] = NULL;
741d8588912SDave May       }
742d8588912SDave May     }
743d8588912SDave May   }
744ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
745841e96a3SJed Brown   /* Give the new MatNest exclusive ownership */
746841e96a3SJed Brown   for (i=0; i<nr*nc; i++) {
7476bf464f9SBarry Smith     ierr = MatDestroy(&b[i]);CHKERRQ(ierr);
748d8588912SDave May   }
749d8588912SDave May   ierr = PetscFree(b);CHKERRQ(ierr);
750d8588912SDave May 
751841e96a3SJed Brown   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
752841e96a3SJed Brown   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
753d8588912SDave May   PetscFunctionReturn(0);
754d8588912SDave May }
755d8588912SDave May 
756d8588912SDave May /* nest api */
757d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
758d8588912SDave May {
759d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
7605fd66863SKarl Rupp 
761d8588912SDave May   PetscFunctionBegin;
762ce94432eSBarry Smith   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
763ce94432eSBarry Smith   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
764d8588912SDave May   *mat = bA->m[idxm][jdxm];
765d8588912SDave May   PetscFunctionReturn(0);
766d8588912SDave May }
767d8588912SDave May 
7689ba0d327SJed Brown /*@
769d8588912SDave May  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
770d8588912SDave May 
771d8588912SDave May  Not collective
772d8588912SDave May 
773d8588912SDave May  Input Parameters:
774629881c0SJed Brown +   A  - nest matrix
775d8588912SDave May .   idxm - index of the matrix within the nest matrix
776629881c0SJed Brown -   jdxm - index of the matrix within the nest matrix
777d8588912SDave May 
778d8588912SDave May  Output Parameter:
779d8588912SDave May .   sub - matrix at index idxm,jdxm within the nest matrix
780d8588912SDave May 
781d8588912SDave May  Level: developer
782d8588912SDave May 
783d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMats()
784d8588912SDave May @*/
7857087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
786d8588912SDave May {
787699a902aSJed Brown   PetscErrorCode ierr;
788d8588912SDave May 
789d8588912SDave May   PetscFunctionBegin;
790699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
791d8588912SDave May   PetscFunctionReturn(0);
792d8588912SDave May }
793d8588912SDave May 
7940782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat)
7950782ca92SJed Brown {
7960782ca92SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
7970782ca92SJed Brown   PetscInt       m,n,M,N,mi,ni,Mi,Ni;
7980782ca92SJed Brown   PetscErrorCode ierr;
7990782ca92SJed Brown 
8000782ca92SJed Brown   PetscFunctionBegin;
801ce94432eSBarry Smith   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
802ce94432eSBarry Smith   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
8030782ca92SJed Brown   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
8040782ca92SJed Brown   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
8050782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr);
8060782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr);
8070782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr);
8080782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr);
809ce94432eSBarry 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);
810ce94432eSBarry 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);
81126fbe8dcSKarl Rupp 
8120782ca92SJed Brown   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
8130782ca92SJed Brown   ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr);
8140782ca92SJed Brown   bA->m[idxm][jdxm] = mat;
8150782ca92SJed Brown   PetscFunctionReturn(0);
8160782ca92SJed Brown }
8170782ca92SJed Brown 
8189ba0d327SJed Brown /*@
8190782ca92SJed Brown  MatNestSetSubMat - Set a single submatrix in the nest matrix.
8200782ca92SJed Brown 
8210782ca92SJed Brown  Logically collective on the submatrix communicator
8220782ca92SJed Brown 
8230782ca92SJed Brown  Input Parameters:
8240782ca92SJed Brown +   A  - nest matrix
8250782ca92SJed Brown .   idxm - index of the matrix within the nest matrix
8260782ca92SJed Brown .   jdxm - index of the matrix within the nest matrix
8270782ca92SJed Brown -   sub - matrix at index idxm,jdxm within the nest matrix
8280782ca92SJed Brown 
8290782ca92SJed Brown  Notes:
8300782ca92SJed Brown  The new submatrix must have the same size and communicator as that block of the nest.
8310782ca92SJed Brown 
8320782ca92SJed Brown  This increments the reference count of the submatrix.
8330782ca92SJed Brown 
8340782ca92SJed Brown  Level: developer
8350782ca92SJed Brown 
836*d5dfb694SBarry Smith .seealso: MatNestSetSubMats(), MatNestGetSubMats()
8370782ca92SJed Brown @*/
8380782ca92SJed Brown PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
8390782ca92SJed Brown {
8400782ca92SJed Brown   PetscErrorCode ierr;
8410782ca92SJed Brown 
8420782ca92SJed Brown   PetscFunctionBegin;
8430782ca92SJed Brown   ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr);
8440782ca92SJed Brown   PetscFunctionReturn(0);
8450782ca92SJed Brown }
8460782ca92SJed Brown 
847d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
848d8588912SDave May {
849d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
8505fd66863SKarl Rupp 
851d8588912SDave May   PetscFunctionBegin;
85226fbe8dcSKarl Rupp   if (M)   *M   = bA->nr;
85326fbe8dcSKarl Rupp   if (N)   *N   = bA->nc;
85426fbe8dcSKarl Rupp   if (mat) *mat = bA->m;
855d8588912SDave May   PetscFunctionReturn(0);
856d8588912SDave May }
857d8588912SDave May 
858d8588912SDave May /*@C
859d8588912SDave May  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
860d8588912SDave May 
861d8588912SDave May  Not collective
862d8588912SDave May 
863d8588912SDave May  Input Parameters:
864629881c0SJed Brown .   A  - nest matrix
865d8588912SDave May 
866d8588912SDave May  Output Parameter:
867629881c0SJed Brown +   M - number of rows in the nest matrix
868d8588912SDave May .   N - number of cols in the nest matrix
869629881c0SJed Brown -   mat - 2d array of matrices
870d8588912SDave May 
871d8588912SDave May  Notes:
872d8588912SDave May 
873d8588912SDave May  The user should not free the array mat.
874d8588912SDave May 
875351962e3SVincent Le Chenadec  In Fortran, this routine has a calling sequence
876351962e3SVincent Le Chenadec $   call MatNestGetSubMats(A, M, N, mat, ierr)
877351962e3SVincent Le Chenadec  where the space allocated for the optional argument mat is assumed large enough (if provided).
878351962e3SVincent Le Chenadec 
879d8588912SDave May  Level: developer
880d8588912SDave May 
881d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMat()
882d8588912SDave May @*/
8837087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
884d8588912SDave May {
885699a902aSJed Brown   PetscErrorCode ierr;
886d8588912SDave May 
887d8588912SDave May   PetscFunctionBegin;
888699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
889d8588912SDave May   PetscFunctionReturn(0);
890d8588912SDave May }
891d8588912SDave May 
8927087cfbeSBarry Smith PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
893d8588912SDave May {
894d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
895d8588912SDave May 
896d8588912SDave May   PetscFunctionBegin;
89726fbe8dcSKarl Rupp   if (M) *M = bA->nr;
89826fbe8dcSKarl Rupp   if (N) *N = bA->nc;
899d8588912SDave May   PetscFunctionReturn(0);
900d8588912SDave May }
901d8588912SDave May 
9029ba0d327SJed Brown /*@
903d8588912SDave May  MatNestGetSize - Returns the size of the nest matrix.
904d8588912SDave May 
905d8588912SDave May  Not collective
906d8588912SDave May 
907d8588912SDave May  Input Parameters:
908d8588912SDave May .   A  - nest matrix
909d8588912SDave May 
910d8588912SDave May  Output Parameter:
911629881c0SJed Brown +   M - number of rows in the nested mat
912629881c0SJed Brown -   N - number of cols in the nested mat
913d8588912SDave May 
914d8588912SDave May  Notes:
915d8588912SDave May 
916d8588912SDave May  Level: developer
917d8588912SDave May 
918d8588912SDave May .seealso: MatNestGetSubMat(), MatNestGetSubMats()
919d8588912SDave May @*/
9207087cfbeSBarry Smith PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
921d8588912SDave May {
922699a902aSJed Brown   PetscErrorCode ierr;
923d8588912SDave May 
924d8588912SDave May   PetscFunctionBegin;
925699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
926d8588912SDave May   PetscFunctionReturn(0);
927d8588912SDave May }
928d8588912SDave May 
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 
9403a4d7b9aSSatish Balay /*@C
941900e7ff2SJed Brown  MatNestGetISs - Returns the index sets partitioning the row and column spaces
942900e7ff2SJed Brown 
943900e7ff2SJed Brown  Not collective
944900e7ff2SJed Brown 
945900e7ff2SJed Brown  Input Parameters:
946900e7ff2SJed Brown .   A  - nest matrix
947900e7ff2SJed Brown 
948900e7ff2SJed Brown  Output Parameter:
949900e7ff2SJed Brown +   rows - array of row index sets
950900e7ff2SJed Brown -   cols - array of column index sets
951900e7ff2SJed Brown 
952900e7ff2SJed Brown  Level: advanced
953900e7ff2SJed Brown 
954900e7ff2SJed Brown  Notes:
955900e7ff2SJed Brown  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
956900e7ff2SJed Brown 
957900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs()
958900e7ff2SJed Brown @*/
959900e7ff2SJed Brown PetscErrorCode  MatNestGetISs(Mat A,IS rows[],IS cols[])
960900e7ff2SJed Brown {
961900e7ff2SJed Brown   PetscErrorCode ierr;
962900e7ff2SJed Brown 
963900e7ff2SJed Brown   PetscFunctionBegin;
964900e7ff2SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
965900e7ff2SJed Brown   ierr = PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
966900e7ff2SJed Brown   PetscFunctionReturn(0);
967900e7ff2SJed Brown }
968900e7ff2SJed Brown 
969f7a08781SBarry Smith static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[])
970900e7ff2SJed Brown {
971900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest*)A->data;
972900e7ff2SJed Brown   PetscInt i;
973900e7ff2SJed Brown 
974900e7ff2SJed Brown   PetscFunctionBegin;
975900e7ff2SJed Brown   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i];
976900e7ff2SJed Brown   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i];
977900e7ff2SJed Brown   PetscFunctionReturn(0);
978900e7ff2SJed Brown }
979900e7ff2SJed Brown 
980900e7ff2SJed Brown /*@C
981900e7ff2SJed Brown  MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces
982900e7ff2SJed Brown 
983900e7ff2SJed Brown  Not collective
984900e7ff2SJed Brown 
985900e7ff2SJed Brown  Input Parameters:
986900e7ff2SJed Brown .   A  - nest matrix
987900e7ff2SJed Brown 
988900e7ff2SJed Brown  Output Parameter:
9890298fd71SBarry Smith +   rows - array of row index sets (or NULL to ignore)
9900298fd71SBarry Smith -   cols - array of column index sets (or NULL to ignore)
991900e7ff2SJed Brown 
992900e7ff2SJed Brown  Level: advanced
993900e7ff2SJed Brown 
994900e7ff2SJed Brown  Notes:
995900e7ff2SJed Brown  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
996900e7ff2SJed Brown 
997900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs()
998900e7ff2SJed Brown @*/
999900e7ff2SJed Brown PetscErrorCode  MatNestGetLocalISs(Mat A,IS rows[],IS cols[])
1000900e7ff2SJed Brown {
1001900e7ff2SJed Brown   PetscErrorCode ierr;
1002900e7ff2SJed Brown 
1003900e7ff2SJed Brown   PetscFunctionBegin;
1004900e7ff2SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1005900e7ff2SJed Brown   ierr = PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
1006900e7ff2SJed Brown   PetscFunctionReturn(0);
1007900e7ff2SJed Brown }
1008900e7ff2SJed Brown 
100919fd82e9SBarry Smith PetscErrorCode  MatNestSetVecType_Nest(Mat A,VecType vtype)
1010207556f9SJed Brown {
1011207556f9SJed Brown   PetscErrorCode ierr;
1012207556f9SJed Brown   PetscBool      flg;
1013207556f9SJed Brown 
1014207556f9SJed Brown   PetscFunctionBegin;
1015207556f9SJed Brown   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
1016207556f9SJed Brown   /* In reality, this only distinguishes VECNEST and "other" */
10172a7a6963SBarry Smith   if (flg) A->ops->getvecs = MatCreateVecs_Nest;
101812b53f24SSatish Balay   else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0;
1019207556f9SJed Brown   PetscFunctionReturn(0);
1020207556f9SJed Brown }
1021207556f9SJed Brown 
1022207556f9SJed Brown /*@C
10232a7a6963SBarry Smith  MatNestSetVecType - Sets the type of Vec returned by MatCreateVecs()
1024207556f9SJed Brown 
1025207556f9SJed Brown  Not collective
1026207556f9SJed Brown 
1027207556f9SJed Brown  Input Parameters:
1028207556f9SJed Brown +  A  - nest matrix
1029207556f9SJed Brown -  vtype - type to use for creating vectors
1030207556f9SJed Brown 
1031207556f9SJed Brown  Notes:
1032207556f9SJed Brown 
1033207556f9SJed Brown  Level: developer
1034207556f9SJed Brown 
10352a7a6963SBarry Smith .seealso: MatCreateVecs()
1036207556f9SJed Brown @*/
103719fd82e9SBarry Smith PetscErrorCode  MatNestSetVecType(Mat A,VecType vtype)
1038207556f9SJed Brown {
1039207556f9SJed Brown   PetscErrorCode ierr;
1040207556f9SJed Brown 
1041207556f9SJed Brown   PetscFunctionBegin;
104219fd82e9SBarry Smith   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));CHKERRQ(ierr);
1043207556f9SJed Brown   PetscFunctionReturn(0);
1044207556f9SJed Brown }
1045207556f9SJed Brown 
1046c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1047d8588912SDave May {
1048c8883902SJed Brown   Mat_Nest       *s = (Mat_Nest*)A->data;
1049c8883902SJed Brown   PetscInt       i,j,m,n,M,N;
1050d8588912SDave May   PetscErrorCode ierr;
1051d8588912SDave May 
1052d8588912SDave May   PetscFunctionBegin;
1053c8883902SJed Brown   s->nr = nr;
1054c8883902SJed Brown   s->nc = nc;
1055d8588912SDave May 
1056c8883902SJed Brown   /* Create space for submatrices */
1057854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&s->m);CHKERRQ(ierr);
1058c8883902SJed Brown   for (i=0; i<nr; i++) {
1059854ce69bSBarry Smith     ierr = PetscMalloc1(nc,&s->m[i]);CHKERRQ(ierr);
1060d8588912SDave May   }
1061c8883902SJed Brown   for (i=0; i<nr; i++) {
1062c8883902SJed Brown     for (j=0; j<nc; j++) {
1063c8883902SJed Brown       s->m[i][j] = a[i*nc+j];
1064c8883902SJed Brown       if (a[i*nc+j]) {
1065c8883902SJed Brown         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
1066d8588912SDave May       }
1067d8588912SDave May     }
1068d8588912SDave May   }
1069d8588912SDave May 
10708188e55aSJed Brown   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
1071d8588912SDave May 
1072854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&s->row_len);CHKERRQ(ierr);
1073854ce69bSBarry Smith   ierr = PetscMalloc1(nc,&s->col_len);CHKERRQ(ierr);
1074c8883902SJed Brown   for (i=0; i<nr; i++) s->row_len[i]=-1;
1075c8883902SJed Brown   for (j=0; j<nc; j++) s->col_len[j]=-1;
1076d8588912SDave May 
10778188e55aSJed Brown   ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr);
1078d8588912SDave May 
1079c8883902SJed Brown   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
1080c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
1081c8883902SJed Brown   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
1082c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
1083c8883902SJed Brown 
1084c8883902SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1085c8883902SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1086c8883902SJed Brown 
10871795a4d1SJed Brown   ierr = PetscCalloc2(nr,&s->left,nc,&s->right);CHKERRQ(ierr);
1088d8588912SDave May   PetscFunctionReturn(0);
1089d8588912SDave May }
1090d8588912SDave May 
1091c8883902SJed Brown /*@
1092c8883902SJed Brown    MatNestSetSubMats - Sets the nested submatrices
1093c8883902SJed Brown 
1094c8883902SJed Brown    Collective on Mat
1095c8883902SJed Brown 
1096c8883902SJed Brown    Input Parameter:
1097c8883902SJed Brown +  N - nested matrix
1098c8883902SJed Brown .  nr - number of nested row blocks
10990298fd71SBarry Smith .  is_row - index sets for each nested row block, or NULL to make contiguous
1100c8883902SJed Brown .  nc - number of nested column blocks
11010298fd71SBarry Smith .  is_col - index sets for each nested column block, or NULL to make contiguous
11020298fd71SBarry Smith -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1103c8883902SJed Brown 
1104c8883902SJed Brown    Level: advanced
1105c8883902SJed Brown 
1106c8883902SJed Brown .seealso: MatCreateNest(), MATNEST
1107c8883902SJed Brown @*/
1108c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1109c8883902SJed Brown {
1110c8883902SJed Brown   PetscErrorCode ierr;
1111c8883902SJed Brown   PetscInt       i;
1112c8883902SJed Brown 
1113c8883902SJed Brown   PetscFunctionBegin;
1114c8883902SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1115ce94432eSBarry Smith   if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
1116c8883902SJed Brown   if (nr && is_row) {
1117c8883902SJed Brown     PetscValidPointer(is_row,3);
1118c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
1119c8883902SJed Brown   }
1120ce94432eSBarry Smith   if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
11211664e352SJed Brown   if (nc && is_col) {
1122c8883902SJed Brown     PetscValidPointer(is_col,5);
11239b30a8f6SBarry Smith     for (i=0; i<nc; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
1124c8883902SJed Brown   }
1125c8883902SJed Brown   if (nr*nc) PetscValidPointer(a,6);
1126c8883902SJed 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);
1127c8883902SJed Brown   PetscFunctionReturn(0);
1128c8883902SJed Brown }
1129d8588912SDave May 
113045b6f7e9SBarry Smith static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog)
113177019fcaSJed Brown {
113277019fcaSJed Brown   PetscErrorCode ierr;
113377019fcaSJed Brown   PetscBool      flg;
113477019fcaSJed Brown   PetscInt       i,j,m,mi,*ix;
113577019fcaSJed Brown 
113677019fcaSJed Brown   PetscFunctionBegin;
113777019fcaSJed Brown   for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) {
113877019fcaSJed Brown     if (islocal[i]) {
113977019fcaSJed Brown       ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
114077019fcaSJed Brown       flg  = PETSC_TRUE;      /* We found a non-trivial entry */
114177019fcaSJed Brown     } else {
114277019fcaSJed Brown       ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
114377019fcaSJed Brown     }
114477019fcaSJed Brown     m += mi;
114577019fcaSJed Brown   }
114677019fcaSJed Brown   if (flg) {
1147785e854fSJed Brown     ierr = PetscMalloc1(m,&ix);CHKERRQ(ierr);
114877019fcaSJed Brown     for (i=0,n=0; i<n; i++) {
11490298fd71SBarry Smith       ISLocalToGlobalMapping smap = NULL;
115077019fcaSJed Brown       VecScatter             scat;
115177019fcaSJed Brown       IS                     isreq;
115277019fcaSJed Brown       Vec                    lvec,gvec;
11533361c9a7SJed Brown       union {char padding[sizeof(PetscScalar)]; PetscInt integer;} *x;
115477019fcaSJed Brown       Mat sub;
115577019fcaSJed Brown 
1156ce94432eSBarry Smith       if (sizeof(*x) != sizeof(PetscScalar)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support when scalars smaller than integers");
115777019fcaSJed Brown       if (colflg) {
115877019fcaSJed Brown         ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
115977019fcaSJed Brown       } else {
116077019fcaSJed Brown         ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
116177019fcaSJed Brown       }
11620298fd71SBarry Smith       if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&smap,NULL);CHKERRQ(ierr);}
116377019fcaSJed Brown       if (islocal[i]) {
116477019fcaSJed Brown         ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
116577019fcaSJed Brown       } else {
116677019fcaSJed Brown         ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
116777019fcaSJed Brown       }
116877019fcaSJed Brown       for (j=0; j<mi; j++) ix[m+j] = j;
116977019fcaSJed Brown       if (smap) {ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr);}
117077019fcaSJed Brown       /*
117177019fcaSJed Brown         Now we need to extract the monolithic global indices that correspond to the given split global indices.
117277019fcaSJed 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.
117377019fcaSJed Brown         The approach here is ugly because it uses VecScatter to move indices.
117477019fcaSJed Brown        */
117577019fcaSJed Brown       ierr = VecCreateSeq(PETSC_COMM_SELF,mi,&lvec);CHKERRQ(ierr);
117677019fcaSJed Brown       ierr = VecCreateMPI(((PetscObject)isglobal[i])->comm,mi,PETSC_DECIDE,&gvec);CHKERRQ(ierr);
117777019fcaSJed Brown       ierr = ISCreateGeneral(((PetscObject)isglobal[i])->comm,mi,ix+m,PETSC_COPY_VALUES,&isreq);CHKERRQ(ierr);
11780298fd71SBarry Smith       ierr = VecScatterCreate(gvec,isreq,lvec,NULL,&scat);CHKERRQ(ierr);
117977019fcaSJed Brown       ierr = VecGetArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
118077019fcaSJed Brown       for (j=0; j<mi; j++) x[j].integer = ix[m+j];
118177019fcaSJed Brown       ierr = VecRestoreArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
118277019fcaSJed Brown       ierr = VecScatterBegin(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
118377019fcaSJed Brown       ierr = VecScatterEnd(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
118477019fcaSJed Brown       ierr = VecGetArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
118577019fcaSJed Brown       for (j=0; j<mi; j++) ix[m+j] = x[j].integer;
118677019fcaSJed Brown       ierr = VecRestoreArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
118777019fcaSJed Brown       ierr = VecDestroy(&lvec);CHKERRQ(ierr);
118877019fcaSJed Brown       ierr = VecDestroy(&gvec);CHKERRQ(ierr);
118977019fcaSJed Brown       ierr = ISDestroy(&isreq);CHKERRQ(ierr);
119077019fcaSJed Brown       ierr = VecScatterDestroy(&scat);CHKERRQ(ierr);
119177019fcaSJed Brown       m   += mi;
119277019fcaSJed Brown     }
1193f0413b6fSBarry Smith     ierr   = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),1,m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
119477019fcaSJed Brown   } else {
11950298fd71SBarry Smith     *ltog  = NULL;
119677019fcaSJed Brown   }
119777019fcaSJed Brown   PetscFunctionReturn(0);
119877019fcaSJed Brown }
119977019fcaSJed Brown 
120077019fcaSJed Brown 
1201d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1202d8588912SDave May /*
1203d8588912SDave May   nprocessors = NP
1204d8588912SDave May   Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1))
1205d8588912SDave May        proc 0: => (g_0,h_0,)
1206d8588912SDave May        proc 1: => (g_1,h_1,)
1207d8588912SDave May        ...
1208d8588912SDave May        proc nprocs-1: => (g_NP-1,h_NP-1,)
1209d8588912SDave May 
1210d8588912SDave May             proc 0:                      proc 1:                    proc nprocs-1:
1211d8588912SDave May     is[0] = (0,1,2,...,nlocal(g_0)-1)  (0,1,...,nlocal(g_1)-1)  (0,1,...,nlocal(g_NP-1))
1212d8588912SDave May 
1213d8588912SDave May             proc 0:
1214d8588912SDave May     is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1)
1215d8588912SDave May             proc 1:
1216d8588912SDave May     is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1)
1217d8588912SDave May 
1218d8588912SDave May             proc NP-1:
1219d8588912SDave May     is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1)
1220d8588912SDave May */
1221841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1222d8588912SDave May {
1223e2d7f03fSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
12248188e55aSJed Brown   PetscInt       i,j,offset,n,nsum,bs;
1225d8588912SDave May   PetscErrorCode ierr;
12260298fd71SBarry Smith   Mat            sub = NULL;
1227d8588912SDave May 
1228d8588912SDave May   PetscFunctionBegin;
1229854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&vs->isglobal.row);CHKERRQ(ierr);
1230854ce69bSBarry Smith   ierr = PetscMalloc1(nc,&vs->isglobal.col);CHKERRQ(ierr);
1231d8588912SDave May   if (is_row) { /* valid IS is passed in */
1232d8588912SDave May     /* refs on is[] are incremeneted */
1233e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1234d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
123526fbe8dcSKarl Rupp 
1236e2d7f03fSJed Brown       vs->isglobal.row[i] = is_row[i];
1237d8588912SDave May     }
12382ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
12398188e55aSJed Brown     nsum = 0;
12408188e55aSJed Brown     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
12418188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1242ce94432eSBarry Smith       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
12430298fd71SBarry Smith       ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
1244ce94432eSBarry Smith       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
12458188e55aSJed Brown       nsum += n;
12468188e55aSJed Brown     }
1247ce94432eSBarry Smith     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
124830bc264bSJed Brown     offset -= nsum;
1249e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1250f349c1fdSJed Brown       ierr    = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
12510298fd71SBarry Smith       ierr    = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
12522ae74bdbSJed Brown       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1253ce94432eSBarry Smith       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
1254e2d7f03fSJed Brown       ierr    = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
12552ae74bdbSJed Brown       offset += n;
1256d8588912SDave May     }
1257d8588912SDave May   }
1258d8588912SDave May 
1259d8588912SDave May   if (is_col) { /* valid IS is passed in */
1260d8588912SDave May     /* refs on is[] are incremeneted */
1261e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1262d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
126326fbe8dcSKarl Rupp 
1264e2d7f03fSJed Brown       vs->isglobal.col[j] = is_col[j];
1265d8588912SDave May     }
12662ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
12672ae74bdbSJed Brown     offset = A->cmap->rstart;
12688188e55aSJed Brown     nsum   = 0;
12698188e55aSJed Brown     for (j=0; j<vs->nc; j++) {
12708188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
1271ce94432eSBarry Smith       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
12720298fd71SBarry Smith       ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
1273ce94432eSBarry Smith       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
12748188e55aSJed Brown       nsum += n;
12758188e55aSJed Brown     }
1276ce94432eSBarry Smith     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
127730bc264bSJed Brown     offset -= nsum;
1278e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1279f349c1fdSJed Brown       ierr    = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
12800298fd71SBarry Smith       ierr    = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
12812ae74bdbSJed Brown       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1282ce94432eSBarry Smith       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
1283e2d7f03fSJed Brown       ierr    = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
12842ae74bdbSJed Brown       offset += n;
1285d8588912SDave May     }
1286d8588912SDave May   }
1287e2d7f03fSJed Brown 
1288e2d7f03fSJed Brown   /* Set up the local ISs */
1289785e854fSJed Brown   ierr = PetscMalloc1(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
1290785e854fSJed Brown   ierr = PetscMalloc1(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
1291e2d7f03fSJed Brown   for (i=0,offset=0; i<vs->nr; i++) {
1292e2d7f03fSJed Brown     IS                     isloc;
12930298fd71SBarry Smith     ISLocalToGlobalMapping rmap = NULL;
1294e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1295e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
12960298fd71SBarry Smith     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,NULL);CHKERRQ(ierr);}
1297207556f9SJed Brown     if (rmap) {
1298e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1299e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
1300e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1301e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1302207556f9SJed Brown     } else {
1303207556f9SJed Brown       nlocal = 0;
13040298fd71SBarry Smith       isloc  = NULL;
1305207556f9SJed Brown     }
1306e2d7f03fSJed Brown     vs->islocal.row[i] = isloc;
1307e2d7f03fSJed Brown     offset            += nlocal;
1308e2d7f03fSJed Brown   }
13098188e55aSJed Brown   for (i=0,offset=0; i<vs->nc; i++) {
1310e2d7f03fSJed Brown     IS                     isloc;
13110298fd71SBarry Smith     ISLocalToGlobalMapping cmap = NULL;
1312e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1313e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
13140298fd71SBarry Smith     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,NULL,&cmap);CHKERRQ(ierr);}
1315207556f9SJed Brown     if (cmap) {
1316e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1317e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(cmap,&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.col[i] = isloc;
1325e2d7f03fSJed Brown     offset            += nlocal;
1326e2d7f03fSJed Brown   }
13270189643fSJed Brown 
132877019fcaSJed Brown   /* Set up the aggregate ISLocalToGlobalMapping */
132977019fcaSJed Brown   {
133045b6f7e9SBarry Smith     ISLocalToGlobalMapping rmap,cmap;
133145b6f7e9SBarry Smith     ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap);CHKERRQ(ierr);
133245b6f7e9SBarry Smith     ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap);CHKERRQ(ierr);
133377019fcaSJed Brown     if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);}
133477019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
133577019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
133677019fcaSJed Brown   }
133777019fcaSJed Brown 
13380189643fSJed Brown #if defined(PETSC_USE_DEBUG)
13390189643fSJed Brown   for (i=0; i<vs->nr; i++) {
13400189643fSJed Brown     for (j=0; j<vs->nc; j++) {
13410189643fSJed Brown       PetscInt m,n,M,N,mi,ni,Mi,Ni;
13420189643fSJed Brown       Mat      B = vs->m[i][j];
13430189643fSJed Brown       if (!B) continue;
13440189643fSJed Brown       ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr);
13450189643fSJed Brown       ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr);
13460189643fSJed Brown       ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr);
13470189643fSJed Brown       ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr);
13480189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr);
13490189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr);
1350ce94432eSBarry 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);
1351ce94432eSBarry 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);
13520189643fSJed Brown     }
13530189643fSJed Brown   }
13540189643fSJed Brown #endif
1355a061e289SJed Brown 
1356a061e289SJed Brown   /* Set A->assembled if all non-null blocks are currently assembled */
1357a061e289SJed Brown   for (i=0; i<vs->nr; i++) {
1358a061e289SJed Brown     for (j=0; j<vs->nc; j++) {
1359a061e289SJed Brown       if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0);
1360a061e289SJed Brown     }
1361a061e289SJed Brown   }
1362a061e289SJed Brown   A->assembled = PETSC_TRUE;
1363d8588912SDave May   PetscFunctionReturn(0);
1364d8588912SDave May }
1365d8588912SDave May 
136645c38901SJed Brown /*@C
1367659c6bb0SJed Brown    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
1368659c6bb0SJed Brown 
1369659c6bb0SJed Brown    Collective on Mat
1370659c6bb0SJed Brown 
1371659c6bb0SJed Brown    Input Parameter:
1372659c6bb0SJed Brown +  comm - Communicator for the new Mat
1373659c6bb0SJed Brown .  nr - number of nested row blocks
13740298fd71SBarry Smith .  is_row - index sets for each nested row block, or NULL to make contiguous
1375659c6bb0SJed Brown .  nc - number of nested column blocks
13760298fd71SBarry Smith .  is_col - index sets for each nested column block, or NULL to make contiguous
13770298fd71SBarry Smith -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1378659c6bb0SJed Brown 
1379659c6bb0SJed Brown    Output Parameter:
1380659c6bb0SJed Brown .  B - new matrix
1381659c6bb0SJed Brown 
1382659c6bb0SJed Brown    Level: advanced
1383659c6bb0SJed Brown 
1384950540a4SJed Brown .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST
1385659c6bb0SJed Brown @*/
13867087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1387d8588912SDave May {
1388d8588912SDave May   Mat            A;
1389d8588912SDave May   PetscErrorCode ierr;
1390d8588912SDave May 
1391d8588912SDave May   PetscFunctionBegin;
1392c8883902SJed Brown   *B   = 0;
1393d8588912SDave May   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1394c8883902SJed Brown   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
139591a28eb3SBarry Smith   A->preallocated = PETSC_TRUE;
1396c8883902SJed Brown   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1397d8588912SDave May   *B   = A;
1398d8588912SDave May   PetscFunctionReturn(0);
1399d8588912SDave May }
1400659c6bb0SJed Brown 
1401629c3df2SDmitry Karpeev #undef __FUNCT__
1402b68353e5Sstefano_zampini #define __FUNCT__ "MatConvert_Nest_SeqAIJ_fast"
1403b68353e5Sstefano_zampini static PetscErrorCode MatConvert_Nest_SeqAIJ_fast(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1404b68353e5Sstefano_zampini {
1405b68353e5Sstefano_zampini   Mat_Nest       *nest = (Mat_Nest*)A->data;
140623875855Sstefano_zampini   Mat            *trans;
1407b68353e5Sstefano_zampini   PetscScalar    **avv;
1408b68353e5Sstefano_zampini   PetscScalar    *vv;
1409b68353e5Sstefano_zampini   PetscInt       **aii,**ajj;
1410b68353e5Sstefano_zampini   PetscInt       *ii,*jj,*ci;
1411b68353e5Sstefano_zampini   PetscInt       nr,nc,nnz,i,j;
1412b68353e5Sstefano_zampini   PetscBool      done;
1413b68353e5Sstefano_zampini   PetscErrorCode ierr;
1414b68353e5Sstefano_zampini 
1415b68353e5Sstefano_zampini   PetscFunctionBegin;
1416b68353e5Sstefano_zampini   ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr);
1417b68353e5Sstefano_zampini   if (reuse == MAT_REUSE_MATRIX) {
1418b68353e5Sstefano_zampini     PetscInt rnr;
1419b68353e5Sstefano_zampini 
1420b68353e5Sstefano_zampini     ierr = MatGetRowIJ(*newmat,0,PETSC_FALSE,PETSC_FALSE,&rnr,(const PetscInt**)&ii,(const PetscInt**)&jj,&done);CHKERRQ(ierr);
1421b68353e5Sstefano_zampini     if (!done) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"MatGetRowIJ");
1422b68353e5Sstefano_zampini     if (rnr != nr) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Cannot reuse matrix, wrong number of rows");
1423b68353e5Sstefano_zampini     ierr = MatSeqAIJGetArray(*newmat,&vv);CHKERRQ(ierr);
1424b68353e5Sstefano_zampini   }
1425b68353e5Sstefano_zampini   /* extract CSR for nested SeqAIJ matrices */
1426b68353e5Sstefano_zampini   nnz  = 0;
142723875855Sstefano_zampini   ierr = PetscCalloc4(nest->nr*nest->nc,&aii,nest->nr*nest->nc,&ajj,nest->nr*nest->nc,&avv,nest->nr*nest->nc,&trans);CHKERRQ(ierr);
1428b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1429b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1430b68353e5Sstefano_zampini       Mat B = nest->m[i][j];
1431b68353e5Sstefano_zampini       if (B) {
1432b68353e5Sstefano_zampini         PetscScalar *naa;
1433b68353e5Sstefano_zampini         PetscInt    *nii,*njj,nnr;
143423875855Sstefano_zampini         PetscBool   istrans;
1435b68353e5Sstefano_zampini 
143623875855Sstefano_zampini         ierr = PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&istrans);CHKERRQ(ierr);
143723875855Sstefano_zampini         if (istrans) {
143823875855Sstefano_zampini           Mat Bt;
143923875855Sstefano_zampini 
144023875855Sstefano_zampini           ierr = MatTransposeGetMat(B,&Bt);CHKERRQ(ierr);
144123875855Sstefano_zampini           ierr = MatTranspose(Bt,MAT_INITIAL_MATRIX,&trans[i*nest->nc+j]);CHKERRQ(ierr);
144223875855Sstefano_zampini           B    = trans[i*nest->nc+j];
144323875855Sstefano_zampini         }
1444b68353e5Sstefano_zampini         ierr = MatGetRowIJ(B,0,PETSC_FALSE,PETSC_FALSE,&nnr,(const PetscInt**)&nii,(const PetscInt**)&njj,&done);CHKERRQ(ierr);
1445b68353e5Sstefano_zampini         if (!done) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"MatGetRowIJ");
1446b68353e5Sstefano_zampini         ierr = MatSeqAIJGetArray(B,&naa);CHKERRQ(ierr);
1447b68353e5Sstefano_zampini         nnz += nii[nnr];
1448b68353e5Sstefano_zampini 
1449b68353e5Sstefano_zampini         aii[i*nest->nc+j] = nii;
1450b68353e5Sstefano_zampini         ajj[i*nest->nc+j] = njj;
1451b68353e5Sstefano_zampini         avv[i*nest->nc+j] = naa;
1452b68353e5Sstefano_zampini       }
1453b68353e5Sstefano_zampini     }
1454b68353e5Sstefano_zampini   }
1455b68353e5Sstefano_zampini   if (reuse != MAT_REUSE_MATRIX) {
1456b68353e5Sstefano_zampini     ierr = PetscMalloc1(nr+1,&ii);CHKERRQ(ierr);
1457b68353e5Sstefano_zampini     ierr = PetscMalloc1(nnz,&jj);CHKERRQ(ierr);
1458b68353e5Sstefano_zampini     ierr = PetscMalloc1(nnz,&vv);CHKERRQ(ierr);
1459b68353e5Sstefano_zampini   } else {
1460b68353e5Sstefano_zampini     if (nnz != ii[nr]) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Cannot reuse matrix, wrong number of nonzeros");
1461b68353e5Sstefano_zampini   }
1462b68353e5Sstefano_zampini 
1463b68353e5Sstefano_zampini   /* new row pointer */
1464b68353e5Sstefano_zampini   ierr = PetscMemzero(ii,(nr+1)*sizeof(PetscInt));CHKERRQ(ierr);
1465b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1466b68353e5Sstefano_zampini     PetscInt       ncr,rst;
1467b68353e5Sstefano_zampini 
1468b68353e5Sstefano_zampini     ierr = ISStrideGetInfo(nest->isglobal.row[i],&rst,NULL);CHKERRQ(ierr);
1469b68353e5Sstefano_zampini     ierr = ISGetLocalSize(nest->isglobal.row[i],&ncr);CHKERRQ(ierr);
1470b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1471b68353e5Sstefano_zampini       if (aii[i*nest->nc+j]) {
1472b68353e5Sstefano_zampini         PetscInt    *nii = aii[i*nest->nc+j];
1473b68353e5Sstefano_zampini         PetscInt    ir;
1474b68353e5Sstefano_zampini 
1475b68353e5Sstefano_zampini         for (ir=rst; ir<ncr+rst; ++ir) {
1476b68353e5Sstefano_zampini           ii[ir+1] += nii[1]-nii[0];
1477b68353e5Sstefano_zampini           nii++;
1478b68353e5Sstefano_zampini         }
1479b68353e5Sstefano_zampini       }
1480b68353e5Sstefano_zampini     }
1481b68353e5Sstefano_zampini   }
1482b68353e5Sstefano_zampini   for (i=0; i<nr; i++) ii[i+1] += ii[i];
1483b68353e5Sstefano_zampini 
1484b68353e5Sstefano_zampini   /* construct CSR for the new matrix */
1485b68353e5Sstefano_zampini   ierr = PetscCalloc1(nr,&ci);CHKERRQ(ierr);
1486b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1487b68353e5Sstefano_zampini     PetscInt       ncr,rst;
1488b68353e5Sstefano_zampini 
1489b68353e5Sstefano_zampini     ierr = ISStrideGetInfo(nest->isglobal.row[i],&rst,NULL);CHKERRQ(ierr);
1490b68353e5Sstefano_zampini     ierr = ISGetLocalSize(nest->isglobal.row[i],&ncr);CHKERRQ(ierr);
1491b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1492b68353e5Sstefano_zampini       if (aii[i*nest->nc+j]) {
1493b68353e5Sstefano_zampini         PetscScalar *nvv = avv[i*nest->nc+j];
1494b68353e5Sstefano_zampini         PetscInt    *nii = aii[i*nest->nc+j];
1495b68353e5Sstefano_zampini         PetscInt    *njj = ajj[i*nest->nc+j];
1496b68353e5Sstefano_zampini         PetscInt    ir,cst;
1497b68353e5Sstefano_zampini 
1498b68353e5Sstefano_zampini         ierr = ISStrideGetInfo(nest->isglobal.col[j],&cst,NULL);CHKERRQ(ierr);
1499b68353e5Sstefano_zampini         for (ir=rst; ir<ncr+rst; ++ir) {
1500b68353e5Sstefano_zampini           PetscInt ij,rsize = nii[1]-nii[0],ist = ii[ir]+ci[ir];
1501b68353e5Sstefano_zampini 
1502b68353e5Sstefano_zampini           for (ij=0;ij<rsize;ij++) {
1503b68353e5Sstefano_zampini             jj[ist+ij] = *njj+cst;
1504b68353e5Sstefano_zampini             vv[ist+ij] = *nvv;
1505b68353e5Sstefano_zampini             njj++;
1506b68353e5Sstefano_zampini             nvv++;
1507b68353e5Sstefano_zampini           }
1508b68353e5Sstefano_zampini           ci[ir] += rsize;
1509b68353e5Sstefano_zampini           nii++;
1510b68353e5Sstefano_zampini         }
1511b68353e5Sstefano_zampini       }
1512b68353e5Sstefano_zampini     }
1513b68353e5Sstefano_zampini   }
1514b68353e5Sstefano_zampini   ierr = PetscFree(ci);CHKERRQ(ierr);
1515b68353e5Sstefano_zampini 
1516b68353e5Sstefano_zampini   /* restore info */
1517b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1518b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1519b68353e5Sstefano_zampini       Mat B = nest->m[i][j];
1520b68353e5Sstefano_zampini       if (B) {
1521b68353e5Sstefano_zampini         PetscInt nnr = 0, k = i*nest->nc+j;
152223875855Sstefano_zampini 
152323875855Sstefano_zampini         B    = (trans[k] ? trans[k] : B);
1524b68353e5Sstefano_zampini         ierr = MatRestoreRowIJ(B,0,PETSC_FALSE,PETSC_FALSE,&nnr,(const PetscInt**)&aii[k],(const PetscInt**)&ajj[k],&done);CHKERRQ(ierr);
1525b68353e5Sstefano_zampini         if (!done) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"MatRestoreRowIJ");
1526b68353e5Sstefano_zampini         ierr = MatSeqAIJRestoreArray(B,&avv[k]);CHKERRQ(ierr);
152723875855Sstefano_zampini         ierr = MatDestroy(&trans[k]);CHKERRQ(ierr);
1528b68353e5Sstefano_zampini       }
1529b68353e5Sstefano_zampini     }
1530b68353e5Sstefano_zampini   }
153123875855Sstefano_zampini   ierr = PetscFree4(aii,ajj,avv,trans);CHKERRQ(ierr);
1532b68353e5Sstefano_zampini 
1533b68353e5Sstefano_zampini   /* finalize newmat */
1534b68353e5Sstefano_zampini   if (reuse == MAT_INITIAL_MATRIX) {
1535b68353e5Sstefano_zampini     ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),nr,nc,ii,jj,vv,newmat);CHKERRQ(ierr);
1536b68353e5Sstefano_zampini   } else if (reuse == MAT_INPLACE_MATRIX) {
1537b68353e5Sstefano_zampini     Mat B;
1538b68353e5Sstefano_zampini 
1539b68353e5Sstefano_zampini     ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),nr,nc,ii,jj,vv,&B);CHKERRQ(ierr);
1540b68353e5Sstefano_zampini     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
1541b68353e5Sstefano_zampini   }
1542b68353e5Sstefano_zampini   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1543b68353e5Sstefano_zampini   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1544b68353e5Sstefano_zampini   {
1545b68353e5Sstefano_zampini     Mat_SeqAIJ *a = (Mat_SeqAIJ*)((*newmat)->data);
1546b68353e5Sstefano_zampini     a->free_a     = PETSC_TRUE;
1547b68353e5Sstefano_zampini     a->free_ij    = PETSC_TRUE;
1548b68353e5Sstefano_zampini   }
1549b68353e5Sstefano_zampini   PetscFunctionReturn(0);
1550b68353e5Sstefano_zampini }
1551b68353e5Sstefano_zampini 
1552b68353e5Sstefano_zampini #undef __FUNCT__
1553629c3df2SDmitry Karpeev #define __FUNCT__ "MatConvert_Nest_AIJ"
1554cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1555629c3df2SDmitry Karpeev {
1556629c3df2SDmitry Karpeev   PetscErrorCode ierr;
1557629c3df2SDmitry Karpeev   Mat_Nest       *nest = (Mat_Nest*)A->data;
155883b1a929SMark Adams   PetscInt       m,n,M,N,i,j,k,*dnnz,*onnz,rstart;
1559649b366bSFande Kong   PetscInt       cstart,cend;
1560b68353e5Sstefano_zampini   PetscMPIInt    size;
1561629c3df2SDmitry Karpeev   Mat            C;
1562629c3df2SDmitry Karpeev 
1563629c3df2SDmitry Karpeev   PetscFunctionBegin;
1564b68353e5Sstefano_zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1565b68353e5Sstefano_zampini   if (size == 1) { /* look for a special case with SeqAIJ matrices and strided-1, contiguous, blocks */
1566b68353e5Sstefano_zampini     PetscInt  nf;
1567b68353e5Sstefano_zampini     PetscBool fast;
1568b68353e5Sstefano_zampini 
1569b68353e5Sstefano_zampini     ierr = PetscStrcmp(newtype,MATAIJ,&fast);CHKERRQ(ierr);
1570b68353e5Sstefano_zampini     if (!fast) {
1571b68353e5Sstefano_zampini       ierr = PetscStrcmp(newtype,MATSEQAIJ,&fast);CHKERRQ(ierr);
1572b68353e5Sstefano_zampini     }
1573b68353e5Sstefano_zampini     for (i=0; i<nest->nr && fast; ++i) {
1574b68353e5Sstefano_zampini       for (j=0; j<nest->nc && fast; ++j) {
1575b68353e5Sstefano_zampini         Mat B = nest->m[i][j];
1576b68353e5Sstefano_zampini         if (B) {
1577b68353e5Sstefano_zampini           ierr = PetscObjectTypeCompare((PetscObject)B,MATSEQAIJ,&fast);CHKERRQ(ierr);
157823875855Sstefano_zampini           if (!fast) {
157923875855Sstefano_zampini             PetscBool istrans;
158023875855Sstefano_zampini 
158123875855Sstefano_zampini             ierr = PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&istrans);CHKERRQ(ierr);
158223875855Sstefano_zampini             if (istrans) {
158323875855Sstefano_zampini               Mat Bt;
158423875855Sstefano_zampini 
158523875855Sstefano_zampini               ierr = MatTransposeGetMat(B,&Bt);CHKERRQ(ierr);
158623875855Sstefano_zampini               ierr = PetscObjectTypeCompare((PetscObject)Bt,MATSEQAIJ,&fast);CHKERRQ(ierr);
158723875855Sstefano_zampini             }
1588b68353e5Sstefano_zampini           }
1589b68353e5Sstefano_zampini         }
1590b68353e5Sstefano_zampini       }
1591b68353e5Sstefano_zampini     }
1592b68353e5Sstefano_zampini     for (i=0, nf=0; i<nest->nr && fast; ++i) {
1593b68353e5Sstefano_zampini       ierr = PetscObjectTypeCompare((PetscObject)nest->isglobal.row[i],ISSTRIDE,&fast);CHKERRQ(ierr);
1594b68353e5Sstefano_zampini       if (fast) {
1595b68353e5Sstefano_zampini         PetscInt f,s;
1596b68353e5Sstefano_zampini 
1597b68353e5Sstefano_zampini         ierr = ISStrideGetInfo(nest->isglobal.row[i],&f,&s);CHKERRQ(ierr);
1598b68353e5Sstefano_zampini         if (f != nf || s != 1) { fast = PETSC_FALSE; }
1599b68353e5Sstefano_zampini         else {
1600b68353e5Sstefano_zampini           ierr = ISGetSize(nest->isglobal.row[i],&f);CHKERRQ(ierr);
1601b68353e5Sstefano_zampini           nf  += f;
1602b68353e5Sstefano_zampini         }
1603b68353e5Sstefano_zampini       }
1604b68353e5Sstefano_zampini     }
1605b68353e5Sstefano_zampini     for (i=0, nf=0; i<nest->nc && fast; ++i) {
1606b68353e5Sstefano_zampini       ierr = PetscObjectTypeCompare((PetscObject)nest->isglobal.col[i],ISSTRIDE,&fast);CHKERRQ(ierr);
1607b68353e5Sstefano_zampini       if (fast) {
1608b68353e5Sstefano_zampini         PetscInt f,s;
1609b68353e5Sstefano_zampini 
1610b68353e5Sstefano_zampini         ierr = ISStrideGetInfo(nest->isglobal.col[i],&f,&s);CHKERRQ(ierr);
1611b68353e5Sstefano_zampini         if (f != nf || s != 1) { fast = PETSC_FALSE; }
1612b68353e5Sstefano_zampini         else {
1613b68353e5Sstefano_zampini           ierr = ISGetSize(nest->isglobal.col[i],&f);CHKERRQ(ierr);
1614b68353e5Sstefano_zampini           nf  += f;
1615b68353e5Sstefano_zampini         }
1616b68353e5Sstefano_zampini       }
1617b68353e5Sstefano_zampini     }
1618b68353e5Sstefano_zampini     if (fast) {
1619b68353e5Sstefano_zampini       ierr = MatConvert_Nest_SeqAIJ_fast(A,newtype,reuse,newmat);CHKERRQ(ierr);
1620b68353e5Sstefano_zampini       PetscFunctionReturn(0);
1621b68353e5Sstefano_zampini     }
1622b68353e5Sstefano_zampini   }
1623629c3df2SDmitry Karpeev   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
1624629c3df2SDmitry Karpeev   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
1625649b366bSFande Kong   ierr = MatGetOwnershipRangeColumn(A,&cstart,&cend);CHKERRQ(ierr);
1626629c3df2SDmitry Karpeev   switch (reuse) {
1627629c3df2SDmitry Karpeev   case MAT_INITIAL_MATRIX:
1628ce94432eSBarry Smith     ierr    = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
1629629c3df2SDmitry Karpeev     ierr    = MatSetType(C,newtype);CHKERRQ(ierr);
1630629c3df2SDmitry Karpeev     ierr    = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr);
1631629c3df2SDmitry Karpeev     *newmat = C;
1632629c3df2SDmitry Karpeev     break;
1633629c3df2SDmitry Karpeev   case MAT_REUSE_MATRIX:
1634629c3df2SDmitry Karpeev     C = *newmat;
1635629c3df2SDmitry Karpeev     break;
1636ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse");
1637629c3df2SDmitry Karpeev   }
1638785e854fSJed Brown   ierr = PetscMalloc1(2*m,&dnnz);CHKERRQ(ierr);
1639629c3df2SDmitry Karpeev   onnz = dnnz + m;
1640629c3df2SDmitry Karpeev   for (k=0; k<m; k++) {
1641629c3df2SDmitry Karpeev     dnnz[k] = 0;
1642629c3df2SDmitry Karpeev     onnz[k] = 0;
1643629c3df2SDmitry Karpeev   }
1644629c3df2SDmitry Karpeev   for (j=0; j<nest->nc; ++j) {
1645629c3df2SDmitry Karpeev     IS             bNis;
1646629c3df2SDmitry Karpeev     PetscInt       bN;
1647629c3df2SDmitry Karpeev     const PetscInt *bNindices;
1648629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
1649629c3df2SDmitry Karpeev     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1650629c3df2SDmitry Karpeev     ierr = ISGetSize(bNis, &bN);CHKERRQ(ierr);
1651629c3df2SDmitry Karpeev     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1652629c3df2SDmitry Karpeev     for (i=0; i<nest->nr; ++i) {
1653629c3df2SDmitry Karpeev       PetscSF        bmsf;
1654649b366bSFande Kong       PetscSFNode    *iremote;
1655629c3df2SDmitry Karpeev       Mat            B;
1656649b366bSFande Kong       PetscInt       bm, *sub_dnnz,*sub_onnz, br;
1657629c3df2SDmitry Karpeev       const PetscInt *bmindices;
1658629c3df2SDmitry Karpeev       B = nest->m[i][j];
1659629c3df2SDmitry Karpeev       if (!B) continue;
1660629c3df2SDmitry Karpeev       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1661629c3df2SDmitry Karpeev       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1662ce94432eSBarry Smith       ierr = PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);CHKERRQ(ierr);
1663649b366bSFande Kong       ierr = PetscMalloc1(bm,&iremote);CHKERRQ(ierr);
1664649b366bSFande Kong       ierr = PetscMalloc1(bm,&sub_dnnz);CHKERRQ(ierr);
1665649b366bSFande Kong       ierr = PetscMalloc1(bm,&sub_onnz);CHKERRQ(ierr);
1666649b366bSFande Kong       for (k = 0; k < bm; ++k){
1667649b366bSFande Kong     	sub_dnnz[k] = 0;
1668649b366bSFande Kong     	sub_onnz[k] = 0;
1669649b366bSFande Kong       }
1670629c3df2SDmitry Karpeev       /*
1671629c3df2SDmitry Karpeev        Locate the owners for all of the locally-owned global row indices for this row block.
1672629c3df2SDmitry Karpeev        These determine the roots of PetscSF used to communicate preallocation data to row owners.
1673629c3df2SDmitry Karpeev        The roots correspond to the dnnz and onnz entries; thus, there are two roots per row.
1674629c3df2SDmitry Karpeev        */
167583b1a929SMark Adams       ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr);
1676629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
1677649b366bSFande Kong         PetscInt       row = bmindices[br], rowowner = 0, brncols, col;
1678629c3df2SDmitry Karpeev         const PetscInt *brcols;
1679a4b3d3acSMatthew G Knepley         PetscInt       rowrel = 0; /* row's relative index on its owner rank */
1680629c3df2SDmitry Karpeev         ierr      = PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);CHKERRQ(ierr);
1681649b366bSFande Kong         /* how many roots  */
1682649b366bSFande Kong         iremote[br].rank = rowowner; iremote[br].index = rowrel;           /* edge from bmdnnz to dnnz */
1683649b366bSFande Kong         /* get nonzero pattern */
168483b1a929SMark Adams         ierr = MatGetRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr);
1685629c3df2SDmitry Karpeev         for (k=0; k<brncols; k++) {
1686629c3df2SDmitry Karpeev           col  = bNindices[brcols[k]];
1687649b366bSFande Kong           if (col>=A->cmap->range[rowowner] && col<A->cmap->range[rowowner+1]) {
1688649b366bSFande Kong             sub_dnnz[br]++;
1689649b366bSFande Kong           } else {
1690649b366bSFande Kong             sub_onnz[br]++;
1691649b366bSFande Kong           }
1692629c3df2SDmitry Karpeev         }
169383b1a929SMark Adams         ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr);
1694629c3df2SDmitry Karpeev       }
1695629c3df2SDmitry Karpeev       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1696629c3df2SDmitry Karpeev       /* bsf will have to take care of disposing of bedges. */
1697649b366bSFande Kong       ierr = PetscSFSetGraph(bmsf,m,bm,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1698649b366bSFande Kong       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr);
1699649b366bSFande Kong       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr);
1700649b366bSFande Kong       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr);
1701649b366bSFande Kong       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr);
1702649b366bSFande Kong       ierr = PetscFree(sub_dnnz);CHKERRQ(ierr);
1703649b366bSFande Kong       ierr = PetscFree(sub_onnz);CHKERRQ(ierr);
1704629c3df2SDmitry Karpeev       ierr = PetscSFDestroy(&bmsf);CHKERRQ(ierr);
1705629c3df2SDmitry Karpeev     }
170622d28d08SBarry Smith     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
1707629c3df2SDmitry Karpeev     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
170865a4a0a3Sstefano_zampini   }
170965a4a0a3Sstefano_zampini   /* Resize preallocation if overestimated */
171065a4a0a3Sstefano_zampini   for (i=0;i<m;i++) {
171165a4a0a3Sstefano_zampini     dnnz[i] = PetscMin(dnnz[i],A->cmap->n);
171265a4a0a3Sstefano_zampini     onnz[i] = PetscMin(onnz[i],A->cmap->N - A->cmap->n);
1713629c3df2SDmitry Karpeev   }
1714629c3df2SDmitry Karpeev   ierr = MatSeqAIJSetPreallocation(C,0,dnnz);CHKERRQ(ierr);
1715629c3df2SDmitry Karpeev   ierr = MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);CHKERRQ(ierr);
1716629c3df2SDmitry Karpeev   ierr = PetscFree(dnnz);CHKERRQ(ierr);
1717629c3df2SDmitry Karpeev 
1718629c3df2SDmitry Karpeev   /* Fill by row */
1719629c3df2SDmitry Karpeev   for (j=0; j<nest->nc; ++j) {
1720629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
1721629c3df2SDmitry Karpeev     IS             bNis;
1722629c3df2SDmitry Karpeev     PetscInt       bN;
1723629c3df2SDmitry Karpeev     const PetscInt *bNindices;
1724629c3df2SDmitry Karpeev     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1725629c3df2SDmitry Karpeev     ierr = ISGetSize(bNis,&bN);CHKERRQ(ierr);
1726629c3df2SDmitry Karpeev     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1727629c3df2SDmitry Karpeev     for (i=0; i<nest->nr; ++i) {
1728629c3df2SDmitry Karpeev       Mat            B;
1729629c3df2SDmitry Karpeev       PetscInt       bm, br;
1730629c3df2SDmitry Karpeev       const PetscInt *bmindices;
1731629c3df2SDmitry Karpeev       B = nest->m[i][j];
1732629c3df2SDmitry Karpeev       if (!B) continue;
1733629c3df2SDmitry Karpeev       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1734629c3df2SDmitry Karpeev       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
173583b1a929SMark Adams       ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr);
1736629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
1737629c3df2SDmitry Karpeev         PetscInt          row = bmindices[br], brncols,  *cols;
1738629c3df2SDmitry Karpeev         const PetscInt    *brcols;
1739629c3df2SDmitry Karpeev         const PetscScalar *brcoldata;
174083b1a929SMark Adams         ierr = MatGetRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
1741785e854fSJed Brown         ierr = PetscMalloc1(brncols,&cols);CHKERRQ(ierr);
174226fbe8dcSKarl Rupp         for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]];
1743629c3df2SDmitry Karpeev         /*
1744629c3df2SDmitry Karpeev           Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match.
1745629c3df2SDmitry Karpeev           Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES.
1746629c3df2SDmitry Karpeev          */
1747a2ea699eSBarry Smith         ierr = MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);CHKERRQ(ierr);
174883b1a929SMark Adams         ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
1749629c3df2SDmitry Karpeev         ierr = PetscFree(cols);CHKERRQ(ierr);
1750629c3df2SDmitry Karpeev       }
1751629c3df2SDmitry Karpeev       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1752629c3df2SDmitry Karpeev     }
1753a2ea699eSBarry Smith     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
1754629c3df2SDmitry Karpeev     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
1755629c3df2SDmitry Karpeev   }
1756629c3df2SDmitry Karpeev   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1757629c3df2SDmitry Karpeev   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1758629c3df2SDmitry Karpeev   PetscFunctionReturn(0);
1759629c3df2SDmitry Karpeev }
1760629c3df2SDmitry Karpeev 
1761659c6bb0SJed Brown /*MC
1762659c6bb0SJed Brown   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
1763659c6bb0SJed Brown 
1764659c6bb0SJed Brown   Level: intermediate
1765659c6bb0SJed Brown 
1766659c6bb0SJed Brown   Notes:
1767659c6bb0SJed Brown   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1768659c6bb0SJed Brown   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1769950540a4SJed Brown   It is usually used with DMComposite and DMCreateMatrix()
1770659c6bb0SJed Brown 
1771659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest()
1772659c6bb0SJed Brown M*/
17738cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A)
1774c8883902SJed Brown {
1775c8883902SJed Brown   Mat_Nest       *s;
1776c8883902SJed Brown   PetscErrorCode ierr;
1777c8883902SJed Brown 
1778c8883902SJed Brown   PetscFunctionBegin;
1779b00a9115SJed Brown   ierr    = PetscNewLog(A,&s);CHKERRQ(ierr);
1780c8883902SJed Brown   A->data = (void*)s;
1781e7c19651SJed Brown 
1782e7c19651SJed Brown   s->nr            = -1;
1783e7c19651SJed Brown   s->nc            = -1;
17840298fd71SBarry Smith   s->m             = NULL;
1785e7c19651SJed Brown   s->splitassembly = PETSC_FALSE;
1786c8883902SJed Brown 
1787c8883902SJed Brown   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
178826fbe8dcSKarl Rupp 
1789c8883902SJed Brown   A->ops->mult                  = MatMult_Nest;
17909194d70fSJed Brown   A->ops->multadd               = MatMultAdd_Nest;
1791c8883902SJed Brown   A->ops->multtranspose         = MatMultTranspose_Nest;
17929194d70fSJed Brown   A->ops->multtransposeadd      = MatMultTransposeAdd_Nest;
1793f8170845SAlex Fikl   A->ops->transpose             = MatTranspose_Nest;
1794c8883902SJed Brown   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1795c8883902SJed Brown   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1796c8883902SJed Brown   A->ops->zeroentries           = MatZeroEntries_Nest;
1797c222c20dSDavid Ham   A->ops->copy                  = MatCopy_Nest;
1798c8883902SJed Brown   A->ops->duplicate             = MatDuplicate_Nest;
1799c8883902SJed Brown   A->ops->getsubmatrix          = MatGetSubMatrix_Nest;
1800c8883902SJed Brown   A->ops->destroy               = MatDestroy_Nest;
1801c8883902SJed Brown   A->ops->view                  = MatView_Nest;
1802c8883902SJed Brown   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1803c8883902SJed Brown   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1804c8883902SJed Brown   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
1805429bac76SJed Brown   A->ops->getdiagonal           = MatGetDiagonal_Nest;
1806429bac76SJed Brown   A->ops->diagonalscale         = MatDiagonalScale_Nest;
1807a061e289SJed Brown   A->ops->scale                 = MatScale_Nest;
1808a061e289SJed Brown   A->ops->shift                 = MatShift_Nest;
180913135bc6SAlex Fikl   A->ops->diagonalset           = MatDiagonalSet_Nest;
1810f8170845SAlex Fikl   A->ops->setrandom             = MatSetRandom_Nest;
1811c8883902SJed Brown 
1812c8883902SJed Brown   A->spptr        = 0;
1813c8883902SJed Brown   A->assembled    = PETSC_FALSE;
1814c8883902SJed Brown 
1815c8883902SJed Brown   /* expose Nest api's */
1816bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",   MatNestGetSubMat_Nest);CHKERRQ(ierr);
1817bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",   MatNestSetSubMat_Nest);CHKERRQ(ierr);
1818bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",  MatNestGetSubMats_Nest);CHKERRQ(ierr);
1819bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",     MatNestGetSize_Nest);CHKERRQ(ierr);
1820bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",      MatNestGetISs_Nest);CHKERRQ(ierr);
1821bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C", MatNestGetLocalISs_Nest);CHKERRQ(ierr);
1822bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",  MatNestSetVecType_Nest);CHKERRQ(ierr);
1823bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",  MatNestSetSubMats_Nest);CHKERRQ(ierr);
182483b1a929SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",MatConvert_Nest_AIJ);CHKERRQ(ierr);
18255e3038f0Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_is_C", MatConvert_Nest_IS);CHKERRQ(ierr);
1826c8883902SJed Brown 
1827c8883902SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
1828c8883902SJed Brown   PetscFunctionReturn(0);
1829c8883902SJed Brown }
1830