xref: /petsc/src/mat/impls/nest/matnest.c (revision 65a4a0a3185e986246d3e0b641268d9edea7e365)
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);
75e3038f0Sstefano_zampini PETSC_INTERN PetscErrorCode MatConvert_Nest_IS(Mat,MatType,MatReuse,Mat*);
8c8883902SJed Brown 
9d8588912SDave May /* private functions */
10d8588912SDave May #undef __FUNCT__
118188e55aSJed Brown #define __FUNCT__ "MatNestGetSizes_Private"
128188e55aSJed Brown static PetscErrorCode MatNestGetSizes_Private(Mat A,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N)
13d8588912SDave May {
14d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
158188e55aSJed Brown   PetscInt       i,j;
16d8588912SDave May   PetscErrorCode ierr;
17d8588912SDave May 
18d8588912SDave May   PetscFunctionBegin;
198188e55aSJed Brown   *m = *n = *M = *N = 0;
208188e55aSJed Brown   for (i=0; i<bA->nr; i++) {  /* rows */
218188e55aSJed Brown     PetscInt sm,sM;
228188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.row[i],&sm);CHKERRQ(ierr);
238188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.row[i],&sM);CHKERRQ(ierr);
248188e55aSJed Brown     *m  += sm;
258188e55aSJed Brown     *M  += sM;
26d8588912SDave May   }
278188e55aSJed Brown   for (j=0; j<bA->nc; j++) {  /* cols */
288188e55aSJed Brown     PetscInt sn,sN;
298188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.col[j],&sn);CHKERRQ(ierr);
308188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.col[j],&sN);CHKERRQ(ierr);
318188e55aSJed Brown     *n  += sn;
328188e55aSJed Brown     *N  += sN;
33d8588912SDave May   }
34d8588912SDave May   PetscFunctionReturn(0);
35d8588912SDave May }
36d8588912SDave May 
37d8588912SDave May /* operations */
38d8588912SDave May #undef __FUNCT__
39d8588912SDave May #define __FUNCT__ "MatMult_Nest"
40207556f9SJed Brown static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y)
41d8588912SDave May {
42d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
43207556f9SJed Brown   Vec            *bx = bA->right,*by = bA->left;
44207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
45d8588912SDave May   PetscErrorCode ierr;
46d8588912SDave May 
47d8588912SDave May   PetscFunctionBegin;
48207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
49207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
50207556f9SJed Brown   for (i=0; i<nr; i++) {
51d8588912SDave May     ierr = VecZeroEntries(by[i]);CHKERRQ(ierr);
52207556f9SJed Brown     for (j=0; j<nc; j++) {
53207556f9SJed Brown       if (!bA->m[i][j]) continue;
54d8588912SDave May       /* y[i] <- y[i] + A[i][j] * x[j] */
55d8588912SDave May       ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr);
56d8588912SDave May     }
57d8588912SDave May   }
58207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
59207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
60d8588912SDave May   PetscFunctionReturn(0);
61d8588912SDave May }
62d8588912SDave May 
63d8588912SDave May #undef __FUNCT__
649194d70fSJed Brown #define __FUNCT__ "MatMultAdd_Nest"
659194d70fSJed Brown static PetscErrorCode MatMultAdd_Nest(Mat A,Vec x,Vec y,Vec z)
669194d70fSJed Brown {
679194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
689194d70fSJed Brown   Vec            *bx = bA->right,*bz = bA->left;
699194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
709194d70fSJed Brown   PetscErrorCode ierr;
719194d70fSJed Brown 
729194d70fSJed Brown   PetscFunctionBegin;
739194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
749194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
759194d70fSJed Brown   for (i=0; i<nr; i++) {
769194d70fSJed Brown     if (y != z) {
779194d70fSJed Brown       Vec by;
789194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
799194d70fSJed Brown       ierr = VecCopy(by,bz[i]);CHKERRQ(ierr);
80336d21e7SJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
819194d70fSJed Brown     }
829194d70fSJed Brown     for (j=0; j<nc; j++) {
839194d70fSJed Brown       if (!bA->m[i][j]) continue;
849194d70fSJed Brown       /* y[i] <- y[i] + A[i][j] * x[j] */
859194d70fSJed Brown       ierr = MatMultAdd(bA->m[i][j],bx[j],bz[i],bz[i]);CHKERRQ(ierr);
869194d70fSJed Brown     }
879194d70fSJed Brown   }
889194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
899194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
909194d70fSJed Brown   PetscFunctionReturn(0);
919194d70fSJed Brown }
929194d70fSJed Brown 
939194d70fSJed Brown #undef __FUNCT__
94d8588912SDave May #define __FUNCT__ "MatMultTranspose_Nest"
95207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y)
96d8588912SDave May {
97d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
98207556f9SJed Brown   Vec            *bx = bA->left,*by = bA->right;
99207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
100d8588912SDave May   PetscErrorCode ierr;
101d8588912SDave May 
102d8588912SDave May   PetscFunctionBegin;
103609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
104609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
105207556f9SJed Brown   for (j=0; j<nc; j++) {
106609e31cbSJed Brown     ierr = VecZeroEntries(by[j]);CHKERRQ(ierr);
107609e31cbSJed Brown     for (i=0; i<nr; i++) {
1086c75ac25SJed Brown       if (!bA->m[i][j]) continue;
109609e31cbSJed Brown       /* y[j] <- y[j] + (A[i][j])^T * x[i] */
110609e31cbSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],by[j],by[j]);CHKERRQ(ierr);
111d8588912SDave May     }
112d8588912SDave May   }
113609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
114609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
115d8588912SDave May   PetscFunctionReturn(0);
116d8588912SDave May }
117d8588912SDave May 
118d8588912SDave May #undef __FUNCT__
1199194d70fSJed Brown #define __FUNCT__ "MatMultTransposeAdd_Nest"
1209194d70fSJed Brown static PetscErrorCode MatMultTransposeAdd_Nest(Mat A,Vec x,Vec y,Vec z)
1219194d70fSJed Brown {
1229194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
1239194d70fSJed Brown   Vec            *bx = bA->left,*bz = bA->right;
1249194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
1259194d70fSJed Brown   PetscErrorCode ierr;
1269194d70fSJed Brown 
1279194d70fSJed Brown   PetscFunctionBegin;
1289194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
1299194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
1309194d70fSJed Brown   for (j=0; j<nc; j++) {
1319194d70fSJed Brown     if (y != z) {
1329194d70fSJed Brown       Vec by;
1339194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
1349194d70fSJed Brown       ierr = VecCopy(by,bz[j]);CHKERRQ(ierr);
1359194d70fSJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
1369194d70fSJed Brown     }
1379194d70fSJed Brown     for (i=0; i<nr; i++) {
1386c75ac25SJed Brown       if (!bA->m[i][j]) continue;
1399194d70fSJed Brown       /* z[j] <- y[j] + (A[i][j])^T * x[i] */
1409194d70fSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],bz[j],bz[j]);CHKERRQ(ierr);
1419194d70fSJed Brown     }
1429194d70fSJed Brown   }
1439194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
1449194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
1459194d70fSJed Brown   PetscFunctionReturn(0);
1469194d70fSJed Brown }
1479194d70fSJed Brown 
1489194d70fSJed Brown #undef __FUNCT__
149f8170845SAlex Fikl #define __FUNCT__ "MatTranspose_Nest"
150f8170845SAlex Fikl static PetscErrorCode MatTranspose_Nest(Mat A,MatReuse reuse,Mat *B)
151f8170845SAlex Fikl {
152f8170845SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data, *bC;
153f8170845SAlex Fikl   Mat            C;
154f8170845SAlex Fikl   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
155f8170845SAlex Fikl   PetscErrorCode ierr;
156f8170845SAlex Fikl 
157f8170845SAlex Fikl   PetscFunctionBegin;
158f8170845SAlex 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");
159f8170845SAlex Fikl 
160f8170845SAlex Fikl   if (reuse == MAT_INITIAL_MATRIX || A == *B) {
161f8170845SAlex Fikl     Mat *subs;
162f8170845SAlex Fikl     IS  *is_row,*is_col;
163f8170845SAlex Fikl 
164f8170845SAlex Fikl     ierr = PetscCalloc1(nr * nc,&subs);CHKERRQ(ierr);
165f8170845SAlex Fikl     ierr = PetscMalloc2(nr,&is_row,nc,&is_col);CHKERRQ(ierr);
166f8170845SAlex Fikl     ierr = MatNestGetISs(A,is_row,is_col);CHKERRQ(ierr);
167ddeb9bd8SAlex Fikl     if (reuse == MAT_REUSE_MATRIX) {
168ddeb9bd8SAlex Fikl       for (i=0; i<nr; i++) {
169ddeb9bd8SAlex Fikl         for (j=0; j<nc; j++) {
170ddeb9bd8SAlex Fikl           subs[i + nr * j] = bA->m[i][j];
171ddeb9bd8SAlex Fikl         }
172ddeb9bd8SAlex Fikl       }
173ddeb9bd8SAlex Fikl     }
174ddeb9bd8SAlex Fikl 
175f8170845SAlex Fikl     ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nc,is_col,nr,is_row,subs,&C);CHKERRQ(ierr);
176f8170845SAlex Fikl     ierr = PetscFree(subs);CHKERRQ(ierr);
1773d994f23SBarry Smith     ierr = PetscFree2(is_row,is_col);CHKERRQ(ierr);
178f8170845SAlex Fikl   } else {
179f8170845SAlex Fikl     C = *B;
180f8170845SAlex Fikl   }
181f8170845SAlex Fikl 
182f8170845SAlex Fikl   bC = (Mat_Nest*)C->data;
183f8170845SAlex Fikl   for (i=0; i<nr; i++) {
184f8170845SAlex Fikl     for (j=0; j<nc; j++) {
185f8170845SAlex Fikl       if (bA->m[i][j]) {
186f8170845SAlex Fikl         ierr = MatTranspose(bA->m[i][j], reuse, &(bC->m[j][i]));CHKERRQ(ierr);
187f8170845SAlex Fikl       } else {
188f8170845SAlex Fikl         bC->m[j][i] = NULL;
189f8170845SAlex Fikl       }
190f8170845SAlex Fikl     }
191f8170845SAlex Fikl   }
192f8170845SAlex Fikl 
193f8170845SAlex Fikl   if (reuse == MAT_INITIAL_MATRIX || A != *B) {
194f8170845SAlex Fikl     *B = C;
195f8170845SAlex Fikl   } else {
196f8170845SAlex Fikl     ierr = MatHeaderMerge(A, &C);CHKERRQ(ierr);
197f8170845SAlex Fikl   }
198f8170845SAlex Fikl   PetscFunctionReturn(0);
199f8170845SAlex Fikl }
200f8170845SAlex Fikl 
201f8170845SAlex Fikl #undef __FUNCT__
202e2d7f03fSJed Brown #define __FUNCT__ "MatNestDestroyISList"
203e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list)
204e2d7f03fSJed Brown {
205e2d7f03fSJed Brown   PetscErrorCode ierr;
206e2d7f03fSJed Brown   IS             *lst = *list;
207e2d7f03fSJed Brown   PetscInt       i;
208e2d7f03fSJed Brown 
209e2d7f03fSJed Brown   PetscFunctionBegin;
210e2d7f03fSJed Brown   if (!lst) PetscFunctionReturn(0);
2116bf464f9SBarry Smith   for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);}
212e2d7f03fSJed Brown   ierr  = PetscFree(lst);CHKERRQ(ierr);
2130298fd71SBarry Smith   *list = NULL;
214e2d7f03fSJed Brown   PetscFunctionReturn(0);
215e2d7f03fSJed Brown }
216e2d7f03fSJed Brown 
217e2d7f03fSJed Brown #undef __FUNCT__
218d8588912SDave May #define __FUNCT__ "MatDestroy_Nest"
219207556f9SJed Brown static PetscErrorCode MatDestroy_Nest(Mat A)
220d8588912SDave May {
221d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
222d8588912SDave May   PetscInt       i,j;
223d8588912SDave May   PetscErrorCode ierr;
224d8588912SDave May 
225d8588912SDave May   PetscFunctionBegin;
226d8588912SDave May   /* release the matrices and the place holders */
227e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr);
228e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr);
229e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
230e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
231d8588912SDave May 
232d8588912SDave May   ierr = PetscFree(vs->row_len);CHKERRQ(ierr);
233d8588912SDave May   ierr = PetscFree(vs->col_len);CHKERRQ(ierr);
234d8588912SDave May 
235207556f9SJed Brown   ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr);
236207556f9SJed Brown 
237d8588912SDave May   /* release the matrices and the place holders */
238d8588912SDave May   if (vs->m) {
239d8588912SDave May     for (i=0; i<vs->nr; i++) {
240d8588912SDave May       for (j=0; j<vs->nc; j++) {
2416bf464f9SBarry Smith         ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr);
242d8588912SDave May       }
243d8588912SDave May       ierr = PetscFree(vs->m[i]);CHKERRQ(ierr);
244d8588912SDave May     }
245d8588912SDave May     ierr = PetscFree(vs->m);CHKERRQ(ierr);
246d8588912SDave May   }
247bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
248d8588912SDave May 
249bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",0);CHKERRQ(ierr);
250bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",0);CHKERRQ(ierr);
251bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",0);CHKERRQ(ierr);
252bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",0);CHKERRQ(ierr);
253bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",0);CHKERRQ(ierr);
254bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",0);CHKERRQ(ierr);
255bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",0);CHKERRQ(ierr);
256bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",0);CHKERRQ(ierr);
2575e3038f0Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",0);CHKERRQ(ierr);
2585e3038f0Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_is_C",0);CHKERRQ(ierr);
259d8588912SDave May   PetscFunctionReturn(0);
260d8588912SDave May }
261d8588912SDave May 
262d8588912SDave May #undef __FUNCT__
263d8588912SDave May #define __FUNCT__ "MatAssemblyBegin_Nest"
264207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
265d8588912SDave May {
266d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
267d8588912SDave May   PetscInt       i,j;
268d8588912SDave May   PetscErrorCode ierr;
269d8588912SDave May 
270d8588912SDave May   PetscFunctionBegin;
271d8588912SDave May   for (i=0; i<vs->nr; i++) {
272d8588912SDave May     for (j=0; j<vs->nc; j++) {
273e7c19651SJed Brown       if (vs->m[i][j]) {
274e7c19651SJed Brown         ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr);
275e7c19651SJed Brown         if (!vs->splitassembly) {
276e7c19651SJed Brown           /* Note: split assembly will fail if the same block appears more than once (even indirectly through a nested
277e7c19651SJed 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
278e7c19651SJed Brown            * already performing an assembly, but the result would by more complicated and appears to offer less
279e7c19651SJed Brown            * potential for diagnostics and correctness checking. Split assembly should be fixed once there is an
280e7c19651SJed Brown            * interface for libraries to make asynchronous progress in "user-defined non-blocking collectives".
281e7c19651SJed Brown            */
282e7c19651SJed Brown           ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr);
283e7c19651SJed Brown         }
284e7c19651SJed Brown       }
285d8588912SDave May     }
286d8588912SDave May   }
287d8588912SDave May   PetscFunctionReturn(0);
288d8588912SDave May }
289d8588912SDave May 
290d8588912SDave May #undef __FUNCT__
291d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest"
292207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
293d8588912SDave May {
294d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
295d8588912SDave May   PetscInt       i,j;
296d8588912SDave May   PetscErrorCode ierr;
297d8588912SDave May 
298d8588912SDave May   PetscFunctionBegin;
299d8588912SDave May   for (i=0; i<vs->nr; i++) {
300d8588912SDave May     for (j=0; j<vs->nc; j++) {
301e7c19651SJed Brown       if (vs->m[i][j]) {
302e7c19651SJed Brown         if (vs->splitassembly) {
303e7c19651SJed Brown           ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr);
304e7c19651SJed Brown         }
305e7c19651SJed Brown       }
306d8588912SDave May     }
307d8588912SDave May   }
308d8588912SDave May   PetscFunctionReturn(0);
309d8588912SDave May }
310d8588912SDave May 
311d8588912SDave May #undef __FUNCT__
312f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatRow"
313f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
314d8588912SDave May {
315207556f9SJed Brown   PetscErrorCode ierr;
316f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
317f349c1fdSJed Brown   PetscInt       j;
318f349c1fdSJed Brown   Mat            sub;
319d8588912SDave May 
320d8588912SDave May   PetscFunctionBegin;
3210298fd71SBarry Smith   sub = (row < vs->nc) ? vs->m[row][row] : (Mat)NULL; /* Prefer to find on the diagonal */
322f349c1fdSJed Brown   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
3234994cf47SJed Brown   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
324f349c1fdSJed Brown   *B = sub;
325f349c1fdSJed Brown   PetscFunctionReturn(0);
326d8588912SDave May }
327d8588912SDave May 
328f349c1fdSJed Brown #undef __FUNCT__
329f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatCol"
330f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
331f349c1fdSJed Brown {
332207556f9SJed Brown   PetscErrorCode ierr;
333f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
334f349c1fdSJed Brown   PetscInt       i;
335f349c1fdSJed Brown   Mat            sub;
336f349c1fdSJed Brown 
337f349c1fdSJed Brown   PetscFunctionBegin;
3380298fd71SBarry Smith   sub = (col < vs->nr) ? vs->m[col][col] : (Mat)NULL; /* Prefer to find on the diagonal */
339f349c1fdSJed Brown   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
3404994cf47SJed Brown   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
341f349c1fdSJed Brown   *B = sub;
342f349c1fdSJed Brown   PetscFunctionReturn(0);
343d8588912SDave May }
344d8588912SDave May 
345f349c1fdSJed Brown #undef __FUNCT__
346f349c1fdSJed Brown #define __FUNCT__ "MatNestFindIS"
347f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
348f349c1fdSJed Brown {
349f349c1fdSJed Brown   PetscErrorCode ierr;
350f349c1fdSJed Brown   PetscInt       i;
351f349c1fdSJed Brown   PetscBool      flg;
352f349c1fdSJed Brown 
353f349c1fdSJed Brown   PetscFunctionBegin;
354f349c1fdSJed Brown   PetscValidPointer(list,3);
355f349c1fdSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,4);
356f349c1fdSJed Brown   PetscValidIntPointer(found,5);
357f349c1fdSJed Brown   *found = -1;
358f349c1fdSJed Brown   for (i=0; i<n; i++) {
359207556f9SJed Brown     if (!list[i]) continue;
360f349c1fdSJed Brown     ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr);
361f349c1fdSJed Brown     if (flg) {
362f349c1fdSJed Brown       *found = i;
363f349c1fdSJed Brown       PetscFunctionReturn(0);
364f349c1fdSJed Brown     }
365f349c1fdSJed Brown   }
366ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Could not find index set");
367f349c1fdSJed Brown   PetscFunctionReturn(0);
368f349c1fdSJed Brown }
369f349c1fdSJed Brown 
370f349c1fdSJed Brown #undef __FUNCT__
3718188e55aSJed Brown #define __FUNCT__ "MatNestGetRow"
3728188e55aSJed Brown /* Get a block row as a new MatNest */
3738188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B)
3748188e55aSJed Brown {
3758188e55aSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
3768188e55aSJed Brown   char           keyname[256];
3778188e55aSJed Brown   PetscErrorCode ierr;
3788188e55aSJed Brown 
3798188e55aSJed Brown   PetscFunctionBegin;
3800298fd71SBarry Smith   *B   = NULL;
3818caf3d72SBarry Smith   ierr = PetscSNPrintf(keyname,sizeof(keyname),"NestRow_%D",row);CHKERRQ(ierr);
3828188e55aSJed Brown   ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr);
3838188e55aSJed Brown   if (*B) PetscFunctionReturn(0);
3848188e55aSJed Brown 
385ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),1,NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr);
38626fbe8dcSKarl Rupp 
3878188e55aSJed Brown   (*B)->assembled = A->assembled;
38826fbe8dcSKarl Rupp 
3898188e55aSJed Brown   ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr);
3908188e55aSJed Brown   ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */
3918188e55aSJed Brown   PetscFunctionReturn(0);
3928188e55aSJed Brown }
3938188e55aSJed Brown 
3948188e55aSJed Brown #undef __FUNCT__
395f349c1fdSJed Brown #define __FUNCT__ "MatNestFindSubMat"
396f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
397f349c1fdSJed Brown {
398f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
3998188e55aSJed Brown   PetscErrorCode ierr;
4006b3a5b13SJed Brown   PetscInt       row,col;
401e072481dSJed Brown   PetscBool      same,isFullCol,isFullColGlobal;
402f349c1fdSJed Brown 
403f349c1fdSJed Brown   PetscFunctionBegin;
4048188e55aSJed Brown   /* Check if full column space. This is a hack */
4058188e55aSJed Brown   isFullCol = PETSC_FALSE;
406251f4c67SDmitry Karpeev   ierr      = PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr);
4078188e55aSJed Brown   if (same) {
40877019fcaSJed Brown     PetscInt n,first,step,i,an,am,afirst,astep;
4098188e55aSJed Brown     ierr      = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
4108188e55aSJed Brown     ierr      = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
41177019fcaSJed Brown     isFullCol = PETSC_TRUE;
41205ce4453SJed Brown     for (i=0,an=A->cmap->rstart; i<vs->nc; i++) {
41377019fcaSJed Brown       ierr = ISStrideGetInfo(is->col[i],&afirst,&astep);CHKERRQ(ierr);
41477019fcaSJed Brown       ierr = ISGetLocalSize(is->col[i],&am);CHKERRQ(ierr);
41577019fcaSJed Brown       if (afirst != an || astep != step) isFullCol = PETSC_FALSE;
41677019fcaSJed Brown       an += am;
41777019fcaSJed Brown     }
41805ce4453SJed Brown     if (an != A->cmap->rstart+n) isFullCol = PETSC_FALSE;
4198188e55aSJed Brown   }
420b2566f29SBarry Smith   ierr = MPIU_Allreduce(&isFullCol,&isFullColGlobal,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)iscol));CHKERRQ(ierr);
4218188e55aSJed Brown 
422e072481dSJed Brown   if (isFullColGlobal) {
4238188e55aSJed Brown     PetscInt row;
4248188e55aSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
4258188e55aSJed Brown     ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr);
4268188e55aSJed Brown   } else {
427f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
428f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr);
429b6480e04SStefano Zampini     if (!vs->m[row][col]) {
430b6480e04SStefano Zampini       PetscInt lr,lc;
431b6480e04SStefano Zampini 
432b6480e04SStefano Zampini       ierr = MatCreate(PetscObjectComm((PetscObject)A),&vs->m[row][col]);CHKERRQ(ierr);
433b6480e04SStefano Zampini       ierr = ISGetLocalSize(vs->isglobal.row[row],&lr);CHKERRQ(ierr);
434b6480e04SStefano Zampini       ierr = ISGetLocalSize(vs->isglobal.col[col],&lc);CHKERRQ(ierr);
435b6480e04SStefano Zampini       ierr = MatSetSizes(vs->m[row][col],lr,lc,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
436b6480e04SStefano Zampini       ierr = MatSetUp(vs->m[row][col]);CHKERRQ(ierr);
437b6480e04SStefano Zampini       ierr = MatAssemblyBegin(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
438b6480e04SStefano Zampini       ierr = MatAssemblyEnd(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
439b6480e04SStefano Zampini     }
440f349c1fdSJed Brown     *B = vs->m[row][col];
4418188e55aSJed Brown   }
442f349c1fdSJed Brown   PetscFunctionReturn(0);
443f349c1fdSJed Brown }
444f349c1fdSJed Brown 
445f349c1fdSJed Brown #undef __FUNCT__
446f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest"
447f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
448f349c1fdSJed Brown {
449f349c1fdSJed Brown   PetscErrorCode ierr;
450f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
451f349c1fdSJed Brown   Mat            sub;
452f349c1fdSJed Brown 
453f349c1fdSJed Brown   PetscFunctionBegin;
454f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr);
455f349c1fdSJed Brown   switch (reuse) {
456f349c1fdSJed Brown   case MAT_INITIAL_MATRIX:
4577874fa86SDave May     if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); }
458f349c1fdSJed Brown     *B = sub;
459f349c1fdSJed Brown     break;
460f349c1fdSJed Brown   case MAT_REUSE_MATRIX:
461ce94432eSBarry Smith     if (sub != *B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
462f349c1fdSJed Brown     break;
463f349c1fdSJed Brown   case MAT_IGNORE_MATRIX:       /* Nothing to do */
464f349c1fdSJed Brown     break;
465511c6705SHong Zhang   case MAT_INPLACE_MATRIX:       /* Nothing to do */
466511c6705SHong Zhang     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX is not supported yet");
467511c6705SHong Zhang     break;
468f349c1fdSJed Brown   }
469f349c1fdSJed Brown   PetscFunctionReturn(0);
470f349c1fdSJed Brown }
471f349c1fdSJed Brown 
472f349c1fdSJed Brown #undef __FUNCT__
473f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest"
474f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
475f349c1fdSJed Brown {
476f349c1fdSJed Brown   PetscErrorCode ierr;
477f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
478f349c1fdSJed Brown   Mat            sub;
479f349c1fdSJed Brown 
480f349c1fdSJed Brown   PetscFunctionBegin;
481f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
482f349c1fdSJed Brown   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
483f349c1fdSJed Brown   if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);}
484f349c1fdSJed Brown   *B = sub;
485d8588912SDave May   PetscFunctionReturn(0);
486d8588912SDave May }
487d8588912SDave May 
488d8588912SDave May #undef __FUNCT__
489d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest"
490207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
491d8588912SDave May {
492d8588912SDave May   PetscErrorCode ierr;
493f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
494f349c1fdSJed Brown   Mat            sub;
495d8588912SDave May 
496d8588912SDave May   PetscFunctionBegin;
497f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
498ce94432eSBarry Smith   if (*B != sub) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
499f349c1fdSJed Brown   if (sub) {
500ce94432eSBarry Smith     if (((PetscObject)sub)->refct <= 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
5016bf464f9SBarry Smith     ierr = MatDestroy(B);CHKERRQ(ierr);
502d8588912SDave May   }
503d8588912SDave May   PetscFunctionReturn(0);
504d8588912SDave May }
505d8588912SDave May 
506d8588912SDave May #undef __FUNCT__
5077874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest"
5087874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
5097874fa86SDave May {
5107874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
5117874fa86SDave May   PetscInt       i;
5127874fa86SDave May   PetscErrorCode ierr;
5137874fa86SDave May 
5147874fa86SDave May   PetscFunctionBegin;
5157874fa86SDave May   for (i=0; i<bA->nr; i++) {
516429bac76SJed Brown     Vec bv;
517429bac76SJed Brown     ierr = VecGetSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
5187874fa86SDave May     if (bA->m[i][i]) {
519429bac76SJed Brown       ierr = MatGetDiagonal(bA->m[i][i],bv);CHKERRQ(ierr);
5207874fa86SDave May     } else {
5215159a857SMatthew G. Knepley       ierr = VecSet(bv,0.0);CHKERRQ(ierr);
5227874fa86SDave May     }
523429bac76SJed Brown     ierr = VecRestoreSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
5247874fa86SDave May   }
5257874fa86SDave May   PetscFunctionReturn(0);
5267874fa86SDave May }
5277874fa86SDave May 
5287874fa86SDave May #undef __FUNCT__
5297874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest"
5307874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
5317874fa86SDave May {
5327874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
533429bac76SJed Brown   Vec            bl,*br;
5347874fa86SDave May   PetscInt       i,j;
5357874fa86SDave May   PetscErrorCode ierr;
5367874fa86SDave May 
5377874fa86SDave May   PetscFunctionBegin;
5383f800ebeSJed Brown   ierr = PetscCalloc1(bA->nc,&br);CHKERRQ(ierr);
5392e6472ebSElliott Sales de Andrade   if (r) {
540429bac76SJed Brown     for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
5412e6472ebSElliott Sales de Andrade   }
5422e6472ebSElliott Sales de Andrade   bl = NULL;
5437874fa86SDave May   for (i=0; i<bA->nr; i++) {
5442e6472ebSElliott Sales de Andrade     if (l) {
545429bac76SJed Brown       ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
5462e6472ebSElliott Sales de Andrade     }
5477874fa86SDave May     for (j=0; j<bA->nc; j++) {
5487874fa86SDave May       if (bA->m[i][j]) {
549429bac76SJed Brown         ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr);
5507874fa86SDave May       }
5517874fa86SDave May     }
5522e6472ebSElliott Sales de Andrade     if (l) {
553a061e289SJed Brown       ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
5547874fa86SDave May     }
5552e6472ebSElliott Sales de Andrade   }
5562e6472ebSElliott Sales de Andrade   if (r) {
557429bac76SJed Brown     for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
5582e6472ebSElliott Sales de Andrade   }
559429bac76SJed Brown   ierr = PetscFree(br);CHKERRQ(ierr);
5607874fa86SDave May   PetscFunctionReturn(0);
5617874fa86SDave May }
5627874fa86SDave May 
5637874fa86SDave May #undef __FUNCT__
564a061e289SJed Brown #define __FUNCT__ "MatScale_Nest"
565a061e289SJed Brown static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a)
566a061e289SJed Brown {
567a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
568a061e289SJed Brown   PetscInt       i,j;
569a061e289SJed Brown   PetscErrorCode ierr;
570a061e289SJed Brown 
571a061e289SJed Brown   PetscFunctionBegin;
572a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
573a061e289SJed Brown     for (j=0; j<bA->nc; j++) {
574a061e289SJed Brown       if (bA->m[i][j]) {
575a061e289SJed Brown         ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr);
576a061e289SJed Brown       }
577a061e289SJed Brown     }
578a061e289SJed Brown   }
579a061e289SJed Brown   PetscFunctionReturn(0);
580a061e289SJed Brown }
581a061e289SJed Brown 
582a061e289SJed Brown #undef __FUNCT__
583a061e289SJed Brown #define __FUNCT__ "MatShift_Nest"
584a061e289SJed Brown static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a)
585a061e289SJed Brown {
586a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
587a061e289SJed Brown   PetscInt       i;
588a061e289SJed Brown   PetscErrorCode ierr;
589a061e289SJed Brown 
590a061e289SJed Brown   PetscFunctionBegin;
591a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
592ce94432eSBarry 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);
593a061e289SJed Brown     ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr);
594a061e289SJed Brown   }
595a061e289SJed Brown   PetscFunctionReturn(0);
596a061e289SJed Brown }
597a061e289SJed Brown 
598a061e289SJed Brown #undef __FUNCT__
59913135bc6SAlex Fikl #define __FUNCT__ "MatDiagonalSet_Nest"
60013135bc6SAlex Fikl static PetscErrorCode MatDiagonalSet_Nest(Mat A,Vec D,InsertMode is)
60113135bc6SAlex Fikl {
60213135bc6SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data;
60313135bc6SAlex Fikl   PetscInt       i;
60413135bc6SAlex Fikl   PetscErrorCode ierr;
60513135bc6SAlex Fikl 
60613135bc6SAlex Fikl   PetscFunctionBegin;
60713135bc6SAlex Fikl   for (i=0; i<bA->nr; i++) {
60813135bc6SAlex Fikl     Vec bv;
60913135bc6SAlex Fikl     ierr = VecGetSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
61013135bc6SAlex Fikl     if (bA->m[i][i]) {
61113135bc6SAlex Fikl       ierr = MatDiagonalSet(bA->m[i][i],bv,is);CHKERRQ(ierr);
61213135bc6SAlex Fikl     }
61313135bc6SAlex Fikl     ierr = VecRestoreSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
61413135bc6SAlex Fikl   }
61513135bc6SAlex Fikl   PetscFunctionReturn(0);
61613135bc6SAlex Fikl }
61713135bc6SAlex Fikl 
61813135bc6SAlex Fikl #undef __FUNCT__
619f8170845SAlex Fikl #define __FUNCT__ "MatSetRandom_Nest"
620f8170845SAlex Fikl static PetscErrorCode MatSetRandom_Nest(Mat A,PetscRandom rctx)
621f8170845SAlex Fikl {
622f8170845SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data;
623f8170845SAlex Fikl   PetscInt       i,j;
624f8170845SAlex Fikl   PetscErrorCode ierr;
625f8170845SAlex Fikl 
626f8170845SAlex Fikl   PetscFunctionBegin;
627f8170845SAlex Fikl   for (i=0; i<bA->nr; i++) {
628f8170845SAlex Fikl     for (j=0; j<bA->nc; j++) {
629f8170845SAlex Fikl       if (bA->m[i][j]) {
630f8170845SAlex Fikl         ierr = MatSetRandom(bA->m[i][j],rctx);CHKERRQ(ierr);
631f8170845SAlex Fikl       }
632f8170845SAlex Fikl     }
633f8170845SAlex Fikl   }
634f8170845SAlex Fikl   PetscFunctionReturn(0);
635f8170845SAlex Fikl }
636f8170845SAlex Fikl 
637f8170845SAlex Fikl #undef __FUNCT__
6382a7a6963SBarry Smith #define __FUNCT__ "MatCreateVecs_Nest"
6392a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left)
640d8588912SDave May {
641d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
642d8588912SDave May   Vec            *L,*R;
643d8588912SDave May   MPI_Comm       comm;
644d8588912SDave May   PetscInt       i,j;
645d8588912SDave May   PetscErrorCode ierr;
646d8588912SDave May 
647d8588912SDave May   PetscFunctionBegin;
648ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
649d8588912SDave May   if (right) {
650d8588912SDave May     /* allocate R */
651854ce69bSBarry Smith     ierr = PetscMalloc1(bA->nc, &R);CHKERRQ(ierr);
652d8588912SDave May     /* Create the right vectors */
653d8588912SDave May     for (j=0; j<bA->nc; j++) {
654d8588912SDave May       for (i=0; i<bA->nr; i++) {
655d8588912SDave May         if (bA->m[i][j]) {
6562a7a6963SBarry Smith           ierr = MatCreateVecs(bA->m[i][j],&R[j],NULL);CHKERRQ(ierr);
657d8588912SDave May           break;
658d8588912SDave May         }
659d8588912SDave May       }
6606c4ed002SBarry Smith       if (i==bA->nr) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
661d8588912SDave May     }
662f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
663d8588912SDave May     /* hand back control to the nest vector */
664d8588912SDave May     for (j=0; j<bA->nc; j++) {
6656bf464f9SBarry Smith       ierr = VecDestroy(&R[j]);CHKERRQ(ierr);
666d8588912SDave May     }
667d8588912SDave May     ierr = PetscFree(R);CHKERRQ(ierr);
668d8588912SDave May   }
669d8588912SDave May 
670d8588912SDave May   if (left) {
671d8588912SDave May     /* allocate L */
672854ce69bSBarry Smith     ierr = PetscMalloc1(bA->nr, &L);CHKERRQ(ierr);
673d8588912SDave May     /* Create the left vectors */
674d8588912SDave May     for (i=0; i<bA->nr; i++) {
675d8588912SDave May       for (j=0; j<bA->nc; j++) {
676d8588912SDave May         if (bA->m[i][j]) {
6772a7a6963SBarry Smith           ierr = MatCreateVecs(bA->m[i][j],NULL,&L[i]);CHKERRQ(ierr);
678d8588912SDave May           break;
679d8588912SDave May         }
680d8588912SDave May       }
6816c4ed002SBarry Smith       if (j==bA->nc) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
682d8588912SDave May     }
683d8588912SDave May 
684f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
685d8588912SDave May     for (i=0; i<bA->nr; i++) {
6866bf464f9SBarry Smith       ierr = VecDestroy(&L[i]);CHKERRQ(ierr);
687d8588912SDave May     }
688d8588912SDave May 
689d8588912SDave May     ierr = PetscFree(L);CHKERRQ(ierr);
690d8588912SDave May   }
691d8588912SDave May   PetscFunctionReturn(0);
692d8588912SDave May }
693d8588912SDave May 
694d8588912SDave May #undef __FUNCT__
695d8588912SDave May #define __FUNCT__ "MatView_Nest"
696207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
697d8588912SDave May {
698d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
699d8588912SDave May   PetscBool      isascii;
700d8588912SDave May   PetscInt       i,j;
701d8588912SDave May   PetscErrorCode ierr;
702d8588912SDave May 
703d8588912SDave May   PetscFunctionBegin;
704251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
705d8588912SDave May   if (isascii) {
706d8588912SDave May 
707d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Matrix object: \n");CHKERRQ(ierr);
708d86155a6SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
709d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer, "type=nest, rows=%D, cols=%D \n",bA->nr,bA->nc);CHKERRQ(ierr);
710d8588912SDave May 
711d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"MatNest structure: \n");CHKERRQ(ierr);
712d8588912SDave May     for (i=0; i<bA->nr; i++) {
713d8588912SDave May       for (j=0; j<bA->nc; j++) {
71419fd82e9SBarry Smith         MatType   type;
715270f95d7SJed Brown         char      name[256] = "",prefix[256] = "";
716d8588912SDave May         PetscInt  NR,NC;
717d8588912SDave May         PetscBool isNest = PETSC_FALSE;
718d8588912SDave May 
719d8588912SDave May         if (!bA->m[i][j]) {
720d86155a6SBarry Smith           CHKERRQ(ierr);PetscViewerASCIIPrintf(viewer, "(%D,%D) : NULL \n",i,j);CHKERRQ(ierr);
721d8588912SDave May           continue;
722d8588912SDave May         }
723d8588912SDave May         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
724d8588912SDave May         ierr = MatGetType(bA->m[i][j], &type);CHKERRQ(ierr);
7258caf3d72SBarry Smith         if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);}
7268caf3d72SBarry Smith         if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);}
727251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
728d8588912SDave May 
729270f95d7SJed Brown         ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr);
730d8588912SDave May 
731d8588912SDave May         if (isNest) {
732270f95d7SJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);  /* push1 */
733d8588912SDave May           ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr);
734270f95d7SJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop1 */
735d8588912SDave May         }
736d8588912SDave May       }
737d8588912SDave May     }
738d86155a6SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop0 */
739d8588912SDave May   }
740d8588912SDave May   PetscFunctionReturn(0);
741d8588912SDave May }
742d8588912SDave May 
743d8588912SDave May #undef __FUNCT__
744d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest"
745207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A)
746d8588912SDave May {
747d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
748d8588912SDave May   PetscInt       i,j;
749d8588912SDave May   PetscErrorCode ierr;
750d8588912SDave May 
751d8588912SDave May   PetscFunctionBegin;
752d8588912SDave May   for (i=0; i<bA->nr; i++) {
753d8588912SDave May     for (j=0; j<bA->nc; j++) {
754d8588912SDave May       if (!bA->m[i][j]) continue;
755d8588912SDave May       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
756d8588912SDave May     }
757d8588912SDave May   }
758d8588912SDave May   PetscFunctionReturn(0);
759d8588912SDave May }
760d8588912SDave May 
761d8588912SDave May #undef __FUNCT__
762c222c20dSDavid Ham #define __FUNCT__ "MatCopy_Nest"
763c222c20dSDavid Ham static PetscErrorCode MatCopy_Nest(Mat A,Mat B,MatStructure str)
764c222c20dSDavid Ham {
765c222c20dSDavid Ham   Mat_Nest       *bA = (Mat_Nest*)A->data,*bB = (Mat_Nest*)B->data;
766c222c20dSDavid Ham   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
767c222c20dSDavid Ham   PetscErrorCode ierr;
768c222c20dSDavid Ham 
769c222c20dSDavid Ham   PetscFunctionBegin;
770c222c20dSDavid 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);
771c222c20dSDavid Ham   for (i=0; i<nr; i++) {
772c222c20dSDavid Ham     for (j=0; j<nc; j++) {
77346a2b97cSJed Brown       if (bA->m[i][j] && bB->m[i][j]) {
774c222c20dSDavid Ham         ierr = MatCopy(bA->m[i][j],bB->m[i][j],str);CHKERRQ(ierr);
77546a2b97cSJed 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);
776c222c20dSDavid Ham     }
777c222c20dSDavid Ham   }
778c222c20dSDavid Ham   PetscFunctionReturn(0);
779c222c20dSDavid Ham }
780c222c20dSDavid Ham 
781c222c20dSDavid Ham #undef __FUNCT__
782d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest"
783207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
784d8588912SDave May {
785d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
786841e96a3SJed Brown   Mat            *b;
787841e96a3SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
788d8588912SDave May   PetscErrorCode ierr;
789d8588912SDave May 
790d8588912SDave May   PetscFunctionBegin;
791785e854fSJed Brown   ierr = PetscMalloc1(nr*nc,&b);CHKERRQ(ierr);
792841e96a3SJed Brown   for (i=0; i<nr; i++) {
793841e96a3SJed Brown     for (j=0; j<nc; j++) {
794841e96a3SJed Brown       if (bA->m[i][j]) {
795841e96a3SJed Brown         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
796841e96a3SJed Brown       } else {
7970298fd71SBarry Smith         b[i*nc+j] = NULL;
798d8588912SDave May       }
799d8588912SDave May     }
800d8588912SDave May   }
801ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
802841e96a3SJed Brown   /* Give the new MatNest exclusive ownership */
803841e96a3SJed Brown   for (i=0; i<nr*nc; i++) {
8046bf464f9SBarry Smith     ierr = MatDestroy(&b[i]);CHKERRQ(ierr);
805d8588912SDave May   }
806d8588912SDave May   ierr = PetscFree(b);CHKERRQ(ierr);
807d8588912SDave May 
808841e96a3SJed Brown   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
809841e96a3SJed Brown   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
810d8588912SDave May   PetscFunctionReturn(0);
811d8588912SDave May }
812d8588912SDave May 
813d8588912SDave May /* nest api */
814d8588912SDave May #undef __FUNCT__
815d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest"
816d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
817d8588912SDave May {
818d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
8195fd66863SKarl Rupp 
820d8588912SDave May   PetscFunctionBegin;
821ce94432eSBarry Smith   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
822ce94432eSBarry Smith   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
823d8588912SDave May   *mat = bA->m[idxm][jdxm];
824d8588912SDave May   PetscFunctionReturn(0);
825d8588912SDave May }
826d8588912SDave May 
827d8588912SDave May #undef __FUNCT__
828d8588912SDave May #define __FUNCT__ "MatNestGetSubMat"
8299ba0d327SJed Brown /*@
830d8588912SDave May  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
831d8588912SDave May 
832d8588912SDave May  Not collective
833d8588912SDave May 
834d8588912SDave May  Input Parameters:
835629881c0SJed Brown +   A  - nest matrix
836d8588912SDave May .   idxm - index of the matrix within the nest matrix
837629881c0SJed Brown -   jdxm - index of the matrix within the nest matrix
838d8588912SDave May 
839d8588912SDave May  Output Parameter:
840d8588912SDave May .   sub - matrix at index idxm,jdxm within the nest matrix
841d8588912SDave May 
842d8588912SDave May  Level: developer
843d8588912SDave May 
844d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMats()
845d8588912SDave May @*/
8467087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
847d8588912SDave May {
848699a902aSJed Brown   PetscErrorCode ierr;
849d8588912SDave May 
850d8588912SDave May   PetscFunctionBegin;
851699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
852d8588912SDave May   PetscFunctionReturn(0);
853d8588912SDave May }
854d8588912SDave May 
855d8588912SDave May #undef __FUNCT__
8560782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat_Nest"
8570782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat)
8580782ca92SJed Brown {
8590782ca92SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
8600782ca92SJed Brown   PetscInt       m,n,M,N,mi,ni,Mi,Ni;
8610782ca92SJed Brown   PetscErrorCode ierr;
8620782ca92SJed Brown 
8630782ca92SJed Brown   PetscFunctionBegin;
864ce94432eSBarry Smith   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
865ce94432eSBarry Smith   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
8660782ca92SJed Brown   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
8670782ca92SJed Brown   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
8680782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr);
8690782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr);
8700782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr);
8710782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr);
872ce94432eSBarry 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);
873ce94432eSBarry 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);
87426fbe8dcSKarl Rupp 
8750782ca92SJed Brown   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
8760782ca92SJed Brown   ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr);
8770782ca92SJed Brown   bA->m[idxm][jdxm] = mat;
8780782ca92SJed Brown   PetscFunctionReturn(0);
8790782ca92SJed Brown }
8800782ca92SJed Brown 
8810782ca92SJed Brown #undef __FUNCT__
8820782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat"
8839ba0d327SJed Brown /*@
8840782ca92SJed Brown  MatNestSetSubMat - Set a single submatrix in the nest matrix.
8850782ca92SJed Brown 
8860782ca92SJed Brown  Logically collective on the submatrix communicator
8870782ca92SJed Brown 
8880782ca92SJed Brown  Input Parameters:
8890782ca92SJed Brown +   A  - nest matrix
8900782ca92SJed Brown .   idxm - index of the matrix within the nest matrix
8910782ca92SJed Brown .   jdxm - index of the matrix within the nest matrix
8920782ca92SJed Brown -   sub - matrix at index idxm,jdxm within the nest matrix
8930782ca92SJed Brown 
8940782ca92SJed Brown  Notes:
8950782ca92SJed Brown  The new submatrix must have the same size and communicator as that block of the nest.
8960782ca92SJed Brown 
8970782ca92SJed Brown  This increments the reference count of the submatrix.
8980782ca92SJed Brown 
8990782ca92SJed Brown  Level: developer
9000782ca92SJed Brown 
9010782ca92SJed Brown .seealso: MatNestSetSubMats(), MatNestGetSubMat()
9020782ca92SJed Brown @*/
9030782ca92SJed Brown PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
9040782ca92SJed Brown {
9050782ca92SJed Brown   PetscErrorCode ierr;
9060782ca92SJed Brown 
9070782ca92SJed Brown   PetscFunctionBegin;
9080782ca92SJed Brown   ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr);
9090782ca92SJed Brown   PetscFunctionReturn(0);
9100782ca92SJed Brown }
9110782ca92SJed Brown 
9120782ca92SJed Brown #undef __FUNCT__
913d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest"
914d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
915d8588912SDave May {
916d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
9175fd66863SKarl Rupp 
918d8588912SDave May   PetscFunctionBegin;
91926fbe8dcSKarl Rupp   if (M)   *M   = bA->nr;
92026fbe8dcSKarl Rupp   if (N)   *N   = bA->nc;
92126fbe8dcSKarl Rupp   if (mat) *mat = bA->m;
922d8588912SDave May   PetscFunctionReturn(0);
923d8588912SDave May }
924d8588912SDave May 
925d8588912SDave May #undef __FUNCT__
926d8588912SDave May #define __FUNCT__ "MatNestGetSubMats"
927d8588912SDave May /*@C
928d8588912SDave May  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
929d8588912SDave May 
930d8588912SDave May  Not collective
931d8588912SDave May 
932d8588912SDave May  Input Parameters:
933629881c0SJed Brown .   A  - nest matrix
934d8588912SDave May 
935d8588912SDave May  Output Parameter:
936629881c0SJed Brown +   M - number of rows in the nest matrix
937d8588912SDave May .   N - number of cols in the nest matrix
938629881c0SJed Brown -   mat - 2d array of matrices
939d8588912SDave May 
940d8588912SDave May  Notes:
941d8588912SDave May 
942d8588912SDave May  The user should not free the array mat.
943d8588912SDave May 
944351962e3SVincent Le Chenadec  In Fortran, this routine has a calling sequence
945351962e3SVincent Le Chenadec $   call MatNestGetSubMats(A, M, N, mat, ierr)
946351962e3SVincent Le Chenadec  where the space allocated for the optional argument mat is assumed large enough (if provided).
947351962e3SVincent Le Chenadec 
948d8588912SDave May  Level: developer
949d8588912SDave May 
950d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMat()
951d8588912SDave May @*/
9527087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
953d8588912SDave May {
954699a902aSJed Brown   PetscErrorCode ierr;
955d8588912SDave May 
956d8588912SDave May   PetscFunctionBegin;
957699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
958d8588912SDave May   PetscFunctionReturn(0);
959d8588912SDave May }
960d8588912SDave May 
961d8588912SDave May #undef __FUNCT__
962d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest"
9637087cfbeSBarry Smith PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
964d8588912SDave May {
965d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
966d8588912SDave May 
967d8588912SDave May   PetscFunctionBegin;
96826fbe8dcSKarl Rupp   if (M) *M = bA->nr;
96926fbe8dcSKarl Rupp   if (N) *N = bA->nc;
970d8588912SDave May   PetscFunctionReturn(0);
971d8588912SDave May }
972d8588912SDave May 
973d8588912SDave May #undef __FUNCT__
974d8588912SDave May #define __FUNCT__ "MatNestGetSize"
9759ba0d327SJed Brown /*@
976d8588912SDave May  MatNestGetSize - Returns the size of the nest matrix.
977d8588912SDave May 
978d8588912SDave May  Not collective
979d8588912SDave May 
980d8588912SDave May  Input Parameters:
981d8588912SDave May .   A  - nest matrix
982d8588912SDave May 
983d8588912SDave May  Output Parameter:
984629881c0SJed Brown +   M - number of rows in the nested mat
985629881c0SJed Brown -   N - number of cols in the nested mat
986d8588912SDave May 
987d8588912SDave May  Notes:
988d8588912SDave May 
989d8588912SDave May  Level: developer
990d8588912SDave May 
991d8588912SDave May .seealso: MatNestGetSubMat(), MatNestGetSubMats()
992d8588912SDave May @*/
9937087cfbeSBarry Smith PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
994d8588912SDave May {
995699a902aSJed Brown   PetscErrorCode ierr;
996d8588912SDave May 
997d8588912SDave May   PetscFunctionBegin;
998699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
999d8588912SDave May   PetscFunctionReturn(0);
1000d8588912SDave May }
1001d8588912SDave May 
1002900e7ff2SJed Brown #undef __FUNCT__
1003900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs_Nest"
1004f7a08781SBarry Smith static PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[])
1005900e7ff2SJed Brown {
1006900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest*)A->data;
1007900e7ff2SJed Brown   PetscInt i;
1008900e7ff2SJed Brown 
1009900e7ff2SJed Brown   PetscFunctionBegin;
1010900e7ff2SJed Brown   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i];
1011900e7ff2SJed Brown   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i];
1012900e7ff2SJed Brown   PetscFunctionReturn(0);
1013900e7ff2SJed Brown }
1014900e7ff2SJed Brown 
1015900e7ff2SJed Brown #undef __FUNCT__
1016900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs"
10173a4d7b9aSSatish Balay /*@C
1018900e7ff2SJed Brown  MatNestGetISs - Returns the index sets partitioning the row and column spaces
1019900e7ff2SJed Brown 
1020900e7ff2SJed Brown  Not collective
1021900e7ff2SJed Brown 
1022900e7ff2SJed Brown  Input Parameters:
1023900e7ff2SJed Brown .   A  - nest matrix
1024900e7ff2SJed Brown 
1025900e7ff2SJed Brown  Output Parameter:
1026900e7ff2SJed Brown +   rows - array of row index sets
1027900e7ff2SJed Brown -   cols - array of column index sets
1028900e7ff2SJed Brown 
1029900e7ff2SJed Brown  Level: advanced
1030900e7ff2SJed Brown 
1031900e7ff2SJed Brown  Notes:
1032900e7ff2SJed Brown  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
1033900e7ff2SJed Brown 
1034900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs()
1035900e7ff2SJed Brown @*/
1036900e7ff2SJed Brown PetscErrorCode  MatNestGetISs(Mat A,IS rows[],IS cols[])
1037900e7ff2SJed Brown {
1038900e7ff2SJed Brown   PetscErrorCode ierr;
1039900e7ff2SJed Brown 
1040900e7ff2SJed Brown   PetscFunctionBegin;
1041900e7ff2SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1042900e7ff2SJed Brown   ierr = PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
1043900e7ff2SJed Brown   PetscFunctionReturn(0);
1044900e7ff2SJed Brown }
1045900e7ff2SJed Brown 
1046900e7ff2SJed Brown #undef __FUNCT__
1047900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs_Nest"
1048f7a08781SBarry Smith static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[])
1049900e7ff2SJed Brown {
1050900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest*)A->data;
1051900e7ff2SJed Brown   PetscInt i;
1052900e7ff2SJed Brown 
1053900e7ff2SJed Brown   PetscFunctionBegin;
1054900e7ff2SJed Brown   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i];
1055900e7ff2SJed Brown   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i];
1056900e7ff2SJed Brown   PetscFunctionReturn(0);
1057900e7ff2SJed Brown }
1058900e7ff2SJed Brown 
1059900e7ff2SJed Brown #undef __FUNCT__
1060900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs"
1061900e7ff2SJed Brown /*@C
1062900e7ff2SJed Brown  MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces
1063900e7ff2SJed Brown 
1064900e7ff2SJed Brown  Not collective
1065900e7ff2SJed Brown 
1066900e7ff2SJed Brown  Input Parameters:
1067900e7ff2SJed Brown .   A  - nest matrix
1068900e7ff2SJed Brown 
1069900e7ff2SJed Brown  Output Parameter:
10700298fd71SBarry Smith +   rows - array of row index sets (or NULL to ignore)
10710298fd71SBarry Smith -   cols - array of column index sets (or NULL to ignore)
1072900e7ff2SJed Brown 
1073900e7ff2SJed Brown  Level: advanced
1074900e7ff2SJed Brown 
1075900e7ff2SJed Brown  Notes:
1076900e7ff2SJed Brown  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
1077900e7ff2SJed Brown 
1078900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs()
1079900e7ff2SJed Brown @*/
1080900e7ff2SJed Brown PetscErrorCode  MatNestGetLocalISs(Mat A,IS rows[],IS cols[])
1081900e7ff2SJed Brown {
1082900e7ff2SJed Brown   PetscErrorCode ierr;
1083900e7ff2SJed Brown 
1084900e7ff2SJed Brown   PetscFunctionBegin;
1085900e7ff2SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1086900e7ff2SJed Brown   ierr = PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
1087900e7ff2SJed Brown   PetscFunctionReturn(0);
1088900e7ff2SJed Brown }
1089900e7ff2SJed Brown 
1090207556f9SJed Brown #undef __FUNCT__
1091207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest"
109219fd82e9SBarry Smith PetscErrorCode  MatNestSetVecType_Nest(Mat A,VecType vtype)
1093207556f9SJed Brown {
1094207556f9SJed Brown   PetscErrorCode ierr;
1095207556f9SJed Brown   PetscBool      flg;
1096207556f9SJed Brown 
1097207556f9SJed Brown   PetscFunctionBegin;
1098207556f9SJed Brown   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
1099207556f9SJed Brown   /* In reality, this only distinguishes VECNEST and "other" */
11002a7a6963SBarry Smith   if (flg) A->ops->getvecs = MatCreateVecs_Nest;
110112b53f24SSatish Balay   else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0;
1102207556f9SJed Brown   PetscFunctionReturn(0);
1103207556f9SJed Brown }
1104207556f9SJed Brown 
1105207556f9SJed Brown #undef __FUNCT__
1106207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType"
1107207556f9SJed Brown /*@C
11082a7a6963SBarry Smith  MatNestSetVecType - Sets the type of Vec returned by MatCreateVecs()
1109207556f9SJed Brown 
1110207556f9SJed Brown  Not collective
1111207556f9SJed Brown 
1112207556f9SJed Brown  Input Parameters:
1113207556f9SJed Brown +  A  - nest matrix
1114207556f9SJed Brown -  vtype - type to use for creating vectors
1115207556f9SJed Brown 
1116207556f9SJed Brown  Notes:
1117207556f9SJed Brown 
1118207556f9SJed Brown  Level: developer
1119207556f9SJed Brown 
11202a7a6963SBarry Smith .seealso: MatCreateVecs()
1121207556f9SJed Brown @*/
112219fd82e9SBarry Smith PetscErrorCode  MatNestSetVecType(Mat A,VecType vtype)
1123207556f9SJed Brown {
1124207556f9SJed Brown   PetscErrorCode ierr;
1125207556f9SJed Brown 
1126207556f9SJed Brown   PetscFunctionBegin;
112719fd82e9SBarry Smith   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));CHKERRQ(ierr);
1128207556f9SJed Brown   PetscFunctionReturn(0);
1129207556f9SJed Brown }
1130207556f9SJed Brown 
1131d8588912SDave May #undef __FUNCT__
1132c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats_Nest"
1133c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1134d8588912SDave May {
1135c8883902SJed Brown   Mat_Nest       *s = (Mat_Nest*)A->data;
1136c8883902SJed Brown   PetscInt       i,j,m,n,M,N;
1137d8588912SDave May   PetscErrorCode ierr;
1138d8588912SDave May 
1139d8588912SDave May   PetscFunctionBegin;
1140c8883902SJed Brown   s->nr = nr;
1141c8883902SJed Brown   s->nc = nc;
1142d8588912SDave May 
1143c8883902SJed Brown   /* Create space for submatrices */
1144854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&s->m);CHKERRQ(ierr);
1145c8883902SJed Brown   for (i=0; i<nr; i++) {
1146854ce69bSBarry Smith     ierr = PetscMalloc1(nc,&s->m[i]);CHKERRQ(ierr);
1147d8588912SDave May   }
1148c8883902SJed Brown   for (i=0; i<nr; i++) {
1149c8883902SJed Brown     for (j=0; j<nc; j++) {
1150c8883902SJed Brown       s->m[i][j] = a[i*nc+j];
1151c8883902SJed Brown       if (a[i*nc+j]) {
1152c8883902SJed Brown         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
1153d8588912SDave May       }
1154d8588912SDave May     }
1155d8588912SDave May   }
1156d8588912SDave May 
11578188e55aSJed Brown   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
1158d8588912SDave May 
1159854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&s->row_len);CHKERRQ(ierr);
1160854ce69bSBarry Smith   ierr = PetscMalloc1(nc,&s->col_len);CHKERRQ(ierr);
1161c8883902SJed Brown   for (i=0; i<nr; i++) s->row_len[i]=-1;
1162c8883902SJed Brown   for (j=0; j<nc; j++) s->col_len[j]=-1;
1163d8588912SDave May 
11648188e55aSJed Brown   ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr);
1165d8588912SDave May 
1166c8883902SJed Brown   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
1167c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
1168c8883902SJed Brown   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
1169c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
1170c8883902SJed Brown 
1171c8883902SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1172c8883902SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1173c8883902SJed Brown 
11741795a4d1SJed Brown   ierr = PetscCalloc2(nr,&s->left,nc,&s->right);CHKERRQ(ierr);
1175d8588912SDave May   PetscFunctionReturn(0);
1176d8588912SDave May }
1177d8588912SDave May 
1178c8883902SJed Brown #undef __FUNCT__
1179c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats"
1180c8883902SJed Brown /*@
1181c8883902SJed Brown    MatNestSetSubMats - Sets the nested submatrices
1182c8883902SJed Brown 
1183c8883902SJed Brown    Collective on Mat
1184c8883902SJed Brown 
1185c8883902SJed Brown    Input Parameter:
1186c8883902SJed Brown +  N - nested matrix
1187c8883902SJed Brown .  nr - number of nested row blocks
11880298fd71SBarry Smith .  is_row - index sets for each nested row block, or NULL to make contiguous
1189c8883902SJed Brown .  nc - number of nested column blocks
11900298fd71SBarry Smith .  is_col - index sets for each nested column block, or NULL to make contiguous
11910298fd71SBarry Smith -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1192c8883902SJed Brown 
1193c8883902SJed Brown    Level: advanced
1194c8883902SJed Brown 
1195c8883902SJed Brown .seealso: MatCreateNest(), MATNEST
1196c8883902SJed Brown @*/
1197c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1198c8883902SJed Brown {
1199c8883902SJed Brown   PetscErrorCode ierr;
1200c8883902SJed Brown   PetscInt       i;
1201c8883902SJed Brown 
1202c8883902SJed Brown   PetscFunctionBegin;
1203c8883902SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1204ce94432eSBarry Smith   if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
1205c8883902SJed Brown   if (nr && is_row) {
1206c8883902SJed Brown     PetscValidPointer(is_row,3);
1207c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
1208c8883902SJed Brown   }
1209ce94432eSBarry Smith   if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
12101664e352SJed Brown   if (nc && is_col) {
1211c8883902SJed Brown     PetscValidPointer(is_col,5);
12129b30a8f6SBarry Smith     for (i=0; i<nc; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
1213c8883902SJed Brown   }
1214c8883902SJed Brown   if (nr*nc) PetscValidPointer(a,6);
1215c8883902SJed 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);
1216c8883902SJed Brown   PetscFunctionReturn(0);
1217c8883902SJed Brown }
1218d8588912SDave May 
121977019fcaSJed Brown #undef __FUNCT__
122077019fcaSJed Brown #define __FUNCT__ "MatNestCreateAggregateL2G_Private"
122145b6f7e9SBarry Smith static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog)
122277019fcaSJed Brown {
122377019fcaSJed Brown   PetscErrorCode ierr;
122477019fcaSJed Brown   PetscBool      flg;
122577019fcaSJed Brown   PetscInt       i,j,m,mi,*ix;
122677019fcaSJed Brown 
122777019fcaSJed Brown   PetscFunctionBegin;
122877019fcaSJed Brown   for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) {
122977019fcaSJed Brown     if (islocal[i]) {
123077019fcaSJed Brown       ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
123177019fcaSJed Brown       flg  = PETSC_TRUE;      /* We found a non-trivial entry */
123277019fcaSJed Brown     } else {
123377019fcaSJed Brown       ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
123477019fcaSJed Brown     }
123577019fcaSJed Brown     m += mi;
123677019fcaSJed Brown   }
123777019fcaSJed Brown   if (flg) {
1238785e854fSJed Brown     ierr = PetscMalloc1(m,&ix);CHKERRQ(ierr);
123977019fcaSJed Brown     for (i=0,n=0; i<n; i++) {
12400298fd71SBarry Smith       ISLocalToGlobalMapping smap = NULL;
124177019fcaSJed Brown       VecScatter             scat;
124277019fcaSJed Brown       IS                     isreq;
124377019fcaSJed Brown       Vec                    lvec,gvec;
12443361c9a7SJed Brown       union {char padding[sizeof(PetscScalar)]; PetscInt integer;} *x;
124577019fcaSJed Brown       Mat sub;
124677019fcaSJed Brown 
1247ce94432eSBarry Smith       if (sizeof(*x) != sizeof(PetscScalar)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support when scalars smaller than integers");
124877019fcaSJed Brown       if (colflg) {
124977019fcaSJed Brown         ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
125077019fcaSJed Brown       } else {
125177019fcaSJed Brown         ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
125277019fcaSJed Brown       }
12530298fd71SBarry Smith       if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&smap,NULL);CHKERRQ(ierr);}
125477019fcaSJed Brown       if (islocal[i]) {
125577019fcaSJed Brown         ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
125677019fcaSJed Brown       } else {
125777019fcaSJed Brown         ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
125877019fcaSJed Brown       }
125977019fcaSJed Brown       for (j=0; j<mi; j++) ix[m+j] = j;
126077019fcaSJed Brown       if (smap) {ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr);}
126177019fcaSJed Brown       /*
126277019fcaSJed Brown         Now we need to extract the monolithic global indices that correspond to the given split global indices.
126377019fcaSJed 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.
126477019fcaSJed Brown         The approach here is ugly because it uses VecScatter to move indices.
126577019fcaSJed Brown        */
126677019fcaSJed Brown       ierr = VecCreateSeq(PETSC_COMM_SELF,mi,&lvec);CHKERRQ(ierr);
126777019fcaSJed Brown       ierr = VecCreateMPI(((PetscObject)isglobal[i])->comm,mi,PETSC_DECIDE,&gvec);CHKERRQ(ierr);
126877019fcaSJed Brown       ierr = ISCreateGeneral(((PetscObject)isglobal[i])->comm,mi,ix+m,PETSC_COPY_VALUES,&isreq);CHKERRQ(ierr);
12690298fd71SBarry Smith       ierr = VecScatterCreate(gvec,isreq,lvec,NULL,&scat);CHKERRQ(ierr);
127077019fcaSJed Brown       ierr = VecGetArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
127177019fcaSJed Brown       for (j=0; j<mi; j++) x[j].integer = ix[m+j];
127277019fcaSJed Brown       ierr = VecRestoreArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
127377019fcaSJed Brown       ierr = VecScatterBegin(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
127477019fcaSJed Brown       ierr = VecScatterEnd(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
127577019fcaSJed Brown       ierr = VecGetArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
127677019fcaSJed Brown       for (j=0; j<mi; j++) ix[m+j] = x[j].integer;
127777019fcaSJed Brown       ierr = VecRestoreArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
127877019fcaSJed Brown       ierr = VecDestroy(&lvec);CHKERRQ(ierr);
127977019fcaSJed Brown       ierr = VecDestroy(&gvec);CHKERRQ(ierr);
128077019fcaSJed Brown       ierr = ISDestroy(&isreq);CHKERRQ(ierr);
128177019fcaSJed Brown       ierr = VecScatterDestroy(&scat);CHKERRQ(ierr);
128277019fcaSJed Brown       m   += mi;
128377019fcaSJed Brown     }
1284f0413b6fSBarry Smith     ierr   = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),1,m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
128577019fcaSJed Brown   } else {
12860298fd71SBarry Smith     *ltog  = NULL;
128777019fcaSJed Brown   }
128877019fcaSJed Brown   PetscFunctionReturn(0);
128977019fcaSJed Brown }
129077019fcaSJed Brown 
129177019fcaSJed Brown 
1292d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1293d8588912SDave May /*
1294d8588912SDave May   nprocessors = NP
1295d8588912SDave May   Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1))
1296d8588912SDave May        proc 0: => (g_0,h_0,)
1297d8588912SDave May        proc 1: => (g_1,h_1,)
1298d8588912SDave May        ...
1299d8588912SDave May        proc nprocs-1: => (g_NP-1,h_NP-1,)
1300d8588912SDave May 
1301d8588912SDave May             proc 0:                      proc 1:                    proc nprocs-1:
1302d8588912SDave May     is[0] = (0,1,2,...,nlocal(g_0)-1)  (0,1,...,nlocal(g_1)-1)  (0,1,...,nlocal(g_NP-1))
1303d8588912SDave May 
1304d8588912SDave May             proc 0:
1305d8588912SDave May     is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1)
1306d8588912SDave May             proc 1:
1307d8588912SDave May     is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1)
1308d8588912SDave May 
1309d8588912SDave May             proc NP-1:
1310d8588912SDave May     is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1)
1311d8588912SDave May */
1312d8588912SDave May #undef __FUNCT__
1313d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private"
1314841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1315d8588912SDave May {
1316e2d7f03fSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
13178188e55aSJed Brown   PetscInt       i,j,offset,n,nsum,bs;
1318d8588912SDave May   PetscErrorCode ierr;
13190298fd71SBarry Smith   Mat            sub = NULL;
1320d8588912SDave May 
1321d8588912SDave May   PetscFunctionBegin;
1322854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&vs->isglobal.row);CHKERRQ(ierr);
1323854ce69bSBarry Smith   ierr = PetscMalloc1(nc,&vs->isglobal.col);CHKERRQ(ierr);
1324d8588912SDave May   if (is_row) { /* valid IS is passed in */
1325d8588912SDave May     /* refs on is[] are incremeneted */
1326e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1327d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
132826fbe8dcSKarl Rupp 
1329e2d7f03fSJed Brown       vs->isglobal.row[i] = is_row[i];
1330d8588912SDave May     }
13312ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
13328188e55aSJed Brown     nsum = 0;
13338188e55aSJed Brown     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
13348188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1335ce94432eSBarry Smith       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
13360298fd71SBarry Smith       ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
1337ce94432eSBarry Smith       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
13388188e55aSJed Brown       nsum += n;
13398188e55aSJed Brown     }
1340ce94432eSBarry Smith     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
134130bc264bSJed Brown     offset -= nsum;
1342e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1343f349c1fdSJed Brown       ierr    = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
13440298fd71SBarry Smith       ierr    = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
13452ae74bdbSJed Brown       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1346ce94432eSBarry Smith       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
1347e2d7f03fSJed Brown       ierr    = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
13482ae74bdbSJed Brown       offset += n;
1349d8588912SDave May     }
1350d8588912SDave May   }
1351d8588912SDave May 
1352d8588912SDave May   if (is_col) { /* valid IS is passed in */
1353d8588912SDave May     /* refs on is[] are incremeneted */
1354e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1355d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
135626fbe8dcSKarl Rupp 
1357e2d7f03fSJed Brown       vs->isglobal.col[j] = is_col[j];
1358d8588912SDave May     }
13592ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
13602ae74bdbSJed Brown     offset = A->cmap->rstart;
13618188e55aSJed Brown     nsum   = 0;
13628188e55aSJed Brown     for (j=0; j<vs->nc; j++) {
13638188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
1364ce94432eSBarry Smith       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
13650298fd71SBarry Smith       ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
1366ce94432eSBarry Smith       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
13678188e55aSJed Brown       nsum += n;
13688188e55aSJed Brown     }
1369ce94432eSBarry Smith     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
137030bc264bSJed Brown     offset -= nsum;
1371e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1372f349c1fdSJed Brown       ierr    = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
13730298fd71SBarry Smith       ierr    = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
13742ae74bdbSJed Brown       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1375ce94432eSBarry Smith       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
1376e2d7f03fSJed Brown       ierr    = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
13772ae74bdbSJed Brown       offset += n;
1378d8588912SDave May     }
1379d8588912SDave May   }
1380e2d7f03fSJed Brown 
1381e2d7f03fSJed Brown   /* Set up the local ISs */
1382785e854fSJed Brown   ierr = PetscMalloc1(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
1383785e854fSJed Brown   ierr = PetscMalloc1(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
1384e2d7f03fSJed Brown   for (i=0,offset=0; i<vs->nr; i++) {
1385e2d7f03fSJed Brown     IS                     isloc;
13860298fd71SBarry Smith     ISLocalToGlobalMapping rmap = NULL;
1387e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1388e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
13890298fd71SBarry Smith     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,NULL);CHKERRQ(ierr);}
1390207556f9SJed Brown     if (rmap) {
1391e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1392e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
1393e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1394e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1395207556f9SJed Brown     } else {
1396207556f9SJed Brown       nlocal = 0;
13970298fd71SBarry Smith       isloc  = NULL;
1398207556f9SJed Brown     }
1399e2d7f03fSJed Brown     vs->islocal.row[i] = isloc;
1400e2d7f03fSJed Brown     offset            += nlocal;
1401e2d7f03fSJed Brown   }
14028188e55aSJed Brown   for (i=0,offset=0; i<vs->nc; i++) {
1403e2d7f03fSJed Brown     IS                     isloc;
14040298fd71SBarry Smith     ISLocalToGlobalMapping cmap = NULL;
1405e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1406e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
14070298fd71SBarry Smith     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,NULL,&cmap);CHKERRQ(ierr);}
1408207556f9SJed Brown     if (cmap) {
1409e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1410e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr);
1411e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1412e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1413207556f9SJed Brown     } else {
1414207556f9SJed Brown       nlocal = 0;
14150298fd71SBarry Smith       isloc  = NULL;
1416207556f9SJed Brown     }
1417e2d7f03fSJed Brown     vs->islocal.col[i] = isloc;
1418e2d7f03fSJed Brown     offset            += nlocal;
1419e2d7f03fSJed Brown   }
14200189643fSJed Brown 
142177019fcaSJed Brown   /* Set up the aggregate ISLocalToGlobalMapping */
142277019fcaSJed Brown   {
142345b6f7e9SBarry Smith     ISLocalToGlobalMapping rmap,cmap;
142445b6f7e9SBarry Smith     ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap);CHKERRQ(ierr);
142545b6f7e9SBarry Smith     ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap);CHKERRQ(ierr);
142677019fcaSJed Brown     if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);}
142777019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
142877019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
142977019fcaSJed Brown   }
143077019fcaSJed Brown 
14310189643fSJed Brown #if defined(PETSC_USE_DEBUG)
14320189643fSJed Brown   for (i=0; i<vs->nr; i++) {
14330189643fSJed Brown     for (j=0; j<vs->nc; j++) {
14340189643fSJed Brown       PetscInt m,n,M,N,mi,ni,Mi,Ni;
14350189643fSJed Brown       Mat      B = vs->m[i][j];
14360189643fSJed Brown       if (!B) continue;
14370189643fSJed Brown       ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr);
14380189643fSJed Brown       ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr);
14390189643fSJed Brown       ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr);
14400189643fSJed Brown       ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr);
14410189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr);
14420189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr);
1443ce94432eSBarry 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);
1444ce94432eSBarry 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);
14450189643fSJed Brown     }
14460189643fSJed Brown   }
14470189643fSJed Brown #endif
1448a061e289SJed Brown 
1449a061e289SJed Brown   /* Set A->assembled if all non-null blocks are currently assembled */
1450a061e289SJed Brown   for (i=0; i<vs->nr; i++) {
1451a061e289SJed Brown     for (j=0; j<vs->nc; j++) {
1452a061e289SJed Brown       if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0);
1453a061e289SJed Brown     }
1454a061e289SJed Brown   }
1455a061e289SJed Brown   A->assembled = PETSC_TRUE;
1456d8588912SDave May   PetscFunctionReturn(0);
1457d8588912SDave May }
1458d8588912SDave May 
1459d8588912SDave May #undef __FUNCT__
1460d8588912SDave May #define __FUNCT__ "MatCreateNest"
146145c38901SJed Brown /*@C
1462659c6bb0SJed Brown    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
1463659c6bb0SJed Brown 
1464659c6bb0SJed Brown    Collective on Mat
1465659c6bb0SJed Brown 
1466659c6bb0SJed Brown    Input Parameter:
1467659c6bb0SJed Brown +  comm - Communicator for the new Mat
1468659c6bb0SJed Brown .  nr - number of nested row blocks
14690298fd71SBarry Smith .  is_row - index sets for each nested row block, or NULL to make contiguous
1470659c6bb0SJed Brown .  nc - number of nested column blocks
14710298fd71SBarry Smith .  is_col - index sets for each nested column block, or NULL to make contiguous
14720298fd71SBarry Smith -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1473659c6bb0SJed Brown 
1474659c6bb0SJed Brown    Output Parameter:
1475659c6bb0SJed Brown .  B - new matrix
1476659c6bb0SJed Brown 
1477659c6bb0SJed Brown    Level: advanced
1478659c6bb0SJed Brown 
1479950540a4SJed Brown .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST
1480659c6bb0SJed Brown @*/
14817087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1482d8588912SDave May {
1483d8588912SDave May   Mat            A;
1484d8588912SDave May   PetscErrorCode ierr;
1485d8588912SDave May 
1486d8588912SDave May   PetscFunctionBegin;
1487c8883902SJed Brown   *B   = 0;
1488d8588912SDave May   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1489c8883902SJed Brown   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
149091a28eb3SBarry Smith   A->preallocated = PETSC_TRUE;
1491c8883902SJed Brown   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1492d8588912SDave May   *B   = A;
1493d8588912SDave May   PetscFunctionReturn(0);
1494d8588912SDave May }
1495659c6bb0SJed Brown 
1496629c3df2SDmitry Karpeev #undef __FUNCT__
1497629c3df2SDmitry Karpeev #define __FUNCT__ "MatConvert_Nest_AIJ"
1498cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1499629c3df2SDmitry Karpeev {
1500629c3df2SDmitry Karpeev   PetscErrorCode ierr;
1501629c3df2SDmitry Karpeev   Mat_Nest       *nest = (Mat_Nest*)A->data;
150283b1a929SMark Adams   PetscInt       m,n,M,N,i,j,k,*dnnz,*onnz,rstart;
1503649b366bSFande Kong   PetscInt       cstart,cend;
1504629c3df2SDmitry Karpeev   Mat            C;
1505629c3df2SDmitry Karpeev 
1506629c3df2SDmitry Karpeev   PetscFunctionBegin;
1507629c3df2SDmitry Karpeev   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
1508629c3df2SDmitry Karpeev   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
1509649b366bSFande Kong   ierr = MatGetOwnershipRangeColumn(A,&cstart,&cend);CHKERRQ(ierr);
1510629c3df2SDmitry Karpeev   switch (reuse) {
1511629c3df2SDmitry Karpeev   case MAT_INITIAL_MATRIX:
1512ce94432eSBarry Smith     ierr    = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
1513629c3df2SDmitry Karpeev     ierr    = MatSetType(C,newtype);CHKERRQ(ierr);
1514629c3df2SDmitry Karpeev     ierr    = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr);
1515629c3df2SDmitry Karpeev     *newmat = C;
1516629c3df2SDmitry Karpeev     break;
1517629c3df2SDmitry Karpeev   case MAT_REUSE_MATRIX:
1518629c3df2SDmitry Karpeev     C = *newmat;
1519629c3df2SDmitry Karpeev     break;
1520ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse");
1521629c3df2SDmitry Karpeev   }
1522785e854fSJed Brown   ierr = PetscMalloc1(2*m,&dnnz);CHKERRQ(ierr);
1523629c3df2SDmitry Karpeev   onnz = dnnz + m;
1524629c3df2SDmitry Karpeev   for (k=0; k<m; k++) {
1525629c3df2SDmitry Karpeev     dnnz[k] = 0;
1526629c3df2SDmitry Karpeev     onnz[k] = 0;
1527629c3df2SDmitry Karpeev   }
1528629c3df2SDmitry Karpeev   for (j=0; j<nest->nc; ++j) {
1529629c3df2SDmitry Karpeev     IS             bNis;
1530629c3df2SDmitry Karpeev     PetscInt       bN;
1531629c3df2SDmitry Karpeev     const PetscInt *bNindices;
1532629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
1533629c3df2SDmitry Karpeev     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1534629c3df2SDmitry Karpeev     ierr = ISGetSize(bNis, &bN);CHKERRQ(ierr);
1535629c3df2SDmitry Karpeev     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1536629c3df2SDmitry Karpeev     for (i=0; i<nest->nr; ++i) {
1537629c3df2SDmitry Karpeev       PetscSF        bmsf;
1538649b366bSFande Kong       PetscSFNode    *iremote;
1539629c3df2SDmitry Karpeev       Mat            B;
1540649b366bSFande Kong       PetscInt       bm, *sub_dnnz,*sub_onnz, br;
1541629c3df2SDmitry Karpeev       const PetscInt *bmindices;
1542629c3df2SDmitry Karpeev       B = nest->m[i][j];
1543629c3df2SDmitry Karpeev       if (!B) continue;
1544629c3df2SDmitry Karpeev       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1545629c3df2SDmitry Karpeev       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1546ce94432eSBarry Smith       ierr = PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);CHKERRQ(ierr);
1547649b366bSFande Kong       ierr = PetscMalloc1(bm,&iremote);CHKERRQ(ierr);
1548649b366bSFande Kong       ierr = PetscMalloc1(bm,&sub_dnnz);CHKERRQ(ierr);
1549649b366bSFande Kong       ierr = PetscMalloc1(bm,&sub_onnz);CHKERRQ(ierr);
1550649b366bSFande Kong       for (k = 0; k < bm; ++k){
1551649b366bSFande Kong     	sub_dnnz[k] = 0;
1552649b366bSFande Kong     	sub_onnz[k] = 0;
1553649b366bSFande Kong       }
1554629c3df2SDmitry Karpeev       /*
1555629c3df2SDmitry Karpeev        Locate the owners for all of the locally-owned global row indices for this row block.
1556629c3df2SDmitry Karpeev        These determine the roots of PetscSF used to communicate preallocation data to row owners.
1557629c3df2SDmitry Karpeev        The roots correspond to the dnnz and onnz entries; thus, there are two roots per row.
1558629c3df2SDmitry Karpeev        */
155983b1a929SMark Adams       ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr);
1560629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
1561649b366bSFande Kong         PetscInt       row = bmindices[br], rowowner = 0, brncols, col;
1562629c3df2SDmitry Karpeev         const PetscInt *brcols;
1563a4b3d3acSMatthew G Knepley         PetscInt       rowrel = 0; /* row's relative index on its owner rank */
1564629c3df2SDmitry Karpeev         ierr      = PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);CHKERRQ(ierr);
1565649b366bSFande Kong         /* how many roots  */
1566649b366bSFande Kong         iremote[br].rank = rowowner; iremote[br].index = rowrel;           /* edge from bmdnnz to dnnz */
1567649b366bSFande Kong         /* get nonzero pattern */
156883b1a929SMark Adams         ierr = MatGetRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr);
1569629c3df2SDmitry Karpeev         for (k=0; k<brncols; k++) {
1570629c3df2SDmitry Karpeev           col  = bNindices[brcols[k]];
1571649b366bSFande Kong           if (col>=A->cmap->range[rowowner] && col<A->cmap->range[rowowner+1]) {
1572649b366bSFande Kong             sub_dnnz[br]++;
1573649b366bSFande Kong           } else {
1574649b366bSFande Kong             sub_onnz[br]++;
1575649b366bSFande Kong           }
1576629c3df2SDmitry Karpeev         }
157783b1a929SMark Adams         ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr);
1578629c3df2SDmitry Karpeev       }
1579629c3df2SDmitry Karpeev       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1580629c3df2SDmitry Karpeev       /* bsf will have to take care of disposing of bedges. */
1581649b366bSFande Kong       ierr = PetscSFSetGraph(bmsf,m,bm,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1582649b366bSFande Kong       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr);
1583649b366bSFande Kong       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr);
1584649b366bSFande Kong       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr);
1585649b366bSFande Kong       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr);
1586649b366bSFande Kong       ierr = PetscFree(sub_dnnz);CHKERRQ(ierr);
1587649b366bSFande Kong       ierr = PetscFree(sub_onnz);CHKERRQ(ierr);
1588629c3df2SDmitry Karpeev       ierr = PetscSFDestroy(&bmsf);CHKERRQ(ierr);
1589629c3df2SDmitry Karpeev     }
159022d28d08SBarry Smith     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
1591629c3df2SDmitry Karpeev     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
1592*65a4a0a3Sstefano_zampini   }
1593*65a4a0a3Sstefano_zampini   /* Resize preallocation if overestimated */
1594*65a4a0a3Sstefano_zampini   for (i=0;i<m;i++) {
1595*65a4a0a3Sstefano_zampini     dnnz[i] = PetscMin(dnnz[i],A->cmap->n);
1596*65a4a0a3Sstefano_zampini     onnz[i] = PetscMin(onnz[i],A->cmap->N - A->cmap->n);
1597629c3df2SDmitry Karpeev   }
1598629c3df2SDmitry Karpeev   ierr = MatSeqAIJSetPreallocation(C,0,dnnz);CHKERRQ(ierr);
1599629c3df2SDmitry Karpeev   ierr = MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);CHKERRQ(ierr);
1600629c3df2SDmitry Karpeev   ierr = PetscFree(dnnz);CHKERRQ(ierr);
1601629c3df2SDmitry Karpeev 
1602629c3df2SDmitry Karpeev   /* Fill by row */
1603629c3df2SDmitry Karpeev   for (j=0; j<nest->nc; ++j) {
1604629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
1605629c3df2SDmitry Karpeev     IS             bNis;
1606629c3df2SDmitry Karpeev     PetscInt       bN;
1607629c3df2SDmitry Karpeev     const PetscInt *bNindices;
1608629c3df2SDmitry Karpeev     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1609629c3df2SDmitry Karpeev     ierr = ISGetSize(bNis,&bN);CHKERRQ(ierr);
1610629c3df2SDmitry Karpeev     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1611629c3df2SDmitry Karpeev     for (i=0; i<nest->nr; ++i) {
1612629c3df2SDmitry Karpeev       Mat            B;
1613629c3df2SDmitry Karpeev       PetscInt       bm, br;
1614629c3df2SDmitry Karpeev       const PetscInt *bmindices;
1615629c3df2SDmitry Karpeev       B = nest->m[i][j];
1616629c3df2SDmitry Karpeev       if (!B) continue;
1617629c3df2SDmitry Karpeev       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1618629c3df2SDmitry Karpeev       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
161983b1a929SMark Adams       ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr);
1620629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
1621629c3df2SDmitry Karpeev         PetscInt          row = bmindices[br], brncols,  *cols;
1622629c3df2SDmitry Karpeev         const PetscInt    *brcols;
1623629c3df2SDmitry Karpeev         const PetscScalar *brcoldata;
162483b1a929SMark Adams         ierr = MatGetRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
1625785e854fSJed Brown         ierr = PetscMalloc1(brncols,&cols);CHKERRQ(ierr);
162626fbe8dcSKarl Rupp         for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]];
1627629c3df2SDmitry Karpeev         /*
1628629c3df2SDmitry Karpeev           Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match.
1629629c3df2SDmitry Karpeev           Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES.
1630629c3df2SDmitry Karpeev          */
1631a2ea699eSBarry Smith         ierr = MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);CHKERRQ(ierr);
163283b1a929SMark Adams         ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
1633629c3df2SDmitry Karpeev         ierr = PetscFree(cols);CHKERRQ(ierr);
1634629c3df2SDmitry Karpeev       }
1635629c3df2SDmitry Karpeev       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1636629c3df2SDmitry Karpeev     }
1637a2ea699eSBarry Smith     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
1638629c3df2SDmitry Karpeev     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
1639629c3df2SDmitry Karpeev   }
1640629c3df2SDmitry Karpeev   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1641629c3df2SDmitry Karpeev   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1642629c3df2SDmitry Karpeev   PetscFunctionReturn(0);
1643629c3df2SDmitry Karpeev }
1644629c3df2SDmitry Karpeev 
1645659c6bb0SJed Brown /*MC
1646659c6bb0SJed Brown   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
1647659c6bb0SJed Brown 
1648659c6bb0SJed Brown   Level: intermediate
1649659c6bb0SJed Brown 
1650659c6bb0SJed Brown   Notes:
1651659c6bb0SJed Brown   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1652659c6bb0SJed Brown   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1653950540a4SJed Brown   It is usually used with DMComposite and DMCreateMatrix()
1654659c6bb0SJed Brown 
1655659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest()
1656659c6bb0SJed Brown M*/
1657c8883902SJed Brown #undef __FUNCT__
1658c8883902SJed Brown #define __FUNCT__ "MatCreate_Nest"
16598cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A)
1660c8883902SJed Brown {
1661c8883902SJed Brown   Mat_Nest       *s;
1662c8883902SJed Brown   PetscErrorCode ierr;
1663c8883902SJed Brown 
1664c8883902SJed Brown   PetscFunctionBegin;
1665b00a9115SJed Brown   ierr    = PetscNewLog(A,&s);CHKERRQ(ierr);
1666c8883902SJed Brown   A->data = (void*)s;
1667e7c19651SJed Brown 
1668e7c19651SJed Brown   s->nr            = -1;
1669e7c19651SJed Brown   s->nc            = -1;
16700298fd71SBarry Smith   s->m             = NULL;
1671e7c19651SJed Brown   s->splitassembly = PETSC_FALSE;
1672c8883902SJed Brown 
1673c8883902SJed Brown   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
167426fbe8dcSKarl Rupp 
1675c8883902SJed Brown   A->ops->mult                  = MatMult_Nest;
16769194d70fSJed Brown   A->ops->multadd               = MatMultAdd_Nest;
1677c8883902SJed Brown   A->ops->multtranspose         = MatMultTranspose_Nest;
16789194d70fSJed Brown   A->ops->multtransposeadd      = MatMultTransposeAdd_Nest;
1679f8170845SAlex Fikl   A->ops->transpose             = MatTranspose_Nest;
1680c8883902SJed Brown   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1681c8883902SJed Brown   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1682c8883902SJed Brown   A->ops->zeroentries           = MatZeroEntries_Nest;
1683c222c20dSDavid Ham   A->ops->copy                  = MatCopy_Nest;
1684c8883902SJed Brown   A->ops->duplicate             = MatDuplicate_Nest;
1685c8883902SJed Brown   A->ops->getsubmatrix          = MatGetSubMatrix_Nest;
1686c8883902SJed Brown   A->ops->destroy               = MatDestroy_Nest;
1687c8883902SJed Brown   A->ops->view                  = MatView_Nest;
1688c8883902SJed Brown   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1689c8883902SJed Brown   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1690c8883902SJed Brown   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
1691429bac76SJed Brown   A->ops->getdiagonal           = MatGetDiagonal_Nest;
1692429bac76SJed Brown   A->ops->diagonalscale         = MatDiagonalScale_Nest;
1693a061e289SJed Brown   A->ops->scale                 = MatScale_Nest;
1694a061e289SJed Brown   A->ops->shift                 = MatShift_Nest;
169513135bc6SAlex Fikl   A->ops->diagonalset           = MatDiagonalSet_Nest;
1696f8170845SAlex Fikl   A->ops->setrandom             = MatSetRandom_Nest;
1697c8883902SJed Brown 
1698c8883902SJed Brown   A->spptr        = 0;
1699c8883902SJed Brown   A->assembled    = PETSC_FALSE;
1700c8883902SJed Brown 
1701c8883902SJed Brown   /* expose Nest api's */
1702bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",   MatNestGetSubMat_Nest);CHKERRQ(ierr);
1703bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",   MatNestSetSubMat_Nest);CHKERRQ(ierr);
1704bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",  MatNestGetSubMats_Nest);CHKERRQ(ierr);
1705bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",     MatNestGetSize_Nest);CHKERRQ(ierr);
1706bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",      MatNestGetISs_Nest);CHKERRQ(ierr);
1707bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C", MatNestGetLocalISs_Nest);CHKERRQ(ierr);
1708bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",  MatNestSetVecType_Nest);CHKERRQ(ierr);
1709bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",  MatNestSetSubMats_Nest);CHKERRQ(ierr);
171083b1a929SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",MatConvert_Nest_AIJ);CHKERRQ(ierr);
17115e3038f0Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_is_C", MatConvert_Nest_IS);CHKERRQ(ierr);
1712c8883902SJed Brown 
1713c8883902SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
1714c8883902SJed Brown   PetscFunctionReturn(0);
1715c8883902SJed Brown }
1716