xref: /petsc/src/mat/impls/nest/matnest.c (revision 77019fcace76d50a517c8cdc493ef0c5ee9b5b6a)
1d8588912SDave May 
2ec9811eeSJed Brown #include "matnestimpl.h" /*I   "petscmat.h"   I*/
3d8588912SDave May 
4c8883902SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]);
57874fa86SDave May static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left);
6c8883902SJed Brown 
7d8588912SDave May /* private functions */
8d8588912SDave May #undef __FUNCT__
98188e55aSJed Brown #define __FUNCT__ "MatNestGetSizes_Private"
108188e55aSJed Brown static PetscErrorCode MatNestGetSizes_Private(Mat A,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N)
11d8588912SDave May {
12d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
138188e55aSJed Brown   PetscInt       i,j;
14d8588912SDave May   PetscErrorCode ierr;
15d8588912SDave May 
16d8588912SDave May   PetscFunctionBegin;
178188e55aSJed Brown   *m = *n = *M = *N = 0;
188188e55aSJed Brown   for (i=0; i<bA->nr; i++) {  /* rows */
198188e55aSJed Brown     PetscInt sm,sM;
208188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.row[i],&sm);CHKERRQ(ierr);
218188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.row[i],&sM);CHKERRQ(ierr);
228188e55aSJed Brown     *m += sm;
238188e55aSJed Brown     *M += sM;
24d8588912SDave May   }
258188e55aSJed Brown   for (j=0; j<bA->nc; j++) {  /* cols */
268188e55aSJed Brown     PetscInt sn,sN;
278188e55aSJed Brown     ierr = ISGetLocalSize(bA->isglobal.col[j],&sn);CHKERRQ(ierr);
288188e55aSJed Brown     ierr = ISGetSize(bA->isglobal.col[j],&sN);CHKERRQ(ierr);
298188e55aSJed Brown     *n += sn;
308188e55aSJed Brown     *N += sN;
31d8588912SDave May   }
32d8588912SDave May   PetscFunctionReturn(0);
33d8588912SDave May }
34d8588912SDave May 
35d8588912SDave May /* operations */
36d8588912SDave May #undef __FUNCT__
37d8588912SDave May #define __FUNCT__ "MatMult_Nest"
38207556f9SJed Brown static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y)
39d8588912SDave May {
40d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
41207556f9SJed Brown   Vec            *bx = bA->right,*by = bA->left;
42207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
43d8588912SDave May   PetscErrorCode ierr;
44d8588912SDave May 
45d8588912SDave May   PetscFunctionBegin;
46207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
47207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
48207556f9SJed Brown   for (i=0; i<nr; i++) {
49d8588912SDave May     ierr = VecZeroEntries(by[i]);CHKERRQ(ierr);
50207556f9SJed Brown     for (j=0; j<nc; j++) {
51207556f9SJed Brown       if (!bA->m[i][j]) continue;
52d8588912SDave May       /* y[i] <- y[i] + A[i][j] * x[j] */
53d8588912SDave May       ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr);
54d8588912SDave May     }
55d8588912SDave May   }
56207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
57207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
58d8588912SDave May   PetscFunctionReturn(0);
59d8588912SDave May }
60d8588912SDave May 
61d8588912SDave May #undef __FUNCT__
629194d70fSJed Brown #define __FUNCT__ "MatMultAdd_Nest"
639194d70fSJed Brown static PetscErrorCode MatMultAdd_Nest(Mat A,Vec x,Vec y,Vec z)
649194d70fSJed Brown {
659194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
669194d70fSJed Brown   Vec            *bx = bA->right,*bz = bA->left;
679194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
689194d70fSJed Brown   PetscErrorCode ierr;
699194d70fSJed Brown 
709194d70fSJed Brown   PetscFunctionBegin;
719194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
729194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
739194d70fSJed Brown   for (i=0; i<nr; i++) {
749194d70fSJed Brown     if (y != z) {
759194d70fSJed Brown       Vec by;
769194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
779194d70fSJed Brown       ierr = VecCopy(by,bz[i]);CHKERRQ(ierr);
78336d21e7SJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by);CHKERRQ(ierr);
799194d70fSJed Brown     }
809194d70fSJed Brown     for (j=0; j<nc; j++) {
819194d70fSJed Brown       if (!bA->m[i][j]) continue;
829194d70fSJed Brown       /* y[i] <- y[i] + A[i][j] * x[j] */
839194d70fSJed Brown       ierr = MatMultAdd(bA->m[i][j],bx[j],bz[i],bz[i]);CHKERRQ(ierr);
849194d70fSJed Brown     }
859194d70fSJed Brown   }
869194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.row[i],&bz[i]);CHKERRQ(ierr);}
879194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
889194d70fSJed Brown   PetscFunctionReturn(0);
899194d70fSJed Brown }
909194d70fSJed Brown 
919194d70fSJed Brown #undef __FUNCT__
92d8588912SDave May #define __FUNCT__ "MatMultTranspose_Nest"
93207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y)
94d8588912SDave May {
95d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
96207556f9SJed Brown   Vec            *bx = bA->left,*by = bA->right;
97207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
98d8588912SDave May   PetscErrorCode ierr;
99d8588912SDave May 
100d8588912SDave May   PetscFunctionBegin;
101609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
102609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
103207556f9SJed Brown   for (j=0; j<nc; j++) {
104609e31cbSJed Brown     ierr = VecZeroEntries(by[j]);CHKERRQ(ierr);
105609e31cbSJed Brown     for (i=0; i<nr; i++) {
106207556f9SJed Brown       if (!bA->m[j][i]) continue;
107609e31cbSJed Brown       /* y[j] <- y[j] + (A[i][j])^T * x[i] */
108609e31cbSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],by[j],by[j]);CHKERRQ(ierr);
109d8588912SDave May     }
110d8588912SDave May   }
111609e31cbSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
112609e31cbSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
113d8588912SDave May   PetscFunctionReturn(0);
114d8588912SDave May }
115d8588912SDave May 
116d8588912SDave May #undef __FUNCT__
1179194d70fSJed Brown #define __FUNCT__ "MatMultTransposeAdd_Nest"
1189194d70fSJed Brown static PetscErrorCode MatMultTransposeAdd_Nest(Mat A,Vec x,Vec y,Vec z)
1199194d70fSJed Brown {
1209194d70fSJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
1219194d70fSJed Brown   Vec            *bx = bA->left,*bz = bA->right;
1229194d70fSJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
1239194d70fSJed Brown   PetscErrorCode ierr;
1249194d70fSJed Brown 
1259194d70fSJed Brown   PetscFunctionBegin;
1269194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
1279194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
1289194d70fSJed Brown   for (j=0; j<nc; j++) {
1299194d70fSJed Brown     if (y != z) {
1309194d70fSJed Brown       Vec by;
1319194d70fSJed Brown       ierr = VecGetSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
1329194d70fSJed Brown       ierr = VecCopy(by,bz[j]);CHKERRQ(ierr);
1339194d70fSJed Brown       ierr = VecRestoreSubVector(y,bA->isglobal.col[j],&by);CHKERRQ(ierr);
1349194d70fSJed Brown     }
1359194d70fSJed Brown     for (i=0; i<nr; i++) {
1369194d70fSJed Brown       if (!bA->m[j][i]) continue;
1379194d70fSJed Brown       /* z[j] <- y[j] + (A[i][j])^T * x[i] */
1389194d70fSJed Brown       ierr = MatMultTransposeAdd(bA->m[i][j],bx[i],bz[j],bz[j]);CHKERRQ(ierr);
1399194d70fSJed Brown     }
1409194d70fSJed Brown   }
1419194d70fSJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
1429194d70fSJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(z,bA->isglobal.col[i],&bz[i]);CHKERRQ(ierr);}
1439194d70fSJed Brown   PetscFunctionReturn(0);
1449194d70fSJed Brown }
1459194d70fSJed Brown 
1469194d70fSJed Brown #undef __FUNCT__
147e2d7f03fSJed Brown #define __FUNCT__ "MatNestDestroyISList"
148e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list)
149e2d7f03fSJed Brown {
150e2d7f03fSJed Brown   PetscErrorCode ierr;
151e2d7f03fSJed Brown   IS             *lst = *list;
152e2d7f03fSJed Brown   PetscInt       i;
153e2d7f03fSJed Brown 
154e2d7f03fSJed Brown   PetscFunctionBegin;
155e2d7f03fSJed Brown   if (!lst) PetscFunctionReturn(0);
1566bf464f9SBarry Smith   for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(&lst[i]);CHKERRQ(ierr);}
157e2d7f03fSJed Brown   ierr = PetscFree(lst);CHKERRQ(ierr);
158e2d7f03fSJed Brown   *list = PETSC_NULL;
159e2d7f03fSJed Brown   PetscFunctionReturn(0);
160e2d7f03fSJed Brown }
161e2d7f03fSJed Brown 
162e2d7f03fSJed Brown #undef __FUNCT__
163d8588912SDave May #define __FUNCT__ "MatDestroy_Nest"
164207556f9SJed Brown static PetscErrorCode MatDestroy_Nest(Mat A)
165d8588912SDave May {
166d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
167d8588912SDave May   PetscInt       i,j;
168d8588912SDave May   PetscErrorCode ierr;
169d8588912SDave May 
170d8588912SDave May   PetscFunctionBegin;
171d8588912SDave May   /* release the matrices and the place holders */
172e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr);
173e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr);
174e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
175e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
176d8588912SDave May 
177d8588912SDave May   ierr = PetscFree(vs->row_len);CHKERRQ(ierr);
178d8588912SDave May   ierr = PetscFree(vs->col_len);CHKERRQ(ierr);
179d8588912SDave May 
180207556f9SJed Brown   ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr);
181207556f9SJed Brown 
182d8588912SDave May   /* release the matrices and the place holders */
183d8588912SDave May   if (vs->m) {
184d8588912SDave May     for (i=0; i<vs->nr; i++) {
185d8588912SDave May       for (j=0; j<vs->nc; j++) {
1866bf464f9SBarry Smith         ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr);
187d8588912SDave May       }
188d8588912SDave May       ierr = PetscFree( vs->m[i] );CHKERRQ(ierr);
189d8588912SDave May     }
190d8588912SDave May     ierr = PetscFree(vs->m);CHKERRQ(ierr);
191d8588912SDave May   }
192bf0cc555SLisandro Dalcin   ierr = PetscFree(A->data);CHKERRQ(ierr);
193d8588912SDave May 
194c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C",   "",0);CHKERRQ(ierr);
1950782ca92SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMat_C",   "",0);CHKERRQ(ierr);
196c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C",  "",0);CHKERRQ(ierr);
197c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C",     "",0);CHKERRQ(ierr);
198c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C",  "",0);CHKERRQ(ierr);
199c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C",   "",0);CHKERRQ(ierr);
200d8588912SDave May   PetscFunctionReturn(0);
201d8588912SDave May }
202d8588912SDave May 
203d8588912SDave May #undef __FUNCT__
204d8588912SDave May #define __FUNCT__ "MatAssemblyBegin_Nest"
205207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
206d8588912SDave May {
207d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
208d8588912SDave May   PetscInt       i,j;
209d8588912SDave May   PetscErrorCode ierr;
210d8588912SDave May 
211d8588912SDave May   PetscFunctionBegin;
212d8588912SDave May   for (i=0; i<vs->nr; i++) {
213d8588912SDave May     for (j=0; j<vs->nc; j++) {
214d8588912SDave May       if (vs->m[i][j]) { ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); }
215d8588912SDave May     }
216d8588912SDave May   }
217d8588912SDave May   PetscFunctionReturn(0);
218d8588912SDave May }
219d8588912SDave May 
220d8588912SDave May #undef __FUNCT__
221d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest"
222207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
223d8588912SDave May {
224d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
225d8588912SDave May   PetscInt       i,j;
226d8588912SDave May   PetscErrorCode ierr;
227d8588912SDave May 
228d8588912SDave May   PetscFunctionBegin;
229d8588912SDave May   for (i=0; i<vs->nr; i++) {
230d8588912SDave May     for (j=0; j<vs->nc; j++) {
231d8588912SDave May       if (vs->m[i][j]) { ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); }
232d8588912SDave May     }
233d8588912SDave May   }
234d8588912SDave May   PetscFunctionReturn(0);
235d8588912SDave May }
236d8588912SDave May 
237d8588912SDave May #undef __FUNCT__
238f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatRow"
239f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
240d8588912SDave May {
241207556f9SJed Brown   PetscErrorCode ierr;
242f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
243f349c1fdSJed Brown   PetscInt       j;
244f349c1fdSJed Brown   Mat            sub;
245d8588912SDave May 
246d8588912SDave May   PetscFunctionBegin;
2478188e55aSJed Brown   sub = (row < vs->nc) ? vs->m[row][row] : PETSC_NULL; /* Prefer to find on the diagonal */
248f349c1fdSJed Brown   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
2498188e55aSJed Brown   if (sub) {ierr = MatPreallocated(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */
250f349c1fdSJed Brown   *B = sub;
251f349c1fdSJed Brown   PetscFunctionReturn(0);
252d8588912SDave May }
253d8588912SDave May 
254f349c1fdSJed Brown #undef __FUNCT__
255f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatCol"
256f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
257f349c1fdSJed Brown {
258207556f9SJed Brown   PetscErrorCode ierr;
259f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
260f349c1fdSJed Brown   PetscInt       i;
261f349c1fdSJed Brown   Mat            sub;
262f349c1fdSJed Brown 
263f349c1fdSJed Brown   PetscFunctionBegin;
2648188e55aSJed Brown   sub = (col < vs->nr) ? vs->m[col][col] : PETSC_NULL; /* Prefer to find on the diagonal */
265f349c1fdSJed Brown   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
2668188e55aSJed Brown   if (sub) {ierr = MatPreallocated(sub);CHKERRQ(ierr);} /* Ensure that the sizes are available */
267f349c1fdSJed Brown   *B = sub;
268f349c1fdSJed Brown   PetscFunctionReturn(0);
269d8588912SDave May }
270d8588912SDave May 
271f349c1fdSJed Brown #undef __FUNCT__
272f349c1fdSJed Brown #define __FUNCT__ "MatNestFindIS"
273f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
274f349c1fdSJed Brown {
275f349c1fdSJed Brown   PetscErrorCode ierr;
276f349c1fdSJed Brown   PetscInt       i;
277f349c1fdSJed Brown   PetscBool      flg;
278f349c1fdSJed Brown 
279f349c1fdSJed Brown   PetscFunctionBegin;
280f349c1fdSJed Brown   PetscValidPointer(list,3);
281f349c1fdSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,4);
282f349c1fdSJed Brown   PetscValidIntPointer(found,5);
283f349c1fdSJed Brown   *found = -1;
284f349c1fdSJed Brown   for (i=0; i<n; i++) {
285207556f9SJed Brown     if (!list[i]) continue;
286f349c1fdSJed Brown     ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr);
287f349c1fdSJed Brown     if (flg) {
288f349c1fdSJed Brown       *found = i;
289f349c1fdSJed Brown       PetscFunctionReturn(0);
290f349c1fdSJed Brown     }
291f349c1fdSJed Brown   }
292f349c1fdSJed Brown   SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Could not find index set");
293f349c1fdSJed Brown   PetscFunctionReturn(0);
294f349c1fdSJed Brown }
295f349c1fdSJed Brown 
296f349c1fdSJed Brown #undef __FUNCT__
2978188e55aSJed Brown #define __FUNCT__ "MatNestGetRow"
2988188e55aSJed Brown /* Get a block row as a new MatNest */
2998188e55aSJed Brown static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B)
3008188e55aSJed Brown {
3018188e55aSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
3028188e55aSJed Brown   char           keyname[256];
3038188e55aSJed Brown   PetscErrorCode ierr;
3048188e55aSJed Brown 
3058188e55aSJed Brown   PetscFunctionBegin;
3068188e55aSJed Brown   *B = PETSC_NULL;
3078188e55aSJed Brown   ierr = PetscSNPrintf(keyname,sizeof keyname,"NestRow_%D",row);CHKERRQ(ierr);
3088188e55aSJed Brown   ierr = PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);CHKERRQ(ierr);
3098188e55aSJed Brown   if (*B) PetscFunctionReturn(0);
3108188e55aSJed Brown 
3118188e55aSJed Brown   ierr = MatCreateNest(((PetscObject)A)->comm,1,PETSC_NULL,vs->nc,vs->isglobal.col,vs->m[row],B);CHKERRQ(ierr);
3128188e55aSJed Brown   (*B)->assembled = A->assembled;
3138188e55aSJed Brown   ierr = PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);CHKERRQ(ierr);
3148188e55aSJed Brown   ierr = PetscObjectDereference((PetscObject)*B);CHKERRQ(ierr); /* Leave the only remaining reference in the composition */
3158188e55aSJed Brown   PetscFunctionReturn(0);
3168188e55aSJed Brown }
3178188e55aSJed Brown 
3188188e55aSJed Brown #undef __FUNCT__
319f349c1fdSJed Brown #define __FUNCT__ "MatNestFindSubMat"
320f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
321f349c1fdSJed Brown {
322f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
3238188e55aSJed Brown   PetscErrorCode ierr;
3246b3a5b13SJed Brown   PetscInt       row,col;
3258188e55aSJed Brown   PetscBool      same,isFullCol;
326f349c1fdSJed Brown 
327f349c1fdSJed Brown   PetscFunctionBegin;
3288188e55aSJed Brown   /* Check if full column space. This is a hack */
3298188e55aSJed Brown   isFullCol = PETSC_FALSE;
3308188e55aSJed Brown   ierr = PetscTypeCompare((PetscObject)iscol,ISSTRIDE,&same);CHKERRQ(ierr);
3318188e55aSJed Brown   if (same) {
332*77019fcaSJed Brown     PetscInt n,first,step,i,an,am,afirst,astep;
3338188e55aSJed Brown     ierr = ISStrideGetInfo(iscol,&first,&step);CHKERRQ(ierr);
3348188e55aSJed Brown     ierr = ISGetLocalSize(iscol,&n);CHKERRQ(ierr);
335*77019fcaSJed Brown     isFullCol = PETSC_TRUE;
336*77019fcaSJed Brown     for (i=0,an=0; i<vs->nc; i++) {
337*77019fcaSJed Brown       ierr = ISStrideGetInfo(is->col[i],&afirst,&astep);CHKERRQ(ierr);
338*77019fcaSJed Brown       ierr = ISGetLocalSize(is->col[i],&am);CHKERRQ(ierr);
339*77019fcaSJed Brown       if (afirst != an || astep != step) isFullCol = PETSC_FALSE;
340*77019fcaSJed Brown       an += am;
341*77019fcaSJed Brown     }
342*77019fcaSJed Brown     if (an != n) isFullCol = PETSC_FALSE;
3438188e55aSJed Brown   }
3448188e55aSJed Brown 
3458188e55aSJed Brown   if (isFullCol) {
3468188e55aSJed Brown     PetscInt row;
3478188e55aSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
3488188e55aSJed Brown     ierr = MatNestGetRow(A,row,B);CHKERRQ(ierr);
3498188e55aSJed Brown   } else {
350f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
351f349c1fdSJed Brown     ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr);
352f349c1fdSJed Brown     *B = vs->m[row][col];
3538188e55aSJed Brown   }
354f349c1fdSJed Brown   PetscFunctionReturn(0);
355f349c1fdSJed Brown }
356f349c1fdSJed Brown 
357f349c1fdSJed Brown #undef __FUNCT__
358f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest"
359f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
360f349c1fdSJed Brown {
361f349c1fdSJed Brown   PetscErrorCode ierr;
362f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
363f349c1fdSJed Brown   Mat            sub;
364f349c1fdSJed Brown 
365f349c1fdSJed Brown   PetscFunctionBegin;
366f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr);
367f349c1fdSJed Brown   switch (reuse) {
368f349c1fdSJed Brown   case MAT_INITIAL_MATRIX:
3697874fa86SDave May     if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); }
370f349c1fdSJed Brown     *B = sub;
371f349c1fdSJed Brown     break;
372f349c1fdSJed Brown   case MAT_REUSE_MATRIX:
373f349c1fdSJed Brown     if (sub != *B) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
374f349c1fdSJed Brown     break;
375f349c1fdSJed Brown   case MAT_IGNORE_MATRIX:       /* Nothing to do */
376f349c1fdSJed Brown     break;
377f349c1fdSJed Brown   }
378f349c1fdSJed Brown   PetscFunctionReturn(0);
379f349c1fdSJed Brown }
380f349c1fdSJed Brown 
381f349c1fdSJed Brown #undef __FUNCT__
382f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest"
383f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
384f349c1fdSJed Brown {
385f349c1fdSJed Brown   PetscErrorCode ierr;
386f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
387f349c1fdSJed Brown   Mat            sub;
388f349c1fdSJed Brown 
389f349c1fdSJed Brown   PetscFunctionBegin;
390f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
391f349c1fdSJed Brown   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
392f349c1fdSJed Brown   if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);}
393f349c1fdSJed Brown   *B = sub;
394d8588912SDave May   PetscFunctionReturn(0);
395d8588912SDave May }
396d8588912SDave May 
397d8588912SDave May #undef __FUNCT__
398d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest"
399207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
400d8588912SDave May {
401d8588912SDave May   PetscErrorCode ierr;
402f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
403f349c1fdSJed Brown   Mat            sub;
404d8588912SDave May 
405d8588912SDave May   PetscFunctionBegin;
406f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
407f349c1fdSJed Brown   if (*B != sub) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
408f349c1fdSJed Brown   if (sub) {
409f349c1fdSJed Brown     if (((PetscObject)sub)->refct <= 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
4106bf464f9SBarry Smith     ierr = MatDestroy(B);CHKERRQ(ierr);
411d8588912SDave May   }
412d8588912SDave May   PetscFunctionReturn(0);
413d8588912SDave May }
414d8588912SDave May 
415d8588912SDave May #undef __FUNCT__
4167874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest"
4177874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
4187874fa86SDave May {
4197874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
4207874fa86SDave May   PetscInt       i;
4217874fa86SDave May   PetscErrorCode ierr;
4227874fa86SDave May 
4237874fa86SDave May   PetscFunctionBegin;
4247874fa86SDave May   for (i=0; i<bA->nr; i++) {
425429bac76SJed Brown     Vec bv;
426429bac76SJed Brown     ierr = VecGetSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
4277874fa86SDave May     if (bA->m[i][i]) {
428429bac76SJed Brown       ierr = MatGetDiagonal(bA->m[i][i],bv);CHKERRQ(ierr);
4297874fa86SDave May     } else {
430429bac76SJed Brown       ierr = VecSet(bv,1.0);CHKERRQ(ierr);
4317874fa86SDave May     }
432429bac76SJed Brown     ierr = VecRestoreSubVector(v,bA->isglobal.row[i],&bv);CHKERRQ(ierr);
4337874fa86SDave May   }
4347874fa86SDave May   PetscFunctionReturn(0);
4357874fa86SDave May }
4367874fa86SDave May 
4377874fa86SDave May #undef __FUNCT__
4387874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest"
4397874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
4407874fa86SDave May {
4417874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
442429bac76SJed Brown   Vec            bl,*br;
4437874fa86SDave May   PetscInt       i,j;
4447874fa86SDave May   PetscErrorCode ierr;
4457874fa86SDave May 
4467874fa86SDave May   PetscFunctionBegin;
447429bac76SJed Brown   ierr = PetscMalloc(bA->nc*sizeof(Vec),&br);CHKERRQ(ierr);
448429bac76SJed Brown   for (j=0; j<bA->nc; j++) {ierr = VecGetSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
4497874fa86SDave May   for (i=0; i<bA->nr; i++) {
450429bac76SJed Brown     ierr = VecGetSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
4517874fa86SDave May     for (j=0; j<bA->nc; j++) {
4527874fa86SDave May       if (bA->m[i][j]) {
453429bac76SJed Brown         ierr = MatDiagonalScale(bA->m[i][j],bl,br[j]);CHKERRQ(ierr);
4547874fa86SDave May       }
4557874fa86SDave May     }
456a061e289SJed Brown     ierr = VecRestoreSubVector(l,bA->isglobal.row[i],&bl);CHKERRQ(ierr);
4577874fa86SDave May   }
458429bac76SJed Brown   for (j=0; j<bA->nc; j++) {ierr = VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);CHKERRQ(ierr);}
459429bac76SJed Brown   ierr = PetscFree(br);CHKERRQ(ierr);
4607874fa86SDave May   PetscFunctionReturn(0);
4617874fa86SDave May }
4627874fa86SDave May 
4637874fa86SDave May #undef __FUNCT__
464a061e289SJed Brown #define __FUNCT__ "MatScale_Nest"
465a061e289SJed Brown static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a)
466a061e289SJed Brown {
467a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
468a061e289SJed Brown   PetscInt       i,j;
469a061e289SJed Brown   PetscErrorCode ierr;
470a061e289SJed Brown 
471a061e289SJed Brown   PetscFunctionBegin;
472a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
473a061e289SJed Brown     for (j=0; j<bA->nc; j++) {
474a061e289SJed Brown       if (bA->m[i][j]) {
475a061e289SJed Brown         ierr = MatScale(bA->m[i][j],a);CHKERRQ(ierr);
476a061e289SJed Brown       }
477a061e289SJed Brown     }
478a061e289SJed Brown   }
479a061e289SJed Brown   PetscFunctionReturn(0);
480a061e289SJed Brown }
481a061e289SJed Brown 
482a061e289SJed Brown #undef __FUNCT__
483a061e289SJed Brown #define __FUNCT__ "MatShift_Nest"
484a061e289SJed Brown static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a)
485a061e289SJed Brown {
486a061e289SJed Brown   Mat_Nest       *bA = (Mat_Nest*)A->data;
487a061e289SJed Brown   PetscInt       i;
488a061e289SJed Brown   PetscErrorCode ierr;
489a061e289SJed Brown 
490a061e289SJed Brown   PetscFunctionBegin;
491a061e289SJed Brown   for (i=0; i<bA->nr; i++) {
492a061e289SJed Brown     if (!bA->m[i][i]) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support for shifting an empty diagonal block, insert a matrix in block (%D,%D)",i,i);
493a061e289SJed Brown     ierr = MatShift(bA->m[i][i],a);CHKERRQ(ierr);
494a061e289SJed Brown   }
495a061e289SJed Brown   PetscFunctionReturn(0);
496a061e289SJed Brown }
497a061e289SJed Brown 
498a061e289SJed Brown #undef __FUNCT__
499d8588912SDave May #define __FUNCT__ "MatGetVecs_Nest"
500207556f9SJed Brown static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left)
501d8588912SDave May {
502d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
503d8588912SDave May   Vec            *L,*R;
504d8588912SDave May   MPI_Comm       comm;
505d8588912SDave May   PetscInt       i,j;
506d8588912SDave May   PetscErrorCode ierr;
507d8588912SDave May 
508d8588912SDave May   PetscFunctionBegin;
509d8588912SDave May   comm = ((PetscObject)A)->comm;
510d8588912SDave May   if (right) {
511d8588912SDave May     /* allocate R */
512d8588912SDave May     ierr = PetscMalloc( sizeof(Vec) * bA->nc, &R );CHKERRQ(ierr);
513d8588912SDave May     /* Create the right vectors */
514d8588912SDave May     for (j=0; j<bA->nc; j++) {
515d8588912SDave May       for (i=0; i<bA->nr; i++) {
516d8588912SDave May         if (bA->m[i][j]) {
517d8588912SDave May           ierr = MatGetVecs(bA->m[i][j],&R[j],PETSC_NULL);CHKERRQ(ierr);
518d8588912SDave May           break;
519d8588912SDave May         }
520d8588912SDave May       }
521d8588912SDave May       if (i==bA->nr) {
522d8588912SDave May         /* have an empty column */
523d8588912SDave May         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
524d8588912SDave May       }
525d8588912SDave May     }
526f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
527d8588912SDave May     /* hand back control to the nest vector */
528d8588912SDave May     for (j=0; j<bA->nc; j++) {
5296bf464f9SBarry Smith       ierr = VecDestroy(&R[j]);CHKERRQ(ierr);
530d8588912SDave May     }
531d8588912SDave May     ierr = PetscFree(R);CHKERRQ(ierr);
532d8588912SDave May   }
533d8588912SDave May 
534d8588912SDave May   if (left) {
535d8588912SDave May     /* allocate L */
536d8588912SDave May     ierr = PetscMalloc( sizeof(Vec) * bA->nr, &L );CHKERRQ(ierr);
537d8588912SDave May     /* Create the left vectors */
538d8588912SDave May     for (i=0; i<bA->nr; i++) {
539d8588912SDave May       for (j=0; j<bA->nc; j++) {
540d8588912SDave May         if (bA->m[i][j]) {
541d8588912SDave May           ierr = MatGetVecs(bA->m[i][j],PETSC_NULL,&L[i]);CHKERRQ(ierr);
542d8588912SDave May           break;
543d8588912SDave May         }
544d8588912SDave May       }
545d8588912SDave May       if (j==bA->nc) {
546d8588912SDave May         /* have an empty row */
547d8588912SDave May         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
548d8588912SDave May       }
549d8588912SDave May     }
550d8588912SDave May 
551f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
552d8588912SDave May     for (i=0; i<bA->nr; i++) {
5536bf464f9SBarry Smith       ierr = VecDestroy(&L[i]);CHKERRQ(ierr);
554d8588912SDave May     }
555d8588912SDave May 
556d8588912SDave May     ierr = PetscFree(L);CHKERRQ(ierr);
557d8588912SDave May   }
558d8588912SDave May   PetscFunctionReturn(0);
559d8588912SDave May }
560d8588912SDave May 
561d8588912SDave May #undef __FUNCT__
562d8588912SDave May #define __FUNCT__ "MatView_Nest"
563207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
564d8588912SDave May {
565d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
566d8588912SDave May   PetscBool      isascii;
567d8588912SDave May   PetscInt       i,j;
568d8588912SDave May   PetscErrorCode ierr;
569d8588912SDave May 
570d8588912SDave May   PetscFunctionBegin;
571d8588912SDave May   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
572d8588912SDave May   if (isascii) {
573d8588912SDave May 
574d8588912SDave May     PetscViewerASCIIPrintf(viewer,"Matrix object: \n" );
575d8588912SDave May     PetscViewerASCIIPushTab( viewer );    /* push0 */
576d8588912SDave May     PetscViewerASCIIPrintf( viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc);
577d8588912SDave May 
578d8588912SDave May     PetscViewerASCIIPrintf(viewer,"MatNest structure: \n" );
579d8588912SDave May     for (i=0; i<bA->nr; i++) {
580d8588912SDave May       for (j=0; j<bA->nc; j++) {
581d8588912SDave May         const MatType type;
582d8588912SDave May         const char *name;
583d8588912SDave May         PetscInt NR,NC;
584d8588912SDave May         PetscBool isNest = PETSC_FALSE;
585d8588912SDave May 
586d8588912SDave May         if (!bA->m[i][j]) {
587d8588912SDave May           PetscViewerASCIIPrintf( viewer, "(%D,%D) : PETSC_NULL \n",i,j);
588d8588912SDave May           continue;
589d8588912SDave May         }
590d8588912SDave May         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
591d8588912SDave May         ierr = MatGetType( bA->m[i][j], &type );CHKERRQ(ierr);
592d8588912SDave May         name = ((PetscObject)bA->m[i][j])->prefix;
593d8588912SDave May         ierr = PetscTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
594d8588912SDave May 
595d8588912SDave May         PetscViewerASCIIPrintf( viewer, "(%D,%D) : name=\"%s\", type=%s, rows=%D, cols=%D \n",i,j,name,type,NR,NC);
596d8588912SDave May 
597d8588912SDave May         if (isNest) {
598d8588912SDave May           PetscViewerASCIIPushTab( viewer );  /* push1 */
599d8588912SDave May           ierr = MatView( bA->m[i][j], viewer );CHKERRQ(ierr);
600d8588912SDave May           PetscViewerASCIIPopTab(viewer);    /* pop1 */
601d8588912SDave May         }
602d8588912SDave May       }
603d8588912SDave May     }
604d8588912SDave May     PetscViewerASCIIPopTab(viewer);    /* pop0 */
605d8588912SDave May   }
606d8588912SDave May   PetscFunctionReturn(0);
607d8588912SDave May }
608d8588912SDave May 
609d8588912SDave May #undef __FUNCT__
610d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest"
611207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A)
612d8588912SDave May {
613d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
614d8588912SDave May   PetscInt       i,j;
615d8588912SDave May   PetscErrorCode ierr;
616d8588912SDave May 
617d8588912SDave May   PetscFunctionBegin;
618d8588912SDave May   for (i=0; i<bA->nr; i++) {
619d8588912SDave May     for (j=0; j<bA->nc; j++) {
620d8588912SDave May       if (!bA->m[i][j]) continue;
621d8588912SDave May       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
622d8588912SDave May     }
623d8588912SDave May   }
624d8588912SDave May   PetscFunctionReturn(0);
625d8588912SDave May }
626d8588912SDave May 
627d8588912SDave May #undef __FUNCT__
628d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest"
629207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
630d8588912SDave May {
631d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
632841e96a3SJed Brown   Mat            *b;
633841e96a3SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
634d8588912SDave May   PetscErrorCode ierr;
635d8588912SDave May 
636d8588912SDave May   PetscFunctionBegin;
637841e96a3SJed Brown   ierr = PetscMalloc(nr*nc*sizeof(Mat),&b);CHKERRQ(ierr);
638841e96a3SJed Brown   for (i=0; i<nr; i++) {
639841e96a3SJed Brown     for (j=0; j<nc; j++) {
640841e96a3SJed Brown       if (bA->m[i][j]) {
641841e96a3SJed Brown         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
642841e96a3SJed Brown       } else {
643841e96a3SJed Brown         b[i*nc+j] = PETSC_NULL;
644d8588912SDave May       }
645d8588912SDave May     }
646d8588912SDave May   }
647f349c1fdSJed Brown   ierr = MatCreateNest(((PetscObject)A)->comm,nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
648841e96a3SJed Brown   /* Give the new MatNest exclusive ownership */
649841e96a3SJed Brown   for (i=0; i<nr*nc; i++) {
6506bf464f9SBarry Smith     ierr = MatDestroy(&b[i]);CHKERRQ(ierr);
651d8588912SDave May   }
652d8588912SDave May   ierr = PetscFree(b);CHKERRQ(ierr);
653d8588912SDave May 
654841e96a3SJed Brown   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
655841e96a3SJed Brown   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
656d8588912SDave May   PetscFunctionReturn(0);
657d8588912SDave May }
658d8588912SDave May 
659d8588912SDave May /* nest api */
660d8588912SDave May EXTERN_C_BEGIN
661d8588912SDave May #undef __FUNCT__
662d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest"
663d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
664d8588912SDave May {
665d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
666d8588912SDave May   PetscFunctionBegin;
667d8588912SDave May   if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
668d8588912SDave May   if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
669d8588912SDave May   *mat = bA->m[idxm][jdxm];
670d8588912SDave May   PetscFunctionReturn(0);
671d8588912SDave May }
672d8588912SDave May EXTERN_C_END
673d8588912SDave May 
674d8588912SDave May #undef __FUNCT__
675d8588912SDave May #define __FUNCT__ "MatNestGetSubMat"
676d8588912SDave May /*@C
677d8588912SDave May  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
678d8588912SDave May 
679d8588912SDave May  Not collective
680d8588912SDave May 
681d8588912SDave May  Input Parameters:
682629881c0SJed Brown +   A  - nest matrix
683d8588912SDave May .   idxm - index of the matrix within the nest matrix
684629881c0SJed Brown -   jdxm - index of the matrix within the nest matrix
685d8588912SDave May 
686d8588912SDave May  Output Parameter:
687d8588912SDave May .   sub - matrix at index idxm,jdxm within the nest matrix
688d8588912SDave May 
689d8588912SDave May  Level: developer
690d8588912SDave May 
691d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMats()
692d8588912SDave May @*/
6937087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
694d8588912SDave May {
695699a902aSJed Brown   PetscErrorCode ierr;
696d8588912SDave May 
697d8588912SDave May   PetscFunctionBegin;
698699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
699d8588912SDave May   PetscFunctionReturn(0);
700d8588912SDave May }
701d8588912SDave May 
702d8588912SDave May EXTERN_C_BEGIN
703d8588912SDave May #undef __FUNCT__
7040782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat_Nest"
7050782ca92SJed Brown PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat)
7060782ca92SJed Brown {
7070782ca92SJed Brown   Mat_Nest *bA = (Mat_Nest*)A->data;
7080782ca92SJed Brown   PetscInt m,n,M,N,mi,ni,Mi,Ni;
7090782ca92SJed Brown   PetscErrorCode ierr;
7100782ca92SJed Brown 
7110782ca92SJed Brown   PetscFunctionBegin;
7120782ca92SJed Brown   if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
7130782ca92SJed Brown   if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
7140782ca92SJed Brown   ierr = MatGetLocalSize(mat,&m,&n);CHKERRQ(ierr);
7150782ca92SJed Brown   ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
7160782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.row[idxm],&mi);CHKERRQ(ierr);
7170782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.row[idxm],&Mi);CHKERRQ(ierr);
7180782ca92SJed Brown   ierr = ISGetLocalSize(bA->isglobal.col[jdxm],&ni);CHKERRQ(ierr);
7190782ca92SJed Brown   ierr = ISGetSize(bA->isglobal.col[jdxm],&Ni);CHKERRQ(ierr);
7200782ca92SJed Brown   if (M != Mi || N != Ni) SETERRQ4(((PetscObject)mat)->comm,PETSC_ERR_ARG_INCOMP,"Submatrix dimension (%D,%D) incompatible with nest block (%D,%D)",M,N,Mi,Ni);
7210782ca92SJed Brown   if (m != mi || n != ni) SETERRQ4(((PetscObject)mat)->comm,PETSC_ERR_ARG_INCOMP,"Submatrix local dimension (%D,%D) incompatible with nest block (%D,%D)",m,n,mi,ni);
7220782ca92SJed Brown   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
7230782ca92SJed Brown   ierr = MatDestroy(&bA->m[idxm][jdxm]);CHKERRQ(ierr);
7240782ca92SJed Brown   bA->m[idxm][jdxm] = mat;
7250782ca92SJed Brown   PetscFunctionReturn(0);
7260782ca92SJed Brown }
7270782ca92SJed Brown EXTERN_C_END
7280782ca92SJed Brown 
7290782ca92SJed Brown #undef __FUNCT__
7300782ca92SJed Brown #define __FUNCT__ "MatNestSetSubMat"
7310782ca92SJed Brown /*@C
7320782ca92SJed Brown  MatNestSetSubMat - Set a single submatrix in the nest matrix.
7330782ca92SJed Brown 
7340782ca92SJed Brown  Logically collective on the submatrix communicator
7350782ca92SJed Brown 
7360782ca92SJed Brown  Input Parameters:
7370782ca92SJed Brown +   A  - nest matrix
7380782ca92SJed Brown .   idxm - index of the matrix within the nest matrix
7390782ca92SJed Brown .   jdxm - index of the matrix within the nest matrix
7400782ca92SJed Brown -   sub - matrix at index idxm,jdxm within the nest matrix
7410782ca92SJed Brown 
7420782ca92SJed Brown  Notes:
7430782ca92SJed Brown  The new submatrix must have the same size and communicator as that block of the nest.
7440782ca92SJed Brown 
7450782ca92SJed Brown  This increments the reference count of the submatrix.
7460782ca92SJed Brown 
7470782ca92SJed Brown  Level: developer
7480782ca92SJed Brown 
7490782ca92SJed Brown .seealso: MatNestSetSubMats(), MatNestGetSubMat()
7500782ca92SJed Brown @*/
7510782ca92SJed Brown PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
7520782ca92SJed Brown {
7530782ca92SJed Brown   PetscErrorCode ierr;
7540782ca92SJed Brown 
7550782ca92SJed Brown   PetscFunctionBegin;
7560782ca92SJed Brown   ierr = PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));CHKERRQ(ierr);
7570782ca92SJed Brown   PetscFunctionReturn(0);
7580782ca92SJed Brown }
7590782ca92SJed Brown 
7600782ca92SJed Brown EXTERN_C_BEGIN
7610782ca92SJed Brown #undef __FUNCT__
762d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest"
763d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
764d8588912SDave May {
765d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
766d8588912SDave May   PetscFunctionBegin;
767d8588912SDave May   if (M)   { *M   = bA->nr; }
768d8588912SDave May   if (N)   { *N   = bA->nc; }
769d8588912SDave May   if (mat) { *mat = bA->m;  }
770d8588912SDave May   PetscFunctionReturn(0);
771d8588912SDave May }
772d8588912SDave May EXTERN_C_END
773d8588912SDave May 
774d8588912SDave May #undef __FUNCT__
775d8588912SDave May #define __FUNCT__ "MatNestGetSubMats"
776d8588912SDave May /*@C
777d8588912SDave May  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
778d8588912SDave May 
779d8588912SDave May  Not collective
780d8588912SDave May 
781d8588912SDave May  Input Parameters:
782629881c0SJed Brown .   A  - nest matrix
783d8588912SDave May 
784d8588912SDave May  Output Parameter:
785629881c0SJed Brown +   M - number of rows in the nest matrix
786d8588912SDave May .   N - number of cols in the nest matrix
787629881c0SJed Brown -   mat - 2d array of matrices
788d8588912SDave May 
789d8588912SDave May  Notes:
790d8588912SDave May 
791d8588912SDave May  The user should not free the array mat.
792d8588912SDave May 
793d8588912SDave May  Level: developer
794d8588912SDave May 
795d8588912SDave May .seealso: MatNestGetSize(), MatNestGetSubMat()
796d8588912SDave May @*/
7977087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
798d8588912SDave May {
799699a902aSJed Brown   PetscErrorCode ierr;
800d8588912SDave May 
801d8588912SDave May   PetscFunctionBegin;
802699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
803d8588912SDave May   PetscFunctionReturn(0);
804d8588912SDave May }
805d8588912SDave May 
806d8588912SDave May EXTERN_C_BEGIN
807d8588912SDave May #undef __FUNCT__
808d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest"
8097087cfbeSBarry Smith PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
810d8588912SDave May {
811d8588912SDave May   Mat_Nest  *bA = (Mat_Nest*)A->data;
812d8588912SDave May 
813d8588912SDave May   PetscFunctionBegin;
814d8588912SDave May   if (M) { *M  = bA->nr; }
815d8588912SDave May   if (N) { *N  = bA->nc; }
816d8588912SDave May   PetscFunctionReturn(0);
817d8588912SDave May }
818d8588912SDave May EXTERN_C_END
819d8588912SDave May 
820d8588912SDave May #undef __FUNCT__
821d8588912SDave May #define __FUNCT__ "MatNestGetSize"
822d8588912SDave May /*@C
823d8588912SDave May  MatNestGetSize - Returns the size of the nest matrix.
824d8588912SDave May 
825d8588912SDave May  Not collective
826d8588912SDave May 
827d8588912SDave May  Input Parameters:
828d8588912SDave May .   A  - nest matrix
829d8588912SDave May 
830d8588912SDave May  Output Parameter:
831629881c0SJed Brown +   M - number of rows in the nested mat
832629881c0SJed Brown -   N - number of cols in the nested mat
833d8588912SDave May 
834d8588912SDave May  Notes:
835d8588912SDave May 
836d8588912SDave May  Level: developer
837d8588912SDave May 
838d8588912SDave May .seealso: MatNestGetSubMat(), MatNestGetSubMats()
839d8588912SDave May @*/
8407087cfbeSBarry Smith PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
841d8588912SDave May {
842699a902aSJed Brown   PetscErrorCode ierr;
843d8588912SDave May 
844d8588912SDave May   PetscFunctionBegin;
845699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
846d8588912SDave May   PetscFunctionReturn(0);
847d8588912SDave May }
848d8588912SDave May 
849207556f9SJed Brown EXTERN_C_BEGIN
850207556f9SJed Brown #undef __FUNCT__
851207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest"
8527087cfbeSBarry Smith PetscErrorCode  MatNestSetVecType_Nest(Mat A,const VecType vtype)
853207556f9SJed Brown {
854207556f9SJed Brown   PetscErrorCode ierr;
855207556f9SJed Brown   PetscBool      flg;
856207556f9SJed Brown 
857207556f9SJed Brown   PetscFunctionBegin;
858207556f9SJed Brown   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
859207556f9SJed Brown   /* In reality, this only distinguishes VECNEST and "other" */
860207556f9SJed Brown   A->ops->getvecs       = flg ? MatGetVecs_Nest       : 0;
861207556f9SJed Brown  PetscFunctionReturn(0);
862207556f9SJed Brown }
863207556f9SJed Brown EXTERN_C_END
864207556f9SJed Brown 
865207556f9SJed Brown #undef __FUNCT__
866207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType"
867207556f9SJed Brown /*@C
868207556f9SJed Brown  MatNestSetVecType - Sets the type of Vec returned by MatGetVecs()
869207556f9SJed Brown 
870207556f9SJed Brown  Not collective
871207556f9SJed Brown 
872207556f9SJed Brown  Input Parameters:
873207556f9SJed Brown +  A  - nest matrix
874207556f9SJed Brown -  vtype - type to use for creating vectors
875207556f9SJed Brown 
876207556f9SJed Brown  Notes:
877207556f9SJed Brown 
878207556f9SJed Brown  Level: developer
879207556f9SJed Brown 
880207556f9SJed Brown .seealso: MatGetVecs()
881207556f9SJed Brown @*/
8827087cfbeSBarry Smith PetscErrorCode  MatNestSetVecType(Mat A,const VecType vtype)
883207556f9SJed Brown {
884207556f9SJed Brown   PetscErrorCode ierr;
885207556f9SJed Brown 
886207556f9SJed Brown   PetscFunctionBegin;
887207556f9SJed Brown   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,const VecType),(A,vtype));CHKERRQ(ierr);
888207556f9SJed Brown   PetscFunctionReturn(0);
889207556f9SJed Brown }
890207556f9SJed Brown 
891c8883902SJed Brown EXTERN_C_BEGIN
892d8588912SDave May #undef __FUNCT__
893c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats_Nest"
894c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
895d8588912SDave May {
896c8883902SJed Brown   Mat_Nest       *s = (Mat_Nest*)A->data;
897c8883902SJed Brown   PetscInt       i,j,m,n,M,N;
898d8588912SDave May   PetscErrorCode ierr;
899d8588912SDave May 
900d8588912SDave May   PetscFunctionBegin;
901c8883902SJed Brown   s->nr = nr;
902c8883902SJed Brown   s->nc = nc;
903d8588912SDave May 
904c8883902SJed Brown   /* Create space for submatrices */
905c8883902SJed Brown   ierr = PetscMalloc(sizeof(Mat*)*nr,&s->m);CHKERRQ(ierr);
906c8883902SJed Brown   for (i=0; i<nr; i++) {
907c8883902SJed Brown     ierr = PetscMalloc(sizeof(Mat)*nc,&s->m[i]);CHKERRQ(ierr);
908d8588912SDave May   }
909c8883902SJed Brown   for (i=0; i<nr; i++) {
910c8883902SJed Brown     for (j=0; j<nc; j++) {
911c8883902SJed Brown       s->m[i][j] = a[i*nc+j];
912c8883902SJed Brown       if (a[i*nc+j]) {
913c8883902SJed Brown         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
914d8588912SDave May       }
915d8588912SDave May     }
916d8588912SDave May   }
917d8588912SDave May 
9188188e55aSJed Brown   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
919d8588912SDave May 
920c8883902SJed Brown   ierr = PetscMalloc(sizeof(PetscInt)*nr,&s->row_len);CHKERRQ(ierr);
921c8883902SJed Brown   ierr = PetscMalloc(sizeof(PetscInt)*nc,&s->col_len);CHKERRQ(ierr);
922c8883902SJed Brown   for (i=0; i<nr; i++) s->row_len[i]=-1;
923c8883902SJed Brown   for (j=0; j<nc; j++) s->col_len[j]=-1;
924d8588912SDave May 
9258188e55aSJed Brown   ierr = MatNestGetSizes_Private(A,&m,&n,&M,&N);CHKERRQ(ierr);
926d8588912SDave May 
927c8883902SJed Brown   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
928c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
929c8883902SJed Brown   ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr);
930c8883902SJed Brown   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
931c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
932c8883902SJed Brown   ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr);
933c8883902SJed Brown 
934c8883902SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
935c8883902SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
936c8883902SJed Brown 
937c8883902SJed Brown   ierr = PetscMalloc2(nr,Vec,&s->left,nc,Vec,&s->right);CHKERRQ(ierr);
938c8883902SJed Brown   ierr = PetscMemzero(s->left,nr*sizeof(Vec));CHKERRQ(ierr);
939c8883902SJed Brown   ierr = PetscMemzero(s->right,nc*sizeof(Vec));CHKERRQ(ierr);
940d8588912SDave May   PetscFunctionReturn(0);
941d8588912SDave May }
942c8883902SJed Brown EXTERN_C_END
943d8588912SDave May 
944c8883902SJed Brown #undef __FUNCT__
945c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats"
946c8883902SJed Brown /*@
947c8883902SJed Brown    MatNestSetSubMats - Sets the nested submatrices
948c8883902SJed Brown 
949c8883902SJed Brown    Collective on Mat
950c8883902SJed Brown 
951c8883902SJed Brown    Input Parameter:
952c8883902SJed Brown +  N - nested matrix
953c8883902SJed Brown .  nr - number of nested row blocks
954c8883902SJed Brown .  is_row - index sets for each nested row block, or PETSC_NULL to make contiguous
955c8883902SJed Brown .  nc - number of nested column blocks
956c8883902SJed Brown .  is_col - index sets for each nested column block, or PETSC_NULL to make contiguous
957c8883902SJed Brown -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL
958c8883902SJed Brown 
959c8883902SJed Brown    Level: advanced
960c8883902SJed Brown 
961c8883902SJed Brown .seealso: MatCreateNest(), MATNEST
962c8883902SJed Brown @*/
963c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
964c8883902SJed Brown {
965c8883902SJed Brown   PetscErrorCode ierr;
966c8883902SJed Brown   PetscInt i;
967c8883902SJed Brown 
968c8883902SJed Brown   PetscFunctionBegin;
969c8883902SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
970c8883902SJed Brown   if (nr < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
971c8883902SJed Brown   if (nr && is_row) {
972c8883902SJed Brown     PetscValidPointer(is_row,3);
973c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
974c8883902SJed Brown   }
975c8883902SJed Brown   if (nc < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
9761664e352SJed Brown   if (nc && is_col) {
977c8883902SJed Brown     PetscValidPointer(is_col,5);
978c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
979c8883902SJed Brown   }
980c8883902SJed Brown   if (nr*nc) PetscValidPointer(a,6);
981c8883902SJed 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);
982c8883902SJed Brown   PetscFunctionReturn(0);
983c8883902SJed Brown }
984d8588912SDave May 
985*77019fcaSJed Brown #undef __FUNCT__
986*77019fcaSJed Brown #define __FUNCT__ "MatNestCreateAggregateL2G_Private"
987*77019fcaSJed Brown static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog,ISLocalToGlobalMapping *ltogb)
988*77019fcaSJed Brown {
989*77019fcaSJed Brown   PetscErrorCode ierr;
990*77019fcaSJed Brown   PetscBool flg;
991*77019fcaSJed Brown   PetscInt i,j,m,mi,*ix;
992*77019fcaSJed Brown 
993*77019fcaSJed Brown   PetscFunctionBegin;
994*77019fcaSJed Brown   for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) {
995*77019fcaSJed Brown     if (islocal[i]) {
996*77019fcaSJed Brown       ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
997*77019fcaSJed Brown       flg = PETSC_TRUE;       /* We found a non-trivial entry */
998*77019fcaSJed Brown     } else {
999*77019fcaSJed Brown       ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
1000*77019fcaSJed Brown     }
1001*77019fcaSJed Brown     m += mi;
1002*77019fcaSJed Brown   }
1003*77019fcaSJed Brown   if (flg) {
1004*77019fcaSJed Brown     ierr = PetscMalloc(m*sizeof(*ix),&ix);CHKERRQ(ierr);
1005*77019fcaSJed Brown     for (i=0,n=0; i<n; i++) {
1006*77019fcaSJed Brown       ISLocalToGlobalMapping smap = PETSC_NULL;
1007*77019fcaSJed Brown       VecScatter scat;
1008*77019fcaSJed Brown       IS isreq;
1009*77019fcaSJed Brown       Vec lvec,gvec;
1010*77019fcaSJed Brown       union {PetscScalar scalar; PetscInt integer;} *x;
1011*77019fcaSJed Brown       Mat sub;
1012*77019fcaSJed Brown 
1013*77019fcaSJed Brown       if (sizeof (*x) != sizeof(PetscScalar)) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"No support when scalars smaller than integers");
1014*77019fcaSJed Brown       if (colflg) {
1015*77019fcaSJed Brown         ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1016*77019fcaSJed Brown       } else {
1017*77019fcaSJed Brown         ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
1018*77019fcaSJed Brown       }
1019*77019fcaSJed Brown       if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&smap,PETSC_NULL);CHKERRQ(ierr);}
1020*77019fcaSJed Brown       if (islocal[i]) {
1021*77019fcaSJed Brown         ierr = ISGetSize(islocal[i],&mi);CHKERRQ(ierr);
1022*77019fcaSJed Brown       } else {
1023*77019fcaSJed Brown         ierr = ISGetSize(isglobal[i],&mi);CHKERRQ(ierr);
1024*77019fcaSJed Brown       }
1025*77019fcaSJed Brown       for (j=0; j<mi; j++) ix[m+j] = j;
1026*77019fcaSJed Brown       if (smap) {ierr = ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);CHKERRQ(ierr);}
1027*77019fcaSJed Brown       /*
1028*77019fcaSJed Brown         Now we need to extract the monolithic global indices that correspond to the given split global indices.
1029*77019fcaSJed 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.
1030*77019fcaSJed Brown         The approach here is ugly because it uses VecScatter to move indices.
1031*77019fcaSJed Brown        */
1032*77019fcaSJed Brown       ierr = VecCreateSeq(PETSC_COMM_SELF,mi,&lvec);CHKERRQ(ierr);
1033*77019fcaSJed Brown       ierr = VecCreateMPI(((PetscObject)isglobal[i])->comm,mi,PETSC_DECIDE,&gvec);CHKERRQ(ierr);
1034*77019fcaSJed Brown       ierr = ISCreateGeneral(((PetscObject)isglobal[i])->comm,mi,ix+m,PETSC_COPY_VALUES,&isreq);CHKERRQ(ierr);
1035*77019fcaSJed Brown       ierr = VecScatterCreate(gvec,isreq,lvec,PETSC_NULL,&scat);CHKERRQ(ierr);
1036*77019fcaSJed Brown       ierr = VecGetArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
1037*77019fcaSJed Brown       for (j=0; j<mi; j++) x[j].integer = ix[m+j];
1038*77019fcaSJed Brown       ierr = VecRestoreArray(gvec,(PetscScalar**)&x);CHKERRQ(ierr);
1039*77019fcaSJed Brown       ierr = VecScatterBegin(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1040*77019fcaSJed Brown       ierr = VecScatterEnd(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1041*77019fcaSJed Brown       ierr = VecGetArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
1042*77019fcaSJed Brown       for (j=0; j<mi; j++) ix[m+j] = x[j].integer;
1043*77019fcaSJed Brown       ierr = VecRestoreArray(lvec,(PetscScalar**)&x);CHKERRQ(ierr);
1044*77019fcaSJed Brown       ierr = VecDestroy(&lvec);CHKERRQ(ierr);
1045*77019fcaSJed Brown       ierr = VecDestroy(&gvec);CHKERRQ(ierr);
1046*77019fcaSJed Brown       ierr = ISDestroy(&isreq);CHKERRQ(ierr);
1047*77019fcaSJed Brown       ierr = VecScatterDestroy(&scat);CHKERRQ(ierr);
1048*77019fcaSJed Brown       m += mi;
1049*77019fcaSJed Brown     }
1050*77019fcaSJed Brown     ierr = ISLocalToGlobalMappingCreate(((PetscObject)A)->comm,m,ix,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
1051*77019fcaSJed Brown     *ltogb = PETSC_NULL;
1052*77019fcaSJed Brown   } else {
1053*77019fcaSJed Brown     *ltog = PETSC_NULL;
1054*77019fcaSJed Brown     *ltogb = PETSC_NULL;
1055*77019fcaSJed Brown   }
1056*77019fcaSJed Brown   PetscFunctionReturn(0);
1057*77019fcaSJed Brown }
1058*77019fcaSJed Brown 
1059*77019fcaSJed Brown 
1060d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1061d8588912SDave May /*
1062d8588912SDave May   nprocessors = NP
1063d8588912SDave May   Nest x^T = ( (g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1) )
1064d8588912SDave May        proc 0: => (g_0,h_0,)
1065d8588912SDave May        proc 1: => (g_1,h_1,)
1066d8588912SDave May        ...
1067d8588912SDave May        proc nprocs-1: => (g_NP-1,h_NP-1,)
1068d8588912SDave May 
1069d8588912SDave May             proc 0:                      proc 1:                    proc nprocs-1:
1070d8588912SDave May     is[0] = ( 0,1,2,...,nlocal(g_0)-1 )  ( 0,1,...,nlocal(g_1)-1 )  ( 0,1,...,nlocal(g_NP-1) )
1071d8588912SDave May 
1072d8588912SDave May             proc 0:
1073d8588912SDave May     is[1] = ( nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1 )
1074d8588912SDave May             proc 1:
1075d8588912SDave May     is[1] = ( nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1 )
1076d8588912SDave May 
1077d8588912SDave May             proc NP-1:
1078d8588912SDave May     is[1] = ( nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1 )
1079d8588912SDave May */
1080d8588912SDave May #undef __FUNCT__
1081d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private"
1082841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1083d8588912SDave May {
1084e2d7f03fSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
10858188e55aSJed Brown   PetscInt       i,j,offset,n,nsum,bs;
1086d8588912SDave May   PetscErrorCode ierr;
10872ae74bdbSJed Brown   Mat            sub;
1088d8588912SDave May 
1089d8588912SDave May   PetscFunctionBegin;
10908188e55aSJed Brown   ierr = PetscMalloc(sizeof(IS)*nr,&vs->isglobal.row);CHKERRQ(ierr);
10918188e55aSJed Brown   ierr = PetscMalloc(sizeof(IS)*nc,&vs->isglobal.col);CHKERRQ(ierr);
1092d8588912SDave May   if (is_row) { /* valid IS is passed in */
1093d8588912SDave May     /* refs on is[] are incremeneted */
1094e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1095d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
1096e2d7f03fSJed Brown       vs->isglobal.row[i] = is_row[i];
1097d8588912SDave May     }
10982ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
10998188e55aSJed Brown     nsum = 0;
11008188e55aSJed Brown     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
11018188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
11028188e55aSJed Brown       if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
11038188e55aSJed Brown       ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr);
11048188e55aSJed Brown       if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
11058188e55aSJed Brown       nsum += n;
11068188e55aSJed Brown     }
11078188e55aSJed Brown     offset = 0;
11083fda65a0SJose Roman #if defined(PETSC_HAVE_MPI_EXSCAN)
11098188e55aSJed Brown     ierr = MPI_Exscan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
11103fda65a0SJose Roman #else
11113fda65a0SJose Roman     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sorry but this code requires MPI_EXSCAN that doesn't exist on your machine's version of MPI, install a MPI2 with PETSc to get this functionality");
11123fda65a0SJose Roman #endif
1113e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1114f349c1fdSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
11152ae74bdbSJed Brown       ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr);
11162ae74bdbSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1117e2d7f03fSJed Brown       ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
1118e2d7f03fSJed Brown       ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
11192ae74bdbSJed Brown       offset += n;
1120d8588912SDave May     }
1121d8588912SDave May   }
1122d8588912SDave May 
1123d8588912SDave May   if (is_col) { /* valid IS is passed in */
1124d8588912SDave May     /* refs on is[] are incremeneted */
1125e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1126d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
1127e2d7f03fSJed Brown       vs->isglobal.col[j] = is_col[j];
1128d8588912SDave May     }
11292ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
11302ae74bdbSJed Brown     offset = A->cmap->rstart;
11318188e55aSJed Brown     nsum = 0;
11328188e55aSJed Brown     for (j=0; j<vs->nc; j++) {
11338188e55aSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
11348188e55aSJed Brown       if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
11358188e55aSJed Brown       ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr);
11368188e55aSJed Brown       if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
11378188e55aSJed Brown       nsum += n;
11388188e55aSJed Brown     }
11398188e55aSJed Brown     offset = 0;
1140d9a7f1c0SMark F. Adams #if defined(PETSC_HAVE_MPI_EXSCAN)
11418188e55aSJed Brown     ierr = MPI_Exscan(&nsum,&offset,1,MPIU_INT,MPI_SUM,((PetscObject)A)->comm);CHKERRQ(ierr);
1142d9a7f1c0SMark F. Adams #else
1143d9a7f1c0SMark F. Adams     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sorry but this code requires MPI_EXSCAN that doesn't exist on your machine's version of MPI, install a MPI2 with PETSc to get this functionality");
1144d9a7f1c0SMark F. Adams #endif
1145e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1146f349c1fdSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
11472ae74bdbSJed Brown       ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr);
11482ae74bdbSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1149e2d7f03fSJed Brown       ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
1150e2d7f03fSJed Brown       ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
11512ae74bdbSJed Brown       offset += n;
1152d8588912SDave May     }
1153d8588912SDave May   }
1154e2d7f03fSJed Brown 
1155e2d7f03fSJed Brown   /* Set up the local ISs */
1156e2d7f03fSJed Brown   ierr = PetscMalloc(vs->nr*sizeof(IS),&vs->islocal.row);CHKERRQ(ierr);
1157e2d7f03fSJed Brown   ierr = PetscMalloc(vs->nc*sizeof(IS),&vs->islocal.col);CHKERRQ(ierr);
1158e2d7f03fSJed Brown   for (i=0,offset=0; i<vs->nr; i++) {
1159e2d7f03fSJed Brown     IS                     isloc;
11608188e55aSJed Brown     ISLocalToGlobalMapping rmap = PETSC_NULL;
1161e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1162e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
11638188e55aSJed Brown     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,&rmap,PETSC_NULL);CHKERRQ(ierr);}
1164207556f9SJed Brown     if (rmap) {
1165e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1166e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
1167e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1168e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1169207556f9SJed Brown     } else {
1170207556f9SJed Brown       nlocal = 0;
1171207556f9SJed Brown       isloc  = PETSC_NULL;
1172207556f9SJed Brown     }
1173e2d7f03fSJed Brown     vs->islocal.row[i] = isloc;
1174e2d7f03fSJed Brown     offset += nlocal;
1175e2d7f03fSJed Brown   }
11768188e55aSJed Brown   for (i=0,offset=0; i<vs->nc; i++) {
1177e2d7f03fSJed Brown     IS                     isloc;
11788188e55aSJed Brown     ISLocalToGlobalMapping cmap = PETSC_NULL;
1179e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1180e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
11818188e55aSJed Brown     if (sub) {ierr = MatGetLocalToGlobalMapping(sub,PETSC_NULL,&cmap);CHKERRQ(ierr);}
1182207556f9SJed Brown     if (cmap) {
1183e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1184e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr);
1185e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1186e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1187207556f9SJed Brown     } else {
1188207556f9SJed Brown       nlocal = 0;
1189207556f9SJed Brown       isloc  = PETSC_NULL;
1190207556f9SJed Brown     }
1191e2d7f03fSJed Brown     vs->islocal.col[i] = isloc;
1192e2d7f03fSJed Brown     offset += nlocal;
1193e2d7f03fSJed Brown   }
11940189643fSJed Brown 
1195*77019fcaSJed Brown   /* Set up the aggregate ISLocalToGlobalMapping */
1196*77019fcaSJed Brown   {
1197*77019fcaSJed Brown     ISLocalToGlobalMapping rmap,rmapb,cmap,cmapb;
1198*77019fcaSJed Brown     ierr = MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap,&rmapb);CHKERRQ(ierr);
1199*77019fcaSJed Brown     ierr = MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap,&cmapb);CHKERRQ(ierr);
1200*77019fcaSJed Brown     if (rmap && cmap) {ierr = MatSetLocalToGlobalMapping(A,rmap,cmap);CHKERRQ(ierr);}
1201*77019fcaSJed Brown     if (rmapb && cmapb) {ierr = MatSetLocalToGlobalMappingBlock(A,rmapb,cmapb);CHKERRQ(ierr);}
1202*77019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&rmap);CHKERRQ(ierr);
1203*77019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&rmapb);CHKERRQ(ierr);
1204*77019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&cmap);CHKERRQ(ierr);
1205*77019fcaSJed Brown     ierr = ISLocalToGlobalMappingDestroy(&cmapb);CHKERRQ(ierr);
1206*77019fcaSJed Brown   }
1207*77019fcaSJed Brown 
12080189643fSJed Brown #if defined(PETSC_USE_DEBUG)
12090189643fSJed Brown   for (i=0; i<vs->nr; i++) {
12100189643fSJed Brown     for (j=0; j<vs->nc; j++) {
12110189643fSJed Brown       PetscInt m,n,M,N,mi,ni,Mi,Ni;
12120189643fSJed Brown       Mat B = vs->m[i][j];
12130189643fSJed Brown       if (!B) continue;
12140189643fSJed Brown       ierr = MatGetSize(B,&M,&N);CHKERRQ(ierr);
12150189643fSJed Brown       ierr = MatGetLocalSize(B,&m,&n);CHKERRQ(ierr);
12160189643fSJed Brown       ierr = ISGetSize(vs->isglobal.row[i],&Mi);CHKERRQ(ierr);
12170189643fSJed Brown       ierr = ISGetSize(vs->isglobal.col[j],&Ni);CHKERRQ(ierr);
12180189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.row[i],&mi);CHKERRQ(ierr);
12190189643fSJed Brown       ierr = ISGetLocalSize(vs->isglobal.col[j],&ni);CHKERRQ(ierr);
12200189643fSJed Brown       if (M != Mi || N != Ni) SETERRQ6(((PetscObject)sub)->comm,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);
12210189643fSJed Brown       if (m != mi || n != ni) SETERRQ6(((PetscObject)sub)->comm,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);
12220189643fSJed Brown     }
12230189643fSJed Brown   }
12240189643fSJed Brown #endif
1225a061e289SJed Brown 
1226a061e289SJed Brown   /* Set A->assembled if all non-null blocks are currently assembled */
1227a061e289SJed Brown   for (i=0; i<vs->nr; i++) {
1228a061e289SJed Brown     for (j=0; j<vs->nc; j++) {
1229a061e289SJed Brown       if (vs->m[i][j] && !vs->m[i][j]->assembled) PetscFunctionReturn(0);
1230a061e289SJed Brown     }
1231a061e289SJed Brown   }
1232a061e289SJed Brown   A->assembled = PETSC_TRUE;
1233d8588912SDave May   PetscFunctionReturn(0);
1234d8588912SDave May }
1235d8588912SDave May 
1236d8588912SDave May #undef __FUNCT__
1237d8588912SDave May #define __FUNCT__ "MatCreateNest"
1238659c6bb0SJed Brown /*@
1239659c6bb0SJed Brown    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
1240659c6bb0SJed Brown 
1241659c6bb0SJed Brown    Collective on Mat
1242659c6bb0SJed Brown 
1243659c6bb0SJed Brown    Input Parameter:
1244659c6bb0SJed Brown +  comm - Communicator for the new Mat
1245659c6bb0SJed Brown .  nr - number of nested row blocks
1246659c6bb0SJed Brown .  is_row - index sets for each nested row block, or PETSC_NULL to make contiguous
1247659c6bb0SJed Brown .  nc - number of nested column blocks
1248659c6bb0SJed Brown .  is_col - index sets for each nested column block, or PETSC_NULL to make contiguous
1249659c6bb0SJed Brown -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL
1250659c6bb0SJed Brown 
1251659c6bb0SJed Brown    Output Parameter:
1252659c6bb0SJed Brown .  B - new matrix
1253659c6bb0SJed Brown 
1254659c6bb0SJed Brown    Level: advanced
1255659c6bb0SJed Brown 
1256659c6bb0SJed Brown .seealso: MatCreate(), VecCreateNest(), DMGetMatrix(), MATNEST
1257659c6bb0SJed Brown @*/
12587087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1259d8588912SDave May {
1260d8588912SDave May   Mat            A;
1261d8588912SDave May   PetscErrorCode ierr;
1262d8588912SDave May 
1263d8588912SDave May   PetscFunctionBegin;
1264c8883902SJed Brown   *B = 0;
1265d8588912SDave May   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1266c8883902SJed Brown   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
1267c8883902SJed Brown   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1268d8588912SDave May   *B = A;
1269d8588912SDave May   PetscFunctionReturn(0);
1270d8588912SDave May }
1271659c6bb0SJed Brown 
1272659c6bb0SJed Brown /*MC
1273659c6bb0SJed Brown   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
1274659c6bb0SJed Brown 
1275659c6bb0SJed Brown   Level: intermediate
1276659c6bb0SJed Brown 
1277659c6bb0SJed Brown   Notes:
1278659c6bb0SJed Brown   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1279659c6bb0SJed Brown   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1280659c6bb0SJed Brown   It is usually used with DMComposite and DMGetMatrix()
1281659c6bb0SJed Brown 
1282659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest()
1283659c6bb0SJed Brown M*/
1284c8883902SJed Brown EXTERN_C_BEGIN
1285c8883902SJed Brown #undef __FUNCT__
1286c8883902SJed Brown #define __FUNCT__ "MatCreate_Nest"
1287c8883902SJed Brown PetscErrorCode MatCreate_Nest(Mat A)
1288c8883902SJed Brown {
1289c8883902SJed Brown   Mat_Nest       *s;
1290c8883902SJed Brown   PetscErrorCode ierr;
1291c8883902SJed Brown 
1292c8883902SJed Brown   PetscFunctionBegin;
1293c8883902SJed Brown   ierr = PetscNewLog(A,Mat_Nest,&s);CHKERRQ(ierr);
1294c8883902SJed Brown   A->data         = (void*)s;
1295c8883902SJed Brown   s->nr = s->nc   = -1;
1296c8883902SJed Brown   s->m            = PETSC_NULL;
1297c8883902SJed Brown 
1298c8883902SJed Brown   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
1299c8883902SJed Brown   A->ops->mult                  = MatMult_Nest;
13009194d70fSJed Brown   A->ops->multadd               = MatMultAdd_Nest;
1301c8883902SJed Brown   A->ops->multtranspose         = MatMultTranspose_Nest;
13029194d70fSJed Brown   A->ops->multtransposeadd      = MatMultTransposeAdd_Nest;
1303c8883902SJed Brown   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1304c8883902SJed Brown   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1305c8883902SJed Brown   A->ops->zeroentries           = MatZeroEntries_Nest;
1306c8883902SJed Brown   A->ops->duplicate             = MatDuplicate_Nest;
1307c8883902SJed Brown   A->ops->getsubmatrix          = MatGetSubMatrix_Nest;
1308c8883902SJed Brown   A->ops->destroy               = MatDestroy_Nest;
1309c8883902SJed Brown   A->ops->view                  = MatView_Nest;
1310c8883902SJed Brown   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1311c8883902SJed Brown   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1312c8883902SJed Brown   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
1313429bac76SJed Brown   A->ops->getdiagonal           = MatGetDiagonal_Nest;
1314429bac76SJed Brown   A->ops->diagonalscale         = MatDiagonalScale_Nest;
1315a061e289SJed Brown   A->ops->scale                 = MatScale_Nest;
1316a061e289SJed Brown   A->ops->shift                 = MatShift_Nest;
1317c8883902SJed Brown 
1318c8883902SJed Brown   A->spptr        = 0;
1319c8883902SJed Brown   A->same_nonzero = PETSC_FALSE;
1320c8883902SJed Brown   A->assembled    = PETSC_FALSE;
1321c8883902SJed Brown 
1322c8883902SJed Brown   /* expose Nest api's */
1323c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C",   "MatNestGetSubMat_Nest",   MatNestGetSubMat_Nest);CHKERRQ(ierr);
13240782ca92SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMat_C",   "MatNestSetSubMat_Nest",   MatNestSetSubMat_Nest);CHKERRQ(ierr);
1325c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C",  "MatNestGetSubMats_Nest",  MatNestGetSubMats_Nest);CHKERRQ(ierr);
1326c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C",     "MatNestGetSize_Nest",     MatNestGetSize_Nest);CHKERRQ(ierr);
1327c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C",  "MatNestSetVecType_Nest",  MatNestSetVecType_Nest);CHKERRQ(ierr);
1328c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C",  "MatNestSetSubMats_Nest",  MatNestSetSubMats_Nest);CHKERRQ(ierr);
1329c8883902SJed Brown 
1330c8883902SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
1331c8883902SJed Brown   PetscFunctionReturn(0);
1332c8883902SJed Brown }
133386e80854SHong Zhang EXTERN_C_END
1334