xref: /petsc/src/mat/impls/nest/matnest.c (revision e108cb999e836831742efd3a44f4221fe84ae696)
1d8588912SDave May 
2aaa7dc30SBarry Smith #include <../src/mat/impls/nest/matnestimpl.h> /*I   "petscmat.h"   I*/
3b68353e5Sstefano_zampini #include <../src/mat/impls/aij/seq/aij.h>
40c312b8eSJed Brown #include <petscsf.h>
5d8588912SDave May 
6c8883902SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]);
72a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left);
85e3038f0Sstefano_zampini PETSC_INTERN PetscErrorCode MatConvert_Nest_IS(Mat,MatType,MatReuse,Mat*);
9c8883902SJed Brown 
10d8588912SDave May /* private functions */
118188e55aSJed Brown static PetscErrorCode MatNestGetSizes_Private(Mat A,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N)
12d8588912SDave May {
13d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
148188e55aSJed Brown   PetscInt       i,j;
15d8588912SDave May   PetscErrorCode ierr;
16d8588912SDave May 
17d8588912SDave May   PetscFunctionBegin;
188188e55aSJed Brown   *m = *n = *M = *N = 0;
198188e55aSJed Brown   for (i=0; i<bA->nr; i++) {  /* rows */
208188e55aSJed Brown     PetscInt sm,sM;
218188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.row[i],&sm);CHKERRQ(ierr);
228188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.row[i],&sM);CHKERRQ(ierr);
238188e55aSJed Brown     *m  += sm;
248188e55aSJed Brown     *M  += sM;
25d8588912SDave May   }
268188e55aSJed Brown   for (j=0; j<bA->nc; j++) {  /* cols */
278188e55aSJed Brown     PetscInt sn,sN;
288188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.col[j],&sn);CHKERRQ(ierr);
298188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.col[j],&sN);CHKERRQ(ierr);
308188e55aSJed Brown     *n  += sn;
318188e55aSJed Brown     *N  += sN;
32d8588912SDave May   }
33d8588912SDave May   PetscFunctionReturn(0);
34d8588912SDave May }
35d8588912SDave May 
36d8588912SDave May /* operations */
37207556f9SJed Brown static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y)
38d8588912SDave May {
39d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
40207556f9SJed Brown   Vec            *bx = bA->right,*by = bA->left;
41207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
42d8588912SDave May   PetscErrorCode ierr;
43d8588912SDave May 
44d8588912SDave May   PetscFunctionBegin;
45207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
46207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
47207556f9SJed Brown   for (i=0; i<nr; i++) {
48d8588912SDave May     ierr = VecZeroEntries(by[i]);CHKERRQ(ierr);
49207556f9SJed Brown     for (j=0; j<nc; j++) {
50207556f9SJed Brown       if (!bA->m[i][j]) continue;
51d8588912SDave May       /* y[i] <- y[i] + A[i][j] * x[j] */
52d8588912SDave May       ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr);
53d8588912SDave May     }
54d8588912SDave May   }
55207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
56207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
57d8588912SDave May   PetscFunctionReturn(0);
58d8588912SDave May }
59d8588912SDave May 
609194d70fSJed Brown static PetscErrorCode MatMultAdd_Nest(Mat A,Vec x,Vec y,Vec z)
619194d70fSJed Brown {
629194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
639194d70fSJed Brown   Vec            *bx = bA->right,*bz = bA->left;
649194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
659194d70fSJed Brown   PetscErrorCode ierr;
669194d70fSJed Brown 
679194d70fSJed Brown   PetscFunctionBegin;
689194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
699194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
709194d70fSJed Brown   for (i=0; i<nr; i++) {
719194d70fSJed Brown     if (y != z) {
729194d70fSJed Brown       Vec by;
739194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
749194d70fSJed Brown       ierr = VecCopy(by,bz[i]);CHKERRQ(ierr);
75336d21e7SJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
769194d70fSJed Brown     }
779194d70fSJed Brown     for (j=0; j<nc; j++) {
789194d70fSJed Brown       if (!bA->m[i][j]) continue;
799194d70fSJed Brown       /* y[i] <- y[i] + A[i][j] * x[j] */
809194d70fSJed Brown       ierr = MatMultAdd(bA->m[i][j],bx[j],bz[i],bz[i]);CHKERRQ(ierr);
819194d70fSJed Brown     }
829194d70fSJed Brown   }
839194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
849194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
859194d70fSJed Brown   PetscFunctionReturn(0);
869194d70fSJed Brown }
879194d70fSJed Brown 
88207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y)
89d8588912SDave May {
90d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
91207556f9SJed Brown   Vec            *bx = bA->left,*by = bA->right;
92207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
93d8588912SDave May   PetscErrorCode ierr;
94d8588912SDave May 
95d8588912SDave May   PetscFunctionBegin;
96609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
97609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
98207556f9SJed Brown   for (j=0; j<nc; j++) {
99609e31cbSJed Brown     ierr = VecZeroEntries(by[j]);CHKERRQ(ierr);
100609e31cbSJed Brown     for (i=0; i<nr; i++) {
1016c75ac25SJed Brown       if (!bA->m[i][j]) continue;
102609e31cbSJed Brown       /* y[j] <- y[j] + (A[i][j])^T * x[i] */
103609e31cbSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],by[j],by[j]);CHKERRQ(ierr);
104d8588912SDave May     }
105d8588912SDave May   }
106609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
107609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
108d8588912SDave May   PetscFunctionReturn(0);
109d8588912SDave May }
110d8588912SDave May 
1119194d70fSJed Brown static PetscErrorCode MatMultTransposeAdd_Nest(Mat A,Vec x,Vec y,Vec z)
1129194d70fSJed Brown {
1139194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
1149194d70fSJed Brown   Vec            *bx = bA->left,*bz = bA->right;
1159194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
1169194d70fSJed Brown   PetscErrorCode ierr;
1179194d70fSJed Brown 
1189194d70fSJed Brown   PetscFunctionBegin;
1199194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
1209194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
1219194d70fSJed Brown   for (j=0; j<nc; j++) {
1229194d70fSJed Brown     if (y != z) {
1239194d70fSJed Brown       Vec by;
1249194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
1259194d70fSJed Brown       ierr = VecCopy(by,bz[j]);CHKERRQ(ierr);
1269194d70fSJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
1279194d70fSJed Brown     }
1289194d70fSJed Brown     for (i=0; i<nr; i++) {
1296c75ac25SJed Brown       if (!bA->m[i][j]) continue;
1309194d70fSJed Brown       /* z[j] <- y[j] + (A[i][j])^T * x[i] */
1319194d70fSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],bz[j],bz[j]);CHKERRQ(ierr);
1329194d70fSJed Brown     }
1339194d70fSJed Brown   }
1349194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
1359194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
1369194d70fSJed Brown   PetscFunctionReturn(0);
1379194d70fSJed Brown }
1389194d70fSJed Brown 
139f8170845SAlex Fikl static PetscErrorCode MatTranspose_Nest(Mat A,MatReuse reuse,Mat *B)
140f8170845SAlex Fikl {
141f8170845SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data, *bC;
142f8170845SAlex Fikl   Mat            C;
143f8170845SAlex Fikl   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
144f8170845SAlex Fikl   PetscErrorCode ierr;
145f8170845SAlex Fikl 
146f8170845SAlex Fikl   PetscFunctionBegin;
147cf37664fSBarry Smith   if (reuse == MAT_INPLACE_MATRIX && nr != nc) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Square nested matrix only for in-place");
148f8170845SAlex Fikl 
149cf37664fSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) {
150f8170845SAlex Fikl     Mat *subs;
151f8170845SAlex Fikl     IS  *is_row,*is_col;
152f8170845SAlex Fikl 
153f8170845SAlex Fikl     ierr = PetscCalloc1(nr * nc,&subs);CHKERRQ(ierr);
154f8170845SAlex Fikl     ierr = PetscMalloc2(nr,&is_row,nc,&is_col);CHKERRQ(ierr);
155f8170845SAlex Fikl     ierr = MatNestGetISs(A,is_row,is_col);CHKERRQ(ierr);
156cf37664fSBarry Smith     if (reuse == MAT_INPLACE_MATRIX) {
157ddeb9bd8SAlex Fikl       for (i=0; i<nr; i++) {
158ddeb9bd8SAlex Fikl         for (j=0; j<nc; j++) {
159ddeb9bd8SAlex Fikl           subs[i + nr * j] = bA->m[i][j];
160ddeb9bd8SAlex Fikl         }
161ddeb9bd8SAlex Fikl       }
162ddeb9bd8SAlex Fikl     }
163ddeb9bd8SAlex Fikl 
164f8170845SAlex Fikl     ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nc,is_col,nr,is_row,subs,&C);CHKERRQ(ierr);
165f8170845SAlex Fikl     ierr = PetscFree(subs);CHKERRQ(ierr);
1663d994f23SBarry Smith     ierr = PetscFree2(is_row,is_col);CHKERRQ(ierr);
167f8170845SAlex Fikl   } else {
168f8170845SAlex Fikl     C = *B;
169f8170845SAlex Fikl   }
170f8170845SAlex Fikl 
171f8170845SAlex Fikl   bC = (Mat_Nest*)C->data;
172f8170845SAlex Fikl   for (i=0; i<nr; i++) {
173f8170845SAlex Fikl     for (j=0; j<nc; j++) {
174f8170845SAlex Fikl       if (bA->m[i][j]) {
175f8170845SAlex Fikl         ierr = MatTranspose(bA->m[i][j], reuse, &(bC->m[j][i]));CHKERRQ(ierr);
176f8170845SAlex Fikl       } else {
177f8170845SAlex Fikl         bC->m[j][i] = NULL;
178f8170845SAlex Fikl       }
179f8170845SAlex Fikl     }
180f8170845SAlex Fikl   }
181f8170845SAlex Fikl 
182cf37664fSBarry Smith   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) {
183f8170845SAlex Fikl     *B = C;
184f8170845SAlex Fikl   } else {
185f8170845SAlex Fikl     ierr = MatHeaderMerge(A, &C);CHKERRQ(ierr);
186f8170845SAlex Fikl   }
187f8170845SAlex Fikl   PetscFunctionReturn(0);
188f8170845SAlex Fikl }
189f8170845SAlex Fikl 
190e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list)
191e2d7f03fSJed Brown {
192e2d7f03fSJed Brown   PetscErrorCode ierr;
193e2d7f03fSJed Brown   IS             *lst = *list;
194e2d7f03fSJed Brown   PetscInt       i;
195e2d7f03fSJed Brown 
196e2d7f03fSJed Brown   PetscFunctionBegin;
197e2d7f03fSJed Brown   if (!lst) PetscFunctionReturn(0);
1986bf464f9SBarry Smith   for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);}
199e2d7f03fSJed Brown   ierr  = PetscFree(lst);CHKERRQ(ierr);
2000298fd71SBarry Smith   *list = NULL;
201e2d7f03fSJed Brown   PetscFunctionReturn(0);
202e2d7f03fSJed Brown }
203e2d7f03fSJed Brown 
204207556f9SJed Brown static PetscErrorCode MatDestroy_Nest(Mat A)
205d8588912SDave May {
206d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
207d8588912SDave May   PetscInt       i,j;
208d8588912SDave May   PetscErrorCode ierr;
209d8588912SDave May 
210d8588912SDave May   PetscFunctionBegin;
211d8588912SDave May   /* release the matrices and the place holders */
212e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr);
213e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr);
214e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
215e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
216d8588912SDave May 
217d8588912SDave May   ierr = PetscFree(vs->row_len);CHKERRQ(ierr);
218d8588912SDave May   ierr = PetscFree(vs->col_len);CHKERRQ(ierr);
219d8588912SDave May 
220207556f9SJed Brown   ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr);
221207556f9SJed Brown 
222d8588912SDave May   /* release the matrices and the place holders */
223d8588912SDave May   if (vs->m) {
224d8588912SDave May     for (i=0; i<vs->nr; i++) {
225d8588912SDave May       for (j=0; j<vs->nc; j++) {
2266bf464f9SBarry Smith         ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr);
227d8588912SDave May       }
228d8588912SDave May       ierr = PetscFree(vs->m[i]);CHKERRQ(ierr);
229d8588912SDave May     }
230d8588912SDave May     ierr = PetscFree(vs->m);CHKERRQ(ierr);
231d8588912SDave May   }
232bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
233d8588912SDave May 
234bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",0);CHKERRQ(ierr);
235bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",0);CHKERRQ(ierr);
236bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",0);CHKERRQ(ierr);
237bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",0);CHKERRQ(ierr);
238bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",0);CHKERRQ(ierr);
239bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",0);CHKERRQ(ierr);
240bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",0);CHKERRQ(ierr);
241bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",0);CHKERRQ(ierr);
2420899c546SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_mpiaij_C",0);CHKERRQ(ierr);
2430899c546SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_seqaij_C",0);CHKERRQ(ierr);
2445e3038f0Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",0);CHKERRQ(ierr);
2455e3038f0Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_is_C",0);CHKERRQ(ierr);
246d8588912SDave May   PetscFunctionReturn(0);
247d8588912SDave May }
248d8588912SDave May 
249207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
250d8588912SDave May {
251d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
252d8588912SDave May   PetscInt       i,j;
253d8588912SDave May   PetscErrorCode ierr;
254d8588912SDave May 
255d8588912SDave May   PetscFunctionBegin;
256d8588912SDave May   for (i=0; i<vs->nr; i++) {
257d8588912SDave May     for (j=0; j<vs->nc; j++) {
258e7c19651SJed Brown       if (vs->m[i][j]) {
259e7c19651SJed Brown         ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr);
260e7c19651SJed Brown         if (!vs->splitassembly) {
261e7c19651SJed Brown           /* Note: split assembly will fail if the same block appears more than once (even indirectly through a nested
262e7c19651SJed 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
263e7c19651SJed Brown            * already performing an assembly, but the result would by more complicated and appears to offer less
264e7c19651SJed Brown            * potential for diagnostics and correctness checking. Split assembly should be fixed once there is an
265e7c19651SJed Brown            * interface for libraries to make asynchronous progress in "user-defined non-blocking collectives".
266e7c19651SJed Brown            */
267e7c19651SJed Brown           ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr);
268e7c19651SJed Brown         }
269e7c19651SJed Brown       }
270d8588912SDave May     }
271d8588912SDave May   }
272d8588912SDave May   PetscFunctionReturn(0);
273d8588912SDave May }
274d8588912SDave May 
275207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
276d8588912SDave May {
277d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
278d8588912SDave May   PetscInt       i,j;
279d8588912SDave May   PetscErrorCode ierr;
280d8588912SDave May 
281d8588912SDave May   PetscFunctionBegin;
282d8588912SDave May   for (i=0; i<vs->nr; i++) {
283d8588912SDave May     for (j=0; j<vs->nc; j++) {
284e7c19651SJed Brown       if (vs->m[i][j]) {
285e7c19651SJed Brown         if (vs->splitassembly) {
286e7c19651SJed Brown           ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr);
287e7c19651SJed Brown         }
288e7c19651SJed Brown       }
289d8588912SDave May     }
290d8588912SDave May   }
291d8588912SDave May   PetscFunctionReturn(0);
292d8588912SDave May }
293d8588912SDave May 
294f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
295d8588912SDave May {
296207556f9SJed Brown   PetscErrorCode ierr;
297f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
298f349c1fdSJed Brown   PetscInt       j;
299f349c1fdSJed Brown   Mat            sub;
300d8588912SDave May 
301d8588912SDave May   PetscFunctionBegin;
3020298fd71SBarry Smith   sub = (row < vs->nc) ? vs->m[row][row] : (Mat)NULL; /* Prefer to find on the diagonal */
303f349c1fdSJed Brown   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
3044994cf47SJed Brown   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
305f349c1fdSJed Brown   *B = sub;
306f349c1fdSJed Brown   PetscFunctionReturn(0);
307d8588912SDave May }
308d8588912SDave May 
309f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
310f349c1fdSJed Brown {
311207556f9SJed Brown   PetscErrorCode ierr;
312f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
313f349c1fdSJed Brown   PetscInt       i;
314f349c1fdSJed Brown   Mat            sub;
315f349c1fdSJed Brown 
316f349c1fdSJed Brown   PetscFunctionBegin;
3170298fd71SBarry Smith   sub = (col < vs->nr) ? vs->m[col][col] : (Mat)NULL; /* Prefer to find on the diagonal */
318f349c1fdSJed Brown   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
3194994cf47SJed Brown   if (sub) {ierr = MatSetUp(sub);CHKERRQ(ierr);}       /* Ensure that the sizes are available */
320f349c1fdSJed Brown   *B = sub;
321f349c1fdSJed Brown   PetscFunctionReturn(0);
322d8588912SDave May }
323d8588912SDave May 
324f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
325f349c1fdSJed Brown {
326f349c1fdSJed Brown   PetscErrorCode ierr;
327f349c1fdSJed Brown   PetscInt       i;
328f349c1fdSJed Brown   PetscBool      flg;
329f349c1fdSJed Brown 
330f349c1fdSJed Brown   PetscFunctionBegin;
331f349c1fdSJed Brown   PetscValidPointer(list,3);
332f349c1fdSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,4);
333f349c1fdSJed Brown   PetscValidIntPointer(found,5);
334f349c1fdSJed Brown   *found = -1;
335f349c1fdSJed Brown   for (i=0; i<n; i++) {
336207556f9SJed Brown     if (!list[i]) continue;
337f349c1fdSJed Brown     ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr);
338f349c1fdSJed Brown     if (flg) {
339f349c1fdSJed Brown       *found = i;
340f349c1fdSJed Brown       PetscFunctionReturn(0);
341f349c1fdSJed Brown     }
342f349c1fdSJed Brown   }
343ce94432eSBarry Smith   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Could not find index set");
344f349c1fdSJed Brown   PetscFunctionReturn(0);
345f349c1fdSJed Brown }
346f349c1fdSJed Brown 
3478188e55aSJed Brown /* Get a block row as a new MatNest */
3488188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B)
3498188e55aSJed Brown {
3508188e55aSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
3518188e55aSJed Brown   char           keyname[256];
3528188e55aSJed Brown   PetscErrorCode ierr;
3538188e55aSJed Brown 
3548188e55aSJed Brown   PetscFunctionBegin;
3550298fd71SBarry Smith   *B   = NULL;
3568caf3d72SBarry Smith   ierr = PetscSNPrintf(keyname,sizeof(keyname),"NestRow_%D",row);CHKERRQ(ierr);
3578188e55aSJed Brown   ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr);
3588188e55aSJed Brown   if (*B) PetscFunctionReturn(0);
3598188e55aSJed Brown 
360ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),1,NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr);
36126fbe8dcSKarl Rupp 
3628188e55aSJed Brown   (*B)->assembled = A->assembled;
36326fbe8dcSKarl Rupp 
3648188e55aSJed Brown   ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr);
3658188e55aSJed Brown   ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */
3668188e55aSJed Brown   PetscFunctionReturn(0);
3678188e55aSJed Brown }
3688188e55aSJed Brown 
369f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
370f349c1fdSJed Brown {
371f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
3728188e55aSJed Brown   PetscErrorCode ierr;
3736b3a5b13SJed Brown   PetscInt       row,col;
374e072481dSJed Brown   PetscBool      same,isFullCol,isFullColGlobal;
375f349c1fdSJed Brown 
376f349c1fdSJed Brown   PetscFunctionBegin;
3778188e55aSJed Brown   /* Check if full column space. This is a hack */
3788188e55aSJed Brown   isFullCol = PETSC_FALSE;
379251f4c67SDmitry Karpeev   ierr      = PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr);
3808188e55aSJed Brown   if (same) {
38177019fcaSJed Brown     PetscInt n,first,step,i,an,am,afirst,astep;
3828188e55aSJed Brown     ierr      = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
3838188e55aSJed Brown     ierr      = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
38477019fcaSJed Brown     isFullCol = PETSC_TRUE;
38505ce4453SJed Brown     for (i=0,an=A->cmap->rstart; i<vs->nc; i++) {
38677019fcaSJed Brown       ierr = ISStrideGetInfo(is->col[i],&afirst,&astep);CHKERRQ(ierr);
38777019fcaSJed Brown       ierr = ISGetLocalSize(is->col[i],&am);CHKERRQ(ierr);
38877019fcaSJed Brown       if (afirst != an || astep != step) isFullCol = PETSC_FALSE;
38977019fcaSJed Brown       an += am;
39077019fcaSJed Brown     }
39105ce4453SJed Brown     if (an != A->cmap->rstart+n) isFullCol = PETSC_FALSE;
3928188e55aSJed Brown   }
393b2566f29SBarry Smith   ierr = MPIU_Allreduce(&isFullCol,&isFullColGlobal,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)iscol));CHKERRQ(ierr);
3948188e55aSJed Brown 
395427230ceSLisandro Dalcin   if (isFullColGlobal && vs->nc > 1) {
3968188e55aSJed Brown     PetscInt row;
3978188e55aSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
3988188e55aSJed Brown     ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr);
3998188e55aSJed Brown   } else {
400f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
401f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr);
402b6480e04SStefano Zampini     if (!vs->m[row][col]) {
403b6480e04SStefano Zampini       PetscInt lr,lc;
404b6480e04SStefano Zampini 
405b6480e04SStefano Zampini       ierr = MatCreate(PetscObjectComm((PetscObject)A),&vs->m[row][col]);CHKERRQ(ierr);
406b6480e04SStefano Zampini       ierr = ISGetLocalSize(vs->isglobal.row[row],&lr);CHKERRQ(ierr);
407b6480e04SStefano Zampini       ierr = ISGetLocalSize(vs->isglobal.col[col],&lc);CHKERRQ(ierr);
408b6480e04SStefano Zampini       ierr = MatSetSizes(vs->m[row][col],lr,lc,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr);
409b6480e04SStefano Zampini       ierr = MatSetUp(vs->m[row][col]);CHKERRQ(ierr);
410b6480e04SStefano Zampini       ierr = MatAssemblyBegin(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
411b6480e04SStefano Zampini       ierr = MatAssemblyEnd(vs->m[row][col],MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
412b6480e04SStefano Zampini     }
413f349c1fdSJed Brown     *B = vs->m[row][col];
4148188e55aSJed Brown   }
415f349c1fdSJed Brown   PetscFunctionReturn(0);
416f349c1fdSJed Brown }
417f349c1fdSJed Brown 
4187dae84e0SHong Zhang static PetscErrorCode MatCreateSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
419f349c1fdSJed Brown {
420f349c1fdSJed Brown   PetscErrorCode ierr;
421f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
422f349c1fdSJed Brown   Mat            sub;
423f349c1fdSJed Brown 
424f349c1fdSJed Brown   PetscFunctionBegin;
425f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr);
426f349c1fdSJed Brown   switch (reuse) {
427f349c1fdSJed Brown   case MAT_INITIAL_MATRIX:
4287874fa86SDave May     if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); }
429f349c1fdSJed Brown     *B = sub;
430f349c1fdSJed Brown     break;
431f349c1fdSJed Brown   case MAT_REUSE_MATRIX:
432ce94432eSBarry Smith     if (sub != *B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
433f349c1fdSJed Brown     break;
434f349c1fdSJed Brown   case MAT_IGNORE_MATRIX:       /* Nothing to do */
435f349c1fdSJed Brown     break;
436511c6705SHong Zhang   case MAT_INPLACE_MATRIX:       /* Nothing to do */
437511c6705SHong Zhang     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX is not supported yet");
438511c6705SHong Zhang     break;
439f349c1fdSJed Brown   }
440f349c1fdSJed Brown   PetscFunctionReturn(0);
441f349c1fdSJed Brown }
442f349c1fdSJed Brown 
443f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
444f349c1fdSJed Brown {
445f349c1fdSJed Brown   PetscErrorCode ierr;
446f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
447f349c1fdSJed Brown   Mat            sub;
448f349c1fdSJed Brown 
449f349c1fdSJed Brown   PetscFunctionBegin;
450f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
451f349c1fdSJed Brown   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
452f349c1fdSJed Brown   if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);}
453f349c1fdSJed Brown   *B = sub;
454d8588912SDave May   PetscFunctionReturn(0);
455d8588912SDave May }
456d8588912SDave May 
457207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
458d8588912SDave May {
459d8588912SDave May   PetscErrorCode ierr;
460f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
461f349c1fdSJed Brown   Mat            sub;
462d8588912SDave May 
463d8588912SDave May   PetscFunctionBegin;
464f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
465ce94432eSBarry Smith   if (*B != sub) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
466f349c1fdSJed Brown   if (sub) {
467ce94432eSBarry Smith     if (((PetscObject)sub)->refct <= 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
4686bf464f9SBarry Smith     ierr = MatDestroy(B);CHKERRQ(ierr);
469d8588912SDave May   }
470d8588912SDave May   PetscFunctionReturn(0);
471d8588912SDave May }
472d8588912SDave May 
4737874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
4747874fa86SDave May {
4757874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
4767874fa86SDave May   PetscInt       i;
4777874fa86SDave May   PetscErrorCode ierr;
4787874fa86SDave May 
4797874fa86SDave May   PetscFunctionBegin;
4807874fa86SDave May   for (i=0; i<bA->nr; i++) {
481429bac76SJed Brown     Vec bv;
482429bac76SJed Brown     ierr = VecGetSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
4837874fa86SDave May     if (bA->m[i][i]) {
484429bac76SJed Brown       ierr = MatGetDiagonal(bA->m[i][i],bv);CHKERRQ(ierr);
4857874fa86SDave May     } else {
4865159a857SMatthew G. Knepley       ierr = VecSet(bv,0.0);CHKERRQ(ierr);
4877874fa86SDave May     }
488429bac76SJed Brown     ierr = VecRestoreSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
4897874fa86SDave May   }
4907874fa86SDave May   PetscFunctionReturn(0);
4917874fa86SDave May }
4927874fa86SDave May 
4937874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
4947874fa86SDave May {
4957874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
496429bac76SJed Brown   Vec            bl,*br;
4977874fa86SDave May   PetscInt       i,j;
4987874fa86SDave May   PetscErrorCode ierr;
4997874fa86SDave May 
5007874fa86SDave May   PetscFunctionBegin;
5013f800ebeSJed Brown   ierr = PetscCalloc1(bA->nc,&br);CHKERRQ(ierr);
5022e6472ebSElliott Sales de Andrade   if (r) {
503429bac76SJed Brown     for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
5042e6472ebSElliott Sales de Andrade   }
5052e6472ebSElliott Sales de Andrade   bl = NULL;
5067874fa86SDave May   for (i=0; i<bA->nr; i++) {
5072e6472ebSElliott Sales de Andrade     if (l) {
508429bac76SJed Brown       ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
5092e6472ebSElliott Sales de Andrade     }
5107874fa86SDave May     for (j=0; j<bA->nc; j++) {
5117874fa86SDave May       if (bA->m[i][j]) {
512429bac76SJed Brown         ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr);
5137874fa86SDave May       }
5147874fa86SDave May     }
5152e6472ebSElliott Sales de Andrade     if (l) {
516a061e289SJed Brown       ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
5177874fa86SDave May     }
5182e6472ebSElliott Sales de Andrade   }
5192e6472ebSElliott Sales de Andrade   if (r) {
520429bac76SJed Brown     for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
5212e6472ebSElliott Sales de Andrade   }
522429bac76SJed Brown   ierr = PetscFree(br);CHKERRQ(ierr);
5237874fa86SDave May   PetscFunctionReturn(0);
5247874fa86SDave May }
5257874fa86SDave May 
526a061e289SJed Brown static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a)
527a061e289SJed Brown {
528a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
529a061e289SJed Brown   PetscInt       i,j;
530a061e289SJed Brown   PetscErrorCode ierr;
531a061e289SJed Brown 
532a061e289SJed Brown   PetscFunctionBegin;
533a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
534a061e289SJed Brown     for (j=0; j<bA->nc; j++) {
535a061e289SJed Brown       if (bA->m[i][j]) {
536a061e289SJed Brown         ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr);
537a061e289SJed Brown       }
538a061e289SJed Brown     }
539a061e289SJed Brown   }
540a061e289SJed Brown   PetscFunctionReturn(0);
541a061e289SJed Brown }
542a061e289SJed Brown 
543a061e289SJed Brown static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a)
544a061e289SJed Brown {
545a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
546a061e289SJed Brown   PetscInt       i;
547a061e289SJed Brown   PetscErrorCode ierr;
548a061e289SJed Brown 
549a061e289SJed Brown   PetscFunctionBegin;
550a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
551ce94432eSBarry 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);
552a061e289SJed Brown     ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr);
553a061e289SJed Brown   }
554a061e289SJed Brown   PetscFunctionReturn(0);
555a061e289SJed Brown }
556a061e289SJed Brown 
55713135bc6SAlex Fikl static PetscErrorCode MatDiagonalSet_Nest(Mat A,Vec D,InsertMode is)
55813135bc6SAlex Fikl {
55913135bc6SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data;
56013135bc6SAlex Fikl   PetscInt       i;
56113135bc6SAlex Fikl   PetscErrorCode ierr;
56213135bc6SAlex Fikl 
56313135bc6SAlex Fikl   PetscFunctionBegin;
56413135bc6SAlex Fikl   for (i=0; i<bA->nr; i++) {
56513135bc6SAlex Fikl     Vec bv;
56613135bc6SAlex Fikl     ierr = VecGetSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
56713135bc6SAlex Fikl     if (bA->m[i][i]) {
56813135bc6SAlex Fikl       ierr = MatDiagonalSet(bA->m[i][i],bv,is);CHKERRQ(ierr);
56913135bc6SAlex Fikl     }
57013135bc6SAlex Fikl     ierr = VecRestoreSubVector(D,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
57113135bc6SAlex Fikl   }
57213135bc6SAlex Fikl   PetscFunctionReturn(0);
57313135bc6SAlex Fikl }
57413135bc6SAlex Fikl 
575f8170845SAlex Fikl static PetscErrorCode MatSetRandom_Nest(Mat A,PetscRandom rctx)
576f8170845SAlex Fikl {
577f8170845SAlex Fikl   Mat_Nest       *bA = (Mat_Nest*)A->data;
578f8170845SAlex Fikl   PetscInt       i,j;
579f8170845SAlex Fikl   PetscErrorCode ierr;
580f8170845SAlex Fikl 
581f8170845SAlex Fikl   PetscFunctionBegin;
582f8170845SAlex Fikl   for (i=0; i<bA->nr; i++) {
583f8170845SAlex Fikl     for (j=0; j<bA->nc; j++) {
584f8170845SAlex Fikl       if (bA->m[i][j]) {
585f8170845SAlex Fikl         ierr = MatSetRandom(bA->m[i][j],rctx);CHKERRQ(ierr);
586f8170845SAlex Fikl       }
587f8170845SAlex Fikl     }
588f8170845SAlex Fikl   }
589f8170845SAlex Fikl   PetscFunctionReturn(0);
590f8170845SAlex Fikl }
591f8170845SAlex Fikl 
5922a7a6963SBarry Smith static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left)
593d8588912SDave May {
594d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
595d8588912SDave May   Vec            *L,*R;
596d8588912SDave May   MPI_Comm       comm;
597d8588912SDave May   PetscInt       i,j;
598d8588912SDave May   PetscErrorCode ierr;
599d8588912SDave May 
600d8588912SDave May   PetscFunctionBegin;
601ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr);
602d8588912SDave May   if (right) {
603d8588912SDave May     /* allocate R */
604854ce69bSBarry Smith     ierr = PetscMalloc1(bA->nc, &R);CHKERRQ(ierr);
605d8588912SDave May     /* Create the right vectors */
606d8588912SDave May     for (j=0; j<bA->nc; j++) {
607d8588912SDave May       for (i=0; i<bA->nr; i++) {
608d8588912SDave May         if (bA->m[i][j]) {
6092a7a6963SBarry Smith           ierr = MatCreateVecs(bA->m[i][j],&R[j],NULL);CHKERRQ(ierr);
610d8588912SDave May           break;
611d8588912SDave May         }
612d8588912SDave May       }
6136c4ed002SBarry Smith       if (i==bA->nr) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
614d8588912SDave May     }
615f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
616d8588912SDave May     /* hand back control to the nest vector */
617d8588912SDave May     for (j=0; j<bA->nc; j++) {
6186bf464f9SBarry Smith       ierr = VecDestroy(&R[j]);CHKERRQ(ierr);
619d8588912SDave May     }
620d8588912SDave May     ierr = PetscFree(R);CHKERRQ(ierr);
621d8588912SDave May   }
622d8588912SDave May 
623d8588912SDave May   if (left) {
624d8588912SDave May     /* allocate L */
625854ce69bSBarry Smith     ierr = PetscMalloc1(bA->nr, &L);CHKERRQ(ierr);
626d8588912SDave May     /* Create the left vectors */
627d8588912SDave May     for (i=0; i<bA->nr; i++) {
628d8588912SDave May       for (j=0; j<bA->nc; j++) {
629d8588912SDave May         if (bA->m[i][j]) {
6302a7a6963SBarry Smith           ierr = MatCreateVecs(bA->m[i][j],NULL,&L[i]);CHKERRQ(ierr);
631d8588912SDave May           break;
632d8588912SDave May         }
633d8588912SDave May       }
6346c4ed002SBarry Smith       if (j==bA->nc) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
635d8588912SDave May     }
636d8588912SDave May 
637f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
638d8588912SDave May     for (i=0; i<bA->nr; i++) {
6396bf464f9SBarry Smith       ierr = VecDestroy(&L[i]);CHKERRQ(ierr);
640d8588912SDave May     }
641d8588912SDave May 
642d8588912SDave May     ierr = PetscFree(L);CHKERRQ(ierr);
643d8588912SDave May   }
644d8588912SDave May   PetscFunctionReturn(0);
645d8588912SDave May }
646d8588912SDave May 
647207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
648d8588912SDave May {
649d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
650d8588912SDave May   PetscBool      isascii;
651d8588912SDave May   PetscInt       i,j;
652d8588912SDave May   PetscErrorCode ierr;
653d8588912SDave May 
654d8588912SDave May   PetscFunctionBegin;
655251f4c67SDmitry Karpeev   ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
656d8588912SDave May   if (isascii) {
657d8588912SDave May 
658d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"Matrix object: \n");CHKERRQ(ierr);
659d86155a6SBarry Smith     ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
660d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer, "type=nest, rows=%D, cols=%D \n",bA->nr,bA->nc);CHKERRQ(ierr);
661d8588912SDave May 
662d86155a6SBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"MatNest structure: \n");CHKERRQ(ierr);
663d8588912SDave May     for (i=0; i<bA->nr; i++) {
664d8588912SDave May       for (j=0; j<bA->nc; j++) {
66519fd82e9SBarry Smith         MatType   type;
666270f95d7SJed Brown         char      name[256] = "",prefix[256] = "";
667d8588912SDave May         PetscInt  NR,NC;
668d8588912SDave May         PetscBool isNest = PETSC_FALSE;
669d8588912SDave May 
670d8588912SDave May         if (!bA->m[i][j]) {
671d86155a6SBarry Smith           CHKERRQ(ierr);PetscViewerASCIIPrintf(viewer, "(%D,%D) : NULL \n",i,j);CHKERRQ(ierr);
672d8588912SDave May           continue;
673d8588912SDave May         }
674d8588912SDave May         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
675d8588912SDave May         ierr = MatGetType(bA->m[i][j], &type);CHKERRQ(ierr);
6768caf3d72SBarry Smith         if (((PetscObject)bA->m[i][j])->name) {ierr = PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);CHKERRQ(ierr);}
6778caf3d72SBarry Smith         if (((PetscObject)bA->m[i][j])->prefix) {ierr = PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);CHKERRQ(ierr);}
678251f4c67SDmitry Karpeev         ierr = PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
679d8588912SDave May 
680270f95d7SJed Brown         ierr = PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);CHKERRQ(ierr);
681d8588912SDave May 
682d8588912SDave May         if (isNest) {
683270f95d7SJed Brown           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);  /* push1 */
684d8588912SDave May           ierr = MatView(bA->m[i][j],viewer);CHKERRQ(ierr);
685270f95d7SJed Brown           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop1 */
686d8588912SDave May         }
687d8588912SDave May       }
688d8588912SDave May     }
689d86155a6SBarry Smith     ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);    /* pop0 */
690d8588912SDave May   }
691d8588912SDave May   PetscFunctionReturn(0);
692d8588912SDave May }
693d8588912SDave May 
694207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A)
695d8588912SDave May {
696d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
697d8588912SDave May   PetscInt       i,j;
698d8588912SDave May   PetscErrorCode ierr;
699d8588912SDave May 
700d8588912SDave May   PetscFunctionBegin;
701d8588912SDave May   for (i=0; i<bA->nr; i++) {
702d8588912SDave May     for (j=0; j<bA->nc; j++) {
703d8588912SDave May       if (!bA->m[i][j]) continue;
704d8588912SDave May       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
705d8588912SDave May     }
706d8588912SDave May   }
707d8588912SDave May   PetscFunctionReturn(0);
708d8588912SDave May }
709d8588912SDave May 
710c222c20dSDavid Ham static PetscErrorCode MatCopy_Nest(Mat A,Mat B,MatStructure str)
711c222c20dSDavid Ham {
712c222c20dSDavid Ham   Mat_Nest       *bA = (Mat_Nest*)A->data,*bB = (Mat_Nest*)B->data;
713c222c20dSDavid Ham   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
714c222c20dSDavid Ham   PetscErrorCode ierr;
715c222c20dSDavid Ham 
716c222c20dSDavid Ham   PetscFunctionBegin;
717c222c20dSDavid 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);
718c222c20dSDavid Ham   for (i=0; i<nr; i++) {
719c222c20dSDavid Ham     for (j=0; j<nc; j++) {
72046a2b97cSJed Brown       if (bA->m[i][j] && bB->m[i][j]) {
721c222c20dSDavid Ham         ierr = MatCopy(bA->m[i][j],bB->m[i][j],str);CHKERRQ(ierr);
72246a2b97cSJed 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);
723c222c20dSDavid Ham     }
724c222c20dSDavid Ham   }
725cdc753b6SBarry Smith   ierr = PetscObjectStateIncrease((PetscObject)B);CHKERRQ(ierr);
726c222c20dSDavid Ham   PetscFunctionReturn(0);
727c222c20dSDavid Ham }
728c222c20dSDavid Ham 
729207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
730d8588912SDave May {
731d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
732841e96a3SJed Brown   Mat            *b;
733841e96a3SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
734d8588912SDave May   PetscErrorCode ierr;
735d8588912SDave May 
736d8588912SDave May   PetscFunctionBegin;
737785e854fSJed Brown   ierr = PetscMalloc1(nr*nc,&b);CHKERRQ(ierr);
738841e96a3SJed Brown   for (i=0; i<nr; i++) {
739841e96a3SJed Brown     for (j=0; j<nc; j++) {
740841e96a3SJed Brown       if (bA->m[i][j]) {
741841e96a3SJed Brown         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
742841e96a3SJed Brown       } else {
7430298fd71SBarry Smith         b[i*nc+j] = NULL;
744d8588912SDave May       }
745d8588912SDave May     }
746d8588912SDave May   }
747ce94432eSBarry Smith   ierr = MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
748841e96a3SJed Brown   /* Give the new MatNest exclusive ownership */
749841e96a3SJed Brown   for (i=0; i<nr*nc; i++) {
7506bf464f9SBarry Smith     ierr = MatDestroy(&b[i]);CHKERRQ(ierr);
751d8588912SDave May   }
752d8588912SDave May   ierr = PetscFree(b);CHKERRQ(ierr);
753d8588912SDave May 
754841e96a3SJed Brown   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
755841e96a3SJed Brown   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
756d8588912SDave May   PetscFunctionReturn(0);
757d8588912SDave May }
758d8588912SDave May 
759d8588912SDave May /* nest api */
760d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
761d8588912SDave May {
762d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
7635fd66863SKarl Rupp 
764d8588912SDave May   PetscFunctionBegin;
765ce94432eSBarry Smith   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
766ce94432eSBarry Smith   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
767d8588912SDave May   *mat = bA->m[idxm][jdxm];
768d8588912SDave May   PetscFunctionReturn(0);
769d8588912SDave May }
770d8588912SDave May 
7719ba0d327SJed Brown /*@
772d8588912SDave May  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
773d8588912SDave May 
774d8588912SDave May  Not collective
775d8588912SDave May 
776d8588912SDave May  Input Parameters:
777629881c0SJed Brown +   A  - nest matrix
778d8588912SDave May .   idxm - index of the matrix within the nest matrix
779629881c0SJed Brown -   jdxm - index of the matrix within the nest matrix
780d8588912SDave May 
781d8588912SDave May  Output Parameter:
782d8588912SDave May .   sub - matrix at index idxm,jdxm within the nest matrix
783d8588912SDave May 
784d8588912SDave May  Level: developer
785d8588912SDave May 
786d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMats()
787d8588912SDave May @*/
7887087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
789d8588912SDave May {
790699a902aSJed Brown   PetscErrorCode ierr;
791d8588912SDave May 
792d8588912SDave May   PetscFunctionBegin;
793699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
794d8588912SDave May   PetscFunctionReturn(0);
795d8588912SDave May }
796d8588912SDave May 
7970782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat)
7980782ca92SJed Brown {
7990782ca92SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
8000782ca92SJed Brown   PetscInt       m,n,M,N,mi,ni,Mi,Ni;
8010782ca92SJed Brown   PetscErrorCode ierr;
8020782ca92SJed Brown 
8030782ca92SJed Brown   PetscFunctionBegin;
804ce94432eSBarry Smith   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
805ce94432eSBarry Smith   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
8060782ca92SJed Brown   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
8070782ca92SJed Brown   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
8080782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr);
8090782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr);
8100782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr);
8110782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr);
812ce94432eSBarry 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);
813ce94432eSBarry 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);
81426fbe8dcSKarl Rupp 
8150782ca92SJed Brown   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
8160782ca92SJed Brown   ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr);
8170782ca92SJed Brown   bA->m[idxm][jdxm] = mat;
8180782ca92SJed Brown   PetscFunctionReturn(0);
8190782ca92SJed Brown }
8200782ca92SJed Brown 
8219ba0d327SJed Brown /*@
8220782ca92SJed Brown  MatNestSetSubMat - Set a single submatrix in the nest matrix.
8230782ca92SJed Brown 
8240782ca92SJed Brown  Logically collective on the submatrix communicator
8250782ca92SJed Brown 
8260782ca92SJed Brown  Input Parameters:
8270782ca92SJed Brown +   A  - nest matrix
8280782ca92SJed Brown .   idxm - index of the matrix within the nest matrix
8290782ca92SJed Brown .   jdxm - index of the matrix within the nest matrix
8300782ca92SJed Brown -   sub - matrix at index idxm,jdxm within the nest matrix
8310782ca92SJed Brown 
8320782ca92SJed Brown  Notes:
8330782ca92SJed Brown  The new submatrix must have the same size and communicator as that block of the nest.
8340782ca92SJed Brown 
8350782ca92SJed Brown  This increments the reference count of the submatrix.
8360782ca92SJed Brown 
8370782ca92SJed Brown  Level: developer
8380782ca92SJed Brown 
839d5dfb694SBarry Smith .seealso: MatNestSetSubMats(), MatNestGetSubMats()
8400782ca92SJed Brown @*/
8410782ca92SJed Brown PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
8420782ca92SJed Brown {
8430782ca92SJed Brown   PetscErrorCode ierr;
8440782ca92SJed Brown 
8450782ca92SJed Brown   PetscFunctionBegin;
8460782ca92SJed Brown   ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr);
8470782ca92SJed Brown   PetscFunctionReturn(0);
8480782ca92SJed Brown }
8490782ca92SJed Brown 
850d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
851d8588912SDave May {
852d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
8535fd66863SKarl Rupp 
854d8588912SDave May   PetscFunctionBegin;
85526fbe8dcSKarl Rupp   if (M)   *M   = bA->nr;
85626fbe8dcSKarl Rupp   if (N)   *N   = bA->nc;
85726fbe8dcSKarl Rupp   if (mat) *mat = bA->m;
858d8588912SDave May   PetscFunctionReturn(0);
859d8588912SDave May }
860d8588912SDave May 
861d8588912SDave May /*@C
862d8588912SDave May  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
863d8588912SDave May 
864d8588912SDave May  Not collective
865d8588912SDave May 
866d8588912SDave May  Input Parameters:
867629881c0SJed Brown .   A  - nest matrix
868d8588912SDave May 
869d8588912SDave May  Output Parameter:
870629881c0SJed Brown +   M - number of rows in the nest matrix
871d8588912SDave May .   N - number of cols in the nest matrix
872629881c0SJed Brown -   mat - 2d array of matrices
873d8588912SDave May 
874d8588912SDave May  Notes:
875d8588912SDave May 
876d8588912SDave May  The user should not free the array mat.
877d8588912SDave May 
878351962e3SVincent Le Chenadec  In Fortran, this routine has a calling sequence
879351962e3SVincent Le Chenadec $   call MatNestGetSubMats(A, M, N, mat, ierr)
880351962e3SVincent Le Chenadec  where the space allocated for the optional argument mat is assumed large enough (if provided).
881351962e3SVincent Le Chenadec 
882d8588912SDave May  Level: developer
883d8588912SDave May 
884d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMat()
885d8588912SDave May @*/
8867087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
887d8588912SDave May {
888699a902aSJed Brown   PetscErrorCode ierr;
889d8588912SDave May 
890d8588912SDave May   PetscFunctionBegin;
891699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
892d8588912SDave May   PetscFunctionReturn(0);
893d8588912SDave May }
894d8588912SDave May 
8957087cfbeSBarry Smith PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
896d8588912SDave May {
897d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
898d8588912SDave May 
899d8588912SDave May   PetscFunctionBegin;
90026fbe8dcSKarl Rupp   if (M) *M = bA->nr;
90126fbe8dcSKarl Rupp   if (N) *N = bA->nc;
902d8588912SDave May   PetscFunctionReturn(0);
903d8588912SDave May }
904d8588912SDave May 
9059ba0d327SJed Brown /*@
906d8588912SDave May  MatNestGetSize - Returns the size of the nest matrix.
907d8588912SDave May 
908d8588912SDave May  Not collective
909d8588912SDave May 
910d8588912SDave May  Input Parameters:
911d8588912SDave May .   A  - nest matrix
912d8588912SDave May 
913d8588912SDave May  Output Parameter:
914629881c0SJed Brown +   M - number of rows in the nested mat
915629881c0SJed Brown -   N - number of cols in the nested mat
916d8588912SDave May 
917d8588912SDave May  Notes:
918d8588912SDave May 
919d8588912SDave May  Level: developer
920d8588912SDave May 
921d8588912SDave May .seealso: MatNestGetSubMat(), MatNestGetSubMats()
922d8588912SDave May @*/
9237087cfbeSBarry Smith PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
924d8588912SDave May {
925699a902aSJed Brown   PetscErrorCode ierr;
926d8588912SDave May 
927d8588912SDave May   PetscFunctionBegin;
928699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
929d8588912SDave May   PetscFunctionReturn(0);
930d8588912SDave May }
931d8588912SDave May 
932f7a08781SBarry Smith static PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[])
933900e7ff2SJed Brown {
934900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest*)A->data;
935900e7ff2SJed Brown   PetscInt i;
936900e7ff2SJed Brown 
937900e7ff2SJed Brown   PetscFunctionBegin;
938900e7ff2SJed Brown   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i];
939900e7ff2SJed Brown   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i];
940900e7ff2SJed Brown   PetscFunctionReturn(0);
941900e7ff2SJed Brown }
942900e7ff2SJed Brown 
9433a4d7b9aSSatish Balay /*@C
944900e7ff2SJed Brown  MatNestGetISs - Returns the index sets partitioning the row and column spaces
945900e7ff2SJed Brown 
946900e7ff2SJed Brown  Not collective
947900e7ff2SJed Brown 
948900e7ff2SJed Brown  Input Parameters:
949900e7ff2SJed Brown .   A  - nest matrix
950900e7ff2SJed Brown 
951900e7ff2SJed Brown  Output Parameter:
952900e7ff2SJed Brown +   rows - array of row index sets
953900e7ff2SJed Brown -   cols - array of column index sets
954900e7ff2SJed Brown 
955900e7ff2SJed Brown  Level: advanced
956900e7ff2SJed Brown 
957900e7ff2SJed Brown  Notes:
958900e7ff2SJed Brown  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
959900e7ff2SJed Brown 
960900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs()
961900e7ff2SJed Brown @*/
962900e7ff2SJed Brown PetscErrorCode  MatNestGetISs(Mat A,IS rows[],IS cols[])
963900e7ff2SJed Brown {
964900e7ff2SJed Brown   PetscErrorCode ierr;
965900e7ff2SJed Brown 
966900e7ff2SJed Brown   PetscFunctionBegin;
967900e7ff2SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
968900e7ff2SJed Brown   ierr = PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
969900e7ff2SJed Brown   PetscFunctionReturn(0);
970900e7ff2SJed Brown }
971900e7ff2SJed Brown 
972f7a08781SBarry Smith static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[])
973900e7ff2SJed Brown {
974900e7ff2SJed Brown   Mat_Nest *vs = (Mat_Nest*)A->data;
975900e7ff2SJed Brown   PetscInt i;
976900e7ff2SJed Brown 
977900e7ff2SJed Brown   PetscFunctionBegin;
978900e7ff2SJed Brown   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i];
979900e7ff2SJed Brown   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i];
980900e7ff2SJed Brown   PetscFunctionReturn(0);
981900e7ff2SJed Brown }
982900e7ff2SJed Brown 
983900e7ff2SJed Brown /*@C
984900e7ff2SJed Brown  MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces
985900e7ff2SJed Brown 
986900e7ff2SJed Brown  Not collective
987900e7ff2SJed Brown 
988900e7ff2SJed Brown  Input Parameters:
989900e7ff2SJed Brown .   A  - nest matrix
990900e7ff2SJed Brown 
991900e7ff2SJed Brown  Output Parameter:
9920298fd71SBarry Smith +   rows - array of row index sets (or NULL to ignore)
9930298fd71SBarry Smith -   cols - array of column index sets (or NULL to ignore)
994900e7ff2SJed Brown 
995900e7ff2SJed Brown  Level: advanced
996900e7ff2SJed Brown 
997900e7ff2SJed Brown  Notes:
998900e7ff2SJed Brown  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.
999900e7ff2SJed Brown 
1000900e7ff2SJed Brown .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs()
1001900e7ff2SJed Brown @*/
1002900e7ff2SJed Brown PetscErrorCode  MatNestGetLocalISs(Mat A,IS rows[],IS cols[])
1003900e7ff2SJed Brown {
1004900e7ff2SJed Brown   PetscErrorCode ierr;
1005900e7ff2SJed Brown 
1006900e7ff2SJed Brown   PetscFunctionBegin;
1007900e7ff2SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1008900e7ff2SJed Brown   ierr = PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));CHKERRQ(ierr);
1009900e7ff2SJed Brown   PetscFunctionReturn(0);
1010900e7ff2SJed Brown }
1011900e7ff2SJed Brown 
101219fd82e9SBarry Smith PetscErrorCode  MatNestSetVecType_Nest(Mat A,VecType vtype)
1013207556f9SJed Brown {
1014207556f9SJed Brown   PetscErrorCode ierr;
1015207556f9SJed Brown   PetscBool      flg;
1016207556f9SJed Brown 
1017207556f9SJed Brown   PetscFunctionBegin;
1018207556f9SJed Brown   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
1019207556f9SJed Brown   /* In reality, this only distinguishes VECNEST and "other" */
10202a7a6963SBarry Smith   if (flg) A->ops->getvecs = MatCreateVecs_Nest;
102112b53f24SSatish Balay   else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0;
1022207556f9SJed Brown   PetscFunctionReturn(0);
1023207556f9SJed Brown }
1024207556f9SJed Brown 
1025207556f9SJed Brown /*@C
10262a7a6963SBarry Smith  MatNestSetVecType - Sets the type of Vec returned by MatCreateVecs()
1027207556f9SJed Brown 
1028207556f9SJed Brown  Not collective
1029207556f9SJed Brown 
1030207556f9SJed Brown  Input Parameters:
1031207556f9SJed Brown +  A  - nest matrix
1032207556f9SJed Brown -  vtype - type to use for creating vectors
1033207556f9SJed Brown 
1034207556f9SJed Brown  Notes:
1035207556f9SJed Brown 
1036207556f9SJed Brown  Level: developer
1037207556f9SJed Brown 
10382a7a6963SBarry Smith .seealso: MatCreateVecs()
1039207556f9SJed Brown @*/
104019fd82e9SBarry Smith PetscErrorCode  MatNestSetVecType(Mat A,VecType vtype)
1041207556f9SJed Brown {
1042207556f9SJed Brown   PetscErrorCode ierr;
1043207556f9SJed Brown 
1044207556f9SJed Brown   PetscFunctionBegin;
104519fd82e9SBarry Smith   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));CHKERRQ(ierr);
1046207556f9SJed Brown   PetscFunctionReturn(0);
1047207556f9SJed Brown }
1048207556f9SJed Brown 
1049c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1050d8588912SDave May {
1051c8883902SJed Brown   Mat_Nest       *s = (Mat_Nest*)A->data;
1052c8883902SJed Brown   PetscInt       i,j,m,n,M,N;
1053d8588912SDave May   PetscErrorCode ierr;
1054d8588912SDave May 
1055d8588912SDave May   PetscFunctionBegin;
1056c8883902SJed Brown   s->nr = nr;
1057c8883902SJed Brown   s->nc = nc;
1058d8588912SDave May 
1059c8883902SJed Brown   /* Create space for submatrices */
1060854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&s->m);CHKERRQ(ierr);
1061c8883902SJed Brown   for (i=0; i<nr; i++) {
1062854ce69bSBarry Smith     ierr = PetscMalloc1(nc,&s->m[i]);CHKERRQ(ierr);
1063d8588912SDave May   }
1064c8883902SJed Brown   for (i=0; i<nr; i++) {
1065c8883902SJed Brown     for (j=0; j<nc; j++) {
1066c8883902SJed Brown       s->m[i][j] = a[i*nc+j];
1067c8883902SJed Brown       if (a[i*nc+j]) {
1068c8883902SJed Brown         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
1069d8588912SDave May       }
1070d8588912SDave May     }
1071d8588912SDave May   }
1072d8588912SDave May 
10738188e55aSJed Brown   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
1074d8588912SDave May 
1075854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&s->row_len);CHKERRQ(ierr);
1076854ce69bSBarry Smith   ierr = PetscMalloc1(nc,&s->col_len);CHKERRQ(ierr);
1077c8883902SJed Brown   for (i=0; i<nr; i++) s->row_len[i]=-1;
1078c8883902SJed Brown   for (j=0; j<nc; j++) s->col_len[j]=-1;
1079d8588912SDave May 
10808188e55aSJed Brown   ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr);
1081d8588912SDave May 
1082c8883902SJed Brown   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
1083c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
1084c8883902SJed Brown   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
1085c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
1086c8883902SJed Brown 
1087c8883902SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1088c8883902SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1089c8883902SJed Brown 
10901795a4d1SJed Brown   ierr = PetscCalloc2(nr,&s->left,nc,&s->right);CHKERRQ(ierr);
1091d8588912SDave May   PetscFunctionReturn(0);
1092d8588912SDave May }
1093d8588912SDave May 
1094c8883902SJed Brown /*@
1095c8883902SJed Brown    MatNestSetSubMats - Sets the nested submatrices
1096c8883902SJed Brown 
1097c8883902SJed Brown    Collective on Mat
1098c8883902SJed Brown 
1099c8883902SJed Brown    Input Parameter:
1100c8883902SJed Brown +  N - nested matrix
1101c8883902SJed Brown .  nr - number of nested row blocks
11020298fd71SBarry Smith .  is_row - index sets for each nested row block, or NULL to make contiguous
1103c8883902SJed Brown .  nc - number of nested column blocks
11040298fd71SBarry Smith .  is_col - index sets for each nested column block, or NULL to make contiguous
11050298fd71SBarry Smith -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1106c8883902SJed Brown 
1107c8883902SJed Brown    Level: advanced
1108c8883902SJed Brown 
1109c8883902SJed Brown .seealso: MatCreateNest(), MATNEST
1110c8883902SJed Brown @*/
1111c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1112c8883902SJed Brown {
1113c8883902SJed Brown   PetscErrorCode ierr;
1114eb6c2100SSatish Balay   PetscInt       i,nr_nc;
1115c8883902SJed Brown 
1116c8883902SJed Brown   PetscFunctionBegin;
1117c8883902SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
1118ce94432eSBarry Smith   if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
1119c8883902SJed Brown   if (nr && is_row) {
1120c8883902SJed Brown     PetscValidPointer(is_row,3);
1121c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
1122c8883902SJed Brown   }
1123ce94432eSBarry Smith   if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
11241664e352SJed Brown   if (nc && is_col) {
1125c8883902SJed Brown     PetscValidPointer(is_col,5);
11269b30a8f6SBarry Smith     for (i=0; i<nc; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
1127c8883902SJed Brown   }
1128eb6c2100SSatish Balay   nr_nc=nr*nc;
1129eb6c2100SSatish Balay   if (nr_nc) PetscValidPointer(a,6);
1130c8883902SJed 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);
1131c8883902SJed Brown   PetscFunctionReturn(0);
1132c8883902SJed Brown }
1133d8588912SDave May 
113445b6f7e9SBarry Smith static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog)
113577019fcaSJed Brown {
113677019fcaSJed Brown   PetscErrorCode ierr;
113777019fcaSJed Brown   PetscBool      flg;
113877019fcaSJed Brown   PetscInt       i,j,m,mi,*ix;
113977019fcaSJed Brown 
114077019fcaSJed Brown   PetscFunctionBegin;
114177019fcaSJed Brown   for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) {
114277019fcaSJed Brown     if (islocal[i]) {
114377019fcaSJed Brown       ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
114477019fcaSJed Brown       flg  = PETSC_TRUE;      /* We found a non-trivial entry */
114577019fcaSJed Brown     } else {
114677019fcaSJed Brown       ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
114777019fcaSJed Brown     }
114877019fcaSJed Brown     m += mi;
114977019fcaSJed Brown   }
115077019fcaSJed Brown   if (flg) {
1151785e854fSJed Brown     ierr = PetscMalloc1(m,&ix);CHKERRQ(ierr);
1152165cd838SBarry Smith     for (i=0,m=0; i<n; i++) {
11530298fd71SBarry Smith       ISLocalToGlobalMapping smap = NULL;
115477019fcaSJed Brown       VecScatter             scat;
115577019fcaSJed Brown       IS                     isreq;
115677019fcaSJed Brown       Vec                    lvec,gvec;
1157*e108cb99SStefano Zampini       Mat                    sub = NULL;
11583361c9a7SJed Brown       union {char padding[sizeof(PetscScalar)]; PetscInt integer;} *x;
115977019fcaSJed Brown 
1160ce94432eSBarry Smith       if (sizeof(*x) != sizeof(PetscScalar)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support when scalars smaller than integers");
1161165cd838SBarry Smith       if (!colflg) {
116277019fcaSJed Brown         ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
116377019fcaSJed Brown       } else {
116477019fcaSJed Brown         ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
116577019fcaSJed Brown       }
1166191fd14bSBarry Smith       if (sub) {
1167191fd14bSBarry Smith         if (!colflg) {
1168191fd14bSBarry Smith           ierr = MatGetLocalToGlobalMapping(sub,&smap,NULL);CHKERRQ(ierr);
1169191fd14bSBarry Smith         } else {
1170191fd14bSBarry Smith           ierr = MatGetLocalToGlobalMapping(sub,NULL,&smap);CHKERRQ(ierr);
1171191fd14bSBarry Smith         }
1172191fd14bSBarry Smith       }
117377019fcaSJed Brown       if (islocal[i]) {
117477019fcaSJed Brown         ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
117577019fcaSJed Brown       } else {
117677019fcaSJed Brown         ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
117777019fcaSJed Brown       }
117877019fcaSJed Brown       for (j=0; j<mi; j++) ix[m+j] = j;
117977019fcaSJed Brown       if (smap) {ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr);}
1180165cd838SBarry Smith 
118177019fcaSJed Brown       /*
118277019fcaSJed Brown         Now we need to extract the monolithic global indices that correspond to the given split global indices.
118377019fcaSJed 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.
118477019fcaSJed Brown         The approach here is ugly because it uses VecScatter to move indices.
118577019fcaSJed Brown        */
118677019fcaSJed Brown       ierr = VecCreateSeq(PETSC_COMM_SELF,mi,&lvec);CHKERRQ(ierr);
118777019fcaSJed Brown       ierr = VecCreateMPI(((PetscObject)isglobal[i])->comm,mi,PETSC_DECIDE,&gvec);CHKERRQ(ierr);
118877019fcaSJed Brown       ierr = ISCreateGeneral(((PetscObject)isglobal[i])->comm,mi,ix+m,PETSC_COPY_VALUES,&isreq);CHKERRQ(ierr);
11890298fd71SBarry Smith       ierr = VecScatterCreate(gvec,isreq,lvec,NULL,&scat);CHKERRQ(ierr);
119077019fcaSJed Brown       ierr = VecGetArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
119177019fcaSJed Brown       for (j=0; j<mi; j++) x[j].integer = ix[m+j];
119277019fcaSJed Brown       ierr = VecRestoreArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
119377019fcaSJed Brown       ierr = VecScatterBegin(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
119477019fcaSJed Brown       ierr = VecScatterEnd(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
119577019fcaSJed Brown       ierr = VecGetArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
119677019fcaSJed Brown       for (j=0; j<mi; j++) ix[m+j] = x[j].integer;
119777019fcaSJed Brown       ierr = VecRestoreArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
119877019fcaSJed Brown       ierr = VecDestroy(&lvec);CHKERRQ(ierr);
119977019fcaSJed Brown       ierr = VecDestroy(&gvec);CHKERRQ(ierr);
120077019fcaSJed Brown       ierr = ISDestroy(&isreq);CHKERRQ(ierr);
120177019fcaSJed Brown       ierr = VecScatterDestroy(&scat);CHKERRQ(ierr);
120277019fcaSJed Brown       m   += mi;
120377019fcaSJed Brown     }
1204f0413b6fSBarry Smith     ierr   = ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),1,m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
120577019fcaSJed Brown   } else {
12060298fd71SBarry Smith     *ltog  = NULL;
120777019fcaSJed Brown   }
120877019fcaSJed Brown   PetscFunctionReturn(0);
120977019fcaSJed Brown }
121077019fcaSJed Brown 
121177019fcaSJed Brown 
1212d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1213d8588912SDave May /*
1214d8588912SDave May   nprocessors = NP
1215d8588912SDave May   Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1))
1216d8588912SDave May        proc 0: => (g_0,h_0,)
1217d8588912SDave May        proc 1: => (g_1,h_1,)
1218d8588912SDave May        ...
1219d8588912SDave May        proc nprocs-1: => (g_NP-1,h_NP-1,)
1220d8588912SDave May 
1221d8588912SDave May             proc 0:                      proc 1:                    proc nprocs-1:
1222d8588912SDave May     is[0] = (0,1,2,...,nlocal(g_0)-1)  (0,1,...,nlocal(g_1)-1)  (0,1,...,nlocal(g_NP-1))
1223d8588912SDave May 
1224d8588912SDave May             proc 0:
1225d8588912SDave May     is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1)
1226d8588912SDave May             proc 1:
1227d8588912SDave May     is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1)
1228d8588912SDave May 
1229d8588912SDave May             proc NP-1:
1230d8588912SDave May     is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1)
1231d8588912SDave May */
1232841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1233d8588912SDave May {
1234e2d7f03fSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
12358188e55aSJed Brown   PetscInt       i,j,offset,n,nsum,bs;
1236d8588912SDave May   PetscErrorCode ierr;
12370298fd71SBarry Smith   Mat            sub = NULL;
1238d8588912SDave May 
1239d8588912SDave May   PetscFunctionBegin;
1240854ce69bSBarry Smith   ierr = PetscMalloc1(nr,&vs->isglobal.row);CHKERRQ(ierr);
1241854ce69bSBarry Smith   ierr = PetscMalloc1(nc,&vs->isglobal.col);CHKERRQ(ierr);
1242d8588912SDave May   if (is_row) { /* valid IS is passed in */
1243d8588912SDave May     /* refs on is[] are incremeneted */
1244e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1245d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
124626fbe8dcSKarl Rupp 
1247e2d7f03fSJed Brown       vs->isglobal.row[i] = is_row[i];
1248d8588912SDave May     }
12492ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
12508188e55aSJed Brown     nsum = 0;
12518188e55aSJed Brown     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
12528188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1253ce94432eSBarry Smith       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
12540298fd71SBarry Smith       ierr = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
1255ce94432eSBarry Smith       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
12568188e55aSJed Brown       nsum += n;
12578188e55aSJed Brown     }
1258ce94432eSBarry Smith     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
125930bc264bSJed Brown     offset -= nsum;
1260e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1261f349c1fdSJed Brown       ierr    = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
12620298fd71SBarry Smith       ierr    = MatGetLocalSize(sub,&n,NULL);CHKERRQ(ierr);
12632ae74bdbSJed Brown       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1264ce94432eSBarry Smith       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
1265e2d7f03fSJed Brown       ierr    = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
12662ae74bdbSJed Brown       offset += n;
1267d8588912SDave May     }
1268d8588912SDave May   }
1269d8588912SDave May 
1270d8588912SDave May   if (is_col) { /* valid IS is passed in */
1271d8588912SDave May     /* refs on is[] are incremeneted */
1272e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1273d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
127426fbe8dcSKarl Rupp 
1275e2d7f03fSJed Brown       vs->isglobal.col[j] = is_col[j];
1276d8588912SDave May     }
12772ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
12782ae74bdbSJed Brown     offset = A->cmap->rstart;
12798188e55aSJed Brown     nsum   = 0;
12808188e55aSJed Brown     for (j=0; j<vs->nc; j++) {
12818188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
1282ce94432eSBarry Smith       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
12830298fd71SBarry Smith       ierr = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
1284ce94432eSBarry Smith       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
12858188e55aSJed Brown       nsum += n;
12868188e55aSJed Brown     }
1287ce94432eSBarry Smith     ierr    = MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));CHKERRQ(ierr);
128830bc264bSJed Brown     offset -= nsum;
1289e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1290f349c1fdSJed Brown       ierr    = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
12910298fd71SBarry Smith       ierr    = MatGetLocalSize(sub,NULL,&n);CHKERRQ(ierr);
12922ae74bdbSJed Brown       ierr    = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1293ce94432eSBarry Smith       ierr    = ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
1294e2d7f03fSJed Brown       ierr    = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
12952ae74bdbSJed Brown       offset += n;
1296d8588912SDave May     }
1297d8588912SDave May   }
1298e2d7f03fSJed Brown 
1299e2d7f03fSJed Brown   /* Set up the local ISs */
1300785e854fSJed Brown   ierr = PetscMalloc1(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
1301785e854fSJed Brown   ierr = PetscMalloc1(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
1302e2d7f03fSJed Brown   for (i=0,offset=0; i<vs->nr; i++) {
1303e2d7f03fSJed Brown     IS                     isloc;
13040298fd71SBarry Smith     ISLocalToGlobalMapping rmap = NULL;
1305e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1306e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
13070298fd71SBarry Smith     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,NULL);CHKERRQ(ierr);}
1308207556f9SJed Brown     if (rmap) {
1309e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1310e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
1311e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1312e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1313207556f9SJed Brown     } else {
1314207556f9SJed Brown       nlocal = 0;
13150298fd71SBarry Smith       isloc  = NULL;
1316207556f9SJed Brown     }
1317e2d7f03fSJed Brown     vs->islocal.row[i] = isloc;
1318e2d7f03fSJed Brown     offset            += nlocal;
1319e2d7f03fSJed Brown   }
13208188e55aSJed Brown   for (i=0,offset=0; i<vs->nc; i++) {
1321e2d7f03fSJed Brown     IS                     isloc;
13220298fd71SBarry Smith     ISLocalToGlobalMapping cmap = NULL;
1323e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1324e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
13250298fd71SBarry Smith     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,NULL,&cmap);CHKERRQ(ierr);}
1326207556f9SJed Brown     if (cmap) {
1327e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1328e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr);
1329e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1330e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1331207556f9SJed Brown     } else {
1332207556f9SJed Brown       nlocal = 0;
13330298fd71SBarry Smith       isloc  = NULL;
1334207556f9SJed Brown     }
1335e2d7f03fSJed Brown     vs->islocal.col[i] = isloc;
1336e2d7f03fSJed Brown     offset            += nlocal;
1337e2d7f03fSJed Brown   }
13380189643fSJed Brown 
133977019fcaSJed Brown   /* Set up the aggregate ISLocalToGlobalMapping */
134077019fcaSJed Brown   {
134145b6f7e9SBarry Smith     ISLocalToGlobalMapping rmap,cmap;
134245b6f7e9SBarry Smith     ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap);CHKERRQ(ierr);
134345b6f7e9SBarry Smith     ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap);CHKERRQ(ierr);
134477019fcaSJed Brown     if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);}
134577019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
134677019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
134777019fcaSJed Brown   }
134877019fcaSJed Brown 
13490189643fSJed Brown #if defined(PETSC_USE_DEBUG)
13500189643fSJed Brown   for (i=0; i<vs->nr; i++) {
13510189643fSJed Brown     for (j=0; j<vs->nc; j++) {
13520189643fSJed Brown       PetscInt m,n,M,N,mi,ni,Mi,Ni;
13530189643fSJed Brown       Mat      B = vs->m[i][j];
13540189643fSJed Brown       if (!B) continue;
13550189643fSJed Brown       ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr);
13560189643fSJed Brown       ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr);
13570189643fSJed Brown       ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr);
13580189643fSJed Brown       ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr);
13590189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr);
13600189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr);
1361ce94432eSBarry 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);
1362ce94432eSBarry 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);
13630189643fSJed Brown     }
13640189643fSJed Brown   }
13650189643fSJed Brown #endif
1366a061e289SJed Brown 
1367a061e289SJed Brown   /* Set A->assembled if all non-null blocks are currently assembled */
1368a061e289SJed Brown   for (i=0; i<vs->nr; i++) {
1369a061e289SJed Brown     for (j=0; j<vs->nc; j++) {
1370a061e289SJed Brown       if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0);
1371a061e289SJed Brown     }
1372a061e289SJed Brown   }
1373a061e289SJed Brown   A->assembled = PETSC_TRUE;
1374d8588912SDave May   PetscFunctionReturn(0);
1375d8588912SDave May }
1376d8588912SDave May 
137745c38901SJed Brown /*@C
1378659c6bb0SJed Brown    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
1379659c6bb0SJed Brown 
1380659c6bb0SJed Brown    Collective on Mat
1381659c6bb0SJed Brown 
1382659c6bb0SJed Brown    Input Parameter:
1383659c6bb0SJed Brown +  comm - Communicator for the new Mat
1384659c6bb0SJed Brown .  nr - number of nested row blocks
13850298fd71SBarry Smith .  is_row - index sets for each nested row block, or NULL to make contiguous
1386659c6bb0SJed Brown .  nc - number of nested column blocks
13870298fd71SBarry Smith .  is_col - index sets for each nested column block, or NULL to make contiguous
13880298fd71SBarry Smith -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL
1389659c6bb0SJed Brown 
1390659c6bb0SJed Brown    Output Parameter:
1391659c6bb0SJed Brown .  B - new matrix
1392659c6bb0SJed Brown 
1393659c6bb0SJed Brown    Level: advanced
1394659c6bb0SJed Brown 
1395950540a4SJed Brown .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST
1396659c6bb0SJed Brown @*/
13977087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1398d8588912SDave May {
1399d8588912SDave May   Mat            A;
1400d8588912SDave May   PetscErrorCode ierr;
1401d8588912SDave May 
1402d8588912SDave May   PetscFunctionBegin;
1403c8883902SJed Brown   *B   = 0;
1404d8588912SDave May   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1405c8883902SJed Brown   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
140691a28eb3SBarry Smith   A->preallocated = PETSC_TRUE;
1407c8883902SJed Brown   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1408d8588912SDave May   *B   = A;
1409d8588912SDave May   PetscFunctionReturn(0);
1410d8588912SDave May }
1411659c6bb0SJed Brown 
1412b68353e5Sstefano_zampini static PetscErrorCode MatConvert_Nest_SeqAIJ_fast(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1413b68353e5Sstefano_zampini {
1414b68353e5Sstefano_zampini   Mat_Nest       *nest = (Mat_Nest*)A->data;
141523875855Sstefano_zampini   Mat            *trans;
1416b68353e5Sstefano_zampini   PetscScalar    **avv;
1417b68353e5Sstefano_zampini   PetscScalar    *vv;
1418b68353e5Sstefano_zampini   PetscInt       **aii,**ajj;
1419b68353e5Sstefano_zampini   PetscInt       *ii,*jj,*ci;
1420b68353e5Sstefano_zampini   PetscInt       nr,nc,nnz,i,j;
1421b68353e5Sstefano_zampini   PetscBool      done;
1422b68353e5Sstefano_zampini   PetscErrorCode ierr;
1423b68353e5Sstefano_zampini 
1424b68353e5Sstefano_zampini   PetscFunctionBegin;
1425b68353e5Sstefano_zampini   ierr = MatGetSize(A,&nr,&nc);CHKERRQ(ierr);
1426b68353e5Sstefano_zampini   if (reuse == MAT_REUSE_MATRIX) {
1427b68353e5Sstefano_zampini     PetscInt rnr;
1428b68353e5Sstefano_zampini 
1429b68353e5Sstefano_zampini     ierr = MatGetRowIJ(*newmat,0,PETSC_FALSE,PETSC_FALSE,&rnr,(const PetscInt**)&ii,(const PetscInt**)&jj,&done);CHKERRQ(ierr);
1430b68353e5Sstefano_zampini     if (!done) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"MatGetRowIJ");
1431b68353e5Sstefano_zampini     if (rnr != nr) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Cannot reuse matrix, wrong number of rows");
1432b68353e5Sstefano_zampini     ierr = MatSeqAIJGetArray(*newmat,&vv);CHKERRQ(ierr);
1433b68353e5Sstefano_zampini   }
1434b68353e5Sstefano_zampini   /* extract CSR for nested SeqAIJ matrices */
1435b68353e5Sstefano_zampini   nnz  = 0;
143623875855Sstefano_zampini   ierr = PetscCalloc4(nest->nr*nest->nc,&aii,nest->nr*nest->nc,&ajj,nest->nr*nest->nc,&avv,nest->nr*nest->nc,&trans);CHKERRQ(ierr);
1437b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1438b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1439b68353e5Sstefano_zampini       Mat B = nest->m[i][j];
1440b68353e5Sstefano_zampini       if (B) {
1441b68353e5Sstefano_zampini         PetscScalar *naa;
1442b68353e5Sstefano_zampini         PetscInt    *nii,*njj,nnr;
144323875855Sstefano_zampini         PetscBool   istrans;
1444b68353e5Sstefano_zampini 
144523875855Sstefano_zampini         ierr = PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&istrans);CHKERRQ(ierr);
144623875855Sstefano_zampini         if (istrans) {
144723875855Sstefano_zampini           Mat Bt;
144823875855Sstefano_zampini 
144923875855Sstefano_zampini           ierr = MatTransposeGetMat(B,&Bt);CHKERRQ(ierr);
145023875855Sstefano_zampini           ierr = MatTranspose(Bt,MAT_INITIAL_MATRIX,&trans[i*nest->nc+j]);CHKERRQ(ierr);
145123875855Sstefano_zampini           B    = trans[i*nest->nc+j];
145223875855Sstefano_zampini         }
1453b68353e5Sstefano_zampini         ierr = MatGetRowIJ(B,0,PETSC_FALSE,PETSC_FALSE,&nnr,(const PetscInt**)&nii,(const PetscInt**)&njj,&done);CHKERRQ(ierr);
1454b68353e5Sstefano_zampini         if (!done) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"MatGetRowIJ");
1455b68353e5Sstefano_zampini         ierr = MatSeqAIJGetArray(B,&naa);CHKERRQ(ierr);
1456b68353e5Sstefano_zampini         nnz += nii[nnr];
1457b68353e5Sstefano_zampini 
1458b68353e5Sstefano_zampini         aii[i*nest->nc+j] = nii;
1459b68353e5Sstefano_zampini         ajj[i*nest->nc+j] = njj;
1460b68353e5Sstefano_zampini         avv[i*nest->nc+j] = naa;
1461b68353e5Sstefano_zampini       }
1462b68353e5Sstefano_zampini     }
1463b68353e5Sstefano_zampini   }
1464b68353e5Sstefano_zampini   if (reuse != MAT_REUSE_MATRIX) {
1465b68353e5Sstefano_zampini     ierr = PetscMalloc1(nr+1,&ii);CHKERRQ(ierr);
1466b68353e5Sstefano_zampini     ierr = PetscMalloc1(nnz,&jj);CHKERRQ(ierr);
1467b68353e5Sstefano_zampini     ierr = PetscMalloc1(nnz,&vv);CHKERRQ(ierr);
1468b68353e5Sstefano_zampini   } else {
1469b68353e5Sstefano_zampini     if (nnz != ii[nr]) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Cannot reuse matrix, wrong number of nonzeros");
1470b68353e5Sstefano_zampini   }
1471b68353e5Sstefano_zampini 
1472b68353e5Sstefano_zampini   /* new row pointer */
1473b68353e5Sstefano_zampini   ierr = PetscMemzero(ii,(nr+1)*sizeof(PetscInt));CHKERRQ(ierr);
1474b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1475b68353e5Sstefano_zampini     PetscInt       ncr,rst;
1476b68353e5Sstefano_zampini 
1477b68353e5Sstefano_zampini     ierr = ISStrideGetInfo(nest->isglobal.row[i],&rst,NULL);CHKERRQ(ierr);
1478b68353e5Sstefano_zampini     ierr = ISGetLocalSize(nest->isglobal.row[i],&ncr);CHKERRQ(ierr);
1479b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1480b68353e5Sstefano_zampini       if (aii[i*nest->nc+j]) {
1481b68353e5Sstefano_zampini         PetscInt    *nii = aii[i*nest->nc+j];
1482b68353e5Sstefano_zampini         PetscInt    ir;
1483b68353e5Sstefano_zampini 
1484b68353e5Sstefano_zampini         for (ir=rst; ir<ncr+rst; ++ir) {
1485b68353e5Sstefano_zampini           ii[ir+1] += nii[1]-nii[0];
1486b68353e5Sstefano_zampini           nii++;
1487b68353e5Sstefano_zampini         }
1488b68353e5Sstefano_zampini       }
1489b68353e5Sstefano_zampini     }
1490b68353e5Sstefano_zampini   }
1491b68353e5Sstefano_zampini   for (i=0; i<nr; i++) ii[i+1] += ii[i];
1492b68353e5Sstefano_zampini 
1493b68353e5Sstefano_zampini   /* construct CSR for the new matrix */
1494b68353e5Sstefano_zampini   ierr = PetscCalloc1(nr,&ci);CHKERRQ(ierr);
1495b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1496b68353e5Sstefano_zampini     PetscInt       ncr,rst;
1497b68353e5Sstefano_zampini 
1498b68353e5Sstefano_zampini     ierr = ISStrideGetInfo(nest->isglobal.row[i],&rst,NULL);CHKERRQ(ierr);
1499b68353e5Sstefano_zampini     ierr = ISGetLocalSize(nest->isglobal.row[i],&ncr);CHKERRQ(ierr);
1500b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1501b68353e5Sstefano_zampini       if (aii[i*nest->nc+j]) {
1502b68353e5Sstefano_zampini         PetscScalar *nvv = avv[i*nest->nc+j];
1503b68353e5Sstefano_zampini         PetscInt    *nii = aii[i*nest->nc+j];
1504b68353e5Sstefano_zampini         PetscInt    *njj = ajj[i*nest->nc+j];
1505b68353e5Sstefano_zampini         PetscInt    ir,cst;
1506b68353e5Sstefano_zampini 
1507b68353e5Sstefano_zampini         ierr = ISStrideGetInfo(nest->isglobal.col[j],&cst,NULL);CHKERRQ(ierr);
1508b68353e5Sstefano_zampini         for (ir=rst; ir<ncr+rst; ++ir) {
1509b68353e5Sstefano_zampini           PetscInt ij,rsize = nii[1]-nii[0],ist = ii[ir]+ci[ir];
1510b68353e5Sstefano_zampini 
1511b68353e5Sstefano_zampini           for (ij=0;ij<rsize;ij++) {
1512b68353e5Sstefano_zampini             jj[ist+ij] = *njj+cst;
1513b68353e5Sstefano_zampini             vv[ist+ij] = *nvv;
1514b68353e5Sstefano_zampini             njj++;
1515b68353e5Sstefano_zampini             nvv++;
1516b68353e5Sstefano_zampini           }
1517b68353e5Sstefano_zampini           ci[ir] += rsize;
1518b68353e5Sstefano_zampini           nii++;
1519b68353e5Sstefano_zampini         }
1520b68353e5Sstefano_zampini       }
1521b68353e5Sstefano_zampini     }
1522b68353e5Sstefano_zampini   }
1523b68353e5Sstefano_zampini   ierr = PetscFree(ci);CHKERRQ(ierr);
1524b68353e5Sstefano_zampini 
1525b68353e5Sstefano_zampini   /* restore info */
1526b68353e5Sstefano_zampini   for (i=0; i<nest->nr; ++i) {
1527b68353e5Sstefano_zampini     for (j=0; j<nest->nc; ++j) {
1528b68353e5Sstefano_zampini       Mat B = nest->m[i][j];
1529b68353e5Sstefano_zampini       if (B) {
1530b68353e5Sstefano_zampini         PetscInt nnr = 0, k = i*nest->nc+j;
153123875855Sstefano_zampini 
153223875855Sstefano_zampini         B    = (trans[k] ? trans[k] : B);
1533b68353e5Sstefano_zampini         ierr = MatRestoreRowIJ(B,0,PETSC_FALSE,PETSC_FALSE,&nnr,(const PetscInt**)&aii[k],(const PetscInt**)&ajj[k],&done);CHKERRQ(ierr);
1534b68353e5Sstefano_zampini         if (!done) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"MatRestoreRowIJ");
1535b68353e5Sstefano_zampini         ierr = MatSeqAIJRestoreArray(B,&avv[k]);CHKERRQ(ierr);
153623875855Sstefano_zampini         ierr = MatDestroy(&trans[k]);CHKERRQ(ierr);
1537b68353e5Sstefano_zampini       }
1538b68353e5Sstefano_zampini     }
1539b68353e5Sstefano_zampini   }
154023875855Sstefano_zampini   ierr = PetscFree4(aii,ajj,avv,trans);CHKERRQ(ierr);
1541b68353e5Sstefano_zampini 
1542b68353e5Sstefano_zampini   /* finalize newmat */
1543b68353e5Sstefano_zampini   if (reuse == MAT_INITIAL_MATRIX) {
1544b68353e5Sstefano_zampini     ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),nr,nc,ii,jj,vv,newmat);CHKERRQ(ierr);
1545b68353e5Sstefano_zampini   } else if (reuse == MAT_INPLACE_MATRIX) {
1546b68353e5Sstefano_zampini     Mat B;
1547b68353e5Sstefano_zampini 
1548b68353e5Sstefano_zampini     ierr = MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),nr,nc,ii,jj,vv,&B);CHKERRQ(ierr);
1549b68353e5Sstefano_zampini     ierr = MatHeaderReplace(A,&B);CHKERRQ(ierr);
1550b68353e5Sstefano_zampini   }
1551b68353e5Sstefano_zampini   ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1552b68353e5Sstefano_zampini   ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1553b68353e5Sstefano_zampini   {
1554b68353e5Sstefano_zampini     Mat_SeqAIJ *a = (Mat_SeqAIJ*)((*newmat)->data);
1555b68353e5Sstefano_zampini     a->free_a     = PETSC_TRUE;
1556b68353e5Sstefano_zampini     a->free_ij    = PETSC_TRUE;
1557b68353e5Sstefano_zampini   }
1558b68353e5Sstefano_zampini   PetscFunctionReturn(0);
1559b68353e5Sstefano_zampini }
1560b68353e5Sstefano_zampini 
1561cc2e6a90SBarry Smith PETSC_INTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1562629c3df2SDmitry Karpeev {
1563629c3df2SDmitry Karpeev   PetscErrorCode ierr;
1564629c3df2SDmitry Karpeev   Mat_Nest       *nest = (Mat_Nest*)A->data;
156583b1a929SMark Adams   PetscInt       m,n,M,N,i,j,k,*dnnz,*onnz,rstart;
1566649b366bSFande Kong   PetscInt       cstart,cend;
1567b68353e5Sstefano_zampini   PetscMPIInt    size;
1568629c3df2SDmitry Karpeev   Mat            C;
1569629c3df2SDmitry Karpeev 
1570629c3df2SDmitry Karpeev   PetscFunctionBegin;
1571b68353e5Sstefano_zampini   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr);
1572b68353e5Sstefano_zampini   if (size == 1) { /* look for a special case with SeqAIJ matrices and strided-1, contiguous, blocks */
1573b68353e5Sstefano_zampini     PetscInt  nf;
1574b68353e5Sstefano_zampini     PetscBool fast;
1575b68353e5Sstefano_zampini 
1576b68353e5Sstefano_zampini     ierr = PetscStrcmp(newtype,MATAIJ,&fast);CHKERRQ(ierr);
1577b68353e5Sstefano_zampini     if (!fast) {
1578b68353e5Sstefano_zampini       ierr = PetscStrcmp(newtype,MATSEQAIJ,&fast);CHKERRQ(ierr);
1579b68353e5Sstefano_zampini     }
1580b68353e5Sstefano_zampini     for (i=0; i<nest->nr && fast; ++i) {
1581b68353e5Sstefano_zampini       for (j=0; j<nest->nc && fast; ++j) {
1582b68353e5Sstefano_zampini         Mat B = nest->m[i][j];
1583b68353e5Sstefano_zampini         if (B) {
1584b68353e5Sstefano_zampini           ierr = PetscObjectTypeCompare((PetscObject)B,MATSEQAIJ,&fast);CHKERRQ(ierr);
158523875855Sstefano_zampini           if (!fast) {
158623875855Sstefano_zampini             PetscBool istrans;
158723875855Sstefano_zampini 
158823875855Sstefano_zampini             ierr = PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&istrans);CHKERRQ(ierr);
158923875855Sstefano_zampini             if (istrans) {
159023875855Sstefano_zampini               Mat Bt;
159123875855Sstefano_zampini 
159223875855Sstefano_zampini               ierr = MatTransposeGetMat(B,&Bt);CHKERRQ(ierr);
159323875855Sstefano_zampini               ierr = PetscObjectTypeCompare((PetscObject)Bt,MATSEQAIJ,&fast);CHKERRQ(ierr);
159423875855Sstefano_zampini             }
1595b68353e5Sstefano_zampini           }
1596b68353e5Sstefano_zampini         }
1597b68353e5Sstefano_zampini       }
1598b68353e5Sstefano_zampini     }
1599b68353e5Sstefano_zampini     for (i=0, nf=0; i<nest->nr && fast; ++i) {
1600b68353e5Sstefano_zampini       ierr = PetscObjectTypeCompare((PetscObject)nest->isglobal.row[i],ISSTRIDE,&fast);CHKERRQ(ierr);
1601b68353e5Sstefano_zampini       if (fast) {
1602b68353e5Sstefano_zampini         PetscInt f,s;
1603b68353e5Sstefano_zampini 
1604b68353e5Sstefano_zampini         ierr = ISStrideGetInfo(nest->isglobal.row[i],&f,&s);CHKERRQ(ierr);
1605b68353e5Sstefano_zampini         if (f != nf || s != 1) { fast = PETSC_FALSE; }
1606b68353e5Sstefano_zampini         else {
1607b68353e5Sstefano_zampini           ierr = ISGetSize(nest->isglobal.row[i],&f);CHKERRQ(ierr);
1608b68353e5Sstefano_zampini           nf  += f;
1609b68353e5Sstefano_zampini         }
1610b68353e5Sstefano_zampini       }
1611b68353e5Sstefano_zampini     }
1612b68353e5Sstefano_zampini     for (i=0, nf=0; i<nest->nc && fast; ++i) {
1613b68353e5Sstefano_zampini       ierr = PetscObjectTypeCompare((PetscObject)nest->isglobal.col[i],ISSTRIDE,&fast);CHKERRQ(ierr);
1614b68353e5Sstefano_zampini       if (fast) {
1615b68353e5Sstefano_zampini         PetscInt f,s;
1616b68353e5Sstefano_zampini 
1617b68353e5Sstefano_zampini         ierr = ISStrideGetInfo(nest->isglobal.col[i],&f,&s);CHKERRQ(ierr);
1618b68353e5Sstefano_zampini         if (f != nf || s != 1) { fast = PETSC_FALSE; }
1619b68353e5Sstefano_zampini         else {
1620b68353e5Sstefano_zampini           ierr = ISGetSize(nest->isglobal.col[i],&f);CHKERRQ(ierr);
1621b68353e5Sstefano_zampini           nf  += f;
1622b68353e5Sstefano_zampini         }
1623b68353e5Sstefano_zampini       }
1624b68353e5Sstefano_zampini     }
1625b68353e5Sstefano_zampini     if (fast) {
1626b68353e5Sstefano_zampini       ierr = MatConvert_Nest_SeqAIJ_fast(A,newtype,reuse,newmat);CHKERRQ(ierr);
1627b68353e5Sstefano_zampini       PetscFunctionReturn(0);
1628b68353e5Sstefano_zampini     }
1629b68353e5Sstefano_zampini   }
1630629c3df2SDmitry Karpeev   ierr = MatGetSize(A,&M,&N);CHKERRQ(ierr);
1631629c3df2SDmitry Karpeev   ierr = MatGetLocalSize(A,&m,&n);CHKERRQ(ierr);
1632649b366bSFande Kong   ierr = MatGetOwnershipRangeColumn(A,&cstart,&cend);CHKERRQ(ierr);
1633629c3df2SDmitry Karpeev   switch (reuse) {
1634629c3df2SDmitry Karpeev   case MAT_INITIAL_MATRIX:
1635ce94432eSBarry Smith     ierr    = MatCreate(PetscObjectComm((PetscObject)A),&C);CHKERRQ(ierr);
1636629c3df2SDmitry Karpeev     ierr    = MatSetType(C,newtype);CHKERRQ(ierr);
1637629c3df2SDmitry Karpeev     ierr    = MatSetSizes(C,m,n,M,N);CHKERRQ(ierr);
1638629c3df2SDmitry Karpeev     *newmat = C;
1639629c3df2SDmitry Karpeev     break;
1640629c3df2SDmitry Karpeev   case MAT_REUSE_MATRIX:
1641629c3df2SDmitry Karpeev     C = *newmat;
1642629c3df2SDmitry Karpeev     break;
1643ce94432eSBarry Smith   default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse");
1644629c3df2SDmitry Karpeev   }
1645785e854fSJed Brown   ierr = PetscMalloc1(2*m,&dnnz);CHKERRQ(ierr);
1646629c3df2SDmitry Karpeev   onnz = dnnz + m;
1647629c3df2SDmitry Karpeev   for (k=0; k<m; k++) {
1648629c3df2SDmitry Karpeev     dnnz[k] = 0;
1649629c3df2SDmitry Karpeev     onnz[k] = 0;
1650629c3df2SDmitry Karpeev   }
1651629c3df2SDmitry Karpeev   for (j=0; j<nest->nc; ++j) {
1652629c3df2SDmitry Karpeev     IS             bNis;
1653629c3df2SDmitry Karpeev     PetscInt       bN;
1654629c3df2SDmitry Karpeev     const PetscInt *bNindices;
1655629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
1656629c3df2SDmitry Karpeev     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1657629c3df2SDmitry Karpeev     ierr = ISGetSize(bNis, &bN);CHKERRQ(ierr);
1658629c3df2SDmitry Karpeev     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1659629c3df2SDmitry Karpeev     for (i=0; i<nest->nr; ++i) {
1660629c3df2SDmitry Karpeev       PetscSF        bmsf;
1661649b366bSFande Kong       PetscSFNode    *iremote;
1662629c3df2SDmitry Karpeev       Mat            B;
1663649b366bSFande Kong       PetscInt       bm, *sub_dnnz,*sub_onnz, br;
1664629c3df2SDmitry Karpeev       const PetscInt *bmindices;
1665629c3df2SDmitry Karpeev       B = nest->m[i][j];
1666629c3df2SDmitry Karpeev       if (!B) continue;
1667629c3df2SDmitry Karpeev       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1668629c3df2SDmitry Karpeev       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1669ce94432eSBarry Smith       ierr = PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);CHKERRQ(ierr);
1670649b366bSFande Kong       ierr = PetscMalloc1(bm,&iremote);CHKERRQ(ierr);
1671649b366bSFande Kong       ierr = PetscMalloc1(bm,&sub_dnnz);CHKERRQ(ierr);
1672649b366bSFande Kong       ierr = PetscMalloc1(bm,&sub_onnz);CHKERRQ(ierr);
1673649b366bSFande Kong       for (k = 0; k < bm; ++k){
1674649b366bSFande Kong     	sub_dnnz[k] = 0;
1675649b366bSFande Kong     	sub_onnz[k] = 0;
1676649b366bSFande Kong       }
1677629c3df2SDmitry Karpeev       /*
1678629c3df2SDmitry Karpeev        Locate the owners for all of the locally-owned global row indices for this row block.
1679629c3df2SDmitry Karpeev        These determine the roots of PetscSF used to communicate preallocation data to row owners.
1680629c3df2SDmitry Karpeev        The roots correspond to the dnnz and onnz entries; thus, there are two roots per row.
1681629c3df2SDmitry Karpeev        */
168283b1a929SMark Adams       ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr);
1683629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
1684649b366bSFande Kong         PetscInt       row = bmindices[br], rowowner = 0, brncols, col;
1685629c3df2SDmitry Karpeev         const PetscInt *brcols;
1686a4b3d3acSMatthew G Knepley         PetscInt       rowrel = 0; /* row's relative index on its owner rank */
1687629c3df2SDmitry Karpeev         ierr      = PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);CHKERRQ(ierr);
1688649b366bSFande Kong         /* how many roots  */
1689649b366bSFande Kong         iremote[br].rank = rowowner; iremote[br].index = rowrel;           /* edge from bmdnnz to dnnz */
1690649b366bSFande Kong         /* get nonzero pattern */
169183b1a929SMark Adams         ierr = MatGetRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr);
1692629c3df2SDmitry Karpeev         for (k=0; k<brncols; k++) {
1693629c3df2SDmitry Karpeev           col  = bNindices[brcols[k]];
1694649b366bSFande Kong           if (col>=A->cmap->range[rowowner] && col<A->cmap->range[rowowner+1]) {
1695649b366bSFande Kong             sub_dnnz[br]++;
1696649b366bSFande Kong           } else {
1697649b366bSFande Kong             sub_onnz[br]++;
1698649b366bSFande Kong           }
1699629c3df2SDmitry Karpeev         }
170083b1a929SMark Adams         ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,NULL);CHKERRQ(ierr);
1701629c3df2SDmitry Karpeev       }
1702629c3df2SDmitry Karpeev       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1703629c3df2SDmitry Karpeev       /* bsf will have to take care of disposing of bedges. */
1704649b366bSFande Kong       ierr = PetscSFSetGraph(bmsf,m,bm,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr);
1705649b366bSFande Kong       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr);
1706649b366bSFande Kong       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);CHKERRQ(ierr);
1707649b366bSFande Kong       ierr = PetscSFReduceBegin(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr);
1708649b366bSFande Kong       ierr = PetscSFReduceEnd(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);CHKERRQ(ierr);
1709649b366bSFande Kong       ierr = PetscFree(sub_dnnz);CHKERRQ(ierr);
1710649b366bSFande Kong       ierr = PetscFree(sub_onnz);CHKERRQ(ierr);
1711629c3df2SDmitry Karpeev       ierr = PetscSFDestroy(&bmsf);CHKERRQ(ierr);
1712629c3df2SDmitry Karpeev     }
171322d28d08SBarry Smith     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
1714629c3df2SDmitry Karpeev     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
171565a4a0a3Sstefano_zampini   }
171665a4a0a3Sstefano_zampini   /* Resize preallocation if overestimated */
171765a4a0a3Sstefano_zampini   for (i=0;i<m;i++) {
171865a4a0a3Sstefano_zampini     dnnz[i] = PetscMin(dnnz[i],A->cmap->n);
171965a4a0a3Sstefano_zampini     onnz[i] = PetscMin(onnz[i],A->cmap->N - A->cmap->n);
1720629c3df2SDmitry Karpeev   }
1721629c3df2SDmitry Karpeev   ierr = MatSeqAIJSetPreallocation(C,0,dnnz);CHKERRQ(ierr);
1722629c3df2SDmitry Karpeev   ierr = MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);CHKERRQ(ierr);
1723629c3df2SDmitry Karpeev   ierr = PetscFree(dnnz);CHKERRQ(ierr);
1724629c3df2SDmitry Karpeev 
1725629c3df2SDmitry Karpeev   /* Fill by row */
1726629c3df2SDmitry Karpeev   for (j=0; j<nest->nc; ++j) {
1727629c3df2SDmitry Karpeev     /* Using global column indices and ISAllGather() is not scalable. */
1728629c3df2SDmitry Karpeev     IS             bNis;
1729629c3df2SDmitry Karpeev     PetscInt       bN;
1730629c3df2SDmitry Karpeev     const PetscInt *bNindices;
1731629c3df2SDmitry Karpeev     ierr = ISAllGather(nest->isglobal.col[j], &bNis);CHKERRQ(ierr);
1732629c3df2SDmitry Karpeev     ierr = ISGetSize(bNis,&bN);CHKERRQ(ierr);
1733629c3df2SDmitry Karpeev     ierr = ISGetIndices(bNis,&bNindices);CHKERRQ(ierr);
1734629c3df2SDmitry Karpeev     for (i=0; i<nest->nr; ++i) {
1735629c3df2SDmitry Karpeev       Mat            B;
1736629c3df2SDmitry Karpeev       PetscInt       bm, br;
1737629c3df2SDmitry Karpeev       const PetscInt *bmindices;
1738629c3df2SDmitry Karpeev       B = nest->m[i][j];
1739629c3df2SDmitry Karpeev       if (!B) continue;
1740629c3df2SDmitry Karpeev       ierr = ISGetLocalSize(nest->isglobal.row[i],&bm);CHKERRQ(ierr);
1741629c3df2SDmitry Karpeev       ierr = ISGetIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
174283b1a929SMark Adams       ierr = MatGetOwnershipRange(B,&rstart,NULL);CHKERRQ(ierr);
1743629c3df2SDmitry Karpeev       for (br = 0; br < bm; ++br) {
1744629c3df2SDmitry Karpeev         PetscInt          row = bmindices[br], brncols,  *cols;
1745629c3df2SDmitry Karpeev         const PetscInt    *brcols;
1746629c3df2SDmitry Karpeev         const PetscScalar *brcoldata;
174783b1a929SMark Adams         ierr = MatGetRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
1748785e854fSJed Brown         ierr = PetscMalloc1(brncols,&cols);CHKERRQ(ierr);
174926fbe8dcSKarl Rupp         for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]];
1750629c3df2SDmitry Karpeev         /*
1751629c3df2SDmitry Karpeev           Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match.
1752629c3df2SDmitry Karpeev           Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES.
1753629c3df2SDmitry Karpeev          */
1754a2ea699eSBarry Smith         ierr = MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);CHKERRQ(ierr);
175583b1a929SMark Adams         ierr = MatRestoreRow(B,br+rstart,&brncols,&brcols,&brcoldata);CHKERRQ(ierr);
1756629c3df2SDmitry Karpeev         ierr = PetscFree(cols);CHKERRQ(ierr);
1757629c3df2SDmitry Karpeev       }
1758629c3df2SDmitry Karpeev       ierr = ISRestoreIndices(nest->isglobal.row[i],&bmindices);CHKERRQ(ierr);
1759629c3df2SDmitry Karpeev     }
1760a2ea699eSBarry Smith     ierr = ISRestoreIndices(bNis,&bNindices);CHKERRQ(ierr);
1761629c3df2SDmitry Karpeev     ierr = ISDestroy(&bNis);CHKERRQ(ierr);
1762629c3df2SDmitry Karpeev   }
1763629c3df2SDmitry Karpeev   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1764629c3df2SDmitry Karpeev   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1765629c3df2SDmitry Karpeev   PetscFunctionReturn(0);
1766629c3df2SDmitry Karpeev }
1767629c3df2SDmitry Karpeev 
17688b7d3b4bSBarry Smith PetscErrorCode  MatHasOperation_Nest(Mat mat,MatOperation op,PetscBool  *has)
17698b7d3b4bSBarry Smith {
17708b7d3b4bSBarry Smith   PetscFunctionBegin;
17718b7d3b4bSBarry Smith   *has = PETSC_FALSE;
17728b7d3b4bSBarry Smith   if (op == MATOP_MULT_TRANSPOSE) {
17738b7d3b4bSBarry Smith     Mat_Nest       *bA = (Mat_Nest*)mat->data;
17748b7d3b4bSBarry Smith     PetscInt       i,j,nr = bA->nr,nc = bA->nc;
17758b7d3b4bSBarry Smith     PetscErrorCode ierr;
17768b7d3b4bSBarry Smith     PetscBool      flg;
17778b7d3b4bSBarry Smith 
17788b7d3b4bSBarry Smith     for (j=0; j<nc; j++) {
17798b7d3b4bSBarry Smith       for (i=0; i<nr; i++) {
17808b7d3b4bSBarry Smith         if (!bA->m[i][j]) continue;
17818b7d3b4bSBarry Smith         ierr = MatHasOperation(bA->m[i][j],MATOP_MULT_TRANSPOSE_ADD,&flg);CHKERRQ(ierr);
17828b7d3b4bSBarry Smith         if (!flg) PetscFunctionReturn(0);
17838b7d3b4bSBarry Smith       }
17848b7d3b4bSBarry Smith     }
17858b7d3b4bSBarry Smith   }
17868b7d3b4bSBarry Smith   if (((void**)mat->ops)[op]) *has =  PETSC_TRUE;
17878b7d3b4bSBarry Smith   PetscFunctionReturn(0);
17888b7d3b4bSBarry Smith }
17898b7d3b4bSBarry Smith 
1790659c6bb0SJed Brown /*MC
1791659c6bb0SJed Brown   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
1792659c6bb0SJed Brown 
1793659c6bb0SJed Brown   Level: intermediate
1794659c6bb0SJed Brown 
1795659c6bb0SJed Brown   Notes:
1796659c6bb0SJed Brown   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1797659c6bb0SJed Brown   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1798950540a4SJed Brown   It is usually used with DMComposite and DMCreateMatrix()
1799659c6bb0SJed Brown 
18008b7d3b4bSBarry Smith   Each of the submatrices lives on the same MPI communicator as the original nest matrix (though they can have zero
18018b7d3b4bSBarry Smith   rows/columns on some processes.) Thus this is not meant for cases where the submatrices live on far fewer processes
18028b7d3b4bSBarry Smith   than the nest matrix.
18038b7d3b4bSBarry Smith 
1804659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest()
1805659c6bb0SJed Brown M*/
18068cc058d9SJed Brown PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A)
1807c8883902SJed Brown {
1808c8883902SJed Brown   Mat_Nest       *s;
1809c8883902SJed Brown   PetscErrorCode ierr;
1810c8883902SJed Brown 
1811c8883902SJed Brown   PetscFunctionBegin;
1812b00a9115SJed Brown   ierr    = PetscNewLog(A,&s);CHKERRQ(ierr);
1813c8883902SJed Brown   A->data = (void*)s;
1814e7c19651SJed Brown 
1815e7c19651SJed Brown   s->nr            = -1;
1816e7c19651SJed Brown   s->nc            = -1;
18170298fd71SBarry Smith   s->m             = NULL;
1818e7c19651SJed Brown   s->splitassembly = PETSC_FALSE;
1819c8883902SJed Brown 
1820c8883902SJed Brown   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
182126fbe8dcSKarl Rupp 
1822c8883902SJed Brown   A->ops->mult                  = MatMult_Nest;
18239194d70fSJed Brown   A->ops->multadd               = MatMultAdd_Nest;
1824c8883902SJed Brown   A->ops->multtranspose         = MatMultTranspose_Nest;
18259194d70fSJed Brown   A->ops->multtransposeadd      = MatMultTransposeAdd_Nest;
1826f8170845SAlex Fikl   A->ops->transpose             = MatTranspose_Nest;
1827c8883902SJed Brown   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1828c8883902SJed Brown   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1829c8883902SJed Brown   A->ops->zeroentries           = MatZeroEntries_Nest;
1830c222c20dSDavid Ham   A->ops->copy                  = MatCopy_Nest;
1831c8883902SJed Brown   A->ops->duplicate             = MatDuplicate_Nest;
18327dae84e0SHong Zhang   A->ops->createsubmatrix       = MatCreateSubMatrix_Nest;
1833c8883902SJed Brown   A->ops->destroy               = MatDestroy_Nest;
1834c8883902SJed Brown   A->ops->view                  = MatView_Nest;
1835c8883902SJed Brown   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1836c8883902SJed Brown   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1837c8883902SJed Brown   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
1838429bac76SJed Brown   A->ops->getdiagonal           = MatGetDiagonal_Nest;
1839429bac76SJed Brown   A->ops->diagonalscale         = MatDiagonalScale_Nest;
1840a061e289SJed Brown   A->ops->scale                 = MatScale_Nest;
1841a061e289SJed Brown   A->ops->shift                 = MatShift_Nest;
184213135bc6SAlex Fikl   A->ops->diagonalset           = MatDiagonalSet_Nest;
1843f8170845SAlex Fikl   A->ops->setrandom             = MatSetRandom_Nest;
18448b7d3b4bSBarry Smith   A->ops->hasoperation          = MatHasOperation_Nest;
1845c8883902SJed Brown 
1846c8883902SJed Brown   A->spptr        = 0;
1847c8883902SJed Brown   A->assembled    = PETSC_FALSE;
1848c8883902SJed Brown 
1849c8883902SJed Brown   /* expose Nest api's */
1850bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",      MatNestGetSubMat_Nest);CHKERRQ(ierr);
1851bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",      MatNestSetSubMat_Nest);CHKERRQ(ierr);
1852bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",     MatNestGetSubMats_Nest);CHKERRQ(ierr);
1853bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",        MatNestGetSize_Nest);CHKERRQ(ierr);
1854bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",         MatNestGetISs_Nest);CHKERRQ(ierr);
1855bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",    MatNestGetLocalISs_Nest);CHKERRQ(ierr);
1856bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",     MatNestSetVecType_Nest);CHKERRQ(ierr);
1857bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",     MatNestSetSubMats_Nest);CHKERRQ(ierr);
18580899c546SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_mpiaij_C",MatConvert_Nest_AIJ);CHKERRQ(ierr);
18590899c546SStefano Zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_seqaij_C",MatConvert_Nest_AIJ);CHKERRQ(ierr);
186083b1a929SMark Adams   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",   MatConvert_Nest_AIJ);CHKERRQ(ierr);
18615e3038f0Sstefano_zampini   ierr = PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_is_C",    MatConvert_Nest_IS);CHKERRQ(ierr);
1862c8883902SJed Brown 
1863c8883902SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
1864c8883902SJed Brown   PetscFunctionReturn(0);
1865c8883902SJed Brown }
1866