xref: /petsc/src/mat/impls/nest/matnest.c (revision d86155a6cb4cbd4961d2c1cfed91c0b6da158021)
1d8588912SDave May 
2aaa7dc30SBarry Smith #include <../src/mat/impls/nest/matnestimpl.h> /*I   "petscmat.h"   I*/
30c312b8eSJed Brown #include <petscsf.h>
4d8588912SDave May 
5c8883902SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]);
62a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left);
7c8883902SJed Brown 
8d8588912SDave May /* private functions */
9d8588912SDave May #undef __FUNCT__
108188e55aSJed Brown #define __FUNCT__ "MatNestGetSizes_Private"
118188e55aSJed Brown static PetscErrorCode MatNestGetSizes_Private(Mat A,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N)
12d8588912SDave May {
13d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
148188e55aSJed Brown   PetscInt       i,j;
15d8588912SDave May   PetscErrorCode ierr;
16d8588912SDave May 
17d8588912SDave May   PetscFunctionBegin;
188188e55aSJed Brown   *m = *n = *M = *N = 0;
198188e55aSJed Brown   for (i=0; i<bA->nr; i++) {  /* rows */
208188e55aSJed Brown     PetscInt sm,sM;
218188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.row[i],&sm);CHKERRQ(ierr);
228188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.row[i],&sM);CHKERRQ(ierr);
238188e55aSJed Brown     *m  += sm;
248188e55aSJed Brown     *M  += sM;
25d8588912SDave May   }
268188e55aSJed Brown   for (j=0; j<bA->nc; j++) {  /* cols */
278188e55aSJed Brown     PetscInt sn,sN;
288188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.col[j],&sn);CHKERRQ(ierr);
298188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.col[j],&sN);CHKERRQ(ierr);
308188e55aSJed Brown     *n  += sn;
318188e55aSJed Brown     *N  += sN;
32d8588912SDave May   }
33d8588912SDave May   PetscFunctionReturn(0);
34d8588912SDave May }
35d8588912SDave May 
36d8588912SDave May /* operations */
37d8588912SDave May #undef __FUNCT__
38d8588912SDave May #define __FUNCT__ "MatMult_Nest"
39207556f9SJed Brown static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y)
40d8588912SDave May {
41d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
42207556f9SJed Brown   Vec            *bx = bA->right,*by = bA->left;
43207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
44d8588912SDave May   PetscErrorCode ierr;
45d8588912SDave May 
46d8588912SDave May   PetscFunctionBegin;
47207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
48207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
49207556f9SJed Brown   for (i=0; i<nr; i++) {
50d8588912SDave May     ierr = VecZeroEntries(by[i]);CHKERRQ(ierr);
51207556f9SJed Brown     for (j=0; j<nc; j++) {
52207556f9SJed Brown       if (!bA->m[i][j]) continue;
53d8588912SDave May       /* y[i] <- y[i] + A[i][j] * x[j] */
54d8588912SDave May       ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr);
55d8588912SDave May     }
56d8588912SDave May   }
57207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
58207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
59d8588912SDave May   PetscFunctionReturn(0);
60d8588912SDave May }
61d8588912SDave May 
62d8588912SDave May #undef __FUNCT__
639194d70fSJed Brown #define __FUNCT__ "MatMultAdd_Nest"
649194d70fSJed Brown static PetscErrorCode MatMultAdd_Nest(Mat A,Vec x,Vec y,Vec z)
659194d70fSJed Brown {
669194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
679194d70fSJed Brown   Vec            *bx = bA->right,*bz = bA->left;
689194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
699194d70fSJed Brown   PetscErrorCode ierr;
709194d70fSJed Brown 
719194d70fSJed Brown   PetscFunctionBegin;
729194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
739194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
749194d70fSJed Brown   for (i=0; i<nr; i++) {
759194d70fSJed Brown     if (y != z) {
769194d70fSJed Brown       Vec by;
779194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
789194d70fSJed Brown       ierr = VecCopy(by,bz[i]);CHKERRQ(ierr);
79336d21e7SJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
809194d70fSJed Brown     }
819194d70fSJed Brown     for (j=0; j<nc; j++) {
829194d70fSJed Brown       if (!bA->m[i][j]) continue;
839194d70fSJed Brown       /* y[i] <- y[i] + A[i][j] * x[j] */
849194d70fSJed Brown       ierr = MatMultAdd(bA->m[i][j],bx[j],bz[i],bz[i]);CHKERRQ(ierr);
859194d70fSJed Brown     }
869194d70fSJed Brown   }
879194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
889194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
899194d70fSJed Brown   PetscFunctionReturn(0);
909194d70fSJed Brown }
919194d70fSJed Brown 
929194d70fSJed Brown #undef __FUNCT__
93d8588912SDave May #define __FUNCT__ "MatMultTranspose_Nest"
94207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y)
95d8588912SDave May {
96d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
97207556f9SJed Brown   Vec            *bx = bA->left,*by = bA->right;
98207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
99d8588912SDave May   PetscErrorCode ierr;
100d8588912SDave May 
101d8588912SDave May   PetscFunctionBegin;
102609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
103609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
104207556f9SJed Brown   for (j=0; j<nc; j++) {
105609e31cbSJed Brown     ierr = VecZeroEntries(by[j]);CHKERRQ(ierr);
106609e31cbSJed Brown     for (i=0; i<nr; i++) {
1076c75ac25SJed Brown       if (!bA->m[i][j]) continue;
108609e31cbSJed Brown       /* y[j] <- y[j] + (A[i][j])^T * x[i] */
109609e31cbSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],by[j],by[j]);CHKERRQ(ierr);
110d8588912SDave May     }
111d8588912SDave May   }
112609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
113609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
114d8588912SDave May   PetscFunctionReturn(0);
115d8588912SDave May }
116d8588912SDave May 
117d8588912SDave May #undef __FUNCT__
1189194d70fSJed Brown #define __FUNCT__ "MatMultTransposeAdd_Nest"
1199194d70fSJed Brown static PetscErrorCode MatMultTransposeAdd_Nest(Mat A,Vec x,Vec y,Vec z)
1209194d70fSJed Brown {
1219194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
1229194d70fSJed Brown   Vec            *bx = bA->left,*bz = bA->right;
1239194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
1249194d70fSJed Brown   PetscErrorCode ierr;
1259194d70fSJed Brown 
1269194d70fSJed Brown   PetscFunctionBegin;
1279194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
1289194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
1299194d70fSJed Brown   for (j=0; j<nc; j++) {
1309194d70fSJed Brown     if (y != z) {
1319194d70fSJed Brown       Vec by;
1329194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
1339194d70fSJed Brown       ierr = VecCopy(by,bz[j]);CHKERRQ(ierr);
1349194d70fSJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
1359194d70fSJed Brown     }
1369194d70fSJed Brown     for (i=0; i<nr; i++) {
1376c75ac25SJed Brown       if (!bA->m[i][j]) continue;
1389194d70fSJed Brown       /* z[j] <- y[j] + (A[i][j])^T * x[i] */
1399194d70fSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],bz[j],bz[j]);CHKERRQ(ierr);
1409194d70fSJed Brown     }
1419194d70fSJed Brown   }
1429194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
1439194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
1449194d70fSJed Brown   PetscFunctionReturn(0);
1459194d70fSJed Brown }
1469194d70fSJed Brown 
1479194d70fSJed Brown #undef __FUNCT__
148f8170845SAlex Fikl #define __FUNCT__ "MatTranspose_Nest"
149f8170845SAlex Fikl static PetscErrorCode MatTranspose_Nest(Mat A,MatReuse reuse,Mat *B)
150f8170845SAlex Fikl {
151f8170845SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data, *bC;
152f8170845SAlex Fikl   Mat            C;
153f8170845SAlex Fikl   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
154f8170845SAlex Fikl   PetscErrorCode ierr;
155f8170845SAlex Fikl 
156f8170845SAlex Fikl   PetscFunctionBegin;
157f8170845SAlex 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");
158f8170845SAlex Fikl 
159f8170845SAlex Fikl   if (reuse == MAT_INITIAL_MATRIX || A == *B) {
160f8170845SAlex Fikl     Mat *subs;
161f8170845SAlex Fikl     IS  *is_row,*is_col;
162f8170845SAlex Fikl 
163f8170845SAlex Fikl     ierr = PetscCalloc1(nr * nc,&subs);CHKERRQ(ierr);
164f8170845SAlex Fikl     ierr = PetscMalloc2(nr,&is_row,nc,&is_col);CHKERRQ(ierr);
165f8170845SAlex Fikl     ierr = MatNestGetISs(A,is_row,is_col);CHKERRQ(ierr);
166ddeb9bd8SAlex Fikl     if (reuse == MAT_REUSE_MATRIX) {
167ddeb9bd8SAlex Fikl       for (i=0; i<nr; i++) {
168ddeb9bd8SAlex Fikl         for (j=0; j<nc; j++) {
169ddeb9bd8SAlex Fikl           subs[i + nr * j] = bA->m[i][j];
170ddeb9bd8SAlex Fikl         }
171ddeb9bd8SAlex Fikl       }
172ddeb9bd8SAlex Fikl     }
173ddeb9bd8SAlex Fikl 
174f8170845SAlex Fikl     ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nc,is_col,nr,is_row,subs,&C);CHKERRQ(ierr);
175f8170845SAlex Fikl     ierr = PetscFree(subs);CHKERRQ(ierr);
1763d994f23SBarry Smith     ierr = PetscFree2(is_row,is_col);CHKERRQ(ierr);
177f8170845SAlex Fikl   } else {
178f8170845SAlex Fikl     C = *B;
179f8170845SAlex Fikl   }
180f8170845SAlex Fikl 
181f8170845SAlex Fikl   bC = (Mat_Nest*)C->data;
182f8170845SAlex Fikl   for (i=0; i<nr; i++) {
183f8170845SAlex Fikl     for (j=0; j<nc; j++) {
184f8170845SAlex Fikl       if (bA->m[i][j]) {
185f8170845SAlex Fikl         ierr = MatTranspose(bA->m[i][j], reuse, &(bC->m[j][i]));CHKERRQ(ierr);
186f8170845SAlex Fikl       } else {
187f8170845SAlex Fikl         bC->m[j][i] = NULL;
188f8170845SAlex Fikl       }
189f8170845SAlex Fikl     }
190f8170845SAlex Fikl   }
191f8170845SAlex Fikl 
192f8170845SAlex Fikl   if (reuse == MAT_INITIAL_MATRIX || A != *B) {
193f8170845SAlex Fikl     *B = C;
194f8170845SAlex Fikl   } else {
195f8170845SAlex Fikl     ierr = MatHeaderMerge(A, &C);CHKERRQ(ierr);
196f8170845SAlex Fikl   }
197f8170845SAlex Fikl   PetscFunctionReturn(0);
198f8170845SAlex Fikl }
199f8170845SAlex Fikl 
200f8170845SAlex Fikl #undef __FUNCT__
201e2d7f03fSJed Brown #define __FUNCT__ "MatNestDestroyISList"
202e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list)
203e2d7f03fSJed Brown {
204e2d7f03fSJed Brown   PetscErrorCode ierr;
205e2d7f03fSJed Brown   IS             *lst = *list;
206e2d7f03fSJed Brown   PetscInt       i;
207e2d7f03fSJed Brown 
208e2d7f03fSJed Brown   PetscFunctionBegin;
209e2d7f03fSJed Brown   if (!lst) PetscFunctionReturn(0);
2106bf464f9SBarry Smith   for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);}
211e2d7f03fSJed Brown   ierr  = PetscFree(lst);CHKERRQ(ierr);
2120298fd71SBarry Smith   *list = NULL;
213e2d7f03fSJed Brown   PetscFunctionReturn(0);
214e2d7f03fSJed Brown }
215e2d7f03fSJed Brown 
216e2d7f03fSJed Brown #undef __FUNCT__
217d8588912SDave May #define __FUNCT__ "MatDestroy_Nest"
218207556f9SJed Brown static PetscErrorCode MatDestroy_Nest(Mat A)
219d8588912SDave May {
220d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
221d8588912SDave May   PetscInt       i,j;
222d8588912SDave May   PetscErrorCode ierr;
223d8588912SDave May 
224d8588912SDave May   PetscFunctionBegin;
225d8588912SDave May   /* release the matrices and the place holders */
226e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr);
227e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr);
228e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
229e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
230d8588912SDave May 
231d8588912SDave May   ierr = PetscFree(vs->row_len);CHKERRQ(ierr);
232d8588912SDave May   ierr = PetscFree(vs->col_len);CHKERRQ(ierr);
233d8588912SDave May 
234207556f9SJed Brown   ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr);
235207556f9SJed Brown 
236d8588912SDave May   /* release the matrices and the place holders */
237d8588912SDave May   if (vs->m) {
238d8588912SDave May     for (i=0; i<vs->nr; i++) {
239d8588912SDave May       for (j=0; j<vs->nc; j++) {
2406bf464f9SBarry Smith         ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr);
241d8588912SDave May       }
242d8588912SDave May       ierr = PetscFree(vs->m[i]);CHKERRQ(ierr);
243d8588912SDave May     }
244d8588912SDave May     ierr = PetscFree(vs->m);CHKERRQ(ierr);
245d8588912SDave May   }
246bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
247d8588912SDave May 
248bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",0);CHKERRQ(ierr);
249bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",0);CHKERRQ(ierr);
250bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",0);CHKERRQ(ierr);
251bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",0);CHKERRQ(ierr);
252bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",0);CHKERRQ(ierr);
253bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",0);CHKERRQ(ierr);
254bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",0);CHKERRQ(ierr);
255bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",0);CHKERRQ(ierr);
256d8588912SDave May   PetscFunctionReturn(0);
257d8588912SDave May }
258d8588912SDave May 
259d8588912SDave May #undef __FUNCT__
260d8588912SDave May #define __FUNCT__ "MatAssemblyBegin_Nest"
261207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
262d8588912SDave May {
263d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
264d8588912SDave May   PetscInt       i,j;
265d8588912SDave May   PetscErrorCode ierr;
266d8588912SDave May 
267d8588912SDave May   PetscFunctionBegin;
268d8588912SDave May   for (i=0; i<vs->nr; i++) {
269d8588912SDave May     for (j=0; j<vs->nc; j++) {
270e7c19651SJed Brown       if (vs->m[i][j]) {
271e7c19651SJed Brown         ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr);
272e7c19651SJed Brown         if (!vs->splitassembly) {
273e7c19651SJed Brown           /* Note: split assembly will fail if the same block appears more than once (even indirectly through a nested
274e7c19651SJed 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
275e7c19651SJed Brown            * already performing an assembly, but the result would by more complicated and appears to offer less
276e7c19651SJed Brown            * potential for diagnostics and correctness checking. Split assembly should be fixed once there is an
277e7c19651SJed Brown            * interface for libraries to make asynchronous progress in "user-defined non-blocking collectives".
278e7c19651SJed Brown            */
279e7c19651SJed Brown           ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr);
280e7c19651SJed Brown         }
281e7c19651SJed Brown       }
282d8588912SDave May     }
283d8588912SDave May   }
284d8588912SDave May   PetscFunctionReturn(0);
285d8588912SDave May }
286d8588912SDave May 
287d8588912SDave May #undef __FUNCT__
288d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest"
289207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
290d8588912SDave May {
291d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
292d8588912SDave May   PetscInt       i,j;
293d8588912SDave May   PetscErrorCode ierr;
294d8588912SDave May 
295d8588912SDave May   PetscFunctionBegin;
296d8588912SDave May   for (i=0; i<vs->nr; i++) {
297d8588912SDave May     for (j=0; j<vs->nc; j++) {
298e7c19651SJed Brown       if (vs->m[i][j]) {
299e7c19651SJed Brown         if (vs->splitassembly) {
300e7c19651SJed Brown           ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr);
301e7c19651SJed Brown         }
302e7c19651SJed Brown       }
303d8588912SDave May     }
304d8588912SDave May   }
305d8588912SDave May   PetscFunctionReturn(0);
306d8588912SDave May }
307d8588912SDave May 
308d8588912SDave May #undef __FUNCT__
309f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatRow"
310f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
311d8588912SDave May {
312207556f9SJed Brown   PetscErrorCode ierr;
313f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
314f349c1fdSJed Brown   PetscInt       j;
315f349c1fdSJed Brown   Mat            sub;
316d8588912SDave May 
317d8588912SDave May   PetscFunctionBegin;
3180298fd71SBarry Smith   sub = (row < vs->nc) ? vs->m[row][row] : (Mat)NULL; /* Prefer to find on the diagonal */
319f349c1fdSJed Brown   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
3204994cf47SJed Brown   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
321f349c1fdSJed Brown   *B = sub;
322f349c1fdSJed Brown   PetscFunctionReturn(0);
323d8588912SDave May }
324d8588912SDave May 
325f349c1fdSJed Brown #undef __FUNCT__
326f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatCol"
327f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
328f349c1fdSJed Brown {
329207556f9SJed Brown   PetscErrorCode ierr;
330f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
331f349c1fdSJed Brown   PetscInt       i;
332f349c1fdSJed Brown   Mat            sub;
333f349c1fdSJed Brown 
334f349c1fdSJed Brown   PetscFunctionBegin;
3350298fd71SBarry Smith   sub = (col < vs->nr) ? vs->m[col][col] : (Mat)NULL; /* Prefer to find on the diagonal */
336f349c1fdSJed Brown   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
3374994cf47SJed Brown   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
338f349c1fdSJed Brown   *B = sub;
339f349c1fdSJed Brown   PetscFunctionReturn(0);
340d8588912SDave May }
341d8588912SDave May 
342f349c1fdSJed Brown #undef __FUNCT__
343f349c1fdSJed Brown #define __FUNCT__ "MatNestFindIS"
344f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
345f349c1fdSJed Brown {
346f349c1fdSJed Brown   PetscErrorCode ierr;
347f349c1fdSJed Brown   PetscInt       i;
348f349c1fdSJed Brown   PetscBool      flg;
349f349c1fdSJed Brown 
350f349c1fdSJed Brown   PetscFunctionBegin;
351f349c1fdSJed Brown   PetscValidPointer(list,3);
352f349c1fdSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,4);
353f349c1fdSJed Brown   PetscValidIntPointer(found,5);
354f349c1fdSJed Brown   *found = -1;
355f349c1fdSJed Brown   for (i=0; i<n; i++) {
356207556f9SJed Brown     if (!list[i]) continue;
357f349c1fdSJed Brown     ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr);
358f349c1fdSJed Brown     if (flg) {
359f349c1fdSJed Brown       *found = i;
360f349c1fdSJed Brown       PetscFunctionReturn(0);
361f349c1fdSJed Brown     }
362f349c1fdSJed Brown   }
363ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Could not find index set");
364f349c1fdSJed Brown   PetscFunctionReturn(0);
365f349c1fdSJed Brown }
366f349c1fdSJed Brown 
367f349c1fdSJed Brown #undef __FUNCT__
3688188e55aSJed Brown #define __FUNCT__ "MatNestGetRow"
3698188e55aSJed Brown /* Get a block row as a new MatNest */
3708188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B)
3718188e55aSJed Brown {
3728188e55aSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
3738188e55aSJed Brown   char           keyname[256];
3748188e55aSJed Brown   PetscErrorCode ierr;
3758188e55aSJed Brown 
3768188e55aSJed Brown   PetscFunctionBegin;
3770298fd71SBarry Smith   *B   = NULL;
3788caf3d72SBarry Smith   ierr = PetscSNPrintf(keyname,sizeof(keyname),"NestRow_%D",row);CHKERRQ(ierr);
3798188e55aSJed Brown   ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr);
3808188e55aSJed Brown   if (*B) PetscFunctionReturn(0);
3818188e55aSJed Brown 
382ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),1,NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr);
38326fbe8dcSKarl Rupp 
3848188e55aSJed Brown   (*B)->assembled = A->assembled;
38526fbe8dcSKarl Rupp 
3868188e55aSJed Brown   ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr);
3878188e55aSJed Brown   ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */
3888188e55aSJed Brown   PetscFunctionReturn(0);
3898188e55aSJed Brown }
3908188e55aSJed Brown 
3918188e55aSJed Brown #undef __FUNCT__
392f349c1fdSJed Brown #define __FUNCT__ "MatNestFindSubMat"
393f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
394f349c1fdSJed Brown {
395f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
3968188e55aSJed Brown   PetscErrorCode ierr;
3976b3a5b13SJed Brown   PetscInt       row,col;
398e072481dSJed Brown   PetscBool      same,isFullCol,isFullColGlobal;
399f349c1fdSJed Brown 
400f349c1fdSJed Brown   PetscFunctionBegin;
4018188e55aSJed Brown   /* Check if full column space. This is a hack */
4028188e55aSJed Brown   isFullCol = PETSC_FALSE;
403251f4c67SDmitry Karpeev   ierr      = PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr);
4048188e55aSJed Brown   if (same) {
40577019fcaSJed Brown     PetscInt n,first,step,i,an,am,afirst,astep;
4068188e55aSJed Brown     ierr      = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
4078188e55aSJed Brown     ierr      = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
40877019fcaSJed Brown     isFullCol = PETSC_TRUE;
40905ce4453SJed Brown     for (i=0,an=A->cmap->rstart; i<vs->nc; i++) {
41077019fcaSJed Brown       ierr = ISStrideGetInfo(is->col[i],&afirst,&astep);CHKERRQ(ierr);
41177019fcaSJed Brown       ierr = ISGetLocalSize(is->col[i],&am);CHKERRQ(ierr);
41277019fcaSJed Brown       if (afirst != an || astep != step) isFullCol = PETSC_FALSE;
41377019fcaSJed Brown       an += am;
41477019fcaSJed Brown     }
41505ce4453SJed Brown     if (an != A->cmap->rstart+n) isFullCol = PETSC_FALSE;
4168188e55aSJed Brown   }
417b2566f29SBarry Smith   ierr = MPIU_Allreduce(&isFullCol,&isFullColGlobal,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)iscol));CHKERRQ(ierr);
4188188e55aSJed Brown 
419e072481dSJed Brown   if (isFullColGlobal) {
4208188e55aSJed Brown     PetscInt row;
4218188e55aSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
4228188e55aSJed Brown     ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr);
4238188e55aSJed Brown   } else {
424f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
425f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr);
426b6480e04SStefano Zampini     if (!vs->m[row][col]) {
427b6480e04SStefano Zampini       PetscInt lr,lc;
428b6480e04SStefano Zampini 
429b6480e04SStefano Zampini       ierr = MatCreate(PetscObjectComm((PetscObject)A),&vs->m[row][col]);CHKERRQ(ierr);
430b6480e04SStefano Zampini       ierr = ISGetLocalSize(vs->isglobal.row[row],&lr);CHKERRQ(ierr);
431b6480e04SStefano Zampini       ierr = ISGetLocalSize(vs->isglobal.col[col],&lc);CHKERRQ(ierr);
432b6480e04SStefano Zampini       ierr = MatSetSizes(vs->m[row][col],lr,lc,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
433b6480e04SStefano Zampini       ierr = MatSetUp(vs->m[row][col]);CHKERRQ(ierr);
434b6480e04SStefano Zampini       ierr = MatAssemblyBegin(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
435b6480e04SStefano Zampini       ierr = MatAssemblyEnd(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
436b6480e04SStefano Zampini     }
437f349c1fdSJed Brown     *B = vs->m[row][col];
4388188e55aSJed Brown   }
439f349c1fdSJed Brown   PetscFunctionReturn(0);
440f349c1fdSJed Brown }
441f349c1fdSJed Brown 
442f349c1fdSJed Brown #undef __FUNCT__
443f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest"
444f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
445f349c1fdSJed Brown {
446f349c1fdSJed Brown   PetscErrorCode ierr;
447f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
448f349c1fdSJed Brown   Mat            sub;
449f349c1fdSJed Brown 
450f349c1fdSJed Brown   PetscFunctionBegin;
451f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr);
452f349c1fdSJed Brown   switch (reuse) {
453f349c1fdSJed Brown   case MAT_INITIAL_MATRIX:
4547874fa86SDave May     if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); }
455f349c1fdSJed Brown     *B = sub;
456f349c1fdSJed Brown     break;
457f349c1fdSJed Brown   case MAT_REUSE_MATRIX:
458ce94432eSBarry Smith     if (sub != *B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
459f349c1fdSJed Brown     break;
460f349c1fdSJed Brown   case MAT_IGNORE_MATRIX:       /* Nothing to do */
461f349c1fdSJed Brown     break;
462511c6705SHong Zhang   case MAT_INPLACE_MATRIX:       /* Nothing to do */
463511c6705SHong Zhang     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX is not supported yet");
464511c6705SHong Zhang     break;
465f349c1fdSJed Brown   }
466f349c1fdSJed Brown   PetscFunctionReturn(0);
467f349c1fdSJed Brown }
468f349c1fdSJed Brown 
469f349c1fdSJed Brown #undef __FUNCT__
470f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest"
471f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
472f349c1fdSJed Brown {
473f349c1fdSJed Brown   PetscErrorCode ierr;
474f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
475f349c1fdSJed Brown   Mat            sub;
476f349c1fdSJed Brown 
477f349c1fdSJed Brown   PetscFunctionBegin;
478f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
479f349c1fdSJed Brown   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
480f349c1fdSJed Brown   if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);}
481f349c1fdSJed Brown   *B = sub;
482d8588912SDave May   PetscFunctionReturn(0);
483d8588912SDave May }
484d8588912SDave May 
485d8588912SDave May #undef __FUNCT__
486d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest"
487207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
488d8588912SDave May {
489d8588912SDave May   PetscErrorCode ierr;
490f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
491f349c1fdSJed Brown   Mat            sub;
492d8588912SDave May 
493d8588912SDave May   PetscFunctionBegin;
494f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
495ce94432eSBarry Smith   if (*B != sub) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
496f349c1fdSJed Brown   if (sub) {
497ce94432eSBarry Smith     if (((PetscObject)sub)->refct <= 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
4986bf464f9SBarry Smith     ierr = MatDestroy(B);CHKERRQ(ierr);
499d8588912SDave May   }
500d8588912SDave May   PetscFunctionReturn(0);
501d8588912SDave May }
502d8588912SDave May 
503d8588912SDave May #undef __FUNCT__
5047874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest"
5057874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
5067874fa86SDave May {
5077874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
5087874fa86SDave May   PetscInt       i;
5097874fa86SDave May   PetscErrorCode ierr;
5107874fa86SDave May 
5117874fa86SDave May   PetscFunctionBegin;
5127874fa86SDave May   for (i=0; i<bA->nr; i++) {
513429bac76SJed Brown     Vec bv;
514429bac76SJed Brown     ierr = VecGetSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
5157874fa86SDave May     if (bA->m[i][i]) {
516429bac76SJed Brown       ierr = MatGetDiagonal(bA->m[i][i],bv);CHKERRQ(ierr);
5177874fa86SDave May     } else {
5185159a857SMatthew G. Knepley       ierr = VecSet(bv,0.0);CHKERRQ(ierr);
5197874fa86SDave May     }
520429bac76SJed Brown     ierr = VecRestoreSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
5217874fa86SDave May   }
5227874fa86SDave May   PetscFunctionReturn(0);
5237874fa86SDave May }
5247874fa86SDave May 
5257874fa86SDave May #undef __FUNCT__
5267874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest"
5277874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
5287874fa86SDave May {
5297874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
530429bac76SJed Brown   Vec            bl,*br;
5317874fa86SDave May   PetscInt       i,j;
5327874fa86SDave May   PetscErrorCode ierr;
5337874fa86SDave May 
5347874fa86SDave May   PetscFunctionBegin;
5353f800ebeSJed Brown   ierr = PetscCalloc1(bA->nc,&br);CHKERRQ(ierr);
5362e6472ebSElliott Sales de Andrade   if (r) {
537429bac76SJed Brown     for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
5382e6472ebSElliott Sales de Andrade   }
5392e6472ebSElliott Sales de Andrade   bl = NULL;
5407874fa86SDave May   for (i=0; i<bA->nr; i++) {
5412e6472ebSElliott Sales de Andrade     if (l) {
542429bac76SJed Brown       ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
5432e6472ebSElliott Sales de Andrade     }
5447874fa86SDave May     for (j=0; j<bA->nc; j++) {
5457874fa86SDave May       if (bA->m[i][j]) {
546429bac76SJed Brown         ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr);
5477874fa86SDave May       }
5487874fa86SDave May     }
5492e6472ebSElliott Sales de Andrade     if (l) {
550a061e289SJed Brown       ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
5517874fa86SDave May     }
5522e6472ebSElliott Sales de Andrade   }
5532e6472ebSElliott Sales de Andrade   if (r) {
554429bac76SJed Brown     for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
5552e6472ebSElliott Sales de Andrade   }
556429bac76SJed Brown   ierr = PetscFree(br);CHKERRQ(ierr);
5577874fa86SDave May   PetscFunctionReturn(0);
5587874fa86SDave May }
5597874fa86SDave May 
5607874fa86SDave May #undef __FUNCT__
561a061e289SJed Brown #define __FUNCT__ "MatScale_Nest"
562a061e289SJed Brown static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a)
563a061e289SJed Brown {
564a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
565a061e289SJed Brown   PetscInt       i,j;
566a061e289SJed Brown   PetscErrorCode ierr;
567a061e289SJed Brown 
568a061e289SJed Brown   PetscFunctionBegin;
569a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
570a061e289SJed Brown     for (j=0; j<bA->nc; j++) {
571a061e289SJed Brown       if (bA->m[i][j]) {
572a061e289SJed Brown         ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr);
573a061e289SJed Brown       }
574a061e289SJed Brown     }
575a061e289SJed Brown   }
576a061e289SJed Brown   PetscFunctionReturn(0);
577a061e289SJed Brown }
578a061e289SJed Brown 
579a061e289SJed Brown #undef __FUNCT__
580a061e289SJed Brown #define __FUNCT__ "MatShift_Nest"
581a061e289SJed Brown static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a)
582a061e289SJed Brown {
583a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
584a061e289SJed Brown   PetscInt       i;
585a061e289SJed Brown   PetscErrorCode ierr;
586a061e289SJed Brown 
587a061e289SJed Brown   PetscFunctionBegin;
588a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
589ce94432eSBarry 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);
590a061e289SJed Brown     ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr);
591a061e289SJed Brown   }
592a061e289SJed Brown   PetscFunctionReturn(0);
593a061e289SJed Brown }
594a061e289SJed Brown 
595a061e289SJed Brown #undef __FUNCT__
59613135bc6SAlex Fikl #define __FUNCT__ "MatDiagonalSet_Nest"
59713135bc6SAlex Fikl static PetscErrorCode MatDiagonalSet_Nest(Mat A,Vec D,InsertMode is)
59813135bc6SAlex Fikl {
59913135bc6SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data;
60013135bc6SAlex Fikl   PetscInt       i;
60113135bc6SAlex Fikl   PetscErrorCode ierr;
60213135bc6SAlex Fikl 
60313135bc6SAlex Fikl   PetscFunctionBegin;
60413135bc6SAlex Fikl   for (i=0; i<bA->nr; i++) {
60513135bc6SAlex Fikl     Vec bv;
60613135bc6SAlex Fikl     ierr = VecGetSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
60713135bc6SAlex Fikl     if (bA->m[i][i]) {
60813135bc6SAlex Fikl       ierr = MatDiagonalSet(bA->m[i][i],bv,is);CHKERRQ(ierr);
60913135bc6SAlex Fikl     }
61013135bc6SAlex Fikl     ierr = VecRestoreSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
61113135bc6SAlex Fikl   }
61213135bc6SAlex Fikl   PetscFunctionReturn(0);
61313135bc6SAlex Fikl }
61413135bc6SAlex Fikl 
61513135bc6SAlex Fikl #undef __FUNCT__
616f8170845SAlex Fikl #define __FUNCT__ "MatSetRandom_Nest"
617f8170845SAlex Fikl static PetscErrorCode MatSetRandom_Nest(Mat A,PetscRandom rctx)
618f8170845SAlex Fikl {
619f8170845SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data;
620f8170845SAlex Fikl   PetscInt       i,j;
621f8170845SAlex Fikl   PetscErrorCode ierr;
622f8170845SAlex Fikl 
623f8170845SAlex Fikl   PetscFunctionBegin;
624f8170845SAlex Fikl   for (i=0; i<bA->nr; i++) {
625f8170845SAlex Fikl     for (j=0; j<bA->nc; j++) {
626f8170845SAlex Fikl       if (bA->m[i][j]) {
627f8170845SAlex Fikl         ierr = MatSetRandom(bA->m[i][j],rctx);CHKERRQ(ierr);
628f8170845SAlex Fikl       }
629f8170845SAlex Fikl     }
630f8170845SAlex Fikl   }
631f8170845SAlex Fikl   PetscFunctionReturn(0);
632f8170845SAlex Fikl }
633f8170845SAlex Fikl 
634f8170845SAlex Fikl #undef __FUNCT__
6352a7a6963SBarry Smith #define __FUNCT__ "MatCreateVecs_Nest"
6362a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left)
637d8588912SDave May {
638d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
639d8588912SDave May   Vec            *L,*R;
640d8588912SDave May   MPI_Comm       comm;
641d8588912SDave May   PetscInt       i,j;
642d8588912SDave May   PetscErrorCode ierr;
643d8588912SDave May 
644d8588912SDave May   PetscFunctionBegin;
645ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
646d8588912SDave May   if (right) {
647d8588912SDave May     /* allocate R */
648854ce69bSBarry Smith     ierr = PetscMalloc1(bA->nc, &R);CHKERRQ(ierr);
649d8588912SDave May     /* Create the right vectors */
650d8588912SDave May     for (j=0; j<bA->nc; j++) {
651d8588912SDave May       for (i=0; i<bA->nr; i++) {
652d8588912SDave May         if (bA->m[i][j]) {
6532a7a6963SBarry Smith           ierr = MatCreateVecs(bA->m[i][j],&R[j],NULL);CHKERRQ(ierr);
654d8588912SDave May           break;
655d8588912SDave May         }
656d8588912SDave May       }
6576c4ed002SBarry Smith       if (i==bA->nr) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
658d8588912SDave May     }
659f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
660d8588912SDave May     /* hand back control to the nest vector */
661d8588912SDave May     for (j=0; j<bA->nc; j++) {
6626bf464f9SBarry Smith       ierr = VecDestroy(&R[j]);CHKERRQ(ierr);
663d8588912SDave May     }
664d8588912SDave May     ierr = PetscFree(R);CHKERRQ(ierr);
665d8588912SDave May   }
666d8588912SDave May 
667d8588912SDave May   if (left) {
668d8588912SDave May     /* allocate L */
669854ce69bSBarry Smith     ierr = PetscMalloc1(bA->nr, &L);CHKERRQ(ierr);
670d8588912SDave May     /* Create the left vectors */
671d8588912SDave May     for (i=0; i<bA->nr; i++) {
672d8588912SDave May       for (j=0; j<bA->nc; j++) {
673d8588912SDave May         if (bA->m[i][j]) {
6742a7a6963SBarry Smith           ierr = MatCreateVecs(bA->m[i][j],NULL,&L[i]);CHKERRQ(ierr);
675d8588912SDave May           break;
676d8588912SDave May         }
677d8588912SDave May       }
6786c4ed002SBarry Smith       if (j==bA->nc) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
679d8588912SDave May     }
680d8588912SDave May 
681f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
682d8588912SDave May     for (i=0; i<bA->nr; i++) {
6836bf464f9SBarry Smith       ierr = VecDestroy(&L[i]);CHKERRQ(ierr);
684d8588912SDave May     }
685d8588912SDave May 
686d8588912SDave May     ierr = PetscFree(L);CHKERRQ(ierr);
687d8588912SDave May   }
688d8588912SDave May   PetscFunctionReturn(0);
689d8588912SDave May }
690d8588912SDave May 
691d8588912SDave May #undef __FUNCT__
692d8588912SDave May #define __FUNCT__ "MatView_Nest"
693207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
694d8588912SDave May {
695d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
696d8588912SDave May   PetscBool      isascii;
697d8588912SDave May   PetscInt       i,j;
698d8588912SDave May   PetscErrorCode ierr;
699d8588912SDave May 
700d8588912SDave May   PetscFunctionBegin;
701251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
702d8588912SDave May   if (isascii) {
703d8588912SDave May 
704*d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Matrix object: \n");CHKERRQ(ierr);
705*d86155a6SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
706*d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer, "type=nest, rows=%D, cols=%D \n",bA->nr,bA->nc);CHKERRQ(ierr);
707d8588912SDave May 
708*d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"MatNest structure: \n");CHKERRQ(ierr);
709d8588912SDave May     for (i=0; i<bA->nr; i++) {
710d8588912SDave May       for (j=0; j<bA->nc; j++) {
71119fd82e9SBarry Smith         MatType   type;
712270f95d7SJed Brown         char      name[256] = "",prefix[256] = "";
713d8588912SDave May         PetscInt  NR,NC;
714d8588912SDave May         PetscBool isNest = PETSC_FALSE;
715d8588912SDave May 
716d8588912SDave May         if (!bA->m[i][j]) {
717*d86155a6SBarry Smith           CHKERRQ(ierr);PetscViewerASCIIPrintf(viewer, "(%D,%D) : NULL \n",i,j);CHKERRQ(ierr);
718d8588912SDave May           continue;
719d8588912SDave May         }
720d8588912SDave May         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
721d8588912SDave May         ierr = MatGetType(bA->m[i][j], &type);CHKERRQ(ierr);
7228caf3d72SBarry Smith         if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);}
7238caf3d72SBarry Smith         if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);}
724251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
725d8588912SDave May 
726270f95d7SJed Brown         ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr);
727d8588912SDave May 
728d8588912SDave May         if (isNest) {
729270f95d7SJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);  /* push1 */
730d8588912SDave May           ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr);
731270f95d7SJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop1 */
732d8588912SDave May         }
733d8588912SDave May       }
734d8588912SDave May     }
735*d86155a6SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop0 */
736d8588912SDave May   }
737d8588912SDave May   PetscFunctionReturn(0);
738d8588912SDave May }
739d8588912SDave May 
740d8588912SDave May #undef __FUNCT__
741d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest"
742207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A)
743d8588912SDave May {
744d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
745d8588912SDave May   PetscInt       i,j;
746d8588912SDave May   PetscErrorCode ierr;
747d8588912SDave May 
748d8588912SDave May   PetscFunctionBegin;
749d8588912SDave May   for (i=0; i<bA->nr; i++) {
750d8588912SDave May     for (j=0; j<bA->nc; j++) {
751d8588912SDave May       if (!bA->m[i][j]) continue;
752d8588912SDave May       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
753d8588912SDave May     }
754d8588912SDave May   }
755d8588912SDave May   PetscFunctionReturn(0);
756d8588912SDave May }
757d8588912SDave May 
758d8588912SDave May #undef __FUNCT__
759c222c20dSDavid Ham #define __FUNCT__ "MatCopy_Nest"
760c222c20dSDavid Ham static PetscErrorCode MatCopy_Nest(Mat A,Mat B,MatStructure str)
761c222c20dSDavid Ham {
762c222c20dSDavid Ham   Mat_Nest       *bA = (Mat_Nest*)A->data,*bB = (Mat_Nest*)B->data;
763c222c20dSDavid Ham   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
764c222c20dSDavid Ham   PetscErrorCode ierr;
765c222c20dSDavid Ham 
766c222c20dSDavid Ham   PetscFunctionBegin;
767c222c20dSDavid 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);
768c222c20dSDavid Ham   for (i=0; i<nr; i++) {
769c222c20dSDavid Ham     for (j=0; j<nc; j++) {
77046a2b97cSJed Brown       if (bA->m[i][j] && bB->m[i][j]) {
771c222c20dSDavid Ham         ierr = MatCopy(bA->m[i][j],bB->m[i][j],str);CHKERRQ(ierr);
77246a2b97cSJed 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);
773c222c20dSDavid Ham     }
774c222c20dSDavid Ham   }
775c222c20dSDavid Ham   PetscFunctionReturn(0);
776c222c20dSDavid Ham }
777c222c20dSDavid Ham 
778c222c20dSDavid Ham #undef __FUNCT__
779d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest"
780207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
781d8588912SDave May {
782d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
783841e96a3SJed Brown   Mat            *b;
784841e96a3SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
785d8588912SDave May   PetscErrorCode ierr;
786d8588912SDave May 
787d8588912SDave May   PetscFunctionBegin;
788785e854fSJed Brown   ierr = PetscMalloc1(nr*nc,&b);CHKERRQ(ierr);
789841e96a3SJed Brown   for (i=0; i<nr; i++) {
790841e96a3SJed Brown     for (j=0; j<nc; j++) {
791841e96a3SJed Brown       if (bA->m[i][j]) {
792841e96a3SJed Brown         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
793841e96a3SJed Brown       } else {
7940298fd71SBarry Smith         b[i*nc+j] = NULL;
795d8588912SDave May       }
796d8588912SDave May     }
797d8588912SDave May   }
798ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
799841e96a3SJed Brown   /* Give the new MatNest exclusive ownership */
800841e96a3SJed Brown   for (i=0; i<nr*nc; i++) {
8016bf464f9SBarry Smith     ierr = MatDestroy(&b[i]);CHKERRQ(ierr);
802d8588912SDave May   }
803d8588912SDave May   ierr = PetscFree(b);CHKERRQ(ierr);
804d8588912SDave May 
805841e96a3SJed Brown   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
806841e96a3SJed Brown   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
807d8588912SDave May   PetscFunctionReturn(0);
808d8588912SDave May }
809d8588912SDave May 
810d8588912SDave May /* nest api */
811d8588912SDave May #undef __FUNCT__
812d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest"
813d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
814d8588912SDave May {
815d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
8165fd66863SKarl Rupp 
817d8588912SDave May   PetscFunctionBegin;
818ce94432eSBarry Smith   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
819ce94432eSBarry Smith   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
820d8588912SDave May   *mat = bA->m[idxm][jdxm];
821d8588912SDave May   PetscFunctionReturn(0);
822d8588912SDave May }
823d8588912SDave May 
824d8588912SDave May #undef __FUNCT__
825d8588912SDave May #define __FUNCT__ "MatNestGetSubMat"
8269ba0d327SJed Brown /*@
827d8588912SDave May  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
828d8588912SDave May 
829d8588912SDave May  Not collective
830d8588912SDave May 
831d8588912SDave May  Input Parameters:
832629881c0SJed Brown +   A  - nest matrix
833d8588912SDave May .   idxm - index of the matrix within the nest matrix
834629881c0SJed Brown -   jdxm - index of the matrix within the nest matrix
835d8588912SDave May 
836d8588912SDave May  Output Parameter:
837d8588912SDave May .   sub - matrix at index idxm,jdxm within the nest matrix
838d8588912SDave May 
839d8588912SDave May  Level: developer
840d8588912SDave May 
841d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMats()
842d8588912SDave May @*/
8437087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
844d8588912SDave May {
845699a902aSJed Brown   PetscErrorCode ierr;
846d8588912SDave May 
847d8588912SDave May   PetscFunctionBegin;
848699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
849d8588912SDave May   PetscFunctionReturn(0);
850d8588912SDave May }
851d8588912SDave May 
852d8588912SDave May #undef __FUNCT__
8530782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat_Nest"
8540782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat)
8550782ca92SJed Brown {
8560782ca92SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
8570782ca92SJed Brown   PetscInt       m,n,M,N,mi,ni,Mi,Ni;
8580782ca92SJed Brown   PetscErrorCode ierr;
8590782ca92SJed Brown 
8600782ca92SJed Brown   PetscFunctionBegin;
861ce94432eSBarry Smith   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
862ce94432eSBarry Smith   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
8630782ca92SJed Brown   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
8640782ca92SJed Brown   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
8650782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr);
8660782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr);
8670782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr);
8680782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr);
869ce94432eSBarry 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);
870ce94432eSBarry 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);
87126fbe8dcSKarl Rupp 
8720782ca92SJed Brown   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
8730782ca92SJed Brown   ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr);
8740782ca92SJed Brown   bA->m[idxm][jdxm] = mat;
8750782ca92SJed Brown   PetscFunctionReturn(0);
8760782ca92SJed Brown }
8770782ca92SJed Brown 
8780782ca92SJed Brown #undef __FUNCT__
8790782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat"
8809ba0d327SJed Brown /*@
8810782ca92SJed Brown  MatNestSetSubMat - Set a single submatrix in the nest matrix.
8820782ca92SJed Brown 
8830782ca92SJed Brown  Logically collective on the submatrix communicator
8840782ca92SJed Brown 
8850782ca92SJed Brown  Input Parameters:
8860782ca92SJed Brown +   A  - nest matrix
8870782ca92SJed Brown .   idxm - index of the matrix within the nest matrix
8880782ca92SJed Brown .   jdxm - index of the matrix within the nest matrix
8890782ca92SJed Brown -   sub - matrix at index idxm,jdxm within the nest matrix
8900782ca92SJed Brown 
8910782ca92SJed Brown  Notes:
8920782ca92SJed Brown  The new submatrix must have the same size and communicator as that block of the nest.
8930782ca92SJed Brown 
8940782ca92SJed Brown  This increments the reference count of the submatrix.
8950782ca92SJed Brown 
8960782ca92SJed Brown  Level: developer
8970782ca92SJed Brown 
8980782ca92SJed Brown .seealso: MatNestSetSubMats(), MatNestGetSubMat()
8990782ca92SJed Brown @*/
9000782ca92SJed Brown PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
9010782ca92SJed Brown {
9020782ca92SJed Brown   PetscErrorCode ierr;
9030782ca92SJed Brown 
9040782ca92SJed Brown   PetscFunctionBegin;
9050782ca92SJed Brown   ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr);
9060782ca92SJed Brown   PetscFunctionReturn(0);
9070782ca92SJed Brown }
9080782ca92SJed Brown 
9090782ca92SJed Brown #undef __FUNCT__
910d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest"
911d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
912d8588912SDave May {
913d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
9145fd66863SKarl Rupp 
915d8588912SDave May   PetscFunctionBegin;
91626fbe8dcSKarl Rupp   if (M)   *M   = bA->nr;
91726fbe8dcSKarl Rupp   if (N)   *N   = bA->nc;
91826fbe8dcSKarl Rupp   if (mat) *mat = bA->m;
919d8588912SDave May   PetscFunctionReturn(0);
920d8588912SDave May }
921d8588912SDave May 
922d8588912SDave May #undef __FUNCT__
923d8588912SDave May #define __FUNCT__ "MatNestGetSubMats"
924d8588912SDave May /*@C
925d8588912SDave May  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
926d8588912SDave May 
927d8588912SDave May  Not collective
928d8588912SDave May 
929d8588912SDave May  Input Parameters:
930629881c0SJed Brown .   A  - nest matrix
931d8588912SDave May 
932d8588912SDave May  Output Parameter:
933629881c0SJed Brown +   M - number of rows in the nest matrix
934d8588912SDave May .   N - number of cols in the nest matrix
935629881c0SJed Brown -   mat - 2d array of matrices
936d8588912SDave May 
937d8588912SDave May  Notes:
938d8588912SDave May 
939d8588912SDave May  The user should not free the array mat.
940d8588912SDave May 
941351962e3SVincent Le Chenadec  In Fortran, this routine has a calling sequence
942351962e3SVincent Le Chenadec $   call MatNestGetSubMats(A, M, N, mat, ierr)
943351962e3SVincent Le Chenadec  where the space allocated for the optional argument mat is assumed large enough (if provided).
944351962e3SVincent Le Chenadec 
945d8588912SDave May  Level: developer
946d8588912SDave May 
947d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMat()
948d8588912SDave May @*/
9497087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
950d8588912SDave May {
951699a902aSJed Brown   PetscErrorCode ierr;
952d8588912SDave May 
953d8588912SDave May   PetscFunctionBegin;
954699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
955d8588912SDave May   PetscFunctionReturn(0);
956d8588912SDave May }
957d8588912SDave May 
958d8588912SDave May #undef __FUNCT__
959d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest"
9607087cfbeSBarry Smith PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
961d8588912SDave May {
962d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
963d8588912SDave May 
964d8588912SDave May   PetscFunctionBegin;
96526fbe8dcSKarl Rupp   if (M) *M = bA->nr;
96626fbe8dcSKarl Rupp   if (N) *N = bA->nc;
967d8588912SDave May   PetscFunctionReturn(0);
968d8588912SDave May }
969d8588912SDave May 
970d8588912SDave May #undef __FUNCT__
971d8588912SDave May #define __FUNCT__ "MatNestGetSize"
9729ba0d327SJed Brown /*@
973d8588912SDave May  MatNestGetSize - Returns the size of the nest matrix.
974d8588912SDave May 
975d8588912SDave May  Not collective
976d8588912SDave May 
977d8588912SDave May  Input Parameters:
978d8588912SDave May .   A  - nest matrix
979d8588912SDave May 
980d8588912SDave May  Output Parameter:
981629881c0SJed Brown +   M - number of rows in the nested mat
982629881c0SJed Brown -   N - number of cols in the nested mat
983d8588912SDave May 
984d8588912SDave May  Notes:
985d8588912SDave May 
986d8588912SDave May  Level: developer
987d8588912SDave May 
988d8588912SDave May .seealso: MatNestGetSubMat(), MatNestGetSubMats()
989d8588912SDave May @*/
9907087cfbeSBarry Smith PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
991d8588912SDave May {
992699a902aSJed Brown   PetscErrorCode ierr;
993d8588912SDave May 
994d8588912SDave May   PetscFunctionBegin;
995699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
996d8588912SDave May   PetscFunctionReturn(0);
997d8588912SDave May }
998d8588912SDave May 
999900e7ff2SJed Brown #undef __FUNCT__
1000900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs_Nest"
1001f7a08781SBarry Smith static PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[])
1002900e7ff2SJed Brown {
1003900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest*)A->data;
1004900e7ff2SJed Brown   PetscInt i;
1005900e7ff2SJed Brown 
1006900e7ff2SJed Brown   PetscFunctionBegin;
1007900e7ff2SJed Brown   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i];
1008900e7ff2SJed Brown   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i];
1009900e7ff2SJed Brown   PetscFunctionReturn(0);
1010900e7ff2SJed Brown }
1011900e7ff2SJed Brown 
1012900e7ff2SJed Brown #undef __FUNCT__
1013900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs"
10143a4d7b9aSSatish Balay /*@C
1015900e7ff2SJed Brown  MatNestGetISs - Returns the index sets partitioning the row and column spaces
1016900e7ff2SJed Brown 
1017900e7ff2SJed Brown  Not collective
1018900e7ff2SJed Brown 
1019900e7ff2SJed Brown  Input Parameters:
1020900e7ff2SJed Brown .   A  - nest matrix
1021900e7ff2SJed Brown 
1022900e7ff2SJed Brown  Output Parameter:
1023900e7ff2SJed Brown +   rows - array of row index sets
1024900e7ff2SJed Brown -   cols - array of column index sets
1025900e7ff2SJed Brown 
1026900e7ff2SJed Brown  Level: advanced
1027900e7ff2SJed Brown 
1028900e7ff2SJed Brown  Notes:
1029900e7ff2SJed Brown  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
1030900e7ff2SJed Brown 
1031900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs()
1032900e7ff2SJed Brown @*/
1033900e7ff2SJed Brown PetscErrorCode  MatNestGetISs(Mat A,IS rows[],IS cols[])
1034900e7ff2SJed Brown {
1035900e7ff2SJed Brown   PetscErrorCode ierr;
1036900e7ff2SJed Brown 
1037900e7ff2SJed Brown   PetscFunctionBegin;
1038900e7ff2SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1039900e7ff2SJed Brown   ierr = PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
1040900e7ff2SJed Brown   PetscFunctionReturn(0);
1041900e7ff2SJed Brown }
1042900e7ff2SJed Brown 
1043900e7ff2SJed Brown #undef __FUNCT__
1044900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs_Nest"
1045f7a08781SBarry Smith static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[])
1046900e7ff2SJed Brown {
1047900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest*)A->data;
1048900e7ff2SJed Brown   PetscInt i;
1049900e7ff2SJed Brown 
1050900e7ff2SJed Brown   PetscFunctionBegin;
1051900e7ff2SJed Brown   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i];
1052900e7ff2SJed Brown   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i];
1053900e7ff2SJed Brown   PetscFunctionReturn(0);
1054900e7ff2SJed Brown }
1055900e7ff2SJed Brown 
1056900e7ff2SJed Brown #undef __FUNCT__
1057900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs"
1058900e7ff2SJed Brown /*@C
1059900e7ff2SJed Brown  MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces
1060900e7ff2SJed Brown 
1061900e7ff2SJed Brown  Not collective
1062900e7ff2SJed Brown 
1063900e7ff2SJed Brown  Input Parameters:
1064900e7ff2SJed Brown .   A  - nest matrix
1065900e7ff2SJed Brown 
1066900e7ff2SJed Brown  Output Parameter:
10670298fd71SBarry Smith +   rows - array of row index sets (or NULL to ignore)
10680298fd71SBarry Smith -   cols - array of column index sets (or NULL to ignore)
1069900e7ff2SJed Brown 
1070900e7ff2SJed Brown  Level: advanced
1071900e7ff2SJed Brown 
1072900e7ff2SJed Brown  Notes:
1073900e7ff2SJed Brown  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
1074900e7ff2SJed Brown 
1075900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs()
1076900e7ff2SJed Brown @*/
1077900e7ff2SJed Brown PetscErrorCode  MatNestGetLocalISs(Mat A,IS rows[],IS cols[])
1078900e7ff2SJed Brown {
1079900e7ff2SJed Brown   PetscErrorCode ierr;
1080900e7ff2SJed Brown 
1081900e7ff2SJed Brown   PetscFunctionBegin;
1082900e7ff2SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1083900e7ff2SJed Brown   ierr = PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
1084900e7ff2SJed Brown   PetscFunctionReturn(0);
1085900e7ff2SJed Brown }
1086900e7ff2SJed Brown 
1087207556f9SJed Brown #undef __FUNCT__
1088207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest"
108919fd82e9SBarry Smith PetscErrorCode  MatNestSetVecType_Nest(Mat A,VecType vtype)
1090207556f9SJed Brown {
1091207556f9SJed Brown   PetscErrorCode ierr;
1092207556f9SJed Brown   PetscBool      flg;
1093207556f9SJed Brown 
1094207556f9SJed Brown   PetscFunctionBegin;
1095207556f9SJed Brown   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
1096207556f9SJed Brown   /* In reality, this only distinguishes VECNEST and "other" */
10972a7a6963SBarry Smith   if (flg) A->ops->getvecs = MatCreateVecs_Nest;
109812b53f24SSatish Balay   else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0;
1099207556f9SJed Brown   PetscFunctionReturn(0);
1100207556f9SJed Brown }
1101207556f9SJed Brown 
1102207556f9SJed Brown #undef __FUNCT__
1103207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType"
1104207556f9SJed Brown /*@C
11052a7a6963SBarry Smith  MatNestSetVecType - Sets the type of Vec returned by MatCreateVecs()
1106207556f9SJed Brown 
1107207556f9SJed Brown  Not collective
1108207556f9SJed Brown 
1109207556f9SJed Brown  Input Parameters:
1110207556f9SJed Brown +  A  - nest matrix
1111207556f9SJed Brown -  vtype - type to use for creating vectors
1112207556f9SJed Brown 
1113207556f9SJed Brown  Notes:
1114207556f9SJed Brown 
1115207556f9SJed Brown  Level: developer
1116207556f9SJed Brown 
11172a7a6963SBarry Smith .seealso: MatCreateVecs()
1118207556f9SJed Brown @*/
111919fd82e9SBarry Smith PetscErrorCode  MatNestSetVecType(Mat A,VecType vtype)
1120207556f9SJed Brown {
1121207556f9SJed Brown   PetscErrorCode ierr;
1122207556f9SJed Brown 
1123207556f9SJed Brown   PetscFunctionBegin;
112419fd82e9SBarry Smith   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));CHKERRQ(ierr);
1125207556f9SJed Brown   PetscFunctionReturn(0);
1126207556f9SJed Brown }
1127207556f9SJed Brown 
1128d8588912SDave May #undef __FUNCT__
1129c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats_Nest"
1130c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1131d8588912SDave May {
1132c8883902SJed Brown   Mat_Nest       *s = (Mat_Nest*)A->data;
1133c8883902SJed Brown   PetscInt       i,j,m,n,M,N;
1134d8588912SDave May   PetscErrorCode ierr;
1135d8588912SDave May 
1136d8588912SDave May   PetscFunctionBegin;
1137c8883902SJed Brown   s->nr = nr;
1138c8883902SJed Brown   s->nc = nc;
1139d8588912SDave May 
1140c8883902SJed Brown   /* Create space for submatrices */
1141854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&s->m);CHKERRQ(ierr);
1142c8883902SJed Brown   for (i=0; i<nr; i++) {
1143854ce69bSBarry Smith     ierr = PetscMalloc1(nc,&s->m[i]);CHKERRQ(ierr);
1144d8588912SDave May   }
1145c8883902SJed Brown   for (i=0; i<nr; i++) {
1146c8883902SJed Brown     for (j=0; j<nc; j++) {
1147c8883902SJed Brown       s->m[i][j] = a[i*nc+j];
1148c8883902SJed Brown       if (a[i*nc+j]) {
1149c8883902SJed Brown         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
1150d8588912SDave May       }
1151d8588912SDave May     }
1152d8588912SDave May   }
1153d8588912SDave May 
11548188e55aSJed Brown   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
1155d8588912SDave May 
1156854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&s->row_len);CHKERRQ(ierr);
1157854ce69bSBarry Smith   ierr = PetscMalloc1(nc,&s->col_len);CHKERRQ(ierr);
1158c8883902SJed Brown   for (i=0; i<nr; i++) s->row_len[i]=-1;
1159c8883902SJed Brown   for (j=0; j<nc; j++) s->col_len[j]=-1;
1160d8588912SDave May 
11618188e55aSJed Brown   ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr);
1162d8588912SDave May 
1163c8883902SJed Brown   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
1164c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
1165c8883902SJed Brown   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
1166c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
1167c8883902SJed Brown 
1168c8883902SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1169c8883902SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1170c8883902SJed Brown 
11711795a4d1SJed Brown   ierr = PetscCalloc2(nr,&s->left,nc,&s->right);CHKERRQ(ierr);
1172d8588912SDave May   PetscFunctionReturn(0);
1173d8588912SDave May }
1174d8588912SDave May 
1175c8883902SJed Brown #undef __FUNCT__
1176c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats"
1177c8883902SJed Brown /*@
1178c8883902SJed Brown    MatNestSetSubMats - Sets the nested submatrices
1179c8883902SJed Brown 
1180c8883902SJed Brown    Collective on Mat
1181c8883902SJed Brown 
1182c8883902SJed Brown    Input Parameter:
1183c8883902SJed Brown +  N - nested matrix
1184c8883902SJed Brown .  nr - number of nested row blocks
11850298fd71SBarry Smith .  is_row - index sets for each nested row block, or NULL to make contiguous
1186c8883902SJed Brown .  nc - number of nested column blocks
11870298fd71SBarry Smith .  is_col - index sets for each nested column block, or NULL to make contiguous
11880298fd71SBarry Smith -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1189c8883902SJed Brown 
1190c8883902SJed Brown    Level: advanced
1191c8883902SJed Brown 
1192c8883902SJed Brown .seealso: MatCreateNest(), MATNEST
1193c8883902SJed Brown @*/
1194c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1195c8883902SJed Brown {
1196c8883902SJed Brown   PetscErrorCode ierr;
1197c8883902SJed Brown   PetscInt       i;
1198c8883902SJed Brown 
1199c8883902SJed Brown   PetscFunctionBegin;
1200c8883902SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1201ce94432eSBarry Smith   if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
1202c8883902SJed Brown   if (nr && is_row) {
1203c8883902SJed Brown     PetscValidPointer(is_row,3);
1204c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
1205c8883902SJed Brown   }
1206ce94432eSBarry Smith   if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
12071664e352SJed Brown   if (nc && is_col) {
1208c8883902SJed Brown     PetscValidPointer(is_col,5);
12099b30a8f6SBarry Smith     for (i=0; i<nc; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
1210c8883902SJed Brown   }
1211c8883902SJed Brown   if (nr*nc) PetscValidPointer(a,6);
1212c8883902SJed 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);
1213c8883902SJed Brown   PetscFunctionReturn(0);
1214c8883902SJed Brown }
1215d8588912SDave May 
121677019fcaSJed Brown #undef __FUNCT__
121777019fcaSJed Brown #define __FUNCT__ "MatNestCreateAggregateL2G_Private"
121845b6f7e9SBarry Smith static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog)
121977019fcaSJed Brown {
122077019fcaSJed Brown   PetscErrorCode ierr;
122177019fcaSJed Brown   PetscBool      flg;
122277019fcaSJed Brown   PetscInt       i,j,m,mi,*ix;
122377019fcaSJed Brown 
122477019fcaSJed Brown   PetscFunctionBegin;
122577019fcaSJed Brown   for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) {
122677019fcaSJed Brown     if (islocal[i]) {
122777019fcaSJed Brown       ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
122877019fcaSJed Brown       flg  = PETSC_TRUE;      /* We found a non-trivial entry */
122977019fcaSJed Brown     } else {
123077019fcaSJed Brown       ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
123177019fcaSJed Brown     }
123277019fcaSJed Brown     m += mi;
123377019fcaSJed Brown   }
123477019fcaSJed Brown   if (flg) {
1235785e854fSJed Brown     ierr = PetscMalloc1(m,&ix);CHKERRQ(ierr);
123677019fcaSJed Brown     for (i=0,n=0; i<n; i++) {
12370298fd71SBarry Smith       ISLocalToGlobalMapping smap = NULL;
123877019fcaSJed Brown       VecScatter             scat;
123977019fcaSJed Brown       IS                     isreq;
124077019fcaSJed Brown       Vec                    lvec,gvec;
12413361c9a7SJed Brown       union {char padding[sizeof(PetscScalar)]; PetscInt integer;} *x;
124277019fcaSJed Brown       Mat sub;
124377019fcaSJed Brown 
1244ce94432eSBarry Smith       if (sizeof(*x) != sizeof(PetscScalar)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support when scalars smaller than integers");
124577019fcaSJed Brown       if (colflg) {
124677019fcaSJed Brown         ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
124777019fcaSJed Brown       } else {
124877019fcaSJed Brown         ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
124977019fcaSJed Brown       }
12500298fd71SBarry Smith       if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&smap,NULL);CHKERRQ(ierr);}
125177019fcaSJed Brown       if (islocal[i]) {
125277019fcaSJed Brown         ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
125377019fcaSJed Brown       } else {
125477019fcaSJed Brown         ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
125577019fcaSJed Brown       }
125677019fcaSJed Brown       for (j=0; j<mi; j++) ix[m+j] = j;
125777019fcaSJed Brown       if (smap) {ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr);}
125877019fcaSJed Brown       /*
125977019fcaSJed Brown         Now we need to extract the monolithic global indices that correspond to the given split global indices.
126077019fcaSJed 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.
126177019fcaSJed Brown         The approach here is ugly because it uses VecScatter to move indices.
126277019fcaSJed Brown        */
126377019fcaSJed Brown       ierr = VecCreateSeq(PETSC_COMM_SELF,mi,&lvec);CHKERRQ(ierr);
126477019fcaSJed Brown       ierr = VecCreateMPI(((PetscObject)isglobal[i])->comm,mi,PETSC_DECIDE,&gvec);CHKERRQ(ierr);
126577019fcaSJed Brown       ierr = ISCreateGeneral(((PetscObject)isglobal[i])->comm,mi,ix+m,PETSC_COPY_VALUES,&isreq);CHKERRQ(ierr);
12660298fd71SBarry Smith       ierr = VecScatterCreate(gvec,isreq,lvec,NULL,&scat);CHKERRQ(ierr);
126777019fcaSJed Brown       ierr = VecGetArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
126877019fcaSJed Brown       for (j=0; j<mi; j++) x[j].integer = ix[m+j];
126977019fcaSJed Brown       ierr = VecRestoreArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
127077019fcaSJed Brown       ierr = VecScatterBegin(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
127177019fcaSJed Brown       ierr = VecScatterEnd(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
127277019fcaSJed Brown       ierr = VecGetArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
127377019fcaSJed Brown       for (j=0; j<mi; j++) ix[m+j] = x[j].integer;
127477019fcaSJed Brown       ierr = VecRestoreArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
127577019fcaSJed Brown       ierr = VecDestroy(&lvec);CHKERRQ(ierr);
127677019fcaSJed Brown       ierr = VecDestroy(&gvec);CHKERRQ(ierr);
127777019fcaSJed Brown       ierr = ISDestroy(&isreq);CHKERRQ(ierr);
127877019fcaSJed Brown       ierr = VecScatterDestroy(&scat);CHKERRQ(ierr);
127977019fcaSJed Brown       m   += mi;
128077019fcaSJed Brown     }
1281f0413b6fSBarry Smith     ierr   = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),1,m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
128277019fcaSJed Brown   } else {
12830298fd71SBarry Smith     *ltog  = NULL;
128477019fcaSJed Brown   }
128577019fcaSJed Brown   PetscFunctionReturn(0);
128677019fcaSJed Brown }
128777019fcaSJed Brown 
128877019fcaSJed Brown 
1289d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1290d8588912SDave May /*
1291d8588912SDave May   nprocessors = NP
1292d8588912SDave May   Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1))
1293d8588912SDave May        proc 0: => (g_0,h_0,)
1294d8588912SDave May        proc 1: => (g_1,h_1,)
1295d8588912SDave May        ...
1296d8588912SDave May        proc nprocs-1: => (g_NP-1,h_NP-1,)
1297d8588912SDave May 
1298d8588912SDave May             proc 0:                      proc 1:                    proc nprocs-1:
1299d8588912SDave May     is[0] = (0,1,2,...,nlocal(g_0)-1)  (0,1,...,nlocal(g_1)-1)  (0,1,...,nlocal(g_NP-1))
1300d8588912SDave May 
1301d8588912SDave May             proc 0:
1302d8588912SDave May     is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1)
1303d8588912SDave May             proc 1:
1304d8588912SDave May     is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1)
1305d8588912SDave May 
1306d8588912SDave May             proc NP-1:
1307d8588912SDave May     is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1)
1308d8588912SDave May */
1309d8588912SDave May #undef __FUNCT__
1310d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private"
1311841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1312d8588912SDave May {
1313e2d7f03fSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
13148188e55aSJed Brown   PetscInt       i,j,offset,n,nsum,bs;
1315d8588912SDave May   PetscErrorCode ierr;
13160298fd71SBarry Smith   Mat            sub = NULL;
1317d8588912SDave May 
1318d8588912SDave May   PetscFunctionBegin;
1319854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&vs->isglobal.row);CHKERRQ(ierr);
1320854ce69bSBarry Smith   ierr = PetscMalloc1(nc,&vs->isglobal.col);CHKERRQ(ierr);
1321d8588912SDave May   if (is_row) { /* valid IS is passed in */
1322d8588912SDave May     /* refs on is[] are incremeneted */
1323e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1324d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
132526fbe8dcSKarl Rupp 
1326e2d7f03fSJed Brown       vs->isglobal.row[i] = is_row[i];
1327d8588912SDave May     }
13282ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
13298188e55aSJed Brown     nsum = 0;
13308188e55aSJed Brown     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
13318188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1332ce94432eSBarry Smith       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
13330298fd71SBarry Smith       ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
1334ce94432eSBarry Smith       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
13358188e55aSJed Brown       nsum += n;
13368188e55aSJed Brown     }
1337ce94432eSBarry Smith     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
133830bc264bSJed Brown     offset -= nsum;
1339e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1340f349c1fdSJed Brown       ierr    = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
13410298fd71SBarry Smith       ierr    = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
13422ae74bdbSJed Brown       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1343ce94432eSBarry Smith       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
1344e2d7f03fSJed Brown       ierr    = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
13452ae74bdbSJed Brown       offset += n;
1346d8588912SDave May     }
1347d8588912SDave May   }
1348d8588912SDave May 
1349d8588912SDave May   if (is_col) { /* valid IS is passed in */
1350d8588912SDave May     /* refs on is[] are incremeneted */
1351e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1352d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
135326fbe8dcSKarl Rupp 
1354e2d7f03fSJed Brown       vs->isglobal.col[j] = is_col[j];
1355d8588912SDave May     }
13562ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
13572ae74bdbSJed Brown     offset = A->cmap->rstart;
13588188e55aSJed Brown     nsum   = 0;
13598188e55aSJed Brown     for (j=0; j<vs->nc; j++) {
13608188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
1361ce94432eSBarry Smith       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
13620298fd71SBarry Smith       ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
1363ce94432eSBarry Smith       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
13648188e55aSJed Brown       nsum += n;
13658188e55aSJed Brown     }
1366ce94432eSBarry Smith     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
136730bc264bSJed Brown     offset -= nsum;
1368e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1369f349c1fdSJed Brown       ierr    = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
13700298fd71SBarry Smith       ierr    = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
13712ae74bdbSJed Brown       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1372ce94432eSBarry Smith       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
1373e2d7f03fSJed Brown       ierr    = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
13742ae74bdbSJed Brown       offset += n;
1375d8588912SDave May     }
1376d8588912SDave May   }
1377e2d7f03fSJed Brown 
1378e2d7f03fSJed Brown   /* Set up the local ISs */
1379785e854fSJed Brown   ierr = PetscMalloc1(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
1380785e854fSJed Brown   ierr = PetscMalloc1(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
1381e2d7f03fSJed Brown   for (i=0,offset=0; i<vs->nr; i++) {
1382e2d7f03fSJed Brown     IS                     isloc;
13830298fd71SBarry Smith     ISLocalToGlobalMapping rmap = NULL;
1384e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1385e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
13860298fd71SBarry Smith     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,NULL);CHKERRQ(ierr);}
1387207556f9SJed Brown     if (rmap) {
1388e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1389e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
1390e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1391e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1392207556f9SJed Brown     } else {
1393207556f9SJed Brown       nlocal = 0;
13940298fd71SBarry Smith       isloc  = NULL;
1395207556f9SJed Brown     }
1396e2d7f03fSJed Brown     vs->islocal.row[i] = isloc;
1397e2d7f03fSJed Brown     offset            += nlocal;
1398e2d7f03fSJed Brown   }
13998188e55aSJed Brown   for (i=0,offset=0; i<vs->nc; i++) {
1400e2d7f03fSJed Brown     IS                     isloc;
14010298fd71SBarry Smith     ISLocalToGlobalMapping cmap = NULL;
1402e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1403e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
14040298fd71SBarry Smith     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,NULL,&cmap);CHKERRQ(ierr);}
1405207556f9SJed Brown     if (cmap) {
1406e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1407e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr);
1408e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1409e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1410207556f9SJed Brown     } else {
1411207556f9SJed Brown       nlocal = 0;
14120298fd71SBarry Smith       isloc  = NULL;
1413207556f9SJed Brown     }
1414e2d7f03fSJed Brown     vs->islocal.col[i] = isloc;
1415e2d7f03fSJed Brown     offset            += nlocal;
1416e2d7f03fSJed Brown   }
14170189643fSJed Brown 
141877019fcaSJed Brown   /* Set up the aggregate ISLocalToGlobalMapping */
141977019fcaSJed Brown   {
142045b6f7e9SBarry Smith     ISLocalToGlobalMapping rmap,cmap;
142145b6f7e9SBarry Smith     ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap);CHKERRQ(ierr);
142245b6f7e9SBarry Smith     ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap);CHKERRQ(ierr);
142377019fcaSJed Brown     if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);}
142477019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
142577019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
142677019fcaSJed Brown   }
142777019fcaSJed Brown 
14280189643fSJed Brown #if defined(PETSC_USE_DEBUG)
14290189643fSJed Brown   for (i=0; i<vs->nr; i++) {
14300189643fSJed Brown     for (j=0; j<vs->nc; j++) {
14310189643fSJed Brown       PetscInt m,n,M,N,mi,ni,Mi,Ni;
14320189643fSJed Brown       Mat      B = vs->m[i][j];
14330189643fSJed Brown       if (!B) continue;
14340189643fSJed Brown       ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr);
14350189643fSJed Brown       ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr);
14360189643fSJed Brown       ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr);
14370189643fSJed Brown       ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr);
14380189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr);
14390189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr);
1440ce94432eSBarry 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);
1441ce94432eSBarry 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);
14420189643fSJed Brown     }
14430189643fSJed Brown   }
14440189643fSJed Brown #endif
1445a061e289SJed Brown 
1446a061e289SJed Brown   /* Set A->assembled if all non-null blocks are currently assembled */
1447a061e289SJed Brown   for (i=0; i<vs->nr; i++) {
1448a061e289SJed Brown     for (j=0; j<vs->nc; j++) {
1449a061e289SJed Brown       if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0);
1450a061e289SJed Brown     }
1451a061e289SJed Brown   }
1452a061e289SJed Brown   A->assembled = PETSC_TRUE;
1453d8588912SDave May   PetscFunctionReturn(0);
1454d8588912SDave May }
1455d8588912SDave May 
1456d8588912SDave May #undef __FUNCT__
1457d8588912SDave May #define __FUNCT__ "MatCreateNest"
145845c38901SJed Brown /*@C
1459659c6bb0SJed Brown    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
1460659c6bb0SJed Brown 
1461659c6bb0SJed Brown    Collective on Mat
1462659c6bb0SJed Brown 
1463659c6bb0SJed Brown    Input Parameter:
1464659c6bb0SJed Brown +  comm - Communicator for the new Mat
1465659c6bb0SJed Brown .  nr - number of nested row blocks
14660298fd71SBarry Smith .  is_row - index sets for each nested row block, or NULL to make contiguous
1467659c6bb0SJed Brown .  nc - number of nested column blocks
14680298fd71SBarry Smith .  is_col - index sets for each nested column block, or NULL to make contiguous
14690298fd71SBarry Smith -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1470659c6bb0SJed Brown 
1471659c6bb0SJed Brown    Output Parameter:
1472659c6bb0SJed Brown .  B - new matrix
1473659c6bb0SJed Brown 
1474659c6bb0SJed Brown    Level: advanced
1475659c6bb0SJed Brown 
1476950540a4SJed Brown .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST
1477659c6bb0SJed Brown @*/
14787087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1479d8588912SDave May {
1480d8588912SDave May   Mat            A;
1481d8588912SDave May   PetscErrorCode ierr;
1482d8588912SDave May 
1483d8588912SDave May   PetscFunctionBegin;
1484c8883902SJed Brown   *B   = 0;
1485d8588912SDave May   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1486c8883902SJed Brown   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
148791a28eb3SBarry Smith   A->preallocated = PETSC_TRUE;
1488c8883902SJed Brown   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1489d8588912SDave May   *B   = A;
1490d8588912SDave May   PetscFunctionReturn(0);
1491d8588912SDave May }
1492659c6bb0SJed Brown 
1493629c3df2SDmitry Karpeev #undef __FUNCT__
1494629c3df2SDmitry Karpeev #define __FUNCT__ "MatConvert_Nest_AIJ"
1495cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1496629c3df2SDmitry Karpeev {
1497629c3df2SDmitry Karpeev   PetscErrorCode ierr;
1498629c3df2SDmitry Karpeev   Mat_Nest       *nest = (Mat_Nest*)A->data;
149983b1a929SMark Adams   PetscInt       m,n,M,N,i,j,k,*dnnz,*onnz,rstart;
1500649b366bSFande Kong   PetscInt       cstart,cend;
1501629c3df2SDmitry Karpeev   Mat            C;
1502629c3df2SDmitry Karpeev 
1503629c3df2SDmitry Karpeev   PetscFunctionBegin;
1504629c3df2SDmitry Karpeev   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
1505629c3df2SDmitry Karpeev   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
1506649b366bSFande Kong   ierr = MatGetOwnershipRangeColumn(A,&cstart,&cend);CHKERRQ(ierr);
1507629c3df2SDmitry Karpeev   switch (reuse) {
1508629c3df2SDmitry Karpeev   case MAT_INITIAL_MATRIX:
1509ce94432eSBarry Smith     ierr    = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
1510629c3df2SDmitry Karpeev     ierr    = MatSetType(C,newtype);CHKERRQ(ierr);
1511629c3df2SDmitry Karpeev     ierr    = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr);
1512629c3df2SDmitry Karpeev     *newmat = C;
1513629c3df2SDmitry Karpeev     break;
1514629c3df2SDmitry Karpeev   case MAT_REUSE_MATRIX:
1515629c3df2SDmitry Karpeev     C = *newmat;
1516629c3df2SDmitry Karpeev     break;
1517ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse");
1518629c3df2SDmitry Karpeev   }
1519785e854fSJed Brown   ierr = PetscMalloc1(2*m,&dnnz);CHKERRQ(ierr);
1520629c3df2SDmitry Karpeev   onnz = dnnz + m;
1521629c3df2SDmitry Karpeev   for (k=0; k<m; k++) {
1522629c3df2SDmitry Karpeev     dnnz[k] = 0;
1523629c3df2SDmitry Karpeev     onnz[k] = 0;
1524629c3df2SDmitry Karpeev   }
1525629c3df2SDmitry Karpeev   for (j=0; j<nest->nc; ++j) {
1526629c3df2SDmitry Karpeev     IS             bNis;
1527629c3df2SDmitry Karpeev     PetscInt       bN;
1528629c3df2SDmitry Karpeev     const PetscInt *bNindices;
1529629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
1530629c3df2SDmitry Karpeev     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1531629c3df2SDmitry Karpeev     ierr = ISGetSize(bNis, &bN);CHKERRQ(ierr);
1532629c3df2SDmitry Karpeev     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1533629c3df2SDmitry Karpeev     for (i=0; i<nest->nr; ++i) {
1534629c3df2SDmitry Karpeev       PetscSF        bmsf;
1535649b366bSFande Kong       PetscSFNode    *iremote;
1536629c3df2SDmitry Karpeev       Mat            B;
1537649b366bSFande Kong       PetscInt       bm, *sub_dnnz,*sub_onnz, br;
1538629c3df2SDmitry Karpeev       const PetscInt *bmindices;
1539629c3df2SDmitry Karpeev       B = nest->m[i][j];
1540629c3df2SDmitry Karpeev       if (!B) continue;
1541629c3df2SDmitry Karpeev       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1542629c3df2SDmitry Karpeev       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1543ce94432eSBarry Smith       ierr = PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);CHKERRQ(ierr);
1544649b366bSFande Kong       ierr = PetscMalloc1(bm,&iremote);CHKERRQ(ierr);
1545649b366bSFande Kong       ierr = PetscMalloc1(bm,&sub_dnnz);CHKERRQ(ierr);
1546649b366bSFande Kong       ierr = PetscMalloc1(bm,&sub_onnz);CHKERRQ(ierr);
1547649b366bSFande Kong       for (k = 0; k < bm; ++k){
1548649b366bSFande Kong     	sub_dnnz[k] = 0;
1549649b366bSFande Kong     	sub_onnz[k] = 0;
1550649b366bSFande Kong       }
1551629c3df2SDmitry Karpeev       /*
1552629c3df2SDmitry Karpeev        Locate the owners for all of the locally-owned global row indices for this row block.
1553629c3df2SDmitry Karpeev        These determine the roots of PetscSF used to communicate preallocation data to row owners.
1554629c3df2SDmitry Karpeev        The roots correspond to the dnnz and onnz entries; thus, there are two roots per row.
1555629c3df2SDmitry Karpeev        */
155683b1a929SMark Adams       ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr);
1557629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
1558649b366bSFande Kong         PetscInt       row = bmindices[br], rowowner = 0, brncols, col;
1559629c3df2SDmitry Karpeev         const PetscInt *brcols;
1560a4b3d3acSMatthew G Knepley         PetscInt       rowrel = 0; /* row's relative index on its owner rank */
1561629c3df2SDmitry Karpeev         ierr      = PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);CHKERRQ(ierr);
1562649b366bSFande Kong         /* how many roots  */
1563649b366bSFande Kong         iremote[br].rank = rowowner; iremote[br].index = rowrel;           /* edge from bmdnnz to dnnz */
1564649b366bSFande Kong         /* get nonzero pattern */
156583b1a929SMark Adams         ierr = MatGetRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr);
1566629c3df2SDmitry Karpeev         for (k=0; k<brncols; k++) {
1567629c3df2SDmitry Karpeev           col  = bNindices[brcols[k]];
1568649b366bSFande Kong           if(col>=A->cmap->range[rowowner] && col<A->cmap->range[rowowner+1]){
1569649b366bSFande Kong         	sub_dnnz[br]++;
1570649b366bSFande Kong           }else{
1571649b366bSFande Kong         	sub_onnz[br]++;
1572649b366bSFande Kong           }
1573629c3df2SDmitry Karpeev         }
157483b1a929SMark Adams         ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr);
1575629c3df2SDmitry Karpeev       }
1576629c3df2SDmitry Karpeev       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1577629c3df2SDmitry Karpeev       /* bsf will have to take care of disposing of bedges. */
1578649b366bSFande Kong       ierr = PetscSFSetGraph(bmsf,m,bm,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1579649b366bSFande Kong       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr);
1580649b366bSFande Kong       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr);
1581649b366bSFande Kong       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr);
1582649b366bSFande Kong       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr);
1583649b366bSFande Kong       ierr = PetscFree(sub_dnnz);CHKERRQ(ierr);
1584649b366bSFande Kong       ierr = PetscFree(sub_onnz);CHKERRQ(ierr);
1585629c3df2SDmitry Karpeev       ierr = PetscSFDestroy(&bmsf);CHKERRQ(ierr);
1586629c3df2SDmitry Karpeev     }
158722d28d08SBarry Smith     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
1588629c3df2SDmitry Karpeev     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
1589629c3df2SDmitry Karpeev   }
1590629c3df2SDmitry Karpeev   ierr = MatSeqAIJSetPreallocation(C,0,dnnz);CHKERRQ(ierr);
1591629c3df2SDmitry Karpeev   ierr = MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);CHKERRQ(ierr);
1592629c3df2SDmitry Karpeev   ierr = PetscFree(dnnz);CHKERRQ(ierr);
1593629c3df2SDmitry Karpeev 
1594629c3df2SDmitry Karpeev   /* Fill by row */
1595629c3df2SDmitry Karpeev   for (j=0; j<nest->nc; ++j) {
1596629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
1597629c3df2SDmitry Karpeev     IS             bNis;
1598629c3df2SDmitry Karpeev     PetscInt       bN;
1599629c3df2SDmitry Karpeev     const PetscInt *bNindices;
1600629c3df2SDmitry Karpeev     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1601629c3df2SDmitry Karpeev     ierr = ISGetSize(bNis,&bN);CHKERRQ(ierr);
1602629c3df2SDmitry Karpeev     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1603629c3df2SDmitry Karpeev     for (i=0; i<nest->nr; ++i) {
1604629c3df2SDmitry Karpeev       Mat            B;
1605629c3df2SDmitry Karpeev       PetscInt       bm, br;
1606629c3df2SDmitry Karpeev       const PetscInt *bmindices;
1607629c3df2SDmitry Karpeev       B = nest->m[i][j];
1608629c3df2SDmitry Karpeev       if (!B) continue;
1609629c3df2SDmitry Karpeev       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1610629c3df2SDmitry Karpeev       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
161183b1a929SMark Adams       ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr);
1612629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
1613629c3df2SDmitry Karpeev         PetscInt          row = bmindices[br], brncols,  *cols;
1614629c3df2SDmitry Karpeev         const PetscInt    *brcols;
1615629c3df2SDmitry Karpeev         const PetscScalar *brcoldata;
161683b1a929SMark Adams         ierr = MatGetRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
1617785e854fSJed Brown         ierr = PetscMalloc1(brncols,&cols);CHKERRQ(ierr);
161826fbe8dcSKarl Rupp         for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]];
1619629c3df2SDmitry Karpeev         /*
1620629c3df2SDmitry Karpeev           Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match.
1621629c3df2SDmitry Karpeev           Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES.
1622629c3df2SDmitry Karpeev          */
1623a2ea699eSBarry Smith         ierr = MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);CHKERRQ(ierr);
162483b1a929SMark Adams         ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
1625629c3df2SDmitry Karpeev         ierr = PetscFree(cols);CHKERRQ(ierr);
1626629c3df2SDmitry Karpeev       }
1627629c3df2SDmitry Karpeev       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1628629c3df2SDmitry Karpeev     }
1629a2ea699eSBarry Smith     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
1630629c3df2SDmitry Karpeev     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
1631629c3df2SDmitry Karpeev   }
1632629c3df2SDmitry Karpeev   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1633629c3df2SDmitry Karpeev   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1634629c3df2SDmitry Karpeev   PetscFunctionReturn(0);
1635629c3df2SDmitry Karpeev }
1636629c3df2SDmitry Karpeev 
1637659c6bb0SJed Brown /*MC
1638659c6bb0SJed Brown   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
1639659c6bb0SJed Brown 
1640659c6bb0SJed Brown   Level: intermediate
1641659c6bb0SJed Brown 
1642659c6bb0SJed Brown   Notes:
1643659c6bb0SJed Brown   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1644659c6bb0SJed Brown   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1645950540a4SJed Brown   It is usually used with DMComposite and DMCreateMatrix()
1646659c6bb0SJed Brown 
1647659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest()
1648659c6bb0SJed Brown M*/
1649c8883902SJed Brown #undef __FUNCT__
1650c8883902SJed Brown #define __FUNCT__ "MatCreate_Nest"
16518cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A)
1652c8883902SJed Brown {
1653c8883902SJed Brown   Mat_Nest       *s;
1654c8883902SJed Brown   PetscErrorCode ierr;
1655c8883902SJed Brown 
1656c8883902SJed Brown   PetscFunctionBegin;
1657b00a9115SJed Brown   ierr    = PetscNewLog(A,&s);CHKERRQ(ierr);
1658c8883902SJed Brown   A->data = (void*)s;
1659e7c19651SJed Brown 
1660e7c19651SJed Brown   s->nr            = -1;
1661e7c19651SJed Brown   s->nc            = -1;
16620298fd71SBarry Smith   s->m             = NULL;
1663e7c19651SJed Brown   s->splitassembly = PETSC_FALSE;
1664c8883902SJed Brown 
1665c8883902SJed Brown   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
166626fbe8dcSKarl Rupp 
1667c8883902SJed Brown   A->ops->mult                  = MatMult_Nest;
16689194d70fSJed Brown   A->ops->multadd               = MatMultAdd_Nest;
1669c8883902SJed Brown   A->ops->multtranspose         = MatMultTranspose_Nest;
16709194d70fSJed Brown   A->ops->multtransposeadd      = MatMultTransposeAdd_Nest;
1671f8170845SAlex Fikl   A->ops->transpose             = MatTranspose_Nest;
1672c8883902SJed Brown   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1673c8883902SJed Brown   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1674c8883902SJed Brown   A->ops->zeroentries           = MatZeroEntries_Nest;
1675c222c20dSDavid Ham   A->ops->copy                  = MatCopy_Nest;
1676c8883902SJed Brown   A->ops->duplicate             = MatDuplicate_Nest;
1677c8883902SJed Brown   A->ops->getsubmatrix          = MatGetSubMatrix_Nest;
1678c8883902SJed Brown   A->ops->destroy               = MatDestroy_Nest;
1679c8883902SJed Brown   A->ops->view                  = MatView_Nest;
1680c8883902SJed Brown   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1681c8883902SJed Brown   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1682c8883902SJed Brown   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
1683429bac76SJed Brown   A->ops->getdiagonal           = MatGetDiagonal_Nest;
1684429bac76SJed Brown   A->ops->diagonalscale         = MatDiagonalScale_Nest;
1685a061e289SJed Brown   A->ops->scale                 = MatScale_Nest;
1686a061e289SJed Brown   A->ops->shift                 = MatShift_Nest;
168713135bc6SAlex Fikl   A->ops->diagonalset           = MatDiagonalSet_Nest;
1688f8170845SAlex Fikl   A->ops->setrandom             = MatSetRandom_Nest;
1689c8883902SJed Brown 
1690c8883902SJed Brown   A->spptr        = 0;
1691c8883902SJed Brown   A->assembled    = PETSC_FALSE;
1692c8883902SJed Brown 
1693c8883902SJed Brown   /* expose Nest api's */
1694bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",   MatNestGetSubMat_Nest);CHKERRQ(ierr);
1695bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",   MatNestSetSubMat_Nest);CHKERRQ(ierr);
1696bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",  MatNestGetSubMats_Nest);CHKERRQ(ierr);
1697bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",     MatNestGetSize_Nest);CHKERRQ(ierr);
1698bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",      MatNestGetISs_Nest);CHKERRQ(ierr);
1699bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C", MatNestGetLocalISs_Nest);CHKERRQ(ierr);
1700bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",  MatNestSetVecType_Nest);CHKERRQ(ierr);
1701bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",  MatNestSetSubMats_Nest);CHKERRQ(ierr);
170283b1a929SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",MatConvert_Nest_AIJ);CHKERRQ(ierr);
1703c8883902SJed Brown 
1704c8883902SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
1705c8883902SJed Brown   PetscFunctionReturn(0);
1706c8883902SJed Brown }
1707