xref: /petsc/src/mat/impls/nest/matnest.c (revision b68353e597724b2b70ae3910f8be994f0f1fe2e3)
1d8588912SDave May 
2aaa7dc30SBarry Smith #include <../src/mat/impls/nest/matnestimpl.h> /*I   "petscmat.h"   I*/
3*b68353e5Sstefano_zampini #include <../src/mat/impls/aij/seq/aij.h>
40c312b8eSJed Brown #include <petscsf.h>
5d8588912SDave May 
6c8883902SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]);
72a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left);
85e3038f0Sstefano_zampini PETSC_INTERN PetscErrorCode MatConvert_Nest_IS(Mat,MatType,MatReuse,Mat*);
9c8883902SJed Brown 
10d8588912SDave May /* private functions */
11d8588912SDave May #undef __FUNCT__
128188e55aSJed Brown #define __FUNCT__ "MatNestGetSizes_Private"
138188e55aSJed Brown static PetscErrorCode MatNestGetSizes_Private(Mat A,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N)
14d8588912SDave May {
15d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
168188e55aSJed Brown   PetscInt       i,j;
17d8588912SDave May   PetscErrorCode ierr;
18d8588912SDave May 
19d8588912SDave May   PetscFunctionBegin;
208188e55aSJed Brown   *m = *n = *M = *N = 0;
218188e55aSJed Brown   for (i=0; i<bA->nr; i++) {  /* rows */
228188e55aSJed Brown     PetscInt sm,sM;
238188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.row[i],&sm);CHKERRQ(ierr);
248188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.row[i],&sM);CHKERRQ(ierr);
258188e55aSJed Brown     *m  += sm;
268188e55aSJed Brown     *M  += sM;
27d8588912SDave May   }
288188e55aSJed Brown   for (j=0; j<bA->nc; j++) {  /* cols */
298188e55aSJed Brown     PetscInt sn,sN;
308188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.col[j],&sn);CHKERRQ(ierr);
318188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.col[j],&sN);CHKERRQ(ierr);
328188e55aSJed Brown     *n  += sn;
338188e55aSJed Brown     *N  += sN;
34d8588912SDave May   }
35d8588912SDave May   PetscFunctionReturn(0);
36d8588912SDave May }
37d8588912SDave May 
38d8588912SDave May /* operations */
39d8588912SDave May #undef __FUNCT__
40d8588912SDave May #define __FUNCT__ "MatMult_Nest"
41207556f9SJed Brown static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y)
42d8588912SDave May {
43d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
44207556f9SJed Brown   Vec            *bx = bA->right,*by = bA->left;
45207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
46d8588912SDave May   PetscErrorCode ierr;
47d8588912SDave May 
48d8588912SDave May   PetscFunctionBegin;
49207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
50207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
51207556f9SJed Brown   for (i=0; i<nr; i++) {
52d8588912SDave May     ierr = VecZeroEntries(by[i]);CHKERRQ(ierr);
53207556f9SJed Brown     for (j=0; j<nc; j++) {
54207556f9SJed Brown       if (!bA->m[i][j]) continue;
55d8588912SDave May       /* y[i] <- y[i] + A[i][j] * x[j] */
56d8588912SDave May       ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr);
57d8588912SDave May     }
58d8588912SDave May   }
59207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
60207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
61d8588912SDave May   PetscFunctionReturn(0);
62d8588912SDave May }
63d8588912SDave May 
64d8588912SDave May #undef __FUNCT__
659194d70fSJed Brown #define __FUNCT__ "MatMultAdd_Nest"
669194d70fSJed Brown static PetscErrorCode MatMultAdd_Nest(Mat A,Vec x,Vec y,Vec z)
679194d70fSJed Brown {
689194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
699194d70fSJed Brown   Vec            *bx = bA->right,*bz = bA->left;
709194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
719194d70fSJed Brown   PetscErrorCode ierr;
729194d70fSJed Brown 
739194d70fSJed Brown   PetscFunctionBegin;
749194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
759194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
769194d70fSJed Brown   for (i=0; i<nr; i++) {
779194d70fSJed Brown     if (y != z) {
789194d70fSJed Brown       Vec by;
799194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
809194d70fSJed Brown       ierr = VecCopy(by,bz[i]);CHKERRQ(ierr);
81336d21e7SJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
829194d70fSJed Brown     }
839194d70fSJed Brown     for (j=0; j<nc; j++) {
849194d70fSJed Brown       if (!bA->m[i][j]) continue;
859194d70fSJed Brown       /* y[i] <- y[i] + A[i][j] * x[j] */
869194d70fSJed Brown       ierr = MatMultAdd(bA->m[i][j],bx[j],bz[i],bz[i]);CHKERRQ(ierr);
879194d70fSJed Brown     }
889194d70fSJed Brown   }
899194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
909194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
919194d70fSJed Brown   PetscFunctionReturn(0);
929194d70fSJed Brown }
939194d70fSJed Brown 
949194d70fSJed Brown #undef __FUNCT__
95d8588912SDave May #define __FUNCT__ "MatMultTranspose_Nest"
96207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y)
97d8588912SDave May {
98d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
99207556f9SJed Brown   Vec            *bx = bA->left,*by = bA->right;
100207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
101d8588912SDave May   PetscErrorCode ierr;
102d8588912SDave May 
103d8588912SDave May   PetscFunctionBegin;
104609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
105609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
106207556f9SJed Brown   for (j=0; j<nc; j++) {
107609e31cbSJed Brown     ierr = VecZeroEntries(by[j]);CHKERRQ(ierr);
108609e31cbSJed Brown     for (i=0; i<nr; i++) {
1096c75ac25SJed Brown       if (!bA->m[i][j]) continue;
110609e31cbSJed Brown       /* y[j] <- y[j] + (A[i][j])^T * x[i] */
111609e31cbSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],by[j],by[j]);CHKERRQ(ierr);
112d8588912SDave May     }
113d8588912SDave May   }
114609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
115609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
116d8588912SDave May   PetscFunctionReturn(0);
117d8588912SDave May }
118d8588912SDave May 
119d8588912SDave May #undef __FUNCT__
1209194d70fSJed Brown #define __FUNCT__ "MatMultTransposeAdd_Nest"
1219194d70fSJed Brown static PetscErrorCode MatMultTransposeAdd_Nest(Mat A,Vec x,Vec y,Vec z)
1229194d70fSJed Brown {
1239194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
1249194d70fSJed Brown   Vec            *bx = bA->left,*bz = bA->right;
1259194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
1269194d70fSJed Brown   PetscErrorCode ierr;
1279194d70fSJed Brown 
1289194d70fSJed Brown   PetscFunctionBegin;
1299194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
1309194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
1319194d70fSJed Brown   for (j=0; j<nc; j++) {
1329194d70fSJed Brown     if (y != z) {
1339194d70fSJed Brown       Vec by;
1349194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
1359194d70fSJed Brown       ierr = VecCopy(by,bz[j]);CHKERRQ(ierr);
1369194d70fSJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
1379194d70fSJed Brown     }
1389194d70fSJed Brown     for (i=0; i<nr; i++) {
1396c75ac25SJed Brown       if (!bA->m[i][j]) continue;
1409194d70fSJed Brown       /* z[j] <- y[j] + (A[i][j])^T * x[i] */
1419194d70fSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],bz[j],bz[j]);CHKERRQ(ierr);
1429194d70fSJed Brown     }
1439194d70fSJed Brown   }
1449194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
1459194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
1469194d70fSJed Brown   PetscFunctionReturn(0);
1479194d70fSJed Brown }
1489194d70fSJed Brown 
1499194d70fSJed Brown #undef __FUNCT__
150f8170845SAlex Fikl #define __FUNCT__ "MatTranspose_Nest"
151f8170845SAlex Fikl static PetscErrorCode MatTranspose_Nest(Mat A,MatReuse reuse,Mat *B)
152f8170845SAlex Fikl {
153f8170845SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data, *bC;
154f8170845SAlex Fikl   Mat            C;
155f8170845SAlex Fikl   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
156f8170845SAlex Fikl   PetscErrorCode ierr;
157f8170845SAlex Fikl 
158f8170845SAlex Fikl   PetscFunctionBegin;
159f8170845SAlex 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");
160f8170845SAlex Fikl 
161f8170845SAlex Fikl   if (reuse == MAT_INITIAL_MATRIX || A == *B) {
162f8170845SAlex Fikl     Mat *subs;
163f8170845SAlex Fikl     IS  *is_row,*is_col;
164f8170845SAlex Fikl 
165f8170845SAlex Fikl     ierr = PetscCalloc1(nr * nc,&subs);CHKERRQ(ierr);
166f8170845SAlex Fikl     ierr = PetscMalloc2(nr,&is_row,nc,&is_col);CHKERRQ(ierr);
167f8170845SAlex Fikl     ierr = MatNestGetISs(A,is_row,is_col);CHKERRQ(ierr);
168ddeb9bd8SAlex Fikl     if (reuse == MAT_REUSE_MATRIX) {
169ddeb9bd8SAlex Fikl       for (i=0; i<nr; i++) {
170ddeb9bd8SAlex Fikl         for (j=0; j<nc; j++) {
171ddeb9bd8SAlex Fikl           subs[i + nr * j] = bA->m[i][j];
172ddeb9bd8SAlex Fikl         }
173ddeb9bd8SAlex Fikl       }
174ddeb9bd8SAlex Fikl     }
175ddeb9bd8SAlex Fikl 
176f8170845SAlex Fikl     ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nc,is_col,nr,is_row,subs,&C);CHKERRQ(ierr);
177f8170845SAlex Fikl     ierr = PetscFree(subs);CHKERRQ(ierr);
1783d994f23SBarry Smith     ierr = PetscFree2(is_row,is_col);CHKERRQ(ierr);
179f8170845SAlex Fikl   } else {
180f8170845SAlex Fikl     C = *B;
181f8170845SAlex Fikl   }
182f8170845SAlex Fikl 
183f8170845SAlex Fikl   bC = (Mat_Nest*)C->data;
184f8170845SAlex Fikl   for (i=0; i<nr; i++) {
185f8170845SAlex Fikl     for (j=0; j<nc; j++) {
186f8170845SAlex Fikl       if (bA->m[i][j]) {
187f8170845SAlex Fikl         ierr = MatTranspose(bA->m[i][j], reuse, &(bC->m[j][i]));CHKERRQ(ierr);
188f8170845SAlex Fikl       } else {
189f8170845SAlex Fikl         bC->m[j][i] = NULL;
190f8170845SAlex Fikl       }
191f8170845SAlex Fikl     }
192f8170845SAlex Fikl   }
193f8170845SAlex Fikl 
194f8170845SAlex Fikl   if (reuse == MAT_INITIAL_MATRIX || A != *B) {
195f8170845SAlex Fikl     *B = C;
196f8170845SAlex Fikl   } else {
197f8170845SAlex Fikl     ierr = MatHeaderMerge(A, &C);CHKERRQ(ierr);
198f8170845SAlex Fikl   }
199f8170845SAlex Fikl   PetscFunctionReturn(0);
200f8170845SAlex Fikl }
201f8170845SAlex Fikl 
202f8170845SAlex Fikl #undef __FUNCT__
203e2d7f03fSJed Brown #define __FUNCT__ "MatNestDestroyISList"
204e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list)
205e2d7f03fSJed Brown {
206e2d7f03fSJed Brown   PetscErrorCode ierr;
207e2d7f03fSJed Brown   IS             *lst = *list;
208e2d7f03fSJed Brown   PetscInt       i;
209e2d7f03fSJed Brown 
210e2d7f03fSJed Brown   PetscFunctionBegin;
211e2d7f03fSJed Brown   if (!lst) PetscFunctionReturn(0);
2126bf464f9SBarry Smith   for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);}
213e2d7f03fSJed Brown   ierr  = PetscFree(lst);CHKERRQ(ierr);
2140298fd71SBarry Smith   *list = NULL;
215e2d7f03fSJed Brown   PetscFunctionReturn(0);
216e2d7f03fSJed Brown }
217e2d7f03fSJed Brown 
218e2d7f03fSJed Brown #undef __FUNCT__
219d8588912SDave May #define __FUNCT__ "MatDestroy_Nest"
220207556f9SJed Brown static PetscErrorCode MatDestroy_Nest(Mat A)
221d8588912SDave May {
222d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
223d8588912SDave May   PetscInt       i,j;
224d8588912SDave May   PetscErrorCode ierr;
225d8588912SDave May 
226d8588912SDave May   PetscFunctionBegin;
227d8588912SDave May   /* release the matrices and the place holders */
228e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr);
229e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr);
230e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
231e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
232d8588912SDave May 
233d8588912SDave May   ierr = PetscFree(vs->row_len);CHKERRQ(ierr);
234d8588912SDave May   ierr = PetscFree(vs->col_len);CHKERRQ(ierr);
235d8588912SDave May 
236207556f9SJed Brown   ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr);
237207556f9SJed Brown 
238d8588912SDave May   /* release the matrices and the place holders */
239d8588912SDave May   if (vs->m) {
240d8588912SDave May     for (i=0; i<vs->nr; i++) {
241d8588912SDave May       for (j=0; j<vs->nc; j++) {
2426bf464f9SBarry Smith         ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr);
243d8588912SDave May       }
244d8588912SDave May       ierr = PetscFree(vs->m[i]);CHKERRQ(ierr);
245d8588912SDave May     }
246d8588912SDave May     ierr = PetscFree(vs->m);CHKERRQ(ierr);
247d8588912SDave May   }
248bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
249d8588912SDave May 
250bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",0);CHKERRQ(ierr);
251bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",0);CHKERRQ(ierr);
252bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",0);CHKERRQ(ierr);
253bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",0);CHKERRQ(ierr);
254bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",0);CHKERRQ(ierr);
255bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",0);CHKERRQ(ierr);
256bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",0);CHKERRQ(ierr);
257bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",0);CHKERRQ(ierr);
2585e3038f0Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",0);CHKERRQ(ierr);
2595e3038f0Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_is_C",0);CHKERRQ(ierr);
260d8588912SDave May   PetscFunctionReturn(0);
261d8588912SDave May }
262d8588912SDave May 
263d8588912SDave May #undef __FUNCT__
264d8588912SDave May #define __FUNCT__ "MatAssemblyBegin_Nest"
265207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
266d8588912SDave May {
267d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
268d8588912SDave May   PetscInt       i,j;
269d8588912SDave May   PetscErrorCode ierr;
270d8588912SDave May 
271d8588912SDave May   PetscFunctionBegin;
272d8588912SDave May   for (i=0; i<vs->nr; i++) {
273d8588912SDave May     for (j=0; j<vs->nc; j++) {
274e7c19651SJed Brown       if (vs->m[i][j]) {
275e7c19651SJed Brown         ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr);
276e7c19651SJed Brown         if (!vs->splitassembly) {
277e7c19651SJed Brown           /* Note: split assembly will fail if the same block appears more than once (even indirectly through a nested
278e7c19651SJed 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
279e7c19651SJed Brown            * already performing an assembly, but the result would by more complicated and appears to offer less
280e7c19651SJed Brown            * potential for diagnostics and correctness checking. Split assembly should be fixed once there is an
281e7c19651SJed Brown            * interface for libraries to make asynchronous progress in "user-defined non-blocking collectives".
282e7c19651SJed Brown            */
283e7c19651SJed Brown           ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr);
284e7c19651SJed Brown         }
285e7c19651SJed Brown       }
286d8588912SDave May     }
287d8588912SDave May   }
288d8588912SDave May   PetscFunctionReturn(0);
289d8588912SDave May }
290d8588912SDave May 
291d8588912SDave May #undef __FUNCT__
292d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest"
293207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
294d8588912SDave May {
295d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
296d8588912SDave May   PetscInt       i,j;
297d8588912SDave May   PetscErrorCode ierr;
298d8588912SDave May 
299d8588912SDave May   PetscFunctionBegin;
300d8588912SDave May   for (i=0; i<vs->nr; i++) {
301d8588912SDave May     for (j=0; j<vs->nc; j++) {
302e7c19651SJed Brown       if (vs->m[i][j]) {
303e7c19651SJed Brown         if (vs->splitassembly) {
304e7c19651SJed Brown           ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr);
305e7c19651SJed Brown         }
306e7c19651SJed Brown       }
307d8588912SDave May     }
308d8588912SDave May   }
309d8588912SDave May   PetscFunctionReturn(0);
310d8588912SDave May }
311d8588912SDave May 
312d8588912SDave May #undef __FUNCT__
313f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatRow"
314f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
315d8588912SDave May {
316207556f9SJed Brown   PetscErrorCode ierr;
317f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
318f349c1fdSJed Brown   PetscInt       j;
319f349c1fdSJed Brown   Mat            sub;
320d8588912SDave May 
321d8588912SDave May   PetscFunctionBegin;
3220298fd71SBarry Smith   sub = (row < vs->nc) ? vs->m[row][row] : (Mat)NULL; /* Prefer to find on the diagonal */
323f349c1fdSJed Brown   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
3244994cf47SJed Brown   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
325f349c1fdSJed Brown   *B = sub;
326f349c1fdSJed Brown   PetscFunctionReturn(0);
327d8588912SDave May }
328d8588912SDave May 
329f349c1fdSJed Brown #undef __FUNCT__
330f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatCol"
331f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
332f349c1fdSJed Brown {
333207556f9SJed Brown   PetscErrorCode ierr;
334f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
335f349c1fdSJed Brown   PetscInt       i;
336f349c1fdSJed Brown   Mat            sub;
337f349c1fdSJed Brown 
338f349c1fdSJed Brown   PetscFunctionBegin;
3390298fd71SBarry Smith   sub = (col < vs->nr) ? vs->m[col][col] : (Mat)NULL; /* Prefer to find on the diagonal */
340f349c1fdSJed Brown   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
3414994cf47SJed Brown   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
342f349c1fdSJed Brown   *B = sub;
343f349c1fdSJed Brown   PetscFunctionReturn(0);
344d8588912SDave May }
345d8588912SDave May 
346f349c1fdSJed Brown #undef __FUNCT__
347f349c1fdSJed Brown #define __FUNCT__ "MatNestFindIS"
348f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
349f349c1fdSJed Brown {
350f349c1fdSJed Brown   PetscErrorCode ierr;
351f349c1fdSJed Brown   PetscInt       i;
352f349c1fdSJed Brown   PetscBool      flg;
353f349c1fdSJed Brown 
354f349c1fdSJed Brown   PetscFunctionBegin;
355f349c1fdSJed Brown   PetscValidPointer(list,3);
356f349c1fdSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,4);
357f349c1fdSJed Brown   PetscValidIntPointer(found,5);
358f349c1fdSJed Brown   *found = -1;
359f349c1fdSJed Brown   for (i=0; i<n; i++) {
360207556f9SJed Brown     if (!list[i]) continue;
361f349c1fdSJed Brown     ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr);
362f349c1fdSJed Brown     if (flg) {
363f349c1fdSJed Brown       *found = i;
364f349c1fdSJed Brown       PetscFunctionReturn(0);
365f349c1fdSJed Brown     }
366f349c1fdSJed Brown   }
367ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Could not find index set");
368f349c1fdSJed Brown   PetscFunctionReturn(0);
369f349c1fdSJed Brown }
370f349c1fdSJed Brown 
371f349c1fdSJed Brown #undef __FUNCT__
3728188e55aSJed Brown #define __FUNCT__ "MatNestGetRow"
3738188e55aSJed Brown /* Get a block row as a new MatNest */
3748188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B)
3758188e55aSJed Brown {
3768188e55aSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
3778188e55aSJed Brown   char           keyname[256];
3788188e55aSJed Brown   PetscErrorCode ierr;
3798188e55aSJed Brown 
3808188e55aSJed Brown   PetscFunctionBegin;
3810298fd71SBarry Smith   *B   = NULL;
3828caf3d72SBarry Smith   ierr = PetscSNPrintf(keyname,sizeof(keyname),"NestRow_%D",row);CHKERRQ(ierr);
3838188e55aSJed Brown   ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr);
3848188e55aSJed Brown   if (*B) PetscFunctionReturn(0);
3858188e55aSJed Brown 
386ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),1,NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr);
38726fbe8dcSKarl Rupp 
3888188e55aSJed Brown   (*B)->assembled = A->assembled;
38926fbe8dcSKarl Rupp 
3908188e55aSJed Brown   ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr);
3918188e55aSJed Brown   ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */
3928188e55aSJed Brown   PetscFunctionReturn(0);
3938188e55aSJed Brown }
3948188e55aSJed Brown 
3958188e55aSJed Brown #undef __FUNCT__
396f349c1fdSJed Brown #define __FUNCT__ "MatNestFindSubMat"
397f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
398f349c1fdSJed Brown {
399f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
4008188e55aSJed Brown   PetscErrorCode ierr;
4016b3a5b13SJed Brown   PetscInt       row,col;
402e072481dSJed Brown   PetscBool      same,isFullCol,isFullColGlobal;
403f349c1fdSJed Brown 
404f349c1fdSJed Brown   PetscFunctionBegin;
4058188e55aSJed Brown   /* Check if full column space. This is a hack */
4068188e55aSJed Brown   isFullCol = PETSC_FALSE;
407251f4c67SDmitry Karpeev   ierr      = PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr);
4088188e55aSJed Brown   if (same) {
40977019fcaSJed Brown     PetscInt n,first,step,i,an,am,afirst,astep;
4108188e55aSJed Brown     ierr      = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
4118188e55aSJed Brown     ierr      = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
41277019fcaSJed Brown     isFullCol = PETSC_TRUE;
41305ce4453SJed Brown     for (i=0,an=A->cmap->rstart; i<vs->nc; i++) {
41477019fcaSJed Brown       ierr = ISStrideGetInfo(is->col[i],&afirst,&astep);CHKERRQ(ierr);
41577019fcaSJed Brown       ierr = ISGetLocalSize(is->col[i],&am);CHKERRQ(ierr);
41677019fcaSJed Brown       if (afirst != an || astep != step) isFullCol = PETSC_FALSE;
41777019fcaSJed Brown       an += am;
41877019fcaSJed Brown     }
41905ce4453SJed Brown     if (an != A->cmap->rstart+n) isFullCol = PETSC_FALSE;
4208188e55aSJed Brown   }
421b2566f29SBarry Smith   ierr = MPIU_Allreduce(&isFullCol,&isFullColGlobal,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)iscol));CHKERRQ(ierr);
4228188e55aSJed Brown 
423e072481dSJed Brown   if (isFullColGlobal) {
4248188e55aSJed Brown     PetscInt row;
4258188e55aSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
4268188e55aSJed Brown     ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr);
4278188e55aSJed Brown   } else {
428f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
429f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr);
430b6480e04SStefano Zampini     if (!vs->m[row][col]) {
431b6480e04SStefano Zampini       PetscInt lr,lc;
432b6480e04SStefano Zampini 
433b6480e04SStefano Zampini       ierr = MatCreate(PetscObjectComm((PetscObject)A),&vs->m[row][col]);CHKERRQ(ierr);
434b6480e04SStefano Zampini       ierr = ISGetLocalSize(vs->isglobal.row[row],&lr);CHKERRQ(ierr);
435b6480e04SStefano Zampini       ierr = ISGetLocalSize(vs->isglobal.col[col],&lc);CHKERRQ(ierr);
436b6480e04SStefano Zampini       ierr = MatSetSizes(vs->m[row][col],lr,lc,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
437b6480e04SStefano Zampini       ierr = MatSetUp(vs->m[row][col]);CHKERRQ(ierr);
438b6480e04SStefano Zampini       ierr = MatAssemblyBegin(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
439b6480e04SStefano Zampini       ierr = MatAssemblyEnd(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
440b6480e04SStefano Zampini     }
441f349c1fdSJed Brown     *B = vs->m[row][col];
4428188e55aSJed Brown   }
443f349c1fdSJed Brown   PetscFunctionReturn(0);
444f349c1fdSJed Brown }
445f349c1fdSJed Brown 
446f349c1fdSJed Brown #undef __FUNCT__
447f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest"
448f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
449f349c1fdSJed Brown {
450f349c1fdSJed Brown   PetscErrorCode ierr;
451f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
452f349c1fdSJed Brown   Mat            sub;
453f349c1fdSJed Brown 
454f349c1fdSJed Brown   PetscFunctionBegin;
455f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr);
456f349c1fdSJed Brown   switch (reuse) {
457f349c1fdSJed Brown   case MAT_INITIAL_MATRIX:
4587874fa86SDave May     if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); }
459f349c1fdSJed Brown     *B = sub;
460f349c1fdSJed Brown     break;
461f349c1fdSJed Brown   case MAT_REUSE_MATRIX:
462ce94432eSBarry Smith     if (sub != *B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
463f349c1fdSJed Brown     break;
464f349c1fdSJed Brown   case MAT_IGNORE_MATRIX:       /* Nothing to do */
465f349c1fdSJed Brown     break;
466511c6705SHong Zhang   case MAT_INPLACE_MATRIX:       /* Nothing to do */
467511c6705SHong Zhang     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX is not supported yet");
468511c6705SHong Zhang     break;
469f349c1fdSJed Brown   }
470f349c1fdSJed Brown   PetscFunctionReturn(0);
471f349c1fdSJed Brown }
472f349c1fdSJed Brown 
473f349c1fdSJed Brown #undef __FUNCT__
474f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest"
475f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
476f349c1fdSJed Brown {
477f349c1fdSJed Brown   PetscErrorCode ierr;
478f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
479f349c1fdSJed Brown   Mat            sub;
480f349c1fdSJed Brown 
481f349c1fdSJed Brown   PetscFunctionBegin;
482f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
483f349c1fdSJed Brown   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
484f349c1fdSJed Brown   if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);}
485f349c1fdSJed Brown   *B = sub;
486d8588912SDave May   PetscFunctionReturn(0);
487d8588912SDave May }
488d8588912SDave May 
489d8588912SDave May #undef __FUNCT__
490d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest"
491207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
492d8588912SDave May {
493d8588912SDave May   PetscErrorCode ierr;
494f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
495f349c1fdSJed Brown   Mat            sub;
496d8588912SDave May 
497d8588912SDave May   PetscFunctionBegin;
498f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
499ce94432eSBarry Smith   if (*B != sub) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
500f349c1fdSJed Brown   if (sub) {
501ce94432eSBarry Smith     if (((PetscObject)sub)->refct <= 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
5026bf464f9SBarry Smith     ierr = MatDestroy(B);CHKERRQ(ierr);
503d8588912SDave May   }
504d8588912SDave May   PetscFunctionReturn(0);
505d8588912SDave May }
506d8588912SDave May 
507d8588912SDave May #undef __FUNCT__
5087874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest"
5097874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
5107874fa86SDave May {
5117874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
5127874fa86SDave May   PetscInt       i;
5137874fa86SDave May   PetscErrorCode ierr;
5147874fa86SDave May 
5157874fa86SDave May   PetscFunctionBegin;
5167874fa86SDave May   for (i=0; i<bA->nr; i++) {
517429bac76SJed Brown     Vec bv;
518429bac76SJed Brown     ierr = VecGetSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
5197874fa86SDave May     if (bA->m[i][i]) {
520429bac76SJed Brown       ierr = MatGetDiagonal(bA->m[i][i],bv);CHKERRQ(ierr);
5217874fa86SDave May     } else {
5225159a857SMatthew G. Knepley       ierr = VecSet(bv,0.0);CHKERRQ(ierr);
5237874fa86SDave May     }
524429bac76SJed Brown     ierr = VecRestoreSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
5257874fa86SDave May   }
5267874fa86SDave May   PetscFunctionReturn(0);
5277874fa86SDave May }
5287874fa86SDave May 
5297874fa86SDave May #undef __FUNCT__
5307874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest"
5317874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
5327874fa86SDave May {
5337874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
534429bac76SJed Brown   Vec            bl,*br;
5357874fa86SDave May   PetscInt       i,j;
5367874fa86SDave May   PetscErrorCode ierr;
5377874fa86SDave May 
5387874fa86SDave May   PetscFunctionBegin;
5393f800ebeSJed Brown   ierr = PetscCalloc1(bA->nc,&br);CHKERRQ(ierr);
5402e6472ebSElliott Sales de Andrade   if (r) {
541429bac76SJed Brown     for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
5422e6472ebSElliott Sales de Andrade   }
5432e6472ebSElliott Sales de Andrade   bl = NULL;
5447874fa86SDave May   for (i=0; i<bA->nr; i++) {
5452e6472ebSElliott Sales de Andrade     if (l) {
546429bac76SJed Brown       ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
5472e6472ebSElliott Sales de Andrade     }
5487874fa86SDave May     for (j=0; j<bA->nc; j++) {
5497874fa86SDave May       if (bA->m[i][j]) {
550429bac76SJed Brown         ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr);
5517874fa86SDave May       }
5527874fa86SDave May     }
5532e6472ebSElliott Sales de Andrade     if (l) {
554a061e289SJed Brown       ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
5557874fa86SDave May     }
5562e6472ebSElliott Sales de Andrade   }
5572e6472ebSElliott Sales de Andrade   if (r) {
558429bac76SJed Brown     for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
5592e6472ebSElliott Sales de Andrade   }
560429bac76SJed Brown   ierr = PetscFree(br);CHKERRQ(ierr);
5617874fa86SDave May   PetscFunctionReturn(0);
5627874fa86SDave May }
5637874fa86SDave May 
5647874fa86SDave May #undef __FUNCT__
565a061e289SJed Brown #define __FUNCT__ "MatScale_Nest"
566a061e289SJed Brown static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a)
567a061e289SJed Brown {
568a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
569a061e289SJed Brown   PetscInt       i,j;
570a061e289SJed Brown   PetscErrorCode ierr;
571a061e289SJed Brown 
572a061e289SJed Brown   PetscFunctionBegin;
573a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
574a061e289SJed Brown     for (j=0; j<bA->nc; j++) {
575a061e289SJed Brown       if (bA->m[i][j]) {
576a061e289SJed Brown         ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr);
577a061e289SJed Brown       }
578a061e289SJed Brown     }
579a061e289SJed Brown   }
580a061e289SJed Brown   PetscFunctionReturn(0);
581a061e289SJed Brown }
582a061e289SJed Brown 
583a061e289SJed Brown #undef __FUNCT__
584a061e289SJed Brown #define __FUNCT__ "MatShift_Nest"
585a061e289SJed Brown static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a)
586a061e289SJed Brown {
587a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
588a061e289SJed Brown   PetscInt       i;
589a061e289SJed Brown   PetscErrorCode ierr;
590a061e289SJed Brown 
591a061e289SJed Brown   PetscFunctionBegin;
592a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
593ce94432eSBarry 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);
594a061e289SJed Brown     ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr);
595a061e289SJed Brown   }
596a061e289SJed Brown   PetscFunctionReturn(0);
597a061e289SJed Brown }
598a061e289SJed Brown 
599a061e289SJed Brown #undef __FUNCT__
60013135bc6SAlex Fikl #define __FUNCT__ "MatDiagonalSet_Nest"
60113135bc6SAlex Fikl static PetscErrorCode MatDiagonalSet_Nest(Mat A,Vec D,InsertMode is)
60213135bc6SAlex Fikl {
60313135bc6SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data;
60413135bc6SAlex Fikl   PetscInt       i;
60513135bc6SAlex Fikl   PetscErrorCode ierr;
60613135bc6SAlex Fikl 
60713135bc6SAlex Fikl   PetscFunctionBegin;
60813135bc6SAlex Fikl   for (i=0; i<bA->nr; i++) {
60913135bc6SAlex Fikl     Vec bv;
61013135bc6SAlex Fikl     ierr = VecGetSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
61113135bc6SAlex Fikl     if (bA->m[i][i]) {
61213135bc6SAlex Fikl       ierr = MatDiagonalSet(bA->m[i][i],bv,is);CHKERRQ(ierr);
61313135bc6SAlex Fikl     }
61413135bc6SAlex Fikl     ierr = VecRestoreSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
61513135bc6SAlex Fikl   }
61613135bc6SAlex Fikl   PetscFunctionReturn(0);
61713135bc6SAlex Fikl }
61813135bc6SAlex Fikl 
61913135bc6SAlex Fikl #undef __FUNCT__
620f8170845SAlex Fikl #define __FUNCT__ "MatSetRandom_Nest"
621f8170845SAlex Fikl static PetscErrorCode MatSetRandom_Nest(Mat A,PetscRandom rctx)
622f8170845SAlex Fikl {
623f8170845SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data;
624f8170845SAlex Fikl   PetscInt       i,j;
625f8170845SAlex Fikl   PetscErrorCode ierr;
626f8170845SAlex Fikl 
627f8170845SAlex Fikl   PetscFunctionBegin;
628f8170845SAlex Fikl   for (i=0; i<bA->nr; i++) {
629f8170845SAlex Fikl     for (j=0; j<bA->nc; j++) {
630f8170845SAlex Fikl       if (bA->m[i][j]) {
631f8170845SAlex Fikl         ierr = MatSetRandom(bA->m[i][j],rctx);CHKERRQ(ierr);
632f8170845SAlex Fikl       }
633f8170845SAlex Fikl     }
634f8170845SAlex Fikl   }
635f8170845SAlex Fikl   PetscFunctionReturn(0);
636f8170845SAlex Fikl }
637f8170845SAlex Fikl 
638f8170845SAlex Fikl #undef __FUNCT__
6392a7a6963SBarry Smith #define __FUNCT__ "MatCreateVecs_Nest"
6402a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left)
641d8588912SDave May {
642d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
643d8588912SDave May   Vec            *L,*R;
644d8588912SDave May   MPI_Comm       comm;
645d8588912SDave May   PetscInt       i,j;
646d8588912SDave May   PetscErrorCode ierr;
647d8588912SDave May 
648d8588912SDave May   PetscFunctionBegin;
649ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
650d8588912SDave May   if (right) {
651d8588912SDave May     /* allocate R */
652854ce69bSBarry Smith     ierr = PetscMalloc1(bA->nc, &R);CHKERRQ(ierr);
653d8588912SDave May     /* Create the right vectors */
654d8588912SDave May     for (j=0; j<bA->nc; j++) {
655d8588912SDave May       for (i=0; i<bA->nr; i++) {
656d8588912SDave May         if (bA->m[i][j]) {
6572a7a6963SBarry Smith           ierr = MatCreateVecs(bA->m[i][j],&R[j],NULL);CHKERRQ(ierr);
658d8588912SDave May           break;
659d8588912SDave May         }
660d8588912SDave May       }
6616c4ed002SBarry Smith       if (i==bA->nr) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
662d8588912SDave May     }
663f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
664d8588912SDave May     /* hand back control to the nest vector */
665d8588912SDave May     for (j=0; j<bA->nc; j++) {
6666bf464f9SBarry Smith       ierr = VecDestroy(&R[j]);CHKERRQ(ierr);
667d8588912SDave May     }
668d8588912SDave May     ierr = PetscFree(R);CHKERRQ(ierr);
669d8588912SDave May   }
670d8588912SDave May 
671d8588912SDave May   if (left) {
672d8588912SDave May     /* allocate L */
673854ce69bSBarry Smith     ierr = PetscMalloc1(bA->nr, &L);CHKERRQ(ierr);
674d8588912SDave May     /* Create the left vectors */
675d8588912SDave May     for (i=0; i<bA->nr; i++) {
676d8588912SDave May       for (j=0; j<bA->nc; j++) {
677d8588912SDave May         if (bA->m[i][j]) {
6782a7a6963SBarry Smith           ierr = MatCreateVecs(bA->m[i][j],NULL,&L[i]);CHKERRQ(ierr);
679d8588912SDave May           break;
680d8588912SDave May         }
681d8588912SDave May       }
6826c4ed002SBarry Smith       if (j==bA->nc) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
683d8588912SDave May     }
684d8588912SDave May 
685f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
686d8588912SDave May     for (i=0; i<bA->nr; i++) {
6876bf464f9SBarry Smith       ierr = VecDestroy(&L[i]);CHKERRQ(ierr);
688d8588912SDave May     }
689d8588912SDave May 
690d8588912SDave May     ierr = PetscFree(L);CHKERRQ(ierr);
691d8588912SDave May   }
692d8588912SDave May   PetscFunctionReturn(0);
693d8588912SDave May }
694d8588912SDave May 
695d8588912SDave May #undef __FUNCT__
696d8588912SDave May #define __FUNCT__ "MatView_Nest"
697207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
698d8588912SDave May {
699d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
700d8588912SDave May   PetscBool      isascii;
701d8588912SDave May   PetscInt       i,j;
702d8588912SDave May   PetscErrorCode ierr;
703d8588912SDave May 
704d8588912SDave May   PetscFunctionBegin;
705251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
706d8588912SDave May   if (isascii) {
707d8588912SDave May 
708d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Matrix object: \n");CHKERRQ(ierr);
709d86155a6SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
710d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer, "type=nest, rows=%D, cols=%D \n",bA->nr,bA->nc);CHKERRQ(ierr);
711d8588912SDave May 
712d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"MatNest structure: \n");CHKERRQ(ierr);
713d8588912SDave May     for (i=0; i<bA->nr; i++) {
714d8588912SDave May       for (j=0; j<bA->nc; j++) {
71519fd82e9SBarry Smith         MatType   type;
716270f95d7SJed Brown         char      name[256] = "",prefix[256] = "";
717d8588912SDave May         PetscInt  NR,NC;
718d8588912SDave May         PetscBool isNest = PETSC_FALSE;
719d8588912SDave May 
720d8588912SDave May         if (!bA->m[i][j]) {
721d86155a6SBarry Smith           CHKERRQ(ierr);PetscViewerASCIIPrintf(viewer, "(%D,%D) : NULL \n",i,j);CHKERRQ(ierr);
722d8588912SDave May           continue;
723d8588912SDave May         }
724d8588912SDave May         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
725d8588912SDave May         ierr = MatGetType(bA->m[i][j], &type);CHKERRQ(ierr);
7268caf3d72SBarry Smith         if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);}
7278caf3d72SBarry Smith         if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);}
728251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
729d8588912SDave May 
730270f95d7SJed Brown         ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr);
731d8588912SDave May 
732d8588912SDave May         if (isNest) {
733270f95d7SJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);  /* push1 */
734d8588912SDave May           ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr);
735270f95d7SJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop1 */
736d8588912SDave May         }
737d8588912SDave May       }
738d8588912SDave May     }
739d86155a6SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop0 */
740d8588912SDave May   }
741d8588912SDave May   PetscFunctionReturn(0);
742d8588912SDave May }
743d8588912SDave May 
744d8588912SDave May #undef __FUNCT__
745d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest"
746207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A)
747d8588912SDave May {
748d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
749d8588912SDave May   PetscInt       i,j;
750d8588912SDave May   PetscErrorCode ierr;
751d8588912SDave May 
752d8588912SDave May   PetscFunctionBegin;
753d8588912SDave May   for (i=0; i<bA->nr; i++) {
754d8588912SDave May     for (j=0; j<bA->nc; j++) {
755d8588912SDave May       if (!bA->m[i][j]) continue;
756d8588912SDave May       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
757d8588912SDave May     }
758d8588912SDave May   }
759d8588912SDave May   PetscFunctionReturn(0);
760d8588912SDave May }
761d8588912SDave May 
762d8588912SDave May #undef __FUNCT__
763c222c20dSDavid Ham #define __FUNCT__ "MatCopy_Nest"
764c222c20dSDavid Ham static PetscErrorCode MatCopy_Nest(Mat A,Mat B,MatStructure str)
765c222c20dSDavid Ham {
766c222c20dSDavid Ham   Mat_Nest       *bA = (Mat_Nest*)A->data,*bB = (Mat_Nest*)B->data;
767c222c20dSDavid Ham   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
768c222c20dSDavid Ham   PetscErrorCode ierr;
769c222c20dSDavid Ham 
770c222c20dSDavid Ham   PetscFunctionBegin;
771c222c20dSDavid 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);
772c222c20dSDavid Ham   for (i=0; i<nr; i++) {
773c222c20dSDavid Ham     for (j=0; j<nc; j++) {
77446a2b97cSJed Brown       if (bA->m[i][j] && bB->m[i][j]) {
775c222c20dSDavid Ham         ierr = MatCopy(bA->m[i][j],bB->m[i][j],str);CHKERRQ(ierr);
77646a2b97cSJed 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);
777c222c20dSDavid Ham     }
778c222c20dSDavid Ham   }
779c222c20dSDavid Ham   PetscFunctionReturn(0);
780c222c20dSDavid Ham }
781c222c20dSDavid Ham 
782c222c20dSDavid Ham #undef __FUNCT__
783d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest"
784207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
785d8588912SDave May {
786d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
787841e96a3SJed Brown   Mat            *b;
788841e96a3SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
789d8588912SDave May   PetscErrorCode ierr;
790d8588912SDave May 
791d8588912SDave May   PetscFunctionBegin;
792785e854fSJed Brown   ierr = PetscMalloc1(nr*nc,&b);CHKERRQ(ierr);
793841e96a3SJed Brown   for (i=0; i<nr; i++) {
794841e96a3SJed Brown     for (j=0; j<nc; j++) {
795841e96a3SJed Brown       if (bA->m[i][j]) {
796841e96a3SJed Brown         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
797841e96a3SJed Brown       } else {
7980298fd71SBarry Smith         b[i*nc+j] = NULL;
799d8588912SDave May       }
800d8588912SDave May     }
801d8588912SDave May   }
802ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
803841e96a3SJed Brown   /* Give the new MatNest exclusive ownership */
804841e96a3SJed Brown   for (i=0; i<nr*nc; i++) {
8056bf464f9SBarry Smith     ierr = MatDestroy(&b[i]);CHKERRQ(ierr);
806d8588912SDave May   }
807d8588912SDave May   ierr = PetscFree(b);CHKERRQ(ierr);
808d8588912SDave May 
809841e96a3SJed Brown   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
810841e96a3SJed Brown   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
811d8588912SDave May   PetscFunctionReturn(0);
812d8588912SDave May }
813d8588912SDave May 
814d8588912SDave May /* nest api */
815d8588912SDave May #undef __FUNCT__
816d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest"
817d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
818d8588912SDave May {
819d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
8205fd66863SKarl Rupp 
821d8588912SDave May   PetscFunctionBegin;
822ce94432eSBarry Smith   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
823ce94432eSBarry Smith   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
824d8588912SDave May   *mat = bA->m[idxm][jdxm];
825d8588912SDave May   PetscFunctionReturn(0);
826d8588912SDave May }
827d8588912SDave May 
828d8588912SDave May #undef __FUNCT__
829d8588912SDave May #define __FUNCT__ "MatNestGetSubMat"
8309ba0d327SJed Brown /*@
831d8588912SDave May  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
832d8588912SDave May 
833d8588912SDave May  Not collective
834d8588912SDave May 
835d8588912SDave May  Input Parameters:
836629881c0SJed Brown +   A  - nest matrix
837d8588912SDave May .   idxm - index of the matrix within the nest matrix
838629881c0SJed Brown -   jdxm - index of the matrix within the nest matrix
839d8588912SDave May 
840d8588912SDave May  Output Parameter:
841d8588912SDave May .   sub - matrix at index idxm,jdxm within the nest matrix
842d8588912SDave May 
843d8588912SDave May  Level: developer
844d8588912SDave May 
845d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMats()
846d8588912SDave May @*/
8477087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
848d8588912SDave May {
849699a902aSJed Brown   PetscErrorCode ierr;
850d8588912SDave May 
851d8588912SDave May   PetscFunctionBegin;
852699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
853d8588912SDave May   PetscFunctionReturn(0);
854d8588912SDave May }
855d8588912SDave May 
856d8588912SDave May #undef __FUNCT__
8570782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat_Nest"
8580782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat)
8590782ca92SJed Brown {
8600782ca92SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
8610782ca92SJed Brown   PetscInt       m,n,M,N,mi,ni,Mi,Ni;
8620782ca92SJed Brown   PetscErrorCode ierr;
8630782ca92SJed Brown 
8640782ca92SJed Brown   PetscFunctionBegin;
865ce94432eSBarry Smith   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
866ce94432eSBarry Smith   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
8670782ca92SJed Brown   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
8680782ca92SJed Brown   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
8690782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr);
8700782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr);
8710782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr);
8720782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr);
873ce94432eSBarry 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);
874ce94432eSBarry 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);
87526fbe8dcSKarl Rupp 
8760782ca92SJed Brown   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
8770782ca92SJed Brown   ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr);
8780782ca92SJed Brown   bA->m[idxm][jdxm] = mat;
8790782ca92SJed Brown   PetscFunctionReturn(0);
8800782ca92SJed Brown }
8810782ca92SJed Brown 
8820782ca92SJed Brown #undef __FUNCT__
8830782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat"
8849ba0d327SJed Brown /*@
8850782ca92SJed Brown  MatNestSetSubMat - Set a single submatrix in the nest matrix.
8860782ca92SJed Brown 
8870782ca92SJed Brown  Logically collective on the submatrix communicator
8880782ca92SJed Brown 
8890782ca92SJed Brown  Input Parameters:
8900782ca92SJed Brown +   A  - nest matrix
8910782ca92SJed Brown .   idxm - index of the matrix within the nest matrix
8920782ca92SJed Brown .   jdxm - index of the matrix within the nest matrix
8930782ca92SJed Brown -   sub - matrix at index idxm,jdxm within the nest matrix
8940782ca92SJed Brown 
8950782ca92SJed Brown  Notes:
8960782ca92SJed Brown  The new submatrix must have the same size and communicator as that block of the nest.
8970782ca92SJed Brown 
8980782ca92SJed Brown  This increments the reference count of the submatrix.
8990782ca92SJed Brown 
9000782ca92SJed Brown  Level: developer
9010782ca92SJed Brown 
9020782ca92SJed Brown .seealso: MatNestSetSubMats(), MatNestGetSubMat()
9030782ca92SJed Brown @*/
9040782ca92SJed Brown PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
9050782ca92SJed Brown {
9060782ca92SJed Brown   PetscErrorCode ierr;
9070782ca92SJed Brown 
9080782ca92SJed Brown   PetscFunctionBegin;
9090782ca92SJed Brown   ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr);
9100782ca92SJed Brown   PetscFunctionReturn(0);
9110782ca92SJed Brown }
9120782ca92SJed Brown 
9130782ca92SJed Brown #undef __FUNCT__
914d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest"
915d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
916d8588912SDave May {
917d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
9185fd66863SKarl Rupp 
919d8588912SDave May   PetscFunctionBegin;
92026fbe8dcSKarl Rupp   if (M)   *M   = bA->nr;
92126fbe8dcSKarl Rupp   if (N)   *N   = bA->nc;
92226fbe8dcSKarl Rupp   if (mat) *mat = bA->m;
923d8588912SDave May   PetscFunctionReturn(0);
924d8588912SDave May }
925d8588912SDave May 
926d8588912SDave May #undef __FUNCT__
927d8588912SDave May #define __FUNCT__ "MatNestGetSubMats"
928d8588912SDave May /*@C
929d8588912SDave May  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
930d8588912SDave May 
931d8588912SDave May  Not collective
932d8588912SDave May 
933d8588912SDave May  Input Parameters:
934629881c0SJed Brown .   A  - nest matrix
935d8588912SDave May 
936d8588912SDave May  Output Parameter:
937629881c0SJed Brown +   M - number of rows in the nest matrix
938d8588912SDave May .   N - number of cols in the nest matrix
939629881c0SJed Brown -   mat - 2d array of matrices
940d8588912SDave May 
941d8588912SDave May  Notes:
942d8588912SDave May 
943d8588912SDave May  The user should not free the array mat.
944d8588912SDave May 
945351962e3SVincent Le Chenadec  In Fortran, this routine has a calling sequence
946351962e3SVincent Le Chenadec $   call MatNestGetSubMats(A, M, N, mat, ierr)
947351962e3SVincent Le Chenadec  where the space allocated for the optional argument mat is assumed large enough (if provided).
948351962e3SVincent Le Chenadec 
949d8588912SDave May  Level: developer
950d8588912SDave May 
951d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMat()
952d8588912SDave May @*/
9537087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
954d8588912SDave May {
955699a902aSJed Brown   PetscErrorCode ierr;
956d8588912SDave May 
957d8588912SDave May   PetscFunctionBegin;
958699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
959d8588912SDave May   PetscFunctionReturn(0);
960d8588912SDave May }
961d8588912SDave May 
962d8588912SDave May #undef __FUNCT__
963d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest"
9647087cfbeSBarry Smith PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
965d8588912SDave May {
966d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
967d8588912SDave May 
968d8588912SDave May   PetscFunctionBegin;
96926fbe8dcSKarl Rupp   if (M) *M = bA->nr;
97026fbe8dcSKarl Rupp   if (N) *N = bA->nc;
971d8588912SDave May   PetscFunctionReturn(0);
972d8588912SDave May }
973d8588912SDave May 
974d8588912SDave May #undef __FUNCT__
975d8588912SDave May #define __FUNCT__ "MatNestGetSize"
9769ba0d327SJed Brown /*@
977d8588912SDave May  MatNestGetSize - Returns the size of the nest matrix.
978d8588912SDave May 
979d8588912SDave May  Not collective
980d8588912SDave May 
981d8588912SDave May  Input Parameters:
982d8588912SDave May .   A  - nest matrix
983d8588912SDave May 
984d8588912SDave May  Output Parameter:
985629881c0SJed Brown +   M - number of rows in the nested mat
986629881c0SJed Brown -   N - number of cols in the nested mat
987d8588912SDave May 
988d8588912SDave May  Notes:
989d8588912SDave May 
990d8588912SDave May  Level: developer
991d8588912SDave May 
992d8588912SDave May .seealso: MatNestGetSubMat(), MatNestGetSubMats()
993d8588912SDave May @*/
9947087cfbeSBarry Smith PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
995d8588912SDave May {
996699a902aSJed Brown   PetscErrorCode ierr;
997d8588912SDave May 
998d8588912SDave May   PetscFunctionBegin;
999699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
1000d8588912SDave May   PetscFunctionReturn(0);
1001d8588912SDave May }
1002d8588912SDave May 
1003900e7ff2SJed Brown #undef __FUNCT__
1004900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs_Nest"
1005f7a08781SBarry Smith static PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[])
1006900e7ff2SJed Brown {
1007900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest*)A->data;
1008900e7ff2SJed Brown   PetscInt i;
1009900e7ff2SJed Brown 
1010900e7ff2SJed Brown   PetscFunctionBegin;
1011900e7ff2SJed Brown   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i];
1012900e7ff2SJed Brown   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i];
1013900e7ff2SJed Brown   PetscFunctionReturn(0);
1014900e7ff2SJed Brown }
1015900e7ff2SJed Brown 
1016900e7ff2SJed Brown #undef __FUNCT__
1017900e7ff2SJed Brown #define __FUNCT__ "MatNestGetISs"
10183a4d7b9aSSatish Balay /*@C
1019900e7ff2SJed Brown  MatNestGetISs - Returns the index sets partitioning the row and column spaces
1020900e7ff2SJed Brown 
1021900e7ff2SJed Brown  Not collective
1022900e7ff2SJed Brown 
1023900e7ff2SJed Brown  Input Parameters:
1024900e7ff2SJed Brown .   A  - nest matrix
1025900e7ff2SJed Brown 
1026900e7ff2SJed Brown  Output Parameter:
1027900e7ff2SJed Brown +   rows - array of row index sets
1028900e7ff2SJed Brown -   cols - array of column index sets
1029900e7ff2SJed Brown 
1030900e7ff2SJed Brown  Level: advanced
1031900e7ff2SJed Brown 
1032900e7ff2SJed Brown  Notes:
1033900e7ff2SJed Brown  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
1034900e7ff2SJed Brown 
1035900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs()
1036900e7ff2SJed Brown @*/
1037900e7ff2SJed Brown PetscErrorCode  MatNestGetISs(Mat A,IS rows[],IS cols[])
1038900e7ff2SJed Brown {
1039900e7ff2SJed Brown   PetscErrorCode ierr;
1040900e7ff2SJed Brown 
1041900e7ff2SJed Brown   PetscFunctionBegin;
1042900e7ff2SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1043900e7ff2SJed Brown   ierr = PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
1044900e7ff2SJed Brown   PetscFunctionReturn(0);
1045900e7ff2SJed Brown }
1046900e7ff2SJed Brown 
1047900e7ff2SJed Brown #undef __FUNCT__
1048900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs_Nest"
1049f7a08781SBarry Smith static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[])
1050900e7ff2SJed Brown {
1051900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest*)A->data;
1052900e7ff2SJed Brown   PetscInt i;
1053900e7ff2SJed Brown 
1054900e7ff2SJed Brown   PetscFunctionBegin;
1055900e7ff2SJed Brown   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i];
1056900e7ff2SJed Brown   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i];
1057900e7ff2SJed Brown   PetscFunctionReturn(0);
1058900e7ff2SJed Brown }
1059900e7ff2SJed Brown 
1060900e7ff2SJed Brown #undef __FUNCT__
1061900e7ff2SJed Brown #define __FUNCT__ "MatNestGetLocalISs"
1062900e7ff2SJed Brown /*@C
1063900e7ff2SJed Brown  MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces
1064900e7ff2SJed Brown 
1065900e7ff2SJed Brown  Not collective
1066900e7ff2SJed Brown 
1067900e7ff2SJed Brown  Input Parameters:
1068900e7ff2SJed Brown .   A  - nest matrix
1069900e7ff2SJed Brown 
1070900e7ff2SJed Brown  Output Parameter:
10710298fd71SBarry Smith +   rows - array of row index sets (or NULL to ignore)
10720298fd71SBarry Smith -   cols - array of column index sets (or NULL to ignore)
1073900e7ff2SJed Brown 
1074900e7ff2SJed Brown  Level: advanced
1075900e7ff2SJed Brown 
1076900e7ff2SJed Brown  Notes:
1077900e7ff2SJed Brown  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
1078900e7ff2SJed Brown 
1079900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs()
1080900e7ff2SJed Brown @*/
1081900e7ff2SJed Brown PetscErrorCode  MatNestGetLocalISs(Mat A,IS rows[],IS cols[])
1082900e7ff2SJed Brown {
1083900e7ff2SJed Brown   PetscErrorCode ierr;
1084900e7ff2SJed Brown 
1085900e7ff2SJed Brown   PetscFunctionBegin;
1086900e7ff2SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1087900e7ff2SJed Brown   ierr = PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
1088900e7ff2SJed Brown   PetscFunctionReturn(0);
1089900e7ff2SJed Brown }
1090900e7ff2SJed Brown 
1091207556f9SJed Brown #undef __FUNCT__
1092207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest"
109319fd82e9SBarry Smith PetscErrorCode  MatNestSetVecType_Nest(Mat A,VecType vtype)
1094207556f9SJed Brown {
1095207556f9SJed Brown   PetscErrorCode ierr;
1096207556f9SJed Brown   PetscBool      flg;
1097207556f9SJed Brown 
1098207556f9SJed Brown   PetscFunctionBegin;
1099207556f9SJed Brown   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
1100207556f9SJed Brown   /* In reality, this only distinguishes VECNEST and "other" */
11012a7a6963SBarry Smith   if (flg) A->ops->getvecs = MatCreateVecs_Nest;
110212b53f24SSatish Balay   else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0;
1103207556f9SJed Brown   PetscFunctionReturn(0);
1104207556f9SJed Brown }
1105207556f9SJed Brown 
1106207556f9SJed Brown #undef __FUNCT__
1107207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType"
1108207556f9SJed Brown /*@C
11092a7a6963SBarry Smith  MatNestSetVecType - Sets the type of Vec returned by MatCreateVecs()
1110207556f9SJed Brown 
1111207556f9SJed Brown  Not collective
1112207556f9SJed Brown 
1113207556f9SJed Brown  Input Parameters:
1114207556f9SJed Brown +  A  - nest matrix
1115207556f9SJed Brown -  vtype - type to use for creating vectors
1116207556f9SJed Brown 
1117207556f9SJed Brown  Notes:
1118207556f9SJed Brown 
1119207556f9SJed Brown  Level: developer
1120207556f9SJed Brown 
11212a7a6963SBarry Smith .seealso: MatCreateVecs()
1122207556f9SJed Brown @*/
112319fd82e9SBarry Smith PetscErrorCode  MatNestSetVecType(Mat A,VecType vtype)
1124207556f9SJed Brown {
1125207556f9SJed Brown   PetscErrorCode ierr;
1126207556f9SJed Brown 
1127207556f9SJed Brown   PetscFunctionBegin;
112819fd82e9SBarry Smith   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));CHKERRQ(ierr);
1129207556f9SJed Brown   PetscFunctionReturn(0);
1130207556f9SJed Brown }
1131207556f9SJed Brown 
1132d8588912SDave May #undef __FUNCT__
1133c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats_Nest"
1134c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1135d8588912SDave May {
1136c8883902SJed Brown   Mat_Nest       *s = (Mat_Nest*)A->data;
1137c8883902SJed Brown   PetscInt       i,j,m,n,M,N;
1138d8588912SDave May   PetscErrorCode ierr;
1139d8588912SDave May 
1140d8588912SDave May   PetscFunctionBegin;
1141c8883902SJed Brown   s->nr = nr;
1142c8883902SJed Brown   s->nc = nc;
1143d8588912SDave May 
1144c8883902SJed Brown   /* Create space for submatrices */
1145854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&s->m);CHKERRQ(ierr);
1146c8883902SJed Brown   for (i=0; i<nr; i++) {
1147854ce69bSBarry Smith     ierr = PetscMalloc1(nc,&s->m[i]);CHKERRQ(ierr);
1148d8588912SDave May   }
1149c8883902SJed Brown   for (i=0; i<nr; i++) {
1150c8883902SJed Brown     for (j=0; j<nc; j++) {
1151c8883902SJed Brown       s->m[i][j] = a[i*nc+j];
1152c8883902SJed Brown       if (a[i*nc+j]) {
1153c8883902SJed Brown         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
1154d8588912SDave May       }
1155d8588912SDave May     }
1156d8588912SDave May   }
1157d8588912SDave May 
11588188e55aSJed Brown   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
1159d8588912SDave May 
1160854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&s->row_len);CHKERRQ(ierr);
1161854ce69bSBarry Smith   ierr = PetscMalloc1(nc,&s->col_len);CHKERRQ(ierr);
1162c8883902SJed Brown   for (i=0; i<nr; i++) s->row_len[i]=-1;
1163c8883902SJed Brown   for (j=0; j<nc; j++) s->col_len[j]=-1;
1164d8588912SDave May 
11658188e55aSJed Brown   ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr);
1166d8588912SDave May 
1167c8883902SJed Brown   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
1168c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
1169c8883902SJed Brown   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
1170c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
1171c8883902SJed Brown 
1172c8883902SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1173c8883902SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1174c8883902SJed Brown 
11751795a4d1SJed Brown   ierr = PetscCalloc2(nr,&s->left,nc,&s->right);CHKERRQ(ierr);
1176d8588912SDave May   PetscFunctionReturn(0);
1177d8588912SDave May }
1178d8588912SDave May 
1179c8883902SJed Brown #undef __FUNCT__
1180c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats"
1181c8883902SJed Brown /*@
1182c8883902SJed Brown    MatNestSetSubMats - Sets the nested submatrices
1183c8883902SJed Brown 
1184c8883902SJed Brown    Collective on Mat
1185c8883902SJed Brown 
1186c8883902SJed Brown    Input Parameter:
1187c8883902SJed Brown +  N - nested matrix
1188c8883902SJed Brown .  nr - number of nested row blocks
11890298fd71SBarry Smith .  is_row - index sets for each nested row block, or NULL to make contiguous
1190c8883902SJed Brown .  nc - number of nested column blocks
11910298fd71SBarry Smith .  is_col - index sets for each nested column block, or NULL to make contiguous
11920298fd71SBarry Smith -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1193c8883902SJed Brown 
1194c8883902SJed Brown    Level: advanced
1195c8883902SJed Brown 
1196c8883902SJed Brown .seealso: MatCreateNest(), MATNEST
1197c8883902SJed Brown @*/
1198c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1199c8883902SJed Brown {
1200c8883902SJed Brown   PetscErrorCode ierr;
1201c8883902SJed Brown   PetscInt       i;
1202c8883902SJed Brown 
1203c8883902SJed Brown   PetscFunctionBegin;
1204c8883902SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1205ce94432eSBarry Smith   if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
1206c8883902SJed Brown   if (nr && is_row) {
1207c8883902SJed Brown     PetscValidPointer(is_row,3);
1208c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
1209c8883902SJed Brown   }
1210ce94432eSBarry Smith   if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
12111664e352SJed Brown   if (nc && is_col) {
1212c8883902SJed Brown     PetscValidPointer(is_col,5);
12139b30a8f6SBarry Smith     for (i=0; i<nc; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
1214c8883902SJed Brown   }
1215c8883902SJed Brown   if (nr*nc) PetscValidPointer(a,6);
1216c8883902SJed 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);
1217c8883902SJed Brown   PetscFunctionReturn(0);
1218c8883902SJed Brown }
1219d8588912SDave May 
122077019fcaSJed Brown #undef __FUNCT__
122177019fcaSJed Brown #define __FUNCT__ "MatNestCreateAggregateL2G_Private"
122245b6f7e9SBarry Smith static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog)
122377019fcaSJed Brown {
122477019fcaSJed Brown   PetscErrorCode ierr;
122577019fcaSJed Brown   PetscBool      flg;
122677019fcaSJed Brown   PetscInt       i,j,m,mi,*ix;
122777019fcaSJed Brown 
122877019fcaSJed Brown   PetscFunctionBegin;
122977019fcaSJed Brown   for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) {
123077019fcaSJed Brown     if (islocal[i]) {
123177019fcaSJed Brown       ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
123277019fcaSJed Brown       flg  = PETSC_TRUE;      /* We found a non-trivial entry */
123377019fcaSJed Brown     } else {
123477019fcaSJed Brown       ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
123577019fcaSJed Brown     }
123677019fcaSJed Brown     m += mi;
123777019fcaSJed Brown   }
123877019fcaSJed Brown   if (flg) {
1239785e854fSJed Brown     ierr = PetscMalloc1(m,&ix);CHKERRQ(ierr);
124077019fcaSJed Brown     for (i=0,n=0; i<n; i++) {
12410298fd71SBarry Smith       ISLocalToGlobalMapping smap = NULL;
124277019fcaSJed Brown       VecScatter             scat;
124377019fcaSJed Brown       IS                     isreq;
124477019fcaSJed Brown       Vec                    lvec,gvec;
12453361c9a7SJed Brown       union {char padding[sizeof(PetscScalar)]; PetscInt integer;} *x;
124677019fcaSJed Brown       Mat sub;
124777019fcaSJed Brown 
1248ce94432eSBarry Smith       if (sizeof(*x) != sizeof(PetscScalar)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support when scalars smaller than integers");
124977019fcaSJed Brown       if (colflg) {
125077019fcaSJed Brown         ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
125177019fcaSJed Brown       } else {
125277019fcaSJed Brown         ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
125377019fcaSJed Brown       }
12540298fd71SBarry Smith       if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&smap,NULL);CHKERRQ(ierr);}
125577019fcaSJed Brown       if (islocal[i]) {
125677019fcaSJed Brown         ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
125777019fcaSJed Brown       } else {
125877019fcaSJed Brown         ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
125977019fcaSJed Brown       }
126077019fcaSJed Brown       for (j=0; j<mi; j++) ix[m+j] = j;
126177019fcaSJed Brown       if (smap) {ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr);}
126277019fcaSJed Brown       /*
126377019fcaSJed Brown         Now we need to extract the monolithic global indices that correspond to the given split global indices.
126477019fcaSJed 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.
126577019fcaSJed Brown         The approach here is ugly because it uses VecScatter to move indices.
126677019fcaSJed Brown        */
126777019fcaSJed Brown       ierr = VecCreateSeq(PETSC_COMM_SELF,mi,&lvec);CHKERRQ(ierr);
126877019fcaSJed Brown       ierr = VecCreateMPI(((PetscObject)isglobal[i])->comm,mi,PETSC_DECIDE,&gvec);CHKERRQ(ierr);
126977019fcaSJed Brown       ierr = ISCreateGeneral(((PetscObject)isglobal[i])->comm,mi,ix+m,PETSC_COPY_VALUES,&isreq);CHKERRQ(ierr);
12700298fd71SBarry Smith       ierr = VecScatterCreate(gvec,isreq,lvec,NULL,&scat);CHKERRQ(ierr);
127177019fcaSJed Brown       ierr = VecGetArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
127277019fcaSJed Brown       for (j=0; j<mi; j++) x[j].integer = ix[m+j];
127377019fcaSJed Brown       ierr = VecRestoreArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
127477019fcaSJed Brown       ierr = VecScatterBegin(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
127577019fcaSJed Brown       ierr = VecScatterEnd(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
127677019fcaSJed Brown       ierr = VecGetArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
127777019fcaSJed Brown       for (j=0; j<mi; j++) ix[m+j] = x[j].integer;
127877019fcaSJed Brown       ierr = VecRestoreArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
127977019fcaSJed Brown       ierr = VecDestroy(&lvec);CHKERRQ(ierr);
128077019fcaSJed Brown       ierr = VecDestroy(&gvec);CHKERRQ(ierr);
128177019fcaSJed Brown       ierr = ISDestroy(&isreq);CHKERRQ(ierr);
128277019fcaSJed Brown       ierr = VecScatterDestroy(&scat);CHKERRQ(ierr);
128377019fcaSJed Brown       m   += mi;
128477019fcaSJed Brown     }
1285f0413b6fSBarry Smith     ierr   = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),1,m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
128677019fcaSJed Brown   } else {
12870298fd71SBarry Smith     *ltog  = NULL;
128877019fcaSJed Brown   }
128977019fcaSJed Brown   PetscFunctionReturn(0);
129077019fcaSJed Brown }
129177019fcaSJed Brown 
129277019fcaSJed Brown 
1293d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1294d8588912SDave May /*
1295d8588912SDave May   nprocessors = NP
1296d8588912SDave May   Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1))
1297d8588912SDave May        proc 0: => (g_0,h_0,)
1298d8588912SDave May        proc 1: => (g_1,h_1,)
1299d8588912SDave May        ...
1300d8588912SDave May        proc nprocs-1: => (g_NP-1,h_NP-1,)
1301d8588912SDave May 
1302d8588912SDave May             proc 0:                      proc 1:                    proc nprocs-1:
1303d8588912SDave May     is[0] = (0,1,2,...,nlocal(g_0)-1)  (0,1,...,nlocal(g_1)-1)  (0,1,...,nlocal(g_NP-1))
1304d8588912SDave May 
1305d8588912SDave May             proc 0:
1306d8588912SDave May     is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1)
1307d8588912SDave May             proc 1:
1308d8588912SDave May     is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1)
1309d8588912SDave May 
1310d8588912SDave May             proc NP-1:
1311d8588912SDave May     is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1)
1312d8588912SDave May */
1313d8588912SDave May #undef __FUNCT__
1314d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private"
1315841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1316d8588912SDave May {
1317e2d7f03fSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
13188188e55aSJed Brown   PetscInt       i,j,offset,n,nsum,bs;
1319d8588912SDave May   PetscErrorCode ierr;
13200298fd71SBarry Smith   Mat            sub = NULL;
1321d8588912SDave May 
1322d8588912SDave May   PetscFunctionBegin;
1323854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&vs->isglobal.row);CHKERRQ(ierr);
1324854ce69bSBarry Smith   ierr = PetscMalloc1(nc,&vs->isglobal.col);CHKERRQ(ierr);
1325d8588912SDave May   if (is_row) { /* valid IS is passed in */
1326d8588912SDave May     /* refs on is[] are incremeneted */
1327e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1328d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
132926fbe8dcSKarl Rupp 
1330e2d7f03fSJed Brown       vs->isglobal.row[i] = is_row[i];
1331d8588912SDave May     }
13322ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
13338188e55aSJed Brown     nsum = 0;
13348188e55aSJed Brown     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
13358188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1336ce94432eSBarry Smith       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
13370298fd71SBarry Smith       ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
1338ce94432eSBarry Smith       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
13398188e55aSJed Brown       nsum += n;
13408188e55aSJed Brown     }
1341ce94432eSBarry Smith     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
134230bc264bSJed Brown     offset -= nsum;
1343e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1344f349c1fdSJed Brown       ierr    = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
13450298fd71SBarry Smith       ierr    = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
13462ae74bdbSJed Brown       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1347ce94432eSBarry Smith       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
1348e2d7f03fSJed Brown       ierr    = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
13492ae74bdbSJed Brown       offset += n;
1350d8588912SDave May     }
1351d8588912SDave May   }
1352d8588912SDave May 
1353d8588912SDave May   if (is_col) { /* valid IS is passed in */
1354d8588912SDave May     /* refs on is[] are incremeneted */
1355e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1356d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
135726fbe8dcSKarl Rupp 
1358e2d7f03fSJed Brown       vs->isglobal.col[j] = is_col[j];
1359d8588912SDave May     }
13602ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
13612ae74bdbSJed Brown     offset = A->cmap->rstart;
13628188e55aSJed Brown     nsum   = 0;
13638188e55aSJed Brown     for (j=0; j<vs->nc; j++) {
13648188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
1365ce94432eSBarry Smith       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
13660298fd71SBarry Smith       ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
1367ce94432eSBarry Smith       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
13688188e55aSJed Brown       nsum += n;
13698188e55aSJed Brown     }
1370ce94432eSBarry Smith     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
137130bc264bSJed Brown     offset -= nsum;
1372e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1373f349c1fdSJed Brown       ierr    = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
13740298fd71SBarry Smith       ierr    = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
13752ae74bdbSJed Brown       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1376ce94432eSBarry Smith       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
1377e2d7f03fSJed Brown       ierr    = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
13782ae74bdbSJed Brown       offset += n;
1379d8588912SDave May     }
1380d8588912SDave May   }
1381e2d7f03fSJed Brown 
1382e2d7f03fSJed Brown   /* Set up the local ISs */
1383785e854fSJed Brown   ierr = PetscMalloc1(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
1384785e854fSJed Brown   ierr = PetscMalloc1(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
1385e2d7f03fSJed Brown   for (i=0,offset=0; i<vs->nr; i++) {
1386e2d7f03fSJed Brown     IS                     isloc;
13870298fd71SBarry Smith     ISLocalToGlobalMapping rmap = NULL;
1388e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1389e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
13900298fd71SBarry Smith     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,NULL);CHKERRQ(ierr);}
1391207556f9SJed Brown     if (rmap) {
1392e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1393e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
1394e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1395e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1396207556f9SJed Brown     } else {
1397207556f9SJed Brown       nlocal = 0;
13980298fd71SBarry Smith       isloc  = NULL;
1399207556f9SJed Brown     }
1400e2d7f03fSJed Brown     vs->islocal.row[i] = isloc;
1401e2d7f03fSJed Brown     offset            += nlocal;
1402e2d7f03fSJed Brown   }
14038188e55aSJed Brown   for (i=0,offset=0; i<vs->nc; i++) {
1404e2d7f03fSJed Brown     IS                     isloc;
14050298fd71SBarry Smith     ISLocalToGlobalMapping cmap = NULL;
1406e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1407e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
14080298fd71SBarry Smith     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,NULL,&cmap);CHKERRQ(ierr);}
1409207556f9SJed Brown     if (cmap) {
1410e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1411e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr);
1412e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1413e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1414207556f9SJed Brown     } else {
1415207556f9SJed Brown       nlocal = 0;
14160298fd71SBarry Smith       isloc  = NULL;
1417207556f9SJed Brown     }
1418e2d7f03fSJed Brown     vs->islocal.col[i] = isloc;
1419e2d7f03fSJed Brown     offset            += nlocal;
1420e2d7f03fSJed Brown   }
14210189643fSJed Brown 
142277019fcaSJed Brown   /* Set up the aggregate ISLocalToGlobalMapping */
142377019fcaSJed Brown   {
142445b6f7e9SBarry Smith     ISLocalToGlobalMapping rmap,cmap;
142545b6f7e9SBarry Smith     ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap);CHKERRQ(ierr);
142645b6f7e9SBarry Smith     ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap);CHKERRQ(ierr);
142777019fcaSJed Brown     if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);}
142877019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
142977019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
143077019fcaSJed Brown   }
143177019fcaSJed Brown 
14320189643fSJed Brown #if defined(PETSC_USE_DEBUG)
14330189643fSJed Brown   for (i=0; i<vs->nr; i++) {
14340189643fSJed Brown     for (j=0; j<vs->nc; j++) {
14350189643fSJed Brown       PetscInt m,n,M,N,mi,ni,Mi,Ni;
14360189643fSJed Brown       Mat      B = vs->m[i][j];
14370189643fSJed Brown       if (!B) continue;
14380189643fSJed Brown       ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr);
14390189643fSJed Brown       ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr);
14400189643fSJed Brown       ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr);
14410189643fSJed Brown       ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr);
14420189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr);
14430189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr);
1444ce94432eSBarry 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);
1445ce94432eSBarry 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);
14460189643fSJed Brown     }
14470189643fSJed Brown   }
14480189643fSJed Brown #endif
1449a061e289SJed Brown 
1450a061e289SJed Brown   /* Set A->assembled if all non-null blocks are currently assembled */
1451a061e289SJed Brown   for (i=0; i<vs->nr; i++) {
1452a061e289SJed Brown     for (j=0; j<vs->nc; j++) {
1453a061e289SJed Brown       if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0);
1454a061e289SJed Brown     }
1455a061e289SJed Brown   }
1456a061e289SJed Brown   A->assembled = PETSC_TRUE;
1457d8588912SDave May   PetscFunctionReturn(0);
1458d8588912SDave May }
1459d8588912SDave May 
1460d8588912SDave May #undef __FUNCT__
1461d8588912SDave May #define __FUNCT__ "MatCreateNest"
146245c38901SJed Brown /*@C
1463659c6bb0SJed Brown    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
1464659c6bb0SJed Brown 
1465659c6bb0SJed Brown    Collective on Mat
1466659c6bb0SJed Brown 
1467659c6bb0SJed Brown    Input Parameter:
1468659c6bb0SJed Brown +  comm - Communicator for the new Mat
1469659c6bb0SJed Brown .  nr - number of nested row blocks
14700298fd71SBarry Smith .  is_row - index sets for each nested row block, or NULL to make contiguous
1471659c6bb0SJed Brown .  nc - number of nested column blocks
14720298fd71SBarry Smith .  is_col - index sets for each nested column block, or NULL to make contiguous
14730298fd71SBarry Smith -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1474659c6bb0SJed Brown 
1475659c6bb0SJed Brown    Output Parameter:
1476659c6bb0SJed Brown .  B - new matrix
1477659c6bb0SJed Brown 
1478659c6bb0SJed Brown    Level: advanced
1479659c6bb0SJed Brown 
1480950540a4SJed Brown .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST
1481659c6bb0SJed Brown @*/
14827087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1483d8588912SDave May {
1484d8588912SDave May   Mat            A;
1485d8588912SDave May   PetscErrorCode ierr;
1486d8588912SDave May 
1487d8588912SDave May   PetscFunctionBegin;
1488c8883902SJed Brown   *B   = 0;
1489d8588912SDave May   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1490c8883902SJed Brown   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
149191a28eb3SBarry Smith   A->preallocated = PETSC_TRUE;
1492c8883902SJed Brown   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1493d8588912SDave May   *B   = A;
1494d8588912SDave May   PetscFunctionReturn(0);
1495d8588912SDave May }
1496659c6bb0SJed Brown 
1497629c3df2SDmitry Karpeev #undef __FUNCT__
1498*b68353e5Sstefano_zampini #define __FUNCT__ "MatConvert_Nest_SeqAIJ_fast"
1499*b68353e5Sstefano_zampini static PetscErrorCode MatConvert_Nest_SeqAIJ_fast(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1500*b68353e5Sstefano_zampini {
1501*b68353e5Sstefano_zampini   Mat_Nest       *nest = (Mat_Nest*)A->data;
1502*b68353e5Sstefano_zampini   PetscScalar    **avv;
1503*b68353e5Sstefano_zampini   PetscScalar    *vv;
1504*b68353e5Sstefano_zampini   PetscInt       **aii,**ajj;
1505*b68353e5Sstefano_zampini   PetscInt       *ii,*jj,*ci;
1506*b68353e5Sstefano_zampini   PetscInt       nr,nc,nnz,i,j;
1507*b68353e5Sstefano_zampini   PetscBool      done;
1508*b68353e5Sstefano_zampini   PetscErrorCode ierr;
1509*b68353e5Sstefano_zampini 
1510*b68353e5Sstefano_zampini   PetscFunctionBegin;
1511*b68353e5Sstefano_zampini   ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr);
1512*b68353e5Sstefano_zampini   if (reuse == MAT_REUSE_MATRIX) {
1513*b68353e5Sstefano_zampini     PetscInt rnr;
1514*b68353e5Sstefano_zampini 
1515*b68353e5Sstefano_zampini     ierr = MatGetRowIJ(*newmat,0,PETSC_FALSE,PETSC_FALSE,&rnr,(const PetscInt**)&ii,(const PetscInt**)&jj,&done);CHKERRQ(ierr);
1516*b68353e5Sstefano_zampini     if (!done) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"MatGetRowIJ");
1517*b68353e5Sstefano_zampini     if (rnr != nr) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Cannot reuse matrix, wrong number of rows");
1518*b68353e5Sstefano_zampini     ierr = MatSeqAIJGetArray(*newmat,&vv);CHKERRQ(ierr);
1519*b68353e5Sstefano_zampini   }
1520*b68353e5Sstefano_zampini   /* extract CSR for nested SeqAIJ matrices */
1521*b68353e5Sstefano_zampini   nnz  = 0;
1522*b68353e5Sstefano_zampini   ierr = PetscCalloc3(nest->nr*nest->nc,&aii,nest->nr*nest->nc,&ajj,nest->nr*nest->nc,&avv);CHKERRQ(ierr);
1523*b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1524*b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1525*b68353e5Sstefano_zampini       Mat B = nest->m[i][j];
1526*b68353e5Sstefano_zampini       if (B) {
1527*b68353e5Sstefano_zampini         PetscScalar *naa;
1528*b68353e5Sstefano_zampini         PetscInt    *nii,*njj,nnr;
1529*b68353e5Sstefano_zampini 
1530*b68353e5Sstefano_zampini         ierr = MatGetRowIJ(B,0,PETSC_FALSE,PETSC_FALSE,&nnr,(const PetscInt**)&nii,(const PetscInt**)&njj,&done);CHKERRQ(ierr);
1531*b68353e5Sstefano_zampini         if (!done) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"MatGetRowIJ");
1532*b68353e5Sstefano_zampini         ierr = MatSeqAIJGetArray(B,&naa);CHKERRQ(ierr);
1533*b68353e5Sstefano_zampini         nnz += nii[nnr];
1534*b68353e5Sstefano_zampini 
1535*b68353e5Sstefano_zampini         aii[i*nest->nc+j] = nii;
1536*b68353e5Sstefano_zampini         ajj[i*nest->nc+j] = njj;
1537*b68353e5Sstefano_zampini         avv[i*nest->nc+j] = naa;
1538*b68353e5Sstefano_zampini       }
1539*b68353e5Sstefano_zampini     }
1540*b68353e5Sstefano_zampini   }
1541*b68353e5Sstefano_zampini   if (reuse != MAT_REUSE_MATRIX) {
1542*b68353e5Sstefano_zampini     ierr = PetscMalloc1(nr+1,&ii);CHKERRQ(ierr);
1543*b68353e5Sstefano_zampini     ierr = PetscMalloc1(nnz,&jj);CHKERRQ(ierr);
1544*b68353e5Sstefano_zampini     ierr = PetscMalloc1(nnz,&vv);CHKERRQ(ierr);
1545*b68353e5Sstefano_zampini   } else {
1546*b68353e5Sstefano_zampini     if (nnz != ii[nr]) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Cannot reuse matrix, wrong number of nonzeros");
1547*b68353e5Sstefano_zampini   }
1548*b68353e5Sstefano_zampini 
1549*b68353e5Sstefano_zampini   /* new row pointer */
1550*b68353e5Sstefano_zampini   ierr = PetscMemzero(ii,(nr+1)*sizeof(PetscInt));CHKERRQ(ierr);
1551*b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1552*b68353e5Sstefano_zampini     PetscInt       ncr,rst;
1553*b68353e5Sstefano_zampini 
1554*b68353e5Sstefano_zampini     ierr = ISStrideGetInfo(nest->isglobal.row[i],&rst,NULL);CHKERRQ(ierr);
1555*b68353e5Sstefano_zampini     ierr = ISGetLocalSize(nest->isglobal.row[i],&ncr);CHKERRQ(ierr);
1556*b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1557*b68353e5Sstefano_zampini       if (aii[i*nest->nc+j]) {
1558*b68353e5Sstefano_zampini         PetscInt    *nii = aii[i*nest->nc+j];
1559*b68353e5Sstefano_zampini         PetscInt    ir;
1560*b68353e5Sstefano_zampini 
1561*b68353e5Sstefano_zampini         for (ir=rst; ir<ncr+rst; ++ir) {
1562*b68353e5Sstefano_zampini           ii[ir+1] += nii[1]-nii[0];
1563*b68353e5Sstefano_zampini           nii++;
1564*b68353e5Sstefano_zampini         }
1565*b68353e5Sstefano_zampini       }
1566*b68353e5Sstefano_zampini     }
1567*b68353e5Sstefano_zampini   }
1568*b68353e5Sstefano_zampini   for (i=0; i<nr; i++) ii[i+1] += ii[i];
1569*b68353e5Sstefano_zampini 
1570*b68353e5Sstefano_zampini   /* construct CSR for the new matrix */
1571*b68353e5Sstefano_zampini   ierr = PetscCalloc1(nr,&ci);CHKERRQ(ierr);
1572*b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1573*b68353e5Sstefano_zampini     PetscInt       ncr,rst;
1574*b68353e5Sstefano_zampini 
1575*b68353e5Sstefano_zampini     ierr = ISStrideGetInfo(nest->isglobal.row[i],&rst,NULL);CHKERRQ(ierr);
1576*b68353e5Sstefano_zampini     ierr = ISGetLocalSize(nest->isglobal.row[i],&ncr);CHKERRQ(ierr);
1577*b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1578*b68353e5Sstefano_zampini       if (aii[i*nest->nc+j]) {
1579*b68353e5Sstefano_zampini         PetscScalar *nvv = avv[i*nest->nc+j];
1580*b68353e5Sstefano_zampini         PetscInt    *nii = aii[i*nest->nc+j];
1581*b68353e5Sstefano_zampini         PetscInt    *njj = ajj[i*nest->nc+j];
1582*b68353e5Sstefano_zampini         PetscInt    ir,cst;
1583*b68353e5Sstefano_zampini 
1584*b68353e5Sstefano_zampini         ierr = ISStrideGetInfo(nest->isglobal.col[j],&cst,NULL);CHKERRQ(ierr);
1585*b68353e5Sstefano_zampini         for (ir=rst; ir<ncr+rst; ++ir) {
1586*b68353e5Sstefano_zampini           PetscInt ij,rsize = nii[1]-nii[0],ist = ii[ir]+ci[ir];
1587*b68353e5Sstefano_zampini 
1588*b68353e5Sstefano_zampini           for (ij=0;ij<rsize;ij++) {
1589*b68353e5Sstefano_zampini             jj[ist+ij] = *njj+cst;
1590*b68353e5Sstefano_zampini             vv[ist+ij] = *nvv;
1591*b68353e5Sstefano_zampini             njj++;
1592*b68353e5Sstefano_zampini             nvv++;
1593*b68353e5Sstefano_zampini           }
1594*b68353e5Sstefano_zampini           ci[ir] += rsize;
1595*b68353e5Sstefano_zampini           nii++;
1596*b68353e5Sstefano_zampini         }
1597*b68353e5Sstefano_zampini       }
1598*b68353e5Sstefano_zampini     }
1599*b68353e5Sstefano_zampini   }
1600*b68353e5Sstefano_zampini   ierr = PetscFree(ci);CHKERRQ(ierr);
1601*b68353e5Sstefano_zampini 
1602*b68353e5Sstefano_zampini   /* restore info */
1603*b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1604*b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1605*b68353e5Sstefano_zampini       Mat B = nest->m[i][j];
1606*b68353e5Sstefano_zampini       if (B) {
1607*b68353e5Sstefano_zampini         PetscInt nnr = 0, k = i*nest->nc+j;
1608*b68353e5Sstefano_zampini         ierr = MatRestoreRowIJ(B,0,PETSC_FALSE,PETSC_FALSE,&nnr,(const PetscInt**)&aii[k],(const PetscInt**)&ajj[k],&done);CHKERRQ(ierr);
1609*b68353e5Sstefano_zampini         if (!done) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"MatRestoreRowIJ");
1610*b68353e5Sstefano_zampini         ierr = MatSeqAIJRestoreArray(B,&avv[k]);CHKERRQ(ierr);
1611*b68353e5Sstefano_zampini       }
1612*b68353e5Sstefano_zampini     }
1613*b68353e5Sstefano_zampini   }
1614*b68353e5Sstefano_zampini   ierr = PetscFree3(aii,ajj,avv);CHKERRQ(ierr);
1615*b68353e5Sstefano_zampini 
1616*b68353e5Sstefano_zampini   /* finalize newmat */
1617*b68353e5Sstefano_zampini   if (reuse == MAT_INITIAL_MATRIX) {
1618*b68353e5Sstefano_zampini     ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),nr,nc,ii,jj,vv,newmat);CHKERRQ(ierr);
1619*b68353e5Sstefano_zampini   } else if (reuse == MAT_INPLACE_MATRIX) {
1620*b68353e5Sstefano_zampini     Mat B;
1621*b68353e5Sstefano_zampini 
1622*b68353e5Sstefano_zampini     ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),nr,nc,ii,jj,vv,&B);CHKERRQ(ierr);
1623*b68353e5Sstefano_zampini     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
1624*b68353e5Sstefano_zampini   }
1625*b68353e5Sstefano_zampini   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1626*b68353e5Sstefano_zampini   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1627*b68353e5Sstefano_zampini   {
1628*b68353e5Sstefano_zampini     Mat_SeqAIJ *a = (Mat_SeqAIJ*)((*newmat)->data);
1629*b68353e5Sstefano_zampini     a->free_a     = PETSC_TRUE;
1630*b68353e5Sstefano_zampini     a->free_ij    = PETSC_TRUE;
1631*b68353e5Sstefano_zampini   }
1632*b68353e5Sstefano_zampini   PetscFunctionReturn(0);
1633*b68353e5Sstefano_zampini }
1634*b68353e5Sstefano_zampini 
1635*b68353e5Sstefano_zampini #undef __FUNCT__
1636629c3df2SDmitry Karpeev #define __FUNCT__ "MatConvert_Nest_AIJ"
1637cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1638629c3df2SDmitry Karpeev {
1639629c3df2SDmitry Karpeev   PetscErrorCode ierr;
1640629c3df2SDmitry Karpeev   Mat_Nest       *nest = (Mat_Nest*)A->data;
164183b1a929SMark Adams   PetscInt       m,n,M,N,i,j,k,*dnnz,*onnz,rstart;
1642649b366bSFande Kong   PetscInt       cstart,cend;
1643*b68353e5Sstefano_zampini   PetscMPIInt    size;
1644629c3df2SDmitry Karpeev   Mat            C;
1645629c3df2SDmitry Karpeev 
1646629c3df2SDmitry Karpeev   PetscFunctionBegin;
1647*b68353e5Sstefano_zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1648*b68353e5Sstefano_zampini   if (size == 1) { /* look for a special case with SeqAIJ matrices and strided-1, contiguous, blocks */
1649*b68353e5Sstefano_zampini     PetscInt  nf;
1650*b68353e5Sstefano_zampini     PetscBool fast;
1651*b68353e5Sstefano_zampini 
1652*b68353e5Sstefano_zampini     ierr = PetscStrcmp(newtype,MATAIJ,&fast);CHKERRQ(ierr);
1653*b68353e5Sstefano_zampini     if (!fast) {
1654*b68353e5Sstefano_zampini       ierr = PetscStrcmp(newtype,MATSEQAIJ,&fast);CHKERRQ(ierr);
1655*b68353e5Sstefano_zampini     }
1656*b68353e5Sstefano_zampini     for (i=0; i<nest->nr && fast; ++i) {
1657*b68353e5Sstefano_zampini       for (j=0; j<nest->nc && fast; ++j) {
1658*b68353e5Sstefano_zampini         Mat B = nest->m[i][j];
1659*b68353e5Sstefano_zampini         if (B) {
1660*b68353e5Sstefano_zampini           ierr = PetscObjectTypeCompare((PetscObject)B,MATAIJ,&fast);CHKERRQ(ierr);
1661*b68353e5Sstefano_zampini           if (!fast) {
1662*b68353e5Sstefano_zampini             ierr = PetscObjectTypeCompare((PetscObject)B,MATSEQAIJ,&fast);CHKERRQ(ierr);
1663*b68353e5Sstefano_zampini           }
1664*b68353e5Sstefano_zampini         }
1665*b68353e5Sstefano_zampini       }
1666*b68353e5Sstefano_zampini     }
1667*b68353e5Sstefano_zampini     for (i=0, nf=0; i<nest->nr && fast; ++i) {
1668*b68353e5Sstefano_zampini       ierr = PetscObjectTypeCompare((PetscObject)nest->isglobal.row[i],ISSTRIDE,&fast);CHKERRQ(ierr);
1669*b68353e5Sstefano_zampini       if (fast) {
1670*b68353e5Sstefano_zampini         PetscInt f,s;
1671*b68353e5Sstefano_zampini 
1672*b68353e5Sstefano_zampini         ierr = ISStrideGetInfo(nest->isglobal.row[i],&f,&s);CHKERRQ(ierr);
1673*b68353e5Sstefano_zampini         if (f != nf || s != 1) { fast = PETSC_FALSE; }
1674*b68353e5Sstefano_zampini         else {
1675*b68353e5Sstefano_zampini           ierr = ISGetSize(nest->isglobal.row[i],&f);CHKERRQ(ierr);
1676*b68353e5Sstefano_zampini           nf  += f;
1677*b68353e5Sstefano_zampini         }
1678*b68353e5Sstefano_zampini       }
1679*b68353e5Sstefano_zampini     }
1680*b68353e5Sstefano_zampini     for (i=0, nf=0; i<nest->nc && fast; ++i) {
1681*b68353e5Sstefano_zampini       ierr = PetscObjectTypeCompare((PetscObject)nest->isglobal.col[i],ISSTRIDE,&fast);CHKERRQ(ierr);
1682*b68353e5Sstefano_zampini       if (fast) {
1683*b68353e5Sstefano_zampini         PetscInt f,s;
1684*b68353e5Sstefano_zampini 
1685*b68353e5Sstefano_zampini         ierr = ISStrideGetInfo(nest->isglobal.col[i],&f,&s);CHKERRQ(ierr);
1686*b68353e5Sstefano_zampini         if (f != nf || s != 1) { fast = PETSC_FALSE; }
1687*b68353e5Sstefano_zampini         else {
1688*b68353e5Sstefano_zampini           ierr = ISGetSize(nest->isglobal.col[i],&f);CHKERRQ(ierr);
1689*b68353e5Sstefano_zampini           nf  += f;
1690*b68353e5Sstefano_zampini         }
1691*b68353e5Sstefano_zampini       }
1692*b68353e5Sstefano_zampini     }
1693*b68353e5Sstefano_zampini     if (fast) {
1694*b68353e5Sstefano_zampini       ierr = MatConvert_Nest_SeqAIJ_fast(A,newtype,reuse,newmat);CHKERRQ(ierr);
1695*b68353e5Sstefano_zampini       PetscFunctionReturn(0);
1696*b68353e5Sstefano_zampini     }
1697*b68353e5Sstefano_zampini   }
1698629c3df2SDmitry Karpeev   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
1699629c3df2SDmitry Karpeev   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
1700649b366bSFande Kong   ierr = MatGetOwnershipRangeColumn(A,&cstart,&cend);CHKERRQ(ierr);
1701629c3df2SDmitry Karpeev   switch (reuse) {
1702629c3df2SDmitry Karpeev   case MAT_INITIAL_MATRIX:
1703ce94432eSBarry Smith     ierr    = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
1704629c3df2SDmitry Karpeev     ierr    = MatSetType(C,newtype);CHKERRQ(ierr);
1705629c3df2SDmitry Karpeev     ierr    = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr);
1706629c3df2SDmitry Karpeev     *newmat = C;
1707629c3df2SDmitry Karpeev     break;
1708629c3df2SDmitry Karpeev   case MAT_REUSE_MATRIX:
1709629c3df2SDmitry Karpeev     C = *newmat;
1710629c3df2SDmitry Karpeev     break;
1711ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse");
1712629c3df2SDmitry Karpeev   }
1713785e854fSJed Brown   ierr = PetscMalloc1(2*m,&dnnz);CHKERRQ(ierr);
1714629c3df2SDmitry Karpeev   onnz = dnnz + m;
1715629c3df2SDmitry Karpeev   for (k=0; k<m; k++) {
1716629c3df2SDmitry Karpeev     dnnz[k] = 0;
1717629c3df2SDmitry Karpeev     onnz[k] = 0;
1718629c3df2SDmitry Karpeev   }
1719629c3df2SDmitry Karpeev   for (j=0; j<nest->nc; ++j) {
1720629c3df2SDmitry Karpeev     IS             bNis;
1721629c3df2SDmitry Karpeev     PetscInt       bN;
1722629c3df2SDmitry Karpeev     const PetscInt *bNindices;
1723629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
1724629c3df2SDmitry Karpeev     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1725629c3df2SDmitry Karpeev     ierr = ISGetSize(bNis, &bN);CHKERRQ(ierr);
1726629c3df2SDmitry Karpeev     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1727629c3df2SDmitry Karpeev     for (i=0; i<nest->nr; ++i) {
1728629c3df2SDmitry Karpeev       PetscSF        bmsf;
1729649b366bSFande Kong       PetscSFNode    *iremote;
1730629c3df2SDmitry Karpeev       Mat            B;
1731649b366bSFande Kong       PetscInt       bm, *sub_dnnz,*sub_onnz, br;
1732629c3df2SDmitry Karpeev       const PetscInt *bmindices;
1733629c3df2SDmitry Karpeev       B = nest->m[i][j];
1734629c3df2SDmitry Karpeev       if (!B) continue;
1735629c3df2SDmitry Karpeev       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1736629c3df2SDmitry Karpeev       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1737ce94432eSBarry Smith       ierr = PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);CHKERRQ(ierr);
1738649b366bSFande Kong       ierr = PetscMalloc1(bm,&iremote);CHKERRQ(ierr);
1739649b366bSFande Kong       ierr = PetscMalloc1(bm,&sub_dnnz);CHKERRQ(ierr);
1740649b366bSFande Kong       ierr = PetscMalloc1(bm,&sub_onnz);CHKERRQ(ierr);
1741649b366bSFande Kong       for (k = 0; k < bm; ++k){
1742649b366bSFande Kong     	sub_dnnz[k] = 0;
1743649b366bSFande Kong     	sub_onnz[k] = 0;
1744649b366bSFande Kong       }
1745629c3df2SDmitry Karpeev       /*
1746629c3df2SDmitry Karpeev        Locate the owners for all of the locally-owned global row indices for this row block.
1747629c3df2SDmitry Karpeev        These determine the roots of PetscSF used to communicate preallocation data to row owners.
1748629c3df2SDmitry Karpeev        The roots correspond to the dnnz and onnz entries; thus, there are two roots per row.
1749629c3df2SDmitry Karpeev        */
175083b1a929SMark Adams       ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr);
1751629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
1752649b366bSFande Kong         PetscInt       row = bmindices[br], rowowner = 0, brncols, col;
1753629c3df2SDmitry Karpeev         const PetscInt *brcols;
1754a4b3d3acSMatthew G Knepley         PetscInt       rowrel = 0; /* row's relative index on its owner rank */
1755629c3df2SDmitry Karpeev         ierr      = PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);CHKERRQ(ierr);
1756649b366bSFande Kong         /* how many roots  */
1757649b366bSFande Kong         iremote[br].rank = rowowner; iremote[br].index = rowrel;           /* edge from bmdnnz to dnnz */
1758649b366bSFande Kong         /* get nonzero pattern */
175983b1a929SMark Adams         ierr = MatGetRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr);
1760629c3df2SDmitry Karpeev         for (k=0; k<brncols; k++) {
1761629c3df2SDmitry Karpeev           col  = bNindices[brcols[k]];
1762649b366bSFande Kong           if (col>=A->cmap->range[rowowner] && col<A->cmap->range[rowowner+1]) {
1763649b366bSFande Kong             sub_dnnz[br]++;
1764649b366bSFande Kong           } else {
1765649b366bSFande Kong             sub_onnz[br]++;
1766649b366bSFande Kong           }
1767629c3df2SDmitry Karpeev         }
176883b1a929SMark Adams         ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr);
1769629c3df2SDmitry Karpeev       }
1770629c3df2SDmitry Karpeev       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1771629c3df2SDmitry Karpeev       /* bsf will have to take care of disposing of bedges. */
1772649b366bSFande Kong       ierr = PetscSFSetGraph(bmsf,m,bm,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1773649b366bSFande Kong       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr);
1774649b366bSFande Kong       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr);
1775649b366bSFande Kong       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr);
1776649b366bSFande Kong       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr);
1777649b366bSFande Kong       ierr = PetscFree(sub_dnnz);CHKERRQ(ierr);
1778649b366bSFande Kong       ierr = PetscFree(sub_onnz);CHKERRQ(ierr);
1779629c3df2SDmitry Karpeev       ierr = PetscSFDestroy(&bmsf);CHKERRQ(ierr);
1780629c3df2SDmitry Karpeev     }
178122d28d08SBarry Smith     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
1782629c3df2SDmitry Karpeev     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
178365a4a0a3Sstefano_zampini   }
178465a4a0a3Sstefano_zampini   /* Resize preallocation if overestimated */
178565a4a0a3Sstefano_zampini   for (i=0;i<m;i++) {
178665a4a0a3Sstefano_zampini     dnnz[i] = PetscMin(dnnz[i],A->cmap->n);
178765a4a0a3Sstefano_zampini     onnz[i] = PetscMin(onnz[i],A->cmap->N - A->cmap->n);
1788629c3df2SDmitry Karpeev   }
1789629c3df2SDmitry Karpeev   ierr = MatSeqAIJSetPreallocation(C,0,dnnz);CHKERRQ(ierr);
1790629c3df2SDmitry Karpeev   ierr = MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);CHKERRQ(ierr);
1791629c3df2SDmitry Karpeev   ierr = PetscFree(dnnz);CHKERRQ(ierr);
1792629c3df2SDmitry Karpeev 
1793629c3df2SDmitry Karpeev   /* Fill by row */
1794629c3df2SDmitry Karpeev   for (j=0; j<nest->nc; ++j) {
1795629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
1796629c3df2SDmitry Karpeev     IS             bNis;
1797629c3df2SDmitry Karpeev     PetscInt       bN;
1798629c3df2SDmitry Karpeev     const PetscInt *bNindices;
1799629c3df2SDmitry Karpeev     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1800629c3df2SDmitry Karpeev     ierr = ISGetSize(bNis,&bN);CHKERRQ(ierr);
1801629c3df2SDmitry Karpeev     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1802629c3df2SDmitry Karpeev     for (i=0; i<nest->nr; ++i) {
1803629c3df2SDmitry Karpeev       Mat            B;
1804629c3df2SDmitry Karpeev       PetscInt       bm, br;
1805629c3df2SDmitry Karpeev       const PetscInt *bmindices;
1806629c3df2SDmitry Karpeev       B = nest->m[i][j];
1807629c3df2SDmitry Karpeev       if (!B) continue;
1808629c3df2SDmitry Karpeev       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1809629c3df2SDmitry Karpeev       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
181083b1a929SMark Adams       ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr);
1811629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
1812629c3df2SDmitry Karpeev         PetscInt          row = bmindices[br], brncols,  *cols;
1813629c3df2SDmitry Karpeev         const PetscInt    *brcols;
1814629c3df2SDmitry Karpeev         const PetscScalar *brcoldata;
181583b1a929SMark Adams         ierr = MatGetRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
1816785e854fSJed Brown         ierr = PetscMalloc1(brncols,&cols);CHKERRQ(ierr);
181726fbe8dcSKarl Rupp         for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]];
1818629c3df2SDmitry Karpeev         /*
1819629c3df2SDmitry Karpeev           Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match.
1820629c3df2SDmitry Karpeev           Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES.
1821629c3df2SDmitry Karpeev          */
1822a2ea699eSBarry Smith         ierr = MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);CHKERRQ(ierr);
182383b1a929SMark Adams         ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
1824629c3df2SDmitry Karpeev         ierr = PetscFree(cols);CHKERRQ(ierr);
1825629c3df2SDmitry Karpeev       }
1826629c3df2SDmitry Karpeev       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1827629c3df2SDmitry Karpeev     }
1828a2ea699eSBarry Smith     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
1829629c3df2SDmitry Karpeev     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
1830629c3df2SDmitry Karpeev   }
1831629c3df2SDmitry Karpeev   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1832629c3df2SDmitry Karpeev   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1833629c3df2SDmitry Karpeev   PetscFunctionReturn(0);
1834629c3df2SDmitry Karpeev }
1835629c3df2SDmitry Karpeev 
1836659c6bb0SJed Brown /*MC
1837659c6bb0SJed Brown   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
1838659c6bb0SJed Brown 
1839659c6bb0SJed Brown   Level: intermediate
1840659c6bb0SJed Brown 
1841659c6bb0SJed Brown   Notes:
1842659c6bb0SJed Brown   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1843659c6bb0SJed Brown   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1844950540a4SJed Brown   It is usually used with DMComposite and DMCreateMatrix()
1845659c6bb0SJed Brown 
1846659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest()
1847659c6bb0SJed Brown M*/
1848c8883902SJed Brown #undef __FUNCT__
1849c8883902SJed Brown #define __FUNCT__ "MatCreate_Nest"
18508cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A)
1851c8883902SJed Brown {
1852c8883902SJed Brown   Mat_Nest       *s;
1853c8883902SJed Brown   PetscErrorCode ierr;
1854c8883902SJed Brown 
1855c8883902SJed Brown   PetscFunctionBegin;
1856b00a9115SJed Brown   ierr    = PetscNewLog(A,&s);CHKERRQ(ierr);
1857c8883902SJed Brown   A->data = (void*)s;
1858e7c19651SJed Brown 
1859e7c19651SJed Brown   s->nr            = -1;
1860e7c19651SJed Brown   s->nc            = -1;
18610298fd71SBarry Smith   s->m             = NULL;
1862e7c19651SJed Brown   s->splitassembly = PETSC_FALSE;
1863c8883902SJed Brown 
1864c8883902SJed Brown   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
186526fbe8dcSKarl Rupp 
1866c8883902SJed Brown   A->ops->mult                  = MatMult_Nest;
18679194d70fSJed Brown   A->ops->multadd               = MatMultAdd_Nest;
1868c8883902SJed Brown   A->ops->multtranspose         = MatMultTranspose_Nest;
18699194d70fSJed Brown   A->ops->multtransposeadd      = MatMultTransposeAdd_Nest;
1870f8170845SAlex Fikl   A->ops->transpose             = MatTranspose_Nest;
1871c8883902SJed Brown   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1872c8883902SJed Brown   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1873c8883902SJed Brown   A->ops->zeroentries           = MatZeroEntries_Nest;
1874c222c20dSDavid Ham   A->ops->copy                  = MatCopy_Nest;
1875c8883902SJed Brown   A->ops->duplicate             = MatDuplicate_Nest;
1876c8883902SJed Brown   A->ops->getsubmatrix          = MatGetSubMatrix_Nest;
1877c8883902SJed Brown   A->ops->destroy               = MatDestroy_Nest;
1878c8883902SJed Brown   A->ops->view                  = MatView_Nest;
1879c8883902SJed Brown   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1880c8883902SJed Brown   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1881c8883902SJed Brown   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
1882429bac76SJed Brown   A->ops->getdiagonal           = MatGetDiagonal_Nest;
1883429bac76SJed Brown   A->ops->diagonalscale         = MatDiagonalScale_Nest;
1884a061e289SJed Brown   A->ops->scale                 = MatScale_Nest;
1885a061e289SJed Brown   A->ops->shift                 = MatShift_Nest;
188613135bc6SAlex Fikl   A->ops->diagonalset           = MatDiagonalSet_Nest;
1887f8170845SAlex Fikl   A->ops->setrandom             = MatSetRandom_Nest;
1888c8883902SJed Brown 
1889c8883902SJed Brown   A->spptr        = 0;
1890c8883902SJed Brown   A->assembled    = PETSC_FALSE;
1891c8883902SJed Brown 
1892c8883902SJed Brown   /* expose Nest api's */
1893bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",   MatNestGetSubMat_Nest);CHKERRQ(ierr);
1894bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",   MatNestSetSubMat_Nest);CHKERRQ(ierr);
1895bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",  MatNestGetSubMats_Nest);CHKERRQ(ierr);
1896bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",     MatNestGetSize_Nest);CHKERRQ(ierr);
1897bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",      MatNestGetISs_Nest);CHKERRQ(ierr);
1898bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C", MatNestGetLocalISs_Nest);CHKERRQ(ierr);
1899bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",  MatNestSetVecType_Nest);CHKERRQ(ierr);
1900bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",  MatNestSetSubMats_Nest);CHKERRQ(ierr);
190183b1a929SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",MatConvert_Nest_AIJ);CHKERRQ(ierr);
19025e3038f0Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_is_C", MatConvert_Nest_IS);CHKERRQ(ierr);
1903c8883902SJed Brown 
1904c8883902SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
1905c8883902SJed Brown   PetscFunctionReturn(0);
1906c8883902SJed Brown }
1907