xref: /petsc/src/mat/impls/nest/matnest.c (revision 86e80854f0b5efa69ecb4ae3900d5e62e4aea4a4)
1d8588912SDave May #define PETSCMAT_DLL
2d8588912SDave May 
3ec9811eeSJed Brown #include "matnestimpl.h" /*I   "petscmat.h"   I*/
4d8588912SDave May 
5c8883902SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]);
6c8883902SJed Brown 
7d8588912SDave May /* private functions */
8d8588912SDave May #undef __FUNCT__
9c8883902SJed Brown #define __FUNCT__ "MatNestGetSize_Private"
10c8883902SJed Brown static PetscErrorCode MatNestGetSize_Private(Mat A,PetscBool globalsize,PetscInt *M,PetscInt *N)
11d8588912SDave May {
12d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
13b9dfa011SHong Zhang   PetscInt       sizeM,sizeN,i,j,index = -1;
14d8588912SDave May   PetscErrorCode ierr;
15d8588912SDave May 
16d8588912SDave May   PetscFunctionBegin;
17d8588912SDave May   /* rows */
18d8588912SDave May   for (i=0; i<bA->nr; i++) {
19d8588912SDave May     for (j=0; j<bA->nc; j++) {
20d8588912SDave May       if (bA->m[i][j]) { index = j; break; }
21d8588912SDave May     }
22d8588912SDave May     if (bA->m[i][index]) {
23d8588912SDave May       if (globalsize) { ierr = MatGetSize(bA->m[i][index],&sizeM,&sizeN);CHKERRQ(ierr); }
24d8588912SDave May       else {            ierr = MatGetLocalSize(bA->m[i][index],&sizeM,&sizeN);CHKERRQ(ierr); }
25d8588912SDave May       *M = *M + sizeM;
26d8588912SDave May     }
27d8588912SDave May   }
28d8588912SDave May 
29d8588912SDave May   /* cols */
30d8588912SDave May   /* now descend recursively */
31d8588912SDave May   for (j=0; j<bA->nc; j++) {
32d8588912SDave May     for (i=0; i<bA->nr; i++) {
33d8588912SDave May       if (bA->m[i][j]) { index = i; break; }
34d8588912SDave May     }
35d8588912SDave May     if (bA->m[index][j]) {
36d8588912SDave May       if (globalsize) { ierr = MatGetSize(bA->m[index][j],&sizeM,&sizeN);CHKERRQ(ierr); }
37d8588912SDave May       else {            ierr = MatGetLocalSize(bA->m[index][j],&sizeM,&sizeN);CHKERRQ(ierr); }
38d8588912SDave May       *N = *N + sizeN;
39d8588912SDave May     }
40d8588912SDave May   }
41d8588912SDave May   PetscFunctionReturn(0);
42d8588912SDave May }
43d8588912SDave May 
44d8588912SDave May #undef __FUNCT__
45d8588912SDave May #define __FUNCT__ "PETSc_MatNest_CheckMatVecCompatibility2"
46207556f9SJed Brown PETSC_UNUSED
47d8588912SDave May static PetscErrorCode PETSc_MatNest_CheckMatVecCompatibility2(Mat A,Vec x,Vec y)
48d8588912SDave May {
49d8588912SDave May   PetscBool      isANest, isxNest, isyNest;
50d8588912SDave May   PetscErrorCode ierr;
51d8588912SDave May 
52d8588912SDave May   PetscFunctionBegin;
53d8588912SDave May   isANest = isxNest = isyNest = PETSC_FALSE;
54d8588912SDave May   ierr = PetscTypeCompare( (PetscObject)A, MATNEST, &isANest );CHKERRQ(ierr);
55d8588912SDave May   ierr = PetscTypeCompare( (PetscObject)x, VECNEST, &isxNest );CHKERRQ(ierr);
56d8588912SDave May   ierr = PetscTypeCompare( (PetscObject)y, VECNEST, &isyNest );CHKERRQ(ierr);
57d8588912SDave May 
58d8588912SDave May   if (isANest && isxNest && isyNest) {
59d8588912SDave May     PetscFunctionReturn(0);
60d8588912SDave May   } else {
61d8588912SDave May     SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_SUP, "Operation only valid for Mat(Nest)-Vec(Nest) combinations");
62d8588912SDave May     PetscFunctionReturn(0);
63d8588912SDave May   }
64d8588912SDave May   PetscFunctionReturn(0);
65d8588912SDave May }
66d8588912SDave May 
67d8588912SDave May /*
68d8588912SDave May  This function is called every time we insert a sub matrix the Nest.
69d8588912SDave May  It ensures that every Nest along row r, and col c has the same dimensions
70d8588912SDave May  as the submat being inserted.
71d8588912SDave May */
72d8588912SDave May #undef __FUNCT__
73d8588912SDave May #define __FUNCT__ "PETSc_MatNest_CheckConsistency"
74063a28d4SJed Brown PETSC_UNUSED
75d8588912SDave May static PetscErrorCode PETSc_MatNest_CheckConsistency(Mat A,Mat submat,PetscInt r,PetscInt c)
76d8588912SDave May {
77d8588912SDave May   Mat_Nest       *b = (Mat_Nest*)A->data;
78d8588912SDave May   PetscInt       i,j;
79d8588912SDave May   PetscInt       nr,nc;
80d8588912SDave May   PetscInt       sM,sN,M,N;
81d8588912SDave May   Mat            mat;
82d8588912SDave May   PetscErrorCode ierr;
83d8588912SDave May 
84d8588912SDave May   PetscFunctionBegin;
85d8588912SDave May   nr = b->nr;
86d8588912SDave May   nc = b->nc;
87d8588912SDave May   ierr = MatGetSize(submat,&sM,&sN);CHKERRQ(ierr);
88d8588912SDave May   for (i=0; i<nr; i++) {
89d8588912SDave May     mat = b->m[i][c];
90d8588912SDave May     if (mat) {
91d8588912SDave May       ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
92d8588912SDave May       /* Check columns match */
93d8588912SDave May       if (sN != N) {
94d8588912SDave May         SETERRQ3(((PetscObject)A)->comm,PETSC_ERR_SUP,"Inserted incompatible submatrix into Nest at (%D,%D). Submatrix must have %D rows",r,c,N );
95d8588912SDave May       }
96d8588912SDave May     }
97d8588912SDave May   }
98d8588912SDave May 
99d8588912SDave May   for (j=0; j<nc; j++) {
100d8588912SDave May     mat = b->m[r][j];
101d8588912SDave May     if (mat) {
102d8588912SDave May       ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
103d8588912SDave May       /* Check rows match */
104d8588912SDave May       if (sM != M) {
105d8588912SDave May         SETERRQ3(((PetscObject)A)->comm,PETSC_ERR_SUP,"Inserted incompatible submatrix into Nest at (%D,%D). Submatrix must have %D cols",r,c,M );
106d8588912SDave May       }
107d8588912SDave May     }
108d8588912SDave May   }
109d8588912SDave May   PetscFunctionReturn(0);
110d8588912SDave May }
111d8588912SDave May 
112d8588912SDave May /*
113d8588912SDave May  Checks if the row/col's contain a complete set of IS's.
114d8588912SDave May  Once they are filled, the offsets are computed. This saves having to update
115d8588912SDave May  the index set entries each time we insert something new.
116d8588912SDave May  */
117d8588912SDave May #undef __FUNCT__
118d8588912SDave May #define __FUNCT__ "PETSc_MatNest_UpdateStructure"
119063a28d4SJed Brown PETSC_UNUSED
120d8588912SDave May static PetscErrorCode PETSc_MatNest_UpdateStructure(Mat A,PetscInt ridx,PetscInt cidx)
121d8588912SDave May {
122d8588912SDave May   Mat_Nest       *b = (Mat_Nest*)A->data;
123d8588912SDave May   PetscInt       m,n,i,j;
124d8588912SDave May   PetscBool      fullrow,fullcol;
125d8588912SDave May   PetscErrorCode ierr;
126d8588912SDave May 
127d8588912SDave May   PetscFunctionBegin;
128d8588912SDave May   ierr = MatGetLocalSize(b->m[ridx][cidx],&m,&n);CHKERRQ(ierr);
129d8588912SDave May   b->row_len[ridx] = m;
130d8588912SDave May   b->col_len[cidx] = n;
131d8588912SDave May 
132d8588912SDave May   fullrow = PETSC_TRUE;
133d8588912SDave May   for (i=0; i<b->nr; i++) {
134d8588912SDave May     if (b->row_len[i]==-1) { fullrow = PETSC_FALSE; break; }
135d8588912SDave May   }
136f349c1fdSJed Brown   if ( (fullrow) && (!b->isglobal.row[0]) ){
137d8588912SDave May     PetscInt cnt = 0;
138d8588912SDave May     for (i=0; i<b->nr; i++) {
139f349c1fdSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,b->row_len[i],cnt,1,&b->isglobal.row[i]);CHKERRQ(ierr);
140d8588912SDave May       cnt = cnt + b->row_len[i];
141d8588912SDave May     }
142d8588912SDave May   }
143d8588912SDave May 
144d8588912SDave May   fullcol = PETSC_TRUE;
145d8588912SDave May   for (j=0; j<b->nc; j++) {
146d8588912SDave May     if (b->col_len[j]==-1) { fullcol = PETSC_FALSE; break; }
147d8588912SDave May   }
148f349c1fdSJed Brown   if( (fullcol) && (!b->isglobal.col[0]) ){
149d8588912SDave May     PetscInt cnt = 0;
150d8588912SDave May     for (j=0; j<b->nc; j++) {
151f349c1fdSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,b->col_len[j],cnt,1,&b->isglobal.col[j]);CHKERRQ(ierr);
152d8588912SDave May       cnt = cnt + b->col_len[j];
153d8588912SDave May     }
154d8588912SDave May   }
155d8588912SDave May   PetscFunctionReturn(0);
156d8588912SDave May }
157d8588912SDave May 
158d8588912SDave May /* operations */
159d8588912SDave May #undef __FUNCT__
160d8588912SDave May #define __FUNCT__ "MatMult_Nest"
161207556f9SJed Brown static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y)
162d8588912SDave May {
163d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
164207556f9SJed Brown   Vec            *bx = bA->right,*by = bA->left;
165207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
166d8588912SDave May   PetscErrorCode ierr;
167d8588912SDave May 
168d8588912SDave May   PetscFunctionBegin;
169207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
170207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
171207556f9SJed Brown   for (i=0; i<nr; i++) {
172d8588912SDave May     ierr = VecZeroEntries(by[i]);CHKERRQ(ierr);
173207556f9SJed Brown     for (j=0; j<nc; j++) {
174207556f9SJed Brown       if (!bA->m[i][j]) continue;
175d8588912SDave May       /* y[i] <- y[i] + A[i][j] * x[j] */
176d8588912SDave May       ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr);
177d8588912SDave May     }
178d8588912SDave May   }
179207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);CHKERRQ(ierr);}
180207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);CHKERRQ(ierr);}
181d8588912SDave May   PetscFunctionReturn(0);
182d8588912SDave May }
183d8588912SDave May 
184d8588912SDave May #undef __FUNCT__
185d8588912SDave May #define __FUNCT__ "MatMultTranspose_Nest"
186207556f9SJed Brown static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y)
187d8588912SDave May {
188d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
189207556f9SJed Brown   Vec            *bx = bA->left,*by = bA->right;
190207556f9SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
191d8588912SDave May   PetscErrorCode ierr;
192d8588912SDave May 
193d8588912SDave May   PetscFunctionBegin;
194d8588912SDave May   if (A->symmetric) {
195d8588912SDave May     ierr = MatMult_Nest(A,x,y);CHKERRQ(ierr);
196d8588912SDave May     PetscFunctionReturn(0);
197d8588912SDave May   }
198207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecGetSubVector(y,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
199207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecGetSubVector(x,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
200207556f9SJed Brown   for (i=0; i<nr; i++) {
201d8588912SDave May     ierr = VecZeroEntries(by[i]);CHKERRQ(ierr);
202207556f9SJed Brown     for (j=0; j<nc; j++) {
203207556f9SJed Brown       if (!bA->m[j][i]) continue;
204d8588912SDave May       /* y[i] <- y[i] + A^T[i][j] * x[j], so we swap i,j in mat[][] */
205d8588912SDave May       ierr = MatMultTransposeAdd(bA->m[j][i],bx[j],by[i],by[i]);CHKERRQ(ierr);
206d8588912SDave May     }
207d8588912SDave May   }
208207556f9SJed Brown   for (i=0; i<nr; i++) {ierr = VecRestoreSubVector(y,bA->isglobal.row[i],&bx[i]);CHKERRQ(ierr);}
209207556f9SJed Brown   for (i=0; i<nc; i++) {ierr = VecRestoreSubVector(x,bA->isglobal.col[i],&by[i]);CHKERRQ(ierr);}
210d8588912SDave May   PetscFunctionReturn(0);
211d8588912SDave May }
212d8588912SDave May 
213d8588912SDave May #undef __FUNCT__
214e2d7f03fSJed Brown #define __FUNCT__ "MatNestDestroyISList"
215e2d7f03fSJed Brown static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list)
216e2d7f03fSJed Brown {
217e2d7f03fSJed Brown   PetscErrorCode ierr;
218e2d7f03fSJed Brown   IS             *lst = *list;
219e2d7f03fSJed Brown   PetscInt       i;
220e2d7f03fSJed Brown 
221e2d7f03fSJed Brown   PetscFunctionBegin;
222e2d7f03fSJed Brown   if (!lst) PetscFunctionReturn(0);
223e2d7f03fSJed Brown   for (i=0; i<n; i++) if (lst[i]) {ierr = ISDestroy(lst[i]);CHKERRQ(ierr);}
224e2d7f03fSJed Brown   ierr = PetscFree(lst);CHKERRQ(ierr);
225e2d7f03fSJed Brown   *list = PETSC_NULL;
226e2d7f03fSJed Brown   PetscFunctionReturn(0);
227e2d7f03fSJed Brown }
228e2d7f03fSJed Brown 
229e2d7f03fSJed Brown #undef __FUNCT__
230d8588912SDave May #define __FUNCT__ "MatDestroy_Nest"
231207556f9SJed Brown static PetscErrorCode MatDestroy_Nest(Mat A)
232d8588912SDave May {
233d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
234d8588912SDave May   PetscInt       i,j;
235d8588912SDave May   PetscErrorCode ierr;
236d8588912SDave May 
237d8588912SDave May   PetscFunctionBegin;
238d8588912SDave May   /* release the matrices and the place holders */
239e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->isglobal.row);CHKERRQ(ierr);
240e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->isglobal.col);CHKERRQ(ierr);
241e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nr,&vs->islocal.row);CHKERRQ(ierr);
242e2d7f03fSJed Brown   ierr = MatNestDestroyISList(vs->nc,&vs->islocal.col);CHKERRQ(ierr);
243d8588912SDave May 
244d8588912SDave May   ierr = PetscFree(vs->row_len);CHKERRQ(ierr);
245d8588912SDave May   ierr = PetscFree(vs->col_len);CHKERRQ(ierr);
246d8588912SDave May 
247207556f9SJed Brown   ierr = PetscFree2(vs->left,vs->right);CHKERRQ(ierr);
248207556f9SJed Brown 
249d8588912SDave May   /* release the matrices and the place holders */
250d8588912SDave May   if (vs->m) {
251d8588912SDave May     for (i=0; i<vs->nr; i++) {
252d8588912SDave May       for (j=0; j<vs->nc; j++) {
253d8588912SDave May 
254d8588912SDave May         if (vs->m[i][j]) {
255d8588912SDave May           ierr = MatDestroy(vs->m[i][j]);CHKERRQ(ierr);
256d8588912SDave May           vs->m[i][j] = PETSC_NULL;
257d8588912SDave May         }
258d8588912SDave May       }
259d8588912SDave May       ierr = PetscFree( vs->m[i] );CHKERRQ(ierr);
260d8588912SDave May       vs->m[i] = PETSC_NULL;
261d8588912SDave May     }
262d8588912SDave May     ierr = PetscFree(vs->m);CHKERRQ(ierr);
263d8588912SDave May     vs->m = PETSC_NULL;
264d8588912SDave May   }
265d8588912SDave May   ierr = PetscFree(vs);CHKERRQ(ierr);
266d8588912SDave May 
267c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C",   "",0);CHKERRQ(ierr);
268c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C",  "",0);CHKERRQ(ierr);
269c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C",     "",0);CHKERRQ(ierr);
270c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C",  "",0);CHKERRQ(ierr);
271c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C",   "",0);CHKERRQ(ierr);
272d8588912SDave May   PetscFunctionReturn(0);
273d8588912SDave May }
274d8588912SDave May 
275d8588912SDave May #undef __FUNCT__
276d8588912SDave May #define __FUNCT__ "MatAssemblyBegin_Nest"
277207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
278d8588912SDave May {
279d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
280d8588912SDave May   PetscInt       i,j;
281d8588912SDave May   PetscErrorCode ierr;
282d8588912SDave May 
283d8588912SDave May   PetscFunctionBegin;
284d8588912SDave May   for (i=0; i<vs->nr; i++) {
285d8588912SDave May     for (j=0; j<vs->nc; j++) {
286d8588912SDave May       if (vs->m[i][j]) { ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); }
287d8588912SDave May     }
288d8588912SDave May   }
289d8588912SDave May   PetscFunctionReturn(0);
290d8588912SDave May }
291d8588912SDave May 
292d8588912SDave May #undef __FUNCT__
293d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest"
294207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
295d8588912SDave May {
296d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
297d8588912SDave May   PetscInt       i,j;
298d8588912SDave May   PetscErrorCode ierr;
299d8588912SDave May 
300d8588912SDave May   PetscFunctionBegin;
301d8588912SDave May   for (i=0; i<vs->nr; i++) {
302d8588912SDave May     for (j=0; j<vs->nc; j++) {
303d8588912SDave May       if (vs->m[i][j]) { ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); }
304d8588912SDave May     }
305d8588912SDave May   }
306d8588912SDave May   PetscFunctionReturn(0);
307d8588912SDave May }
308d8588912SDave May 
309d8588912SDave May #undef __FUNCT__
310f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatRow"
311f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
312d8588912SDave May {
313207556f9SJed Brown   PetscErrorCode ierr;
314f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
315f349c1fdSJed Brown   PetscInt       j;
316f349c1fdSJed Brown   Mat            sub;
317d8588912SDave May 
318d8588912SDave May   PetscFunctionBegin;
319f349c1fdSJed Brown   sub = vs->m[row][row];        /* Prefer to find on the diagonal */
320f349c1fdSJed Brown   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
321f349c1fdSJed Brown   if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",row);
322207556f9SJed Brown   ierr = MatPreallocated(sub);CHKERRQ(ierr); /* Ensure that the sizes are available */
323f349c1fdSJed Brown   *B = sub;
324f349c1fdSJed Brown   PetscFunctionReturn(0);
325d8588912SDave May }
326d8588912SDave May 
327f349c1fdSJed Brown #undef __FUNCT__
328f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatCol"
329f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
330f349c1fdSJed Brown {
331207556f9SJed Brown   PetscErrorCode ierr;
332f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
333f349c1fdSJed Brown   PetscInt       i;
334f349c1fdSJed Brown   Mat            sub;
335f349c1fdSJed Brown 
336f349c1fdSJed Brown   PetscFunctionBegin;
337f349c1fdSJed Brown   sub = vs->m[col][col];        /* Prefer to find on the diagonal */
338f349c1fdSJed Brown   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
339f349c1fdSJed Brown   if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",col);
340207556f9SJed Brown   ierr = MatPreallocated(sub);CHKERRQ(ierr); /* Ensure that the sizes are available */
341f349c1fdSJed Brown   *B = sub;
342f349c1fdSJed Brown   PetscFunctionReturn(0);
343d8588912SDave May }
344d8588912SDave May 
345f349c1fdSJed Brown #undef __FUNCT__
346f349c1fdSJed Brown #define __FUNCT__ "MatNestFindIS"
347f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
348f349c1fdSJed Brown {
349f349c1fdSJed Brown   PetscErrorCode ierr;
350f349c1fdSJed Brown   PetscInt       i;
351f349c1fdSJed Brown   PetscBool      flg;
352f349c1fdSJed Brown 
353f349c1fdSJed Brown   PetscFunctionBegin;
354f349c1fdSJed Brown   PetscValidPointer(list,3);
355f349c1fdSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,4);
356f349c1fdSJed Brown   PetscValidIntPointer(found,5);
357f349c1fdSJed Brown   *found = -1;
358f349c1fdSJed Brown   for (i=0; i<n; i++) {
359207556f9SJed Brown     if (!list[i]) continue;
360f349c1fdSJed Brown     ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr);
361f349c1fdSJed Brown     if (flg) {
362f349c1fdSJed Brown       *found = i;
363f349c1fdSJed Brown       PetscFunctionReturn(0);
364f349c1fdSJed Brown     }
365f349c1fdSJed Brown   }
366f349c1fdSJed Brown   SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Could not find index set");
367f349c1fdSJed Brown   PetscFunctionReturn(0);
368f349c1fdSJed Brown }
369f349c1fdSJed Brown 
370f349c1fdSJed Brown #undef __FUNCT__
371f349c1fdSJed Brown #define __FUNCT__ "MatNestFindSubMat"
372f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
373f349c1fdSJed Brown {
374f349c1fdSJed Brown   PetscErrorCode ierr;
375f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
376f349c1fdSJed Brown   PetscInt       row,col;
377f349c1fdSJed Brown 
378f349c1fdSJed Brown   PetscFunctionBegin;
379f349c1fdSJed Brown   ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
380f349c1fdSJed Brown   ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr);
381f349c1fdSJed Brown   *B = vs->m[row][col];
382f349c1fdSJed Brown   PetscFunctionReturn(0);
383f349c1fdSJed Brown }
384f349c1fdSJed Brown 
385f349c1fdSJed Brown #undef __FUNCT__
386f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest"
387f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
388f349c1fdSJed Brown {
389f349c1fdSJed Brown   PetscErrorCode ierr;
390f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
391f349c1fdSJed Brown   Mat            sub;
392f349c1fdSJed Brown 
393f349c1fdSJed Brown   PetscFunctionBegin;
394f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr);
395f349c1fdSJed Brown   switch (reuse) {
396f349c1fdSJed Brown   case MAT_INITIAL_MATRIX:
397f349c1fdSJed Brown     ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);
398f349c1fdSJed Brown     *B = sub;
399f349c1fdSJed Brown     break;
400f349c1fdSJed Brown   case MAT_REUSE_MATRIX:
401f349c1fdSJed Brown     if (sub != *B) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
402f349c1fdSJed Brown     break;
403f349c1fdSJed Brown   case MAT_IGNORE_MATRIX:       /* Nothing to do */
404f349c1fdSJed Brown     break;
405f349c1fdSJed Brown   }
406f349c1fdSJed Brown   PetscFunctionReturn(0);
407f349c1fdSJed Brown }
408f349c1fdSJed Brown 
409f349c1fdSJed Brown #undef __FUNCT__
410f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest"
411f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
412f349c1fdSJed Brown {
413f349c1fdSJed Brown   PetscErrorCode ierr;
414f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
415f349c1fdSJed Brown   Mat            sub;
416f349c1fdSJed Brown 
417f349c1fdSJed Brown   PetscFunctionBegin;
418f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
419f349c1fdSJed Brown   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
420f349c1fdSJed Brown   if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);}
421f349c1fdSJed Brown   *B = sub;
422d8588912SDave May   PetscFunctionReturn(0);
423d8588912SDave May }
424d8588912SDave May 
425d8588912SDave May #undef __FUNCT__
426d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest"
427207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
428d8588912SDave May {
429d8588912SDave May   PetscErrorCode ierr;
430f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
431f349c1fdSJed Brown   Mat            sub;
432d8588912SDave May 
433d8588912SDave May   PetscFunctionBegin;
434f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
435f349c1fdSJed Brown   if (*B != sub) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
436f349c1fdSJed Brown   if (sub) {
437f349c1fdSJed Brown     if (((PetscObject)sub)->refct <= 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
438f349c1fdSJed Brown     ierr = MatDestroy(*B);CHKERRQ(ierr);
439d8588912SDave May   }
440f349c1fdSJed Brown   *B = PETSC_NULL;
441d8588912SDave May   PetscFunctionReturn(0);
442d8588912SDave May }
443d8588912SDave May 
444d8588912SDave May #undef __FUNCT__
445d8588912SDave May #define __FUNCT__ "MatGetVecs_Nest"
446207556f9SJed Brown static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left)
447d8588912SDave May {
448d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
449d8588912SDave May   Vec            *L,*R;
450d8588912SDave May   MPI_Comm       comm;
451d8588912SDave May   PetscInt       i,j;
452d8588912SDave May   PetscErrorCode ierr;
453d8588912SDave May 
454d8588912SDave May   PetscFunctionBegin;
455d8588912SDave May   comm = ((PetscObject)A)->comm;
456d8588912SDave May   if (right) {
457d8588912SDave May     /* allocate R */
458d8588912SDave May     ierr = PetscMalloc( sizeof(Vec) * bA->nc, &R );CHKERRQ(ierr);
459d8588912SDave May     /* Create the right vectors */
460d8588912SDave May     for (j=0; j<bA->nc; j++) {
461d8588912SDave May       for (i=0; i<bA->nr; i++) {
462d8588912SDave May         if (bA->m[i][j]) {
463d8588912SDave May           ierr = MatGetVecs(bA->m[i][j],&R[j],PETSC_NULL);CHKERRQ(ierr);
464d8588912SDave May           break;
465d8588912SDave May         }
466d8588912SDave May       }
467d8588912SDave May       if (i==bA->nr) {
468d8588912SDave May         /* have an empty column */
469d8588912SDave May         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
470d8588912SDave May       }
471d8588912SDave May     }
472f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
473d8588912SDave May     /* hand back control to the nest vector */
474d8588912SDave May     for (j=0; j<bA->nc; j++) {
475d8588912SDave May       ierr = VecDestroy(R[j]);CHKERRQ(ierr);
476d8588912SDave May     }
477d8588912SDave May     ierr = PetscFree(R);CHKERRQ(ierr);
478d8588912SDave May   }
479d8588912SDave May 
480d8588912SDave May   if (left) {
481d8588912SDave May     /* allocate L */
482d8588912SDave May     ierr = PetscMalloc( sizeof(Vec) * bA->nr, &L );CHKERRQ(ierr);
483d8588912SDave May     /* Create the left vectors */
484d8588912SDave May     for (i=0; i<bA->nr; i++) {
485d8588912SDave May       for (j=0; j<bA->nc; j++) {
486d8588912SDave May         if (bA->m[i][j]) {
487d8588912SDave May           ierr = MatGetVecs(bA->m[i][j],PETSC_NULL,&L[i]);CHKERRQ(ierr);
488d8588912SDave May           break;
489d8588912SDave May         }
490d8588912SDave May       }
491d8588912SDave May       if (j==bA->nc) {
492d8588912SDave May         /* have an empty row */
493d8588912SDave May         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
494d8588912SDave May       }
495d8588912SDave May     }
496d8588912SDave May 
497f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
498d8588912SDave May     for (i=0; i<bA->nr; i++) {
499d8588912SDave May       ierr = VecDestroy(L[i]);CHKERRQ(ierr);
500d8588912SDave May     }
501d8588912SDave May 
502d8588912SDave May     ierr = PetscFree(L);CHKERRQ(ierr);
503d8588912SDave May   }
504d8588912SDave May   PetscFunctionReturn(0);
505d8588912SDave May }
506d8588912SDave May 
507d8588912SDave May #undef __FUNCT__
508d8588912SDave May #define __FUNCT__ "MatView_Nest"
509207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
510d8588912SDave May {
511d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
512d8588912SDave May   PetscBool      isascii;
513d8588912SDave May   PetscInt       i,j;
514d8588912SDave May   PetscErrorCode ierr;
515d8588912SDave May 
516d8588912SDave May   PetscFunctionBegin;
517d8588912SDave May   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
518d8588912SDave May   if (isascii) {
519d8588912SDave May 
520d8588912SDave May     PetscViewerASCIIPrintf(viewer,"Matrix object: \n" );
521d8588912SDave May     PetscViewerASCIIPushTab( viewer );    /* push0 */
522d8588912SDave May     PetscViewerASCIIPrintf( viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc);
523d8588912SDave May 
524d8588912SDave May     PetscViewerASCIIPrintf(viewer,"MatNest structure: \n" );
525d8588912SDave May     for (i=0; i<bA->nr; i++) {
526d8588912SDave May       for (j=0; j<bA->nc; j++) {
527d8588912SDave May         const MatType type;
528d8588912SDave May         const char *name;
529d8588912SDave May         PetscInt NR,NC;
530d8588912SDave May         PetscBool isNest = PETSC_FALSE;
531d8588912SDave May 
532d8588912SDave May         if (!bA->m[i][j]) {
533d8588912SDave May           PetscViewerASCIIPrintf( viewer, "(%D,%D) : PETSC_NULL \n",i,j);
534d8588912SDave May           continue;
535d8588912SDave May         }
536d8588912SDave May         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
537d8588912SDave May         ierr = MatGetType( bA->m[i][j], &type );CHKERRQ(ierr);
538d8588912SDave May         name = ((PetscObject)bA->m[i][j])->prefix;
539d8588912SDave May         ierr = PetscTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
540d8588912SDave May 
541d8588912SDave May         PetscViewerASCIIPrintf( viewer, "(%D,%D) : name=\"%s\", type=%s, rows=%D, cols=%D \n",i,j,name,type,NR,NC);
542d8588912SDave May 
543d8588912SDave May         if (isNest) {
544d8588912SDave May           PetscViewerASCIIPushTab( viewer );  /* push1 */
545d8588912SDave May           ierr = MatView( bA->m[i][j], viewer );CHKERRQ(ierr);
546d8588912SDave May           PetscViewerASCIIPopTab(viewer);    /* pop1 */
547d8588912SDave May         }
548d8588912SDave May       }
549d8588912SDave May     }
550d8588912SDave May     PetscViewerASCIIPopTab(viewer);    /* pop0 */
551d8588912SDave May   }
552d8588912SDave May   PetscFunctionReturn(0);
553d8588912SDave May }
554d8588912SDave May 
555d8588912SDave May #undef __FUNCT__
556d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest"
557207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A)
558d8588912SDave May {
559d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
560d8588912SDave May   PetscInt       i,j;
561d8588912SDave May   PetscErrorCode ierr;
562d8588912SDave May 
563d8588912SDave May   PetscFunctionBegin;
564d8588912SDave May   for (i=0; i<bA->nr; i++) {
565d8588912SDave May     for (j=0; j<bA->nc; j++) {
566d8588912SDave May       if (!bA->m[i][j]) continue;
567d8588912SDave May       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
568d8588912SDave May     }
569d8588912SDave May   }
570d8588912SDave May   PetscFunctionReturn(0);
571d8588912SDave May }
572d8588912SDave May 
573d8588912SDave May #undef __FUNCT__
574d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest"
575207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
576d8588912SDave May {
577d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
578841e96a3SJed Brown   Mat            *b;
579841e96a3SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
580d8588912SDave May   PetscErrorCode ierr;
581d8588912SDave May 
582d8588912SDave May   PetscFunctionBegin;
583841e96a3SJed Brown   ierr = PetscMalloc(nr*nc*sizeof(Mat),&b);CHKERRQ(ierr);
584841e96a3SJed Brown   for (i=0; i<nr; i++) {
585841e96a3SJed Brown     for (j=0; j<nc; j++) {
586841e96a3SJed Brown       if (bA->m[i][j]) {
587841e96a3SJed Brown         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
588841e96a3SJed Brown       } else {
589841e96a3SJed Brown         b[i*nc+j] = PETSC_NULL;
590d8588912SDave May       }
591d8588912SDave May     }
592d8588912SDave May   }
593f349c1fdSJed Brown   ierr = MatCreateNest(((PetscObject)A)->comm,nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
594841e96a3SJed Brown   /* Give the new MatNest exclusive ownership */
595841e96a3SJed Brown   for (i=0; i<nr*nc; i++) {
596841e96a3SJed Brown     if (b[i]) {ierr = MatDestroy(b[i]);CHKERRQ(ierr);}
597d8588912SDave May   }
598d8588912SDave May   ierr = PetscFree(b);CHKERRQ(ierr);
599d8588912SDave May 
600841e96a3SJed Brown   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
601841e96a3SJed Brown   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
602d8588912SDave May   PetscFunctionReturn(0);
603d8588912SDave May }
604d8588912SDave May 
605d8588912SDave May /* nest api */
606d8588912SDave May EXTERN_C_BEGIN
607d8588912SDave May #undef __FUNCT__
608d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest"
609d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
610d8588912SDave May {
611d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
612d8588912SDave May   PetscFunctionBegin;
613d8588912SDave May   if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
614d8588912SDave May   if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
615d8588912SDave May   *mat = bA->m[idxm][jdxm];
616d8588912SDave May   PetscFunctionReturn(0);
617d8588912SDave May }
618d8588912SDave May EXTERN_C_END
619d8588912SDave May 
620d8588912SDave May #undef __FUNCT__
621d8588912SDave May #define __FUNCT__ "MatNestGetSubMat"
622d8588912SDave May /*@C
623d8588912SDave May  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
624d8588912SDave May 
625d8588912SDave May  Not collective
626d8588912SDave May 
627d8588912SDave May  Input Parameters:
628d8588912SDave May  .  A  - nest matrix
629d8588912SDave May  .  idxm - index of the matrix within the nest matrix
630d8588912SDave May  .  jdxm - index of the matrix within the nest matrix
631d8588912SDave May 
632d8588912SDave May  Output Parameter:
633d8588912SDave May  .  sub - matrix at index idxm,jdxm within the nest matrix
634d8588912SDave May 
635d8588912SDave May  Notes:
636d8588912SDave May 
637d8588912SDave May  Level: developer
638d8588912SDave May 
639d8588912SDave May  .seealso: MatNestGetSize(), MatNestGetSubMats()
640d8588912SDave May @*/
6417087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
642d8588912SDave May {
643699a902aSJed Brown   PetscErrorCode ierr;
644d8588912SDave May 
645d8588912SDave May   PetscFunctionBegin;
646699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
647d8588912SDave May   PetscFunctionReturn(0);
648d8588912SDave May }
649d8588912SDave May 
650d8588912SDave May EXTERN_C_BEGIN
651d8588912SDave May #undef __FUNCT__
652d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest"
653d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
654d8588912SDave May {
655d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
656d8588912SDave May   PetscFunctionBegin;
657d8588912SDave May   if (M)   { *M   = bA->nr; }
658d8588912SDave May   if (N)   { *N   = bA->nc; }
659d8588912SDave May   if (mat) { *mat = bA->m;  }
660d8588912SDave May   PetscFunctionReturn(0);
661d8588912SDave May }
662d8588912SDave May EXTERN_C_END
663d8588912SDave May 
664d8588912SDave May #undef __FUNCT__
665d8588912SDave May #define __FUNCT__ "MatNestGetSubMats"
666d8588912SDave May /*@C
667d8588912SDave May  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
668d8588912SDave May 
669d8588912SDave May  Not collective
670d8588912SDave May 
671d8588912SDave May  Input Parameters:
672d8588912SDave May  .  A  - nest matri
673d8588912SDave May 
674d8588912SDave May  Output Parameter:
675d8588912SDave May  .  M - number of rows in the nest matrix
676d8588912SDave May  .  N - number of cols in the nest matrix
677d8588912SDave May  .  mat - 2d array of matrices
678d8588912SDave May 
679d8588912SDave May  Notes:
680d8588912SDave May 
681d8588912SDave May  The user should not free the array mat.
682d8588912SDave May 
683d8588912SDave May  Level: developer
684d8588912SDave May 
685d8588912SDave May  .seealso: MatNestGetSize(), MatNestGetSubMat()
686d8588912SDave May @*/
6877087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
688d8588912SDave May {
689699a902aSJed Brown   PetscErrorCode ierr;
690d8588912SDave May 
691d8588912SDave May   PetscFunctionBegin;
692699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
693d8588912SDave May   PetscFunctionReturn(0);
694d8588912SDave May }
695d8588912SDave May 
696d8588912SDave May EXTERN_C_BEGIN
697d8588912SDave May #undef __FUNCT__
698d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest"
6997087cfbeSBarry Smith PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
700d8588912SDave May {
701d8588912SDave May   Mat_Nest  *bA = (Mat_Nest*)A->data;
702d8588912SDave May 
703d8588912SDave May   PetscFunctionBegin;
704d8588912SDave May   if (M) { *M  = bA->nr; }
705d8588912SDave May   if (N) { *N  = bA->nc; }
706d8588912SDave May   PetscFunctionReturn(0);
707d8588912SDave May }
708d8588912SDave May EXTERN_C_END
709d8588912SDave May 
710d8588912SDave May #undef __FUNCT__
711d8588912SDave May #define __FUNCT__ "MatNestGetSize"
712d8588912SDave May /*@C
713d8588912SDave May  MatNestGetSize - Returns the size of the nest matrix.
714d8588912SDave May 
715d8588912SDave May  Not collective
716d8588912SDave May 
717d8588912SDave May  Input Parameters:
718d8588912SDave May  .  A  - nest matrix
719d8588912SDave May 
720d8588912SDave May  Output Parameter:
721d8588912SDave May  .  M - number of rows in the nested mat
722d8588912SDave May  .  N - number of cols in the nested mat
723d8588912SDave May 
724d8588912SDave May  Notes:
725d8588912SDave May 
726d8588912SDave May  Level: developer
727d8588912SDave May 
728d8588912SDave May  .seealso: MatNestGetSubMat(), MatNestGetSubMats()
729d8588912SDave May @*/
7307087cfbeSBarry Smith PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
731d8588912SDave May {
732699a902aSJed Brown   PetscErrorCode ierr;
733d8588912SDave May 
734d8588912SDave May   PetscFunctionBegin;
735699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
736d8588912SDave May   PetscFunctionReturn(0);
737d8588912SDave May }
738d8588912SDave May 
739207556f9SJed Brown EXTERN_C_BEGIN
740207556f9SJed Brown #undef __FUNCT__
741207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest"
7427087cfbeSBarry Smith PetscErrorCode  MatNestSetVecType_Nest(Mat A,const VecType vtype)
743207556f9SJed Brown {
744207556f9SJed Brown   PetscErrorCode ierr;
745207556f9SJed Brown   PetscBool      flg;
746207556f9SJed Brown 
747207556f9SJed Brown   PetscFunctionBegin;
748207556f9SJed Brown   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
749207556f9SJed Brown   /* In reality, this only distinguishes VECNEST and "other" */
750207556f9SJed Brown   A->ops->getvecs = flg ? MatGetVecs_Nest : 0;
751207556f9SJed Brown   PetscFunctionReturn(0);
752207556f9SJed Brown }
753207556f9SJed Brown EXTERN_C_END
754207556f9SJed Brown 
755207556f9SJed Brown #undef __FUNCT__
756207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType"
757207556f9SJed Brown /*@C
758207556f9SJed Brown  MatNestSetVecType - Sets the type of Vec returned by MatGetVecs()
759207556f9SJed Brown 
760207556f9SJed Brown  Not collective
761207556f9SJed Brown 
762207556f9SJed Brown  Input Parameters:
763207556f9SJed Brown +  A  - nest matrix
764207556f9SJed Brown -  vtype - type to use for creating vectors
765207556f9SJed Brown 
766207556f9SJed Brown  Notes:
767207556f9SJed Brown 
768207556f9SJed Brown  Level: developer
769207556f9SJed Brown 
770207556f9SJed Brown  .seealso: MatGetVecs()
771207556f9SJed Brown @*/
7727087cfbeSBarry Smith PetscErrorCode  MatNestSetVecType(Mat A,const VecType vtype)
773207556f9SJed Brown {
774207556f9SJed Brown   PetscErrorCode ierr;
775207556f9SJed Brown 
776207556f9SJed Brown   PetscFunctionBegin;
777207556f9SJed Brown   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,const VecType),(A,vtype));CHKERRQ(ierr);
778207556f9SJed Brown   PetscFunctionReturn(0);
779207556f9SJed Brown }
780207556f9SJed Brown 
781c8883902SJed Brown EXTERN_C_BEGIN
782d8588912SDave May #undef __FUNCT__
783c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats_Nest"
784c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
785d8588912SDave May {
786c8883902SJed Brown   Mat_Nest       *s = (Mat_Nest*)A->data;
787c8883902SJed Brown   PetscInt       i,j,m,n,M,N;
788d8588912SDave May   PetscErrorCode ierr;
789d8588912SDave May 
790d8588912SDave May   PetscFunctionBegin;
791c8883902SJed Brown   s->nr = nr;
792c8883902SJed Brown   s->nc = nc;
793d8588912SDave May 
794c8883902SJed Brown   /* Create space for submatrices */
795c8883902SJed Brown   ierr = PetscMalloc(sizeof(Mat*)*nr,&s->m);CHKERRQ(ierr);
796c8883902SJed Brown   for (i=0; i<nr; i++) {
797c8883902SJed Brown     ierr = PetscMalloc(sizeof(Mat)*nc,&s->m[i]);CHKERRQ(ierr);
798d8588912SDave May   }
799c8883902SJed Brown   for (i=0; i<nr; i++) {
800c8883902SJed Brown     for (j=0; j<nc; j++) {
801c8883902SJed Brown       s->m[i][j] = a[i*nc+j];
802c8883902SJed Brown       if (a[i*nc+j]) {
803c8883902SJed Brown         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
804d8588912SDave May       }
805d8588912SDave May     }
806d8588912SDave May   }
807d8588912SDave May 
808c8883902SJed Brown   ierr = PetscMalloc(sizeof(IS)*nr,&s->isglobal.row);CHKERRQ(ierr);
809c8883902SJed Brown   ierr = PetscMalloc(sizeof(IS)*nc,&s->isglobal.col);CHKERRQ(ierr);
810d8588912SDave May 
811c8883902SJed Brown   ierr = PetscMalloc(sizeof(PetscInt)*nr,&s->row_len);CHKERRQ(ierr);
812c8883902SJed Brown   ierr = PetscMalloc(sizeof(PetscInt)*nc,&s->col_len);CHKERRQ(ierr);
813c8883902SJed Brown   for (i=0; i<nr; i++) s->row_len[i]=-1;
814c8883902SJed Brown   for (j=0; j<nc; j++) s->col_len[j]=-1;
815d8588912SDave May 
816c8883902SJed Brown   m = n = 0;
817c8883902SJed Brown   M = N = 0;
818c8883902SJed Brown   ierr = MatNestGetSize_Private(A,PETSC_TRUE,&M,&N);CHKERRQ(ierr);
819c8883902SJed Brown   ierr = MatNestGetSize_Private(A,PETSC_FALSE,&m,&n);CHKERRQ(ierr);
820d8588912SDave May 
821c8883902SJed Brown   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
822c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
823c8883902SJed Brown   ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr);
824c8883902SJed Brown   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
825c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
826c8883902SJed Brown   ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr);
827c8883902SJed Brown 
828c8883902SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
829c8883902SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
830c8883902SJed Brown 
831c8883902SJed Brown   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
832c8883902SJed Brown   ierr = PetscMalloc2(nr,Vec,&s->left,nc,Vec,&s->right);CHKERRQ(ierr);
833c8883902SJed Brown   ierr = PetscMemzero(s->left,nr*sizeof(Vec));CHKERRQ(ierr);
834c8883902SJed Brown   ierr = PetscMemzero(s->right,nc*sizeof(Vec));CHKERRQ(ierr);
835d8588912SDave May   PetscFunctionReturn(0);
836d8588912SDave May }
837c8883902SJed Brown EXTERN_C_END
838d8588912SDave May 
839c8883902SJed Brown #undef __FUNCT__
840c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats"
841c8883902SJed Brown /*@
842c8883902SJed Brown    MatNestSetSubMats - Sets the nested submatrices
843c8883902SJed Brown 
844c8883902SJed Brown    Collective on Mat
845c8883902SJed Brown 
846c8883902SJed Brown    Input Parameter:
847c8883902SJed Brown +  N - nested matrix
848c8883902SJed Brown .  nr - number of nested row blocks
849c8883902SJed Brown .  is_row - index sets for each nested row block, or PETSC_NULL to make contiguous
850c8883902SJed Brown .  nc - number of nested column blocks
851c8883902SJed Brown .  is_col - index sets for each nested column block, or PETSC_NULL to make contiguous
852c8883902SJed Brown -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL
853c8883902SJed Brown 
854c8883902SJed Brown    Level: advanced
855c8883902SJed Brown 
856c8883902SJed Brown .seealso: MatCreateNest(), MATNEST
857c8883902SJed Brown @*/
858c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
859c8883902SJed Brown {
860c8883902SJed Brown   PetscErrorCode ierr;
861c8883902SJed Brown   PetscInt i;
862c8883902SJed Brown 
863c8883902SJed Brown   PetscFunctionBegin;
864c8883902SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
865c8883902SJed Brown   if (nr < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
866c8883902SJed Brown   if (nr && is_row) {
867c8883902SJed Brown     PetscValidPointer(is_row,3);
868c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
869c8883902SJed Brown   }
870c8883902SJed Brown   if (nc < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
871c8883902SJed Brown   if (nc && is_row) {
872c8883902SJed Brown     PetscValidPointer(is_col,5);
873c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
874c8883902SJed Brown   }
875c8883902SJed Brown   if (nr*nc) PetscValidPointer(a,6);
876c8883902SJed 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);
877c8883902SJed Brown   PetscFunctionReturn(0);
878c8883902SJed Brown }
879d8588912SDave May 
880d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
881d8588912SDave May /*
882d8588912SDave May   nprocessors = NP
883d8588912SDave May   Nest x^T = ( (g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1) )
884d8588912SDave May        proc 0: => (g_0,h_0,)
885d8588912SDave May        proc 1: => (g_1,h_1,)
886d8588912SDave May        ...
887d8588912SDave May        proc nprocs-1: => (g_NP-1,h_NP-1,)
888d8588912SDave May 
889d8588912SDave May             proc 0:                      proc 1:                    proc nprocs-1:
890d8588912SDave May     is[0] = ( 0,1,2,...,nlocal(g_0)-1 )  ( 0,1,...,nlocal(g_1)-1 )  ( 0,1,...,nlocal(g_NP-1) )
891d8588912SDave May 
892d8588912SDave May             proc 0:
893d8588912SDave May     is[1] = ( nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1 )
894d8588912SDave May             proc 1:
895d8588912SDave May     is[1] = ( nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1 )
896d8588912SDave May 
897d8588912SDave May             proc NP-1:
898d8588912SDave May     is[1] = ( nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1 )
899d8588912SDave May */
900d8588912SDave May #undef __FUNCT__
901d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private"
902841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
903d8588912SDave May {
904e2d7f03fSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
9052ae74bdbSJed Brown   PetscInt       i,j,offset,n,bs;
906d8588912SDave May   PetscErrorCode ierr;
9072ae74bdbSJed Brown   Mat            sub;
908d8588912SDave May 
909d8588912SDave May   PetscFunctionBegin;
910d8588912SDave May   if (is_row) { /* valid IS is passed in */
911d8588912SDave May     /* refs on is[] are incremeneted */
912e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
913d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
914e2d7f03fSJed Brown       vs->isglobal.row[i] = is_row[i];
915d8588912SDave May     }
9162ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
9172ae74bdbSJed Brown     offset = A->rmap->rstart;
918e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
919f349c1fdSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
9202ae74bdbSJed Brown       ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr);
921207556f9SJed Brown       if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
9222ae74bdbSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
923e2d7f03fSJed Brown       ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
924e2d7f03fSJed Brown       ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
9252ae74bdbSJed Brown       offset += n;
926d8588912SDave May     }
927d8588912SDave May   }
928d8588912SDave May 
929d8588912SDave May   if (is_col) { /* valid IS is passed in */
930d8588912SDave May     /* refs on is[] are incremeneted */
931e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
932d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
933e2d7f03fSJed Brown       vs->isglobal.col[j] = is_col[j];
934d8588912SDave May     }
9352ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
9362ae74bdbSJed Brown     offset = A->cmap->rstart;
937e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
938f349c1fdSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
9392ae74bdbSJed Brown       ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr);
940207556f9SJed Brown       if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
9412ae74bdbSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
942e2d7f03fSJed Brown       ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
943e2d7f03fSJed Brown       ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
9442ae74bdbSJed Brown       offset += n;
945d8588912SDave May     }
946d8588912SDave May   }
947e2d7f03fSJed Brown 
948e2d7f03fSJed Brown   /* Set up the local ISs */
949e2d7f03fSJed Brown   ierr = PetscMalloc(vs->nr*sizeof(IS),&vs->islocal.row);CHKERRQ(ierr);
950e2d7f03fSJed Brown   ierr = PetscMalloc(vs->nc*sizeof(IS),&vs->islocal.col);CHKERRQ(ierr);
951e2d7f03fSJed Brown   for (i=0,offset=0; i<vs->nr; i++) {
952e2d7f03fSJed Brown     IS                     isloc;
953e2d7f03fSJed Brown     ISLocalToGlobalMapping rmap;
954e2d7f03fSJed Brown     PetscInt               nlocal,bs;
955e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
956e2d7f03fSJed Brown     ierr = MatGetLocalToGlobalMapping(sub,&rmap,PETSC_NULL);CHKERRQ(ierr);
957207556f9SJed Brown     if (rmap) {
958e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
959e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
960e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
961e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
962207556f9SJed Brown     } else {
963207556f9SJed Brown       nlocal = 0;
964207556f9SJed Brown       isloc  = PETSC_NULL;
965207556f9SJed Brown     }
966e2d7f03fSJed Brown     vs->islocal.row[i] = isloc;
967e2d7f03fSJed Brown     offset += nlocal;
968e2d7f03fSJed Brown   }
969e2d7f03fSJed Brown   for (i=0,offset=0; i<vs->nr; i++) {
970e2d7f03fSJed Brown     IS                     isloc;
971e2d7f03fSJed Brown     ISLocalToGlobalMapping cmap;
972e2d7f03fSJed Brown     PetscInt               nlocal,bs;
973e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
974e2d7f03fSJed Brown     ierr = MatGetLocalToGlobalMapping(sub,PETSC_NULL,&cmap);CHKERRQ(ierr);
975207556f9SJed Brown     if (cmap) {
976e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
977e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr);
978e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
979e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
980207556f9SJed Brown     } else {
981207556f9SJed Brown       nlocal = 0;
982207556f9SJed Brown       isloc  = PETSC_NULL;
983207556f9SJed Brown     }
984e2d7f03fSJed Brown     vs->islocal.col[i] = isloc;
985e2d7f03fSJed Brown     offset += nlocal;
986e2d7f03fSJed Brown   }
987d8588912SDave May   PetscFunctionReturn(0);
988d8588912SDave May }
989d8588912SDave May 
990d8588912SDave May #undef __FUNCT__
991d8588912SDave May #define __FUNCT__ "MatCreateNest"
992659c6bb0SJed Brown /*@
993659c6bb0SJed Brown    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
994659c6bb0SJed Brown 
995659c6bb0SJed Brown    Collective on Mat
996659c6bb0SJed Brown 
997659c6bb0SJed Brown    Input Parameter:
998659c6bb0SJed Brown +  comm - Communicator for the new Mat
999659c6bb0SJed Brown .  nr - number of nested row blocks
1000659c6bb0SJed Brown .  is_row - index sets for each nested row block, or PETSC_NULL to make contiguous
1001659c6bb0SJed Brown .  nc - number of nested column blocks
1002659c6bb0SJed Brown .  is_col - index sets for each nested column block, or PETSC_NULL to make contiguous
1003659c6bb0SJed Brown -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL
1004659c6bb0SJed Brown 
1005659c6bb0SJed Brown    Output Parameter:
1006659c6bb0SJed Brown .  B - new matrix
1007659c6bb0SJed Brown 
1008659c6bb0SJed Brown    Level: advanced
1009659c6bb0SJed Brown 
1010659c6bb0SJed Brown .seealso: MatCreate(), VecCreateNest(), DMGetMatrix(), MATNEST
1011659c6bb0SJed Brown @*/
10127087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1013d8588912SDave May {
1014d8588912SDave May   Mat            A;
1015d8588912SDave May   PetscErrorCode ierr;
1016d8588912SDave May 
1017d8588912SDave May   PetscFunctionBegin;
1018c8883902SJed Brown   *B = 0;
1019d8588912SDave May   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1020c8883902SJed Brown   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
1021c8883902SJed Brown   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1022d8588912SDave May   *B = A;
1023d8588912SDave May   PetscFunctionReturn(0);
1024d8588912SDave May }
1025659c6bb0SJed Brown 
1026659c6bb0SJed Brown /*MC
1027659c6bb0SJed Brown   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
1028659c6bb0SJed Brown 
1029659c6bb0SJed Brown   Level: intermediate
1030659c6bb0SJed Brown 
1031659c6bb0SJed Brown   Notes:
1032659c6bb0SJed Brown   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1033659c6bb0SJed Brown   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1034659c6bb0SJed Brown   It is usually used with DMComposite and DMGetMatrix()
1035659c6bb0SJed Brown 
1036659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest()
1037659c6bb0SJed Brown M*/
1038c8883902SJed Brown EXTERN_C_BEGIN
1039c8883902SJed Brown #undef __FUNCT__
1040c8883902SJed Brown #define __FUNCT__ "MatCreate_Nest"
1041c8883902SJed Brown PetscErrorCode MatCreate_Nest(Mat A)
1042c8883902SJed Brown {
1043c8883902SJed Brown   Mat_Nest       *s;
1044c8883902SJed Brown   PetscErrorCode ierr;
1045c8883902SJed Brown 
1046c8883902SJed Brown   PetscFunctionBegin;
1047c8883902SJed Brown   ierr = PetscNewLog(A,Mat_Nest,&s);CHKERRQ(ierr);
1048c8883902SJed Brown   A->data         = (void*)s;
1049c8883902SJed Brown   s->nr = s->nc   = -1;
1050c8883902SJed Brown   s->m            = PETSC_NULL;
1051c8883902SJed Brown 
1052c8883902SJed Brown   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
1053c8883902SJed Brown   A->ops->mult                  = MatMult_Nest;
1054c8883902SJed Brown   A->ops->multadd               = 0;
1055c8883902SJed Brown   A->ops->multtranspose         = MatMultTranspose_Nest;
1056c8883902SJed Brown   A->ops->multtransposeadd      = 0;
1057c8883902SJed Brown   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1058c8883902SJed Brown   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1059c8883902SJed Brown   A->ops->zeroentries           = MatZeroEntries_Nest;
1060c8883902SJed Brown   A->ops->duplicate             = MatDuplicate_Nest;
1061c8883902SJed Brown   A->ops->getsubmatrix          = MatGetSubMatrix_Nest;
1062c8883902SJed Brown   A->ops->destroy               = MatDestroy_Nest;
1063c8883902SJed Brown   A->ops->view                  = MatView_Nest;
1064c8883902SJed Brown   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1065c8883902SJed Brown   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1066c8883902SJed Brown   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
1067c8883902SJed Brown 
1068c8883902SJed Brown   A->spptr        = 0;
1069c8883902SJed Brown   A->same_nonzero = PETSC_FALSE;
1070c8883902SJed Brown   A->assembled    = PETSC_FALSE;
1071c8883902SJed Brown 
1072c8883902SJed Brown   /* expose Nest api's */
1073c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C",   "MatNestGetSubMat_Nest",   MatNestGetSubMat_Nest);CHKERRQ(ierr);
1074c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C",  "MatNestGetSubMats_Nest",  MatNestGetSubMats_Nest);CHKERRQ(ierr);
1075c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C",     "MatNestGetSize_Nest",     MatNestGetSize_Nest);CHKERRQ(ierr);
1076c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C",  "MatNestSetVecType_Nest",  MatNestSetVecType_Nest);CHKERRQ(ierr);
1077c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C",  "MatNestSetSubMats_Nest",  MatNestSetSubMats_Nest);CHKERRQ(ierr);
1078c8883902SJed Brown 
1079c8883902SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
1080c8883902SJed Brown   PetscFunctionReturn(0);
1081c8883902SJed Brown }
1082*86e80854SHong Zhang EXTERN_C_END
1083