xref: /petsc/src/mat/impls/nest/matnest.c (revision a0769a71603f0d6f68a6e87fcac4dba2599b7624)
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__
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);
2236bf464f9SBarry Smith   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++) {
2536bf464f9SBarry Smith         ierr = MatDestroy(&vs->m[i][j]);CHKERRQ(ierr);
254d8588912SDave May       }
255d8588912SDave May       ierr = PetscFree( vs->m[i] );CHKERRQ(ierr);
256d8588912SDave May     }
257d8588912SDave May     ierr = PetscFree(vs->m);CHKERRQ(ierr);
258d8588912SDave May   }
259d8588912SDave May   ierr = PetscFree(vs);CHKERRQ(ierr);
260d8588912SDave May 
261c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C",   "",0);CHKERRQ(ierr);
262c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C",  "",0);CHKERRQ(ierr);
263c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C",     "",0);CHKERRQ(ierr);
264c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C",  "",0);CHKERRQ(ierr);
265c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C",   "",0);CHKERRQ(ierr);
266d8588912SDave May   PetscFunctionReturn(0);
267d8588912SDave May }
268d8588912SDave May 
269d8588912SDave May #undef __FUNCT__
270d8588912SDave May #define __FUNCT__ "MatAssemblyBegin_Nest"
271207556f9SJed Brown static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
272d8588912SDave May {
273d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
274d8588912SDave May   PetscInt       i,j;
275d8588912SDave May   PetscErrorCode ierr;
276d8588912SDave May 
277d8588912SDave May   PetscFunctionBegin;
278d8588912SDave May   for (i=0; i<vs->nr; i++) {
279d8588912SDave May     for (j=0; j<vs->nc; j++) {
280d8588912SDave May       if (vs->m[i][j]) { ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); }
281d8588912SDave May     }
282d8588912SDave May   }
283d8588912SDave May   PetscFunctionReturn(0);
284d8588912SDave May }
285d8588912SDave May 
286d8588912SDave May #undef __FUNCT__
287d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest"
288207556f9SJed Brown static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
289d8588912SDave May {
290d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
291d8588912SDave May   PetscInt       i,j;
292d8588912SDave May   PetscErrorCode ierr;
293d8588912SDave May 
294d8588912SDave May   PetscFunctionBegin;
295d8588912SDave May   for (i=0; i<vs->nr; i++) {
296d8588912SDave May     for (j=0; j<vs->nc; j++) {
297d8588912SDave May       if (vs->m[i][j]) { ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); }
298d8588912SDave May     }
299d8588912SDave May   }
300d8588912SDave May   PetscFunctionReturn(0);
301d8588912SDave May }
302d8588912SDave May 
303d8588912SDave May #undef __FUNCT__
304f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatRow"
305f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
306d8588912SDave May {
307207556f9SJed Brown   PetscErrorCode ierr;
308f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
309f349c1fdSJed Brown   PetscInt       j;
310f349c1fdSJed Brown   Mat            sub;
311d8588912SDave May 
312d8588912SDave May   PetscFunctionBegin;
313f349c1fdSJed Brown   sub = vs->m[row][row];        /* Prefer to find on the diagonal */
314f349c1fdSJed Brown   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
315f349c1fdSJed Brown   if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",row);
316207556f9SJed Brown   ierr = MatPreallocated(sub);CHKERRQ(ierr); /* Ensure that the sizes are available */
317f349c1fdSJed Brown   *B = sub;
318f349c1fdSJed Brown   PetscFunctionReturn(0);
319d8588912SDave May }
320d8588912SDave May 
321f349c1fdSJed Brown #undef __FUNCT__
322f349c1fdSJed Brown #define __FUNCT__ "MatNestFindNonzeroSubMatCol"
323f349c1fdSJed Brown static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
324f349c1fdSJed Brown {
325207556f9SJed Brown   PetscErrorCode ierr;
326f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
327f349c1fdSJed Brown   PetscInt       i;
328f349c1fdSJed Brown   Mat            sub;
329f349c1fdSJed Brown 
330f349c1fdSJed Brown   PetscFunctionBegin;
331f349c1fdSJed Brown   sub = vs->m[col][col];        /* Prefer to find on the diagonal */
332f349c1fdSJed Brown   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
333f349c1fdSJed Brown   if (!sub) SETERRQ1(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",col);
334207556f9SJed Brown   ierr = MatPreallocated(sub);CHKERRQ(ierr); /* Ensure that the sizes are available */
335f349c1fdSJed Brown   *B = sub;
336f349c1fdSJed Brown   PetscFunctionReturn(0);
337d8588912SDave May }
338d8588912SDave May 
339f349c1fdSJed Brown #undef __FUNCT__
340f349c1fdSJed Brown #define __FUNCT__ "MatNestFindIS"
341f349c1fdSJed Brown static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
342f349c1fdSJed Brown {
343f349c1fdSJed Brown   PetscErrorCode ierr;
344f349c1fdSJed Brown   PetscInt       i;
345f349c1fdSJed Brown   PetscBool      flg;
346f349c1fdSJed Brown 
347f349c1fdSJed Brown   PetscFunctionBegin;
348f349c1fdSJed Brown   PetscValidPointer(list,3);
349f349c1fdSJed Brown   PetscValidHeaderSpecific(is,IS_CLASSID,4);
350f349c1fdSJed Brown   PetscValidIntPointer(found,5);
351f349c1fdSJed Brown   *found = -1;
352f349c1fdSJed Brown   for (i=0; i<n; i++) {
353207556f9SJed Brown     if (!list[i]) continue;
354f349c1fdSJed Brown     ierr = ISEqual(list[i],is,&flg);CHKERRQ(ierr);
355f349c1fdSJed Brown     if (flg) {
356f349c1fdSJed Brown       *found = i;
357f349c1fdSJed Brown       PetscFunctionReturn(0);
358f349c1fdSJed Brown     }
359f349c1fdSJed Brown   }
360f349c1fdSJed Brown   SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_INCOMP,"Could not find index set");
361f349c1fdSJed Brown   PetscFunctionReturn(0);
362f349c1fdSJed Brown }
363f349c1fdSJed Brown 
364f349c1fdSJed Brown #undef __FUNCT__
365f349c1fdSJed Brown #define __FUNCT__ "MatNestFindSubMat"
366f349c1fdSJed Brown static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
367f349c1fdSJed Brown {
368f349c1fdSJed Brown   PetscErrorCode ierr;
369f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
370f349c1fdSJed Brown   PetscInt       row,col;
371f349c1fdSJed Brown 
372f349c1fdSJed Brown   PetscFunctionBegin;
373f349c1fdSJed Brown   ierr = MatNestFindIS(A,vs->nr,is->row,isrow,&row);CHKERRQ(ierr);
374f349c1fdSJed Brown   ierr = MatNestFindIS(A,vs->nc,is->col,iscol,&col);CHKERRQ(ierr);
375f349c1fdSJed Brown   *B = vs->m[row][col];
376f349c1fdSJed Brown   PetscFunctionReturn(0);
377f349c1fdSJed Brown }
378f349c1fdSJed Brown 
379f349c1fdSJed Brown #undef __FUNCT__
380f349c1fdSJed Brown #define __FUNCT__ "MatGetSubMatrix_Nest"
381f349c1fdSJed Brown static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
382f349c1fdSJed Brown {
383f349c1fdSJed Brown   PetscErrorCode ierr;
384f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
385f349c1fdSJed Brown   Mat            sub;
386f349c1fdSJed Brown 
387f349c1fdSJed Brown   PetscFunctionBegin;
388f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);CHKERRQ(ierr);
389f349c1fdSJed Brown   switch (reuse) {
390f349c1fdSJed Brown   case MAT_INITIAL_MATRIX:
3917874fa86SDave May     if (sub) { ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr); }
392f349c1fdSJed Brown     *B = sub;
393f349c1fdSJed Brown     break;
394f349c1fdSJed Brown   case MAT_REUSE_MATRIX:
395f349c1fdSJed Brown     if (sub != *B) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
396f349c1fdSJed Brown     break;
397f349c1fdSJed Brown   case MAT_IGNORE_MATRIX:       /* Nothing to do */
398f349c1fdSJed Brown     break;
399f349c1fdSJed Brown   }
400f349c1fdSJed Brown   PetscFunctionReturn(0);
401f349c1fdSJed Brown }
402f349c1fdSJed Brown 
403f349c1fdSJed Brown #undef __FUNCT__
404f349c1fdSJed Brown #define __FUNCT__ "MatGetLocalSubMatrix_Nest"
405f349c1fdSJed Brown PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
406f349c1fdSJed Brown {
407f349c1fdSJed Brown   PetscErrorCode ierr;
408f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
409f349c1fdSJed Brown   Mat            sub;
410f349c1fdSJed Brown 
411f349c1fdSJed Brown   PetscFunctionBegin;
412f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
413f349c1fdSJed Brown   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
414f349c1fdSJed Brown   if (sub) {ierr = PetscObjectReference((PetscObject)sub);CHKERRQ(ierr);}
415f349c1fdSJed Brown   *B = sub;
416d8588912SDave May   PetscFunctionReturn(0);
417d8588912SDave May }
418d8588912SDave May 
419d8588912SDave May #undef __FUNCT__
420d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest"
421207556f9SJed Brown static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
422d8588912SDave May {
423d8588912SDave May   PetscErrorCode ierr;
424f349c1fdSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
425f349c1fdSJed Brown   Mat            sub;
426d8588912SDave May 
427d8588912SDave May   PetscFunctionBegin;
428f349c1fdSJed Brown   ierr = MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);CHKERRQ(ierr);
429f349c1fdSJed Brown   if (*B != sub) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
430f349c1fdSJed Brown   if (sub) {
431f349c1fdSJed Brown     if (((PetscObject)sub)->refct <= 1) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
4326bf464f9SBarry Smith     ierr = MatDestroy(B);CHKERRQ(ierr);
433d8588912SDave May   }
434d8588912SDave May   PetscFunctionReturn(0);
435d8588912SDave May }
436d8588912SDave May 
437d8588912SDave May #undef __FUNCT__
4387874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest"
4397874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
4407874fa86SDave May {
4417874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
4427874fa86SDave May   Vec            *bdiag;
4437874fa86SDave May   PetscInt       i;
4447874fa86SDave May   PetscErrorCode ierr;
4457874fa86SDave May 
4467874fa86SDave May   PetscFunctionBegin;
447*a0769a71SSatish Balay   /*  ierr = MatGetVecs_Nest(A,&diag,PETSC_NULL);CHKERRQ(ierr); */
448*a0769a71SSatish Balay   /*  ierr = VecNestGetSubVecs(diag,PETSC_NULL,&bdiag);CHKERRQ(ierr); */
4497874fa86SDave May   ierr = VecNestGetSubVecs(v,PETSC_NULL,&bdiag);CHKERRQ(ierr);
4507874fa86SDave May   for (i=0; i<bA->nr; i++) {
4517874fa86SDave May     if (bA->m[i][i]) {
4527874fa86SDave May       ierr = MatGetDiagonal(bA->m[i][i],bdiag[i]);CHKERRQ(ierr);
4537874fa86SDave May     } else {
4547874fa86SDave May       ierr = VecSet(bdiag[i],1.0);CHKERRQ(ierr);
4557874fa86SDave May     }
4567874fa86SDave May   }
4577874fa86SDave May   PetscFunctionReturn(0);
4587874fa86SDave May }
4597874fa86SDave May 
4607874fa86SDave May #undef __FUNCT__
4617874fa86SDave May #define __FUNCT__ "MatGetDiagonal_Nest2"
4627874fa86SDave May static PetscErrorCode MatGetDiagonal_Nest2(Mat A,Vec v)
4637874fa86SDave May {
4647874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
4657874fa86SDave May   Vec            diag,*bdiag;
4667874fa86SDave May   VecScatter     *vscat;
4677874fa86SDave May   PetscInt       i;
4687874fa86SDave May   PetscErrorCode ierr;
4697874fa86SDave May 
4707874fa86SDave May   PetscFunctionBegin;
4717874fa86SDave May   ierr = MatGetVecs_Nest(A,&diag,PETSC_NULL);CHKERRQ(ierr);
4727874fa86SDave May   ierr = MatGetDiagonal_Nest(A,diag);CHKERRQ(ierr);
4737874fa86SDave May 
4747874fa86SDave May   /* scatter diag into v */
4757874fa86SDave May   ierr = VecNestGetSubVecs(diag,PETSC_NULL,&bdiag);CHKERRQ(ierr);
4767874fa86SDave May   ierr = PetscMalloc( sizeof(VecScatter) * bA->nr, &vscat );CHKERRQ(ierr);
4777874fa86SDave May   for (i=0; i<bA->nr; i++) {
4787874fa86SDave May     ierr = VecScatterCreate(v,bA->isglobal.row[i], bdiag[i],PETSC_NULL,&vscat[i]);CHKERRQ(ierr);
4797874fa86SDave May     ierr = VecScatterBegin(vscat[i],bdiag[i],v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4807874fa86SDave May   }
4817874fa86SDave May   for (i=0; i<bA->nr; i++) {
4827874fa86SDave May     ierr = VecScatterEnd(vscat[i],bdiag[i],v,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr);
4837874fa86SDave May   }
4847874fa86SDave May   for (i=0; i<bA->nr; i++) {
4856bf464f9SBarry Smith     ierr = VecScatterDestroy(&vscat[i]);CHKERRQ(ierr);
4867874fa86SDave May   }
4877874fa86SDave May   ierr = PetscFree(vscat);CHKERRQ(ierr);
4886bf464f9SBarry Smith   ierr = VecDestroy(&diag);CHKERRQ(ierr);
4897874fa86SDave May   PetscFunctionReturn(0);
4907874fa86SDave May }
4917874fa86SDave May 
4927874fa86SDave May #undef __FUNCT__
4937874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest"
4947874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
4957874fa86SDave May {
4967874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
4977874fa86SDave May   Vec            *bl,*br;
4987874fa86SDave May   PetscInt       i,j;
4997874fa86SDave May   PetscErrorCode ierr;
5007874fa86SDave May 
5017874fa86SDave May   PetscFunctionBegin;
5027874fa86SDave May   ierr = VecNestGetSubVecs(l,PETSC_NULL,&bl);CHKERRQ(ierr);
5037874fa86SDave May   ierr = VecNestGetSubVecs(r,PETSC_NULL,&br);CHKERRQ(ierr);
5047874fa86SDave May   for (i=0; i<bA->nr; i++) {
5057874fa86SDave May     for (j=0; j<bA->nc; j++) {
5067874fa86SDave May       if (bA->m[i][j]) {
5077874fa86SDave May         ierr = MatDiagonalScale(bA->m[i][j],bl[i],br[j]);CHKERRQ(ierr);
5087874fa86SDave May       }
5097874fa86SDave May     }
5107874fa86SDave May   }
5117874fa86SDave May   PetscFunctionReturn(0);
5127874fa86SDave May }
5137874fa86SDave May 
5147874fa86SDave May #undef __FUNCT__
5157874fa86SDave May #define __FUNCT__ "MatDiagonalScale_Nest2"
5167874fa86SDave May static PetscErrorCode MatDiagonalScale_Nest2(Mat A,Vec l,Vec r)
5177874fa86SDave May {
5187874fa86SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
5197874fa86SDave May   Vec            bl,br,*ble,*bre;
5207874fa86SDave May   VecScatter     *vscatl,*vscatr;
5217874fa86SDave May   PetscInt       i;
5227874fa86SDave May   PetscErrorCode ierr;
5237874fa86SDave May 
5247874fa86SDave May   PetscFunctionBegin;
5257874fa86SDave May   /* scatter l,r into bl,br */
5267874fa86SDave May   ierr = MatGetVecs_Nest(A,&bl,&br);CHKERRQ(ierr);
5277874fa86SDave May   ierr = VecNestGetSubVecs(bl,PETSC_NULL,&ble);CHKERRQ(ierr);
5287874fa86SDave May   ierr = VecNestGetSubVecs(br,PETSC_NULL,&bre);CHKERRQ(ierr);
5297874fa86SDave May 
5307874fa86SDave May   /* row */
5317874fa86SDave May   ierr = PetscMalloc( sizeof(VecScatter) * bA->nr, &vscatl );CHKERRQ(ierr);
5327874fa86SDave May   for (i=0; i<bA->nr; i++) {
5337874fa86SDave May     ierr = VecScatterCreate(l,bA->isglobal.row[i], ble[i],PETSC_NULL,&vscatl[i]);CHKERRQ(ierr);
5347874fa86SDave May     ierr = VecScatterBegin(vscatl[i],l,ble[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5357874fa86SDave May   }
5367874fa86SDave May   for (i=0; i<bA->nr; i++) {
5377874fa86SDave May     ierr = VecScatterEnd(vscatl[i],l,ble[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5387874fa86SDave May   }
5397874fa86SDave May   /* col */
5407874fa86SDave May   ierr = PetscMalloc( sizeof(VecScatter) * bA->nc, &vscatr );CHKERRQ(ierr);
5417874fa86SDave May   for (i=0; i<bA->nc; i++) {
5427874fa86SDave May     ierr = VecScatterCreate(r,bA->isglobal.col[i], bre[i],PETSC_NULL,&vscatr[i]);CHKERRQ(ierr);
5437874fa86SDave May     ierr = VecScatterBegin(vscatr[i],l,bre[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5447874fa86SDave May   }
5457874fa86SDave May   for (i=0; i<bA->nc; i++) {
5467874fa86SDave May     ierr = VecScatterEnd(vscatr[i],r,bre[i],INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
5477874fa86SDave May   }
5487874fa86SDave May 
5497874fa86SDave May   ierr = MatDiagonalScale_Nest(A,bl,br);CHKERRQ(ierr);
5507874fa86SDave May 
5517874fa86SDave May   for (i=0; i<bA->nr; i++) {
5526bf464f9SBarry Smith     ierr = VecScatterDestroy(&vscatl[i]);CHKERRQ(ierr);
5537874fa86SDave May   }
5547874fa86SDave May   for (i=0; i<bA->nc; i++) {
5556bf464f9SBarry Smith     ierr = VecScatterDestroy(&vscatr[i]);CHKERRQ(ierr);
5567874fa86SDave May   }
5577874fa86SDave May   ierr = PetscFree(vscatl);CHKERRQ(ierr);
5587874fa86SDave May   ierr = PetscFree(vscatr);CHKERRQ(ierr);
5596bf464f9SBarry Smith   ierr = VecDestroy(&bl);CHKERRQ(ierr);
5606bf464f9SBarry Smith   ierr = VecDestroy(&br);CHKERRQ(ierr);
5617874fa86SDave May   PetscFunctionReturn(0);
5627874fa86SDave May }
5637874fa86SDave May 
5647874fa86SDave May #undef __FUNCT__
565d8588912SDave May #define __FUNCT__ "MatGetVecs_Nest"
566207556f9SJed Brown static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left)
567d8588912SDave May {
568d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
569d8588912SDave May   Vec            *L,*R;
570d8588912SDave May   MPI_Comm       comm;
571d8588912SDave May   PetscInt       i,j;
572d8588912SDave May   PetscErrorCode ierr;
573d8588912SDave May 
574d8588912SDave May   PetscFunctionBegin;
575d8588912SDave May   comm = ((PetscObject)A)->comm;
576d8588912SDave May   if (right) {
577d8588912SDave May     /* allocate R */
578d8588912SDave May     ierr = PetscMalloc( sizeof(Vec) * bA->nc, &R );CHKERRQ(ierr);
579d8588912SDave May     /* Create the right vectors */
580d8588912SDave May     for (j=0; j<bA->nc; j++) {
581d8588912SDave May       for (i=0; i<bA->nr; i++) {
582d8588912SDave May         if (bA->m[i][j]) {
583d8588912SDave May           ierr = MatGetVecs(bA->m[i][j],&R[j],PETSC_NULL);CHKERRQ(ierr);
584d8588912SDave May           break;
585d8588912SDave May         }
586d8588912SDave May       }
587d8588912SDave May       if (i==bA->nr) {
588d8588912SDave May         /* have an empty column */
589d8588912SDave May         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
590d8588912SDave May       }
591d8588912SDave May     }
592f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);CHKERRQ(ierr);
593d8588912SDave May     /* hand back control to the nest vector */
594d8588912SDave May     for (j=0; j<bA->nc; j++) {
5956bf464f9SBarry Smith       ierr = VecDestroy(&R[j]);CHKERRQ(ierr);
596d8588912SDave May     }
597d8588912SDave May     ierr = PetscFree(R);CHKERRQ(ierr);
598d8588912SDave May   }
599d8588912SDave May 
600d8588912SDave May   if (left) {
601d8588912SDave May     /* allocate L */
602d8588912SDave May     ierr = PetscMalloc( sizeof(Vec) * bA->nr, &L );CHKERRQ(ierr);
603d8588912SDave May     /* Create the left vectors */
604d8588912SDave May     for (i=0; i<bA->nr; i++) {
605d8588912SDave May       for (j=0; j<bA->nc; j++) {
606d8588912SDave May         if (bA->m[i][j]) {
607d8588912SDave May           ierr = MatGetVecs(bA->m[i][j],PETSC_NULL,&L[i]);CHKERRQ(ierr);
608d8588912SDave May           break;
609d8588912SDave May         }
610d8588912SDave May       }
611d8588912SDave May       if (j==bA->nc) {
612d8588912SDave May         /* have an empty row */
613d8588912SDave May         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
614d8588912SDave May       }
615d8588912SDave May     }
616d8588912SDave May 
617f349c1fdSJed Brown     ierr = VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);CHKERRQ(ierr);
618d8588912SDave May     for (i=0; i<bA->nr; i++) {
6196bf464f9SBarry Smith       ierr = VecDestroy(&L[i]);CHKERRQ(ierr);
620d8588912SDave May     }
621d8588912SDave May 
622d8588912SDave May     ierr = PetscFree(L);CHKERRQ(ierr);
623d8588912SDave May   }
624d8588912SDave May   PetscFunctionReturn(0);
625d8588912SDave May }
626d8588912SDave May 
627d8588912SDave May #undef __FUNCT__
628d8588912SDave May #define __FUNCT__ "MatView_Nest"
629207556f9SJed Brown static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
630d8588912SDave May {
631d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
632d8588912SDave May   PetscBool      isascii;
633d8588912SDave May   PetscInt       i,j;
634d8588912SDave May   PetscErrorCode ierr;
635d8588912SDave May 
636d8588912SDave May   PetscFunctionBegin;
637d8588912SDave May   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
638d8588912SDave May   if (isascii) {
639d8588912SDave May 
640d8588912SDave May     PetscViewerASCIIPrintf(viewer,"Matrix object: \n" );
641d8588912SDave May     PetscViewerASCIIPushTab( viewer );    /* push0 */
642d8588912SDave May     PetscViewerASCIIPrintf( viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc);
643d8588912SDave May 
644d8588912SDave May     PetscViewerASCIIPrintf(viewer,"MatNest structure: \n" );
645d8588912SDave May     for (i=0; i<bA->nr; i++) {
646d8588912SDave May       for (j=0; j<bA->nc; j++) {
647d8588912SDave May         const MatType type;
648d8588912SDave May         const char *name;
649d8588912SDave May         PetscInt NR,NC;
650d8588912SDave May         PetscBool isNest = PETSC_FALSE;
651d8588912SDave May 
652d8588912SDave May         if (!bA->m[i][j]) {
653d8588912SDave May           PetscViewerASCIIPrintf( viewer, "(%D,%D) : PETSC_NULL \n",i,j);
654d8588912SDave May           continue;
655d8588912SDave May         }
656d8588912SDave May         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
657d8588912SDave May         ierr = MatGetType( bA->m[i][j], &type );CHKERRQ(ierr);
658d8588912SDave May         name = ((PetscObject)bA->m[i][j])->prefix;
659d8588912SDave May         ierr = PetscTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
660d8588912SDave May 
661d8588912SDave May         PetscViewerASCIIPrintf( viewer, "(%D,%D) : name=\"%s\", type=%s, rows=%D, cols=%D \n",i,j,name,type,NR,NC);
662d8588912SDave May 
663d8588912SDave May         if (isNest) {
664d8588912SDave May           PetscViewerASCIIPushTab( viewer );  /* push1 */
665d8588912SDave May           ierr = MatView( bA->m[i][j], viewer );CHKERRQ(ierr);
666d8588912SDave May           PetscViewerASCIIPopTab(viewer);    /* pop1 */
667d8588912SDave May         }
668d8588912SDave May       }
669d8588912SDave May     }
670d8588912SDave May     PetscViewerASCIIPopTab(viewer);    /* pop0 */
671d8588912SDave May   }
672d8588912SDave May   PetscFunctionReturn(0);
673d8588912SDave May }
674d8588912SDave May 
675d8588912SDave May #undef __FUNCT__
676d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest"
677207556f9SJed Brown static PetscErrorCode MatZeroEntries_Nest(Mat A)
678d8588912SDave May {
679d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
680d8588912SDave May   PetscInt       i,j;
681d8588912SDave May   PetscErrorCode ierr;
682d8588912SDave May 
683d8588912SDave May   PetscFunctionBegin;
684d8588912SDave May   for (i=0; i<bA->nr; i++) {
685d8588912SDave May     for (j=0; j<bA->nc; j++) {
686d8588912SDave May       if (!bA->m[i][j]) continue;
687d8588912SDave May       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
688d8588912SDave May     }
689d8588912SDave May   }
690d8588912SDave May   PetscFunctionReturn(0);
691d8588912SDave May }
692d8588912SDave May 
693d8588912SDave May #undef __FUNCT__
694d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest"
695207556f9SJed Brown static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
696d8588912SDave May {
697d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
698841e96a3SJed Brown   Mat            *b;
699841e96a3SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
700d8588912SDave May   PetscErrorCode ierr;
701d8588912SDave May 
702d8588912SDave May   PetscFunctionBegin;
703841e96a3SJed Brown   ierr = PetscMalloc(nr*nc*sizeof(Mat),&b);CHKERRQ(ierr);
704841e96a3SJed Brown   for (i=0; i<nr; i++) {
705841e96a3SJed Brown     for (j=0; j<nc; j++) {
706841e96a3SJed Brown       if (bA->m[i][j]) {
707841e96a3SJed Brown         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
708841e96a3SJed Brown       } else {
709841e96a3SJed Brown         b[i*nc+j] = PETSC_NULL;
710d8588912SDave May       }
711d8588912SDave May     }
712d8588912SDave May   }
713f349c1fdSJed Brown   ierr = MatCreateNest(((PetscObject)A)->comm,nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);CHKERRQ(ierr);
714841e96a3SJed Brown   /* Give the new MatNest exclusive ownership */
715841e96a3SJed Brown   for (i=0; i<nr*nc; i++) {
7166bf464f9SBarry Smith     ierr = MatDestroy(&b[i]);CHKERRQ(ierr);
717d8588912SDave May   }
718d8588912SDave May   ierr = PetscFree(b);CHKERRQ(ierr);
719d8588912SDave May 
720841e96a3SJed Brown   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
721841e96a3SJed Brown   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
722d8588912SDave May   PetscFunctionReturn(0);
723d8588912SDave May }
724d8588912SDave May 
725d8588912SDave May /* nest api */
726d8588912SDave May EXTERN_C_BEGIN
727d8588912SDave May #undef __FUNCT__
728d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest"
729d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
730d8588912SDave May {
731d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
732d8588912SDave May   PetscFunctionBegin;
733d8588912SDave May   if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
734d8588912SDave May   if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
735d8588912SDave May   *mat = bA->m[idxm][jdxm];
736d8588912SDave May   PetscFunctionReturn(0);
737d8588912SDave May }
738d8588912SDave May EXTERN_C_END
739d8588912SDave May 
740d8588912SDave May #undef __FUNCT__
741d8588912SDave May #define __FUNCT__ "MatNestGetSubMat"
742d8588912SDave May /*@C
743d8588912SDave May  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
744d8588912SDave May 
745d8588912SDave May  Not collective
746d8588912SDave May 
747d8588912SDave May  Input Parameters:
748d8588912SDave May  .  A  - nest matrix
749d8588912SDave May  .  idxm - index of the matrix within the nest matrix
750d8588912SDave May  .  jdxm - index of the matrix within the nest matrix
751d8588912SDave May 
752d8588912SDave May  Output Parameter:
753d8588912SDave May  .  sub - matrix at index idxm,jdxm within the nest matrix
754d8588912SDave May 
755d8588912SDave May  Notes:
756d8588912SDave May 
757d8588912SDave May  Level: developer
758d8588912SDave May 
759d8588912SDave May  .seealso: MatNestGetSize(), MatNestGetSubMats()
760d8588912SDave May @*/
7617087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
762d8588912SDave May {
763699a902aSJed Brown   PetscErrorCode ierr;
764d8588912SDave May 
765d8588912SDave May   PetscFunctionBegin;
766699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));CHKERRQ(ierr);
767d8588912SDave May   PetscFunctionReturn(0);
768d8588912SDave May }
769d8588912SDave May 
770d8588912SDave May EXTERN_C_BEGIN
771d8588912SDave May #undef __FUNCT__
772d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest"
773d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
774d8588912SDave May {
775d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
776d8588912SDave May   PetscFunctionBegin;
777d8588912SDave May   if (M)   { *M   = bA->nr; }
778d8588912SDave May   if (N)   { *N   = bA->nc; }
779d8588912SDave May   if (mat) { *mat = bA->m;  }
780d8588912SDave May   PetscFunctionReturn(0);
781d8588912SDave May }
782d8588912SDave May EXTERN_C_END
783d8588912SDave May 
784d8588912SDave May #undef __FUNCT__
785d8588912SDave May #define __FUNCT__ "MatNestGetSubMats"
786d8588912SDave May /*@C
787d8588912SDave May  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
788d8588912SDave May 
789d8588912SDave May  Not collective
790d8588912SDave May 
791d8588912SDave May  Input Parameters:
792d8588912SDave May  .  A  - nest matri
793d8588912SDave May 
794d8588912SDave May  Output Parameter:
795d8588912SDave May  .  M - number of rows in the nest matrix
796d8588912SDave May  .  N - number of cols in the nest matrix
797d8588912SDave May  .  mat - 2d array of matrices
798d8588912SDave May 
799d8588912SDave May  Notes:
800d8588912SDave May 
801d8588912SDave May  The user should not free the array mat.
802d8588912SDave May 
803d8588912SDave May  Level: developer
804d8588912SDave May 
805d8588912SDave May  .seealso: MatNestGetSize(), MatNestGetSubMat()
806d8588912SDave May @*/
8077087cfbeSBarry Smith PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
808d8588912SDave May {
809699a902aSJed Brown   PetscErrorCode ierr;
810d8588912SDave May 
811d8588912SDave May   PetscFunctionBegin;
812699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));CHKERRQ(ierr);
813d8588912SDave May   PetscFunctionReturn(0);
814d8588912SDave May }
815d8588912SDave May 
816d8588912SDave May EXTERN_C_BEGIN
817d8588912SDave May #undef __FUNCT__
818d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest"
8197087cfbeSBarry Smith PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
820d8588912SDave May {
821d8588912SDave May   Mat_Nest  *bA = (Mat_Nest*)A->data;
822d8588912SDave May 
823d8588912SDave May   PetscFunctionBegin;
824d8588912SDave May   if (M) { *M  = bA->nr; }
825d8588912SDave May   if (N) { *N  = bA->nc; }
826d8588912SDave May   PetscFunctionReturn(0);
827d8588912SDave May }
828d8588912SDave May EXTERN_C_END
829d8588912SDave May 
830d8588912SDave May #undef __FUNCT__
831d8588912SDave May #define __FUNCT__ "MatNestGetSize"
832d8588912SDave May /*@C
833d8588912SDave May  MatNestGetSize - Returns the size of the nest matrix.
834d8588912SDave May 
835d8588912SDave May  Not collective
836d8588912SDave May 
837d8588912SDave May  Input Parameters:
838d8588912SDave May  .  A  - nest matrix
839d8588912SDave May 
840d8588912SDave May  Output Parameter:
841d8588912SDave May  .  M - number of rows in the nested mat
842d8588912SDave May  .  N - number of cols in the nested mat
843d8588912SDave May 
844d8588912SDave May  Notes:
845d8588912SDave May 
846d8588912SDave May  Level: developer
847d8588912SDave May 
848d8588912SDave May  .seealso: MatNestGetSubMat(), MatNestGetSubMats()
849d8588912SDave May @*/
8507087cfbeSBarry Smith PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
851d8588912SDave May {
852699a902aSJed Brown   PetscErrorCode ierr;
853d8588912SDave May 
854d8588912SDave May   PetscFunctionBegin;
855699a902aSJed Brown   ierr = PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));CHKERRQ(ierr);
856d8588912SDave May   PetscFunctionReturn(0);
857d8588912SDave May }
858d8588912SDave May 
859207556f9SJed Brown EXTERN_C_BEGIN
860207556f9SJed Brown #undef __FUNCT__
861207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType_Nest"
8627087cfbeSBarry Smith PetscErrorCode  MatNestSetVecType_Nest(Mat A,const VecType vtype)
863207556f9SJed Brown {
864207556f9SJed Brown   PetscErrorCode ierr;
865207556f9SJed Brown   PetscBool      flg;
866207556f9SJed Brown 
867207556f9SJed Brown   PetscFunctionBegin;
868207556f9SJed Brown   ierr = PetscStrcmp(vtype,VECNEST,&flg);CHKERRQ(ierr);
869207556f9SJed Brown   /* In reality, this only distinguishes VECNEST and "other" */
870207556f9SJed Brown   A->ops->getvecs       = flg ? MatGetVecs_Nest       : 0;
8717874fa86SDave May   A->ops->getdiagonal   = flg ? MatGetDiagonal_Nest   : 0;
8727874fa86SDave May   A->ops->diagonalscale = flg ? MatDiagonalScale_Nest : 0;
8737874fa86SDave May 
874207556f9SJed Brown  PetscFunctionReturn(0);
875207556f9SJed Brown }
876207556f9SJed Brown EXTERN_C_END
877207556f9SJed Brown 
878207556f9SJed Brown #undef __FUNCT__
879207556f9SJed Brown #define __FUNCT__ "MatNestSetVecType"
880207556f9SJed Brown /*@C
881207556f9SJed Brown  MatNestSetVecType - Sets the type of Vec returned by MatGetVecs()
882207556f9SJed Brown 
883207556f9SJed Brown  Not collective
884207556f9SJed Brown 
885207556f9SJed Brown  Input Parameters:
886207556f9SJed Brown +  A  - nest matrix
887207556f9SJed Brown -  vtype - type to use for creating vectors
888207556f9SJed Brown 
889207556f9SJed Brown  Notes:
890207556f9SJed Brown 
891207556f9SJed Brown  Level: developer
892207556f9SJed Brown 
893207556f9SJed Brown  .seealso: MatGetVecs()
894207556f9SJed Brown @*/
8957087cfbeSBarry Smith PetscErrorCode  MatNestSetVecType(Mat A,const VecType vtype)
896207556f9SJed Brown {
897207556f9SJed Brown   PetscErrorCode ierr;
898207556f9SJed Brown 
899207556f9SJed Brown   PetscFunctionBegin;
900207556f9SJed Brown   ierr = PetscTryMethod(A,"MatNestSetVecType_C",(Mat,const VecType),(A,vtype));CHKERRQ(ierr);
901207556f9SJed Brown   PetscFunctionReturn(0);
902207556f9SJed Brown }
903207556f9SJed Brown 
904c8883902SJed Brown EXTERN_C_BEGIN
905d8588912SDave May #undef __FUNCT__
906c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats_Nest"
907c8883902SJed Brown PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
908d8588912SDave May {
909c8883902SJed Brown   Mat_Nest       *s = (Mat_Nest*)A->data;
910c8883902SJed Brown   PetscInt       i,j,m,n,M,N;
911d8588912SDave May   PetscErrorCode ierr;
912d8588912SDave May 
913d8588912SDave May   PetscFunctionBegin;
914c8883902SJed Brown   s->nr = nr;
915c8883902SJed Brown   s->nc = nc;
916d8588912SDave May 
917c8883902SJed Brown   /* Create space for submatrices */
918c8883902SJed Brown   ierr = PetscMalloc(sizeof(Mat*)*nr,&s->m);CHKERRQ(ierr);
919c8883902SJed Brown   for (i=0; i<nr; i++) {
920c8883902SJed Brown     ierr = PetscMalloc(sizeof(Mat)*nc,&s->m[i]);CHKERRQ(ierr);
921d8588912SDave May   }
922c8883902SJed Brown   for (i=0; i<nr; i++) {
923c8883902SJed Brown     for (j=0; j<nc; j++) {
924c8883902SJed Brown       s->m[i][j] = a[i*nc+j];
925c8883902SJed Brown       if (a[i*nc+j]) {
926c8883902SJed Brown         ierr = PetscObjectReference((PetscObject)a[i*nc+j]);CHKERRQ(ierr);
927d8588912SDave May       }
928d8588912SDave May     }
929d8588912SDave May   }
930d8588912SDave May 
931c8883902SJed Brown   ierr = PetscMalloc(sizeof(IS)*nr,&s->isglobal.row);CHKERRQ(ierr);
932c8883902SJed Brown   ierr = PetscMalloc(sizeof(IS)*nc,&s->isglobal.col);CHKERRQ(ierr);
933d8588912SDave May 
934c8883902SJed Brown   ierr = PetscMalloc(sizeof(PetscInt)*nr,&s->row_len);CHKERRQ(ierr);
935c8883902SJed Brown   ierr = PetscMalloc(sizeof(PetscInt)*nc,&s->col_len);CHKERRQ(ierr);
936c8883902SJed Brown   for (i=0; i<nr; i++) s->row_len[i]=-1;
937c8883902SJed Brown   for (j=0; j<nc; j++) s->col_len[j]=-1;
938d8588912SDave May 
939c8883902SJed Brown   m = n = 0;
940c8883902SJed Brown   M = N = 0;
941c8883902SJed Brown   ierr = MatNestGetSize_Private(A,PETSC_TRUE,&M,&N);CHKERRQ(ierr);
942c8883902SJed Brown   ierr = MatNestGetSize_Private(A,PETSC_FALSE,&m,&n);CHKERRQ(ierr);
943d8588912SDave May 
944c8883902SJed Brown   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
945c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
946c8883902SJed Brown   ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr);
947c8883902SJed Brown   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
948c8883902SJed Brown   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
949c8883902SJed Brown   ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr);
950c8883902SJed Brown 
951c8883902SJed Brown   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
952c8883902SJed Brown   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
953c8883902SJed Brown 
954c8883902SJed Brown   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
955c8883902SJed Brown   ierr = PetscMalloc2(nr,Vec,&s->left,nc,Vec,&s->right);CHKERRQ(ierr);
956c8883902SJed Brown   ierr = PetscMemzero(s->left,nr*sizeof(Vec));CHKERRQ(ierr);
957c8883902SJed Brown   ierr = PetscMemzero(s->right,nc*sizeof(Vec));CHKERRQ(ierr);
958d8588912SDave May   PetscFunctionReturn(0);
959d8588912SDave May }
960c8883902SJed Brown EXTERN_C_END
961d8588912SDave May 
962c8883902SJed Brown #undef __FUNCT__
963c8883902SJed Brown #define __FUNCT__ "MatNestSetSubMats"
964c8883902SJed Brown /*@
965c8883902SJed Brown    MatNestSetSubMats - Sets the nested submatrices
966c8883902SJed Brown 
967c8883902SJed Brown    Collective on Mat
968c8883902SJed Brown 
969c8883902SJed Brown    Input Parameter:
970c8883902SJed Brown +  N - nested matrix
971c8883902SJed Brown .  nr - number of nested row blocks
972c8883902SJed Brown .  is_row - index sets for each nested row block, or PETSC_NULL to make contiguous
973c8883902SJed Brown .  nc - number of nested column blocks
974c8883902SJed Brown .  is_col - index sets for each nested column block, or PETSC_NULL to make contiguous
975c8883902SJed Brown -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL
976c8883902SJed Brown 
977c8883902SJed Brown    Level: advanced
978c8883902SJed Brown 
979c8883902SJed Brown .seealso: MatCreateNest(), MATNEST
980c8883902SJed Brown @*/
981c8883902SJed Brown PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
982c8883902SJed Brown {
983c8883902SJed Brown   PetscErrorCode ierr;
984c8883902SJed Brown   PetscInt i;
985c8883902SJed Brown 
986c8883902SJed Brown   PetscFunctionBegin;
987c8883902SJed Brown   PetscValidHeaderSpecific(A,MAT_CLASSID,1);
988c8883902SJed Brown   if (nr < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
989c8883902SJed Brown   if (nr && is_row) {
990c8883902SJed Brown     PetscValidPointer(is_row,3);
991c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_row[i],IS_CLASSID,3);
992c8883902SJed Brown   }
993c8883902SJed Brown   if (nc < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
994c8883902SJed Brown   if (nc && is_row) {
995c8883902SJed Brown     PetscValidPointer(is_col,5);
996c8883902SJed Brown     for (i=0; i<nr; i++) PetscValidHeaderSpecific(is_col[i],IS_CLASSID,5);
997c8883902SJed Brown   }
998c8883902SJed Brown   if (nr*nc) PetscValidPointer(a,6);
999c8883902SJed 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);
1000c8883902SJed Brown   PetscFunctionReturn(0);
1001c8883902SJed Brown }
1002d8588912SDave May 
1003d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1004d8588912SDave May /*
1005d8588912SDave May   nprocessors = NP
1006d8588912SDave May   Nest x^T = ( (g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1) )
1007d8588912SDave May        proc 0: => (g_0,h_0,)
1008d8588912SDave May        proc 1: => (g_1,h_1,)
1009d8588912SDave May        ...
1010d8588912SDave May        proc nprocs-1: => (g_NP-1,h_NP-1,)
1011d8588912SDave May 
1012d8588912SDave May             proc 0:                      proc 1:                    proc nprocs-1:
1013d8588912SDave May     is[0] = ( 0,1,2,...,nlocal(g_0)-1 )  ( 0,1,...,nlocal(g_1)-1 )  ( 0,1,...,nlocal(g_NP-1) )
1014d8588912SDave May 
1015d8588912SDave May             proc 0:
1016d8588912SDave May     is[1] = ( nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1 )
1017d8588912SDave May             proc 1:
1018d8588912SDave May     is[1] = ( nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1 )
1019d8588912SDave May 
1020d8588912SDave May             proc NP-1:
1021d8588912SDave May     is[1] = ( nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1 )
1022d8588912SDave May */
1023d8588912SDave May #undef __FUNCT__
1024d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private"
1025841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1026d8588912SDave May {
1027e2d7f03fSJed Brown   Mat_Nest       *vs = (Mat_Nest*)A->data;
10282ae74bdbSJed Brown   PetscInt       i,j,offset,n,bs;
1029d8588912SDave May   PetscErrorCode ierr;
10302ae74bdbSJed Brown   Mat            sub;
1031d8588912SDave May 
1032d8588912SDave May   PetscFunctionBegin;
1033d8588912SDave May   if (is_row) { /* valid IS is passed in */
1034d8588912SDave May     /* refs on is[] are incremeneted */
1035e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1036d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
1037e2d7f03fSJed Brown       vs->isglobal.row[i] = is_row[i];
1038d8588912SDave May     }
10392ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
10402ae74bdbSJed Brown     offset = A->rmap->rstart;
1041e2d7f03fSJed Brown     for (i=0; i<vs->nr; i++) {
1042f349c1fdSJed Brown       ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
10432ae74bdbSJed Brown       ierr = MatGetLocalSize(sub,&n,PETSC_NULL);CHKERRQ(ierr);
1044207556f9SJed Brown       if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
10452ae74bdbSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1046e2d7f03fSJed Brown       ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.row[i]);CHKERRQ(ierr);
1047e2d7f03fSJed Brown       ierr = ISSetBlockSize(vs->isglobal.row[i],bs);CHKERRQ(ierr);
10482ae74bdbSJed Brown       offset += n;
1049d8588912SDave May     }
1050d8588912SDave May   }
1051d8588912SDave May 
1052d8588912SDave May   if (is_col) { /* valid IS is passed in */
1053d8588912SDave May     /* refs on is[] are incremeneted */
1054e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1055d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
1056e2d7f03fSJed Brown       vs->isglobal.col[j] = is_col[j];
1057d8588912SDave May     }
10582ae74bdbSJed Brown   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
10592ae74bdbSJed Brown     offset = A->cmap->rstart;
1060e2d7f03fSJed Brown     for (j=0; j<vs->nc; j++) {
1061f349c1fdSJed Brown       ierr = MatNestFindNonzeroSubMatCol(A,j,&sub);CHKERRQ(ierr);
10622ae74bdbSJed Brown       ierr = MatGetLocalSize(sub,PETSC_NULL,&n);CHKERRQ(ierr);
1063207556f9SJed Brown       if (n < 0) SETERRQ(((PetscObject)A)->comm,PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
10642ae74bdbSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1065e2d7f03fSJed Brown       ierr = ISCreateStride(((PetscObject)sub)->comm,n,offset,1,&vs->isglobal.col[j]);CHKERRQ(ierr);
1066e2d7f03fSJed Brown       ierr = ISSetBlockSize(vs->isglobal.col[j],bs);CHKERRQ(ierr);
10672ae74bdbSJed Brown       offset += n;
1068d8588912SDave May     }
1069d8588912SDave May   }
1070e2d7f03fSJed Brown 
1071e2d7f03fSJed Brown   /* Set up the local ISs */
1072e2d7f03fSJed Brown   ierr = PetscMalloc(vs->nr*sizeof(IS),&vs->islocal.row);CHKERRQ(ierr);
1073e2d7f03fSJed Brown   ierr = PetscMalloc(vs->nc*sizeof(IS),&vs->islocal.col);CHKERRQ(ierr);
1074e2d7f03fSJed Brown   for (i=0,offset=0; i<vs->nr; i++) {
1075e2d7f03fSJed Brown     IS                     isloc;
1076e2d7f03fSJed Brown     ISLocalToGlobalMapping rmap;
1077e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1078e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatRow(A,i,&sub);CHKERRQ(ierr);
1079e2d7f03fSJed Brown     ierr = MatGetLocalToGlobalMapping(sub,&rmap,PETSC_NULL);CHKERRQ(ierr);
1080207556f9SJed Brown     if (rmap) {
1081e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1082e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(rmap,&nlocal);CHKERRQ(ierr);
1083e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1084e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1085207556f9SJed Brown     } else {
1086207556f9SJed Brown       nlocal = 0;
1087207556f9SJed Brown       isloc  = PETSC_NULL;
1088207556f9SJed Brown     }
1089e2d7f03fSJed Brown     vs->islocal.row[i] = isloc;
1090e2d7f03fSJed Brown     offset += nlocal;
1091e2d7f03fSJed Brown   }
1092e2d7f03fSJed Brown   for (i=0,offset=0; i<vs->nr; i++) {
1093e2d7f03fSJed Brown     IS                     isloc;
1094e2d7f03fSJed Brown     ISLocalToGlobalMapping cmap;
1095e2d7f03fSJed Brown     PetscInt               nlocal,bs;
1096e2d7f03fSJed Brown     ierr = MatNestFindNonzeroSubMatCol(A,i,&sub);CHKERRQ(ierr);
1097e2d7f03fSJed Brown     ierr = MatGetLocalToGlobalMapping(sub,PETSC_NULL,&cmap);CHKERRQ(ierr);
1098207556f9SJed Brown     if (cmap) {
1099e2d7f03fSJed Brown       ierr = MatGetBlockSize(sub,&bs);CHKERRQ(ierr);
1100e2d7f03fSJed Brown       ierr = ISLocalToGlobalMappingGetSize(cmap,&nlocal);CHKERRQ(ierr);
1101e2d7f03fSJed Brown       ierr = ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);CHKERRQ(ierr);
1102e2d7f03fSJed Brown       ierr = ISSetBlockSize(isloc,bs);CHKERRQ(ierr);
1103207556f9SJed Brown     } else {
1104207556f9SJed Brown       nlocal = 0;
1105207556f9SJed Brown       isloc  = PETSC_NULL;
1106207556f9SJed Brown     }
1107e2d7f03fSJed Brown     vs->islocal.col[i] = isloc;
1108e2d7f03fSJed Brown     offset += nlocal;
1109e2d7f03fSJed Brown   }
1110d8588912SDave May   PetscFunctionReturn(0);
1111d8588912SDave May }
1112d8588912SDave May 
1113d8588912SDave May #undef __FUNCT__
1114d8588912SDave May #define __FUNCT__ "MatCreateNest"
1115659c6bb0SJed Brown /*@
1116659c6bb0SJed Brown    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately
1117659c6bb0SJed Brown 
1118659c6bb0SJed Brown    Collective on Mat
1119659c6bb0SJed Brown 
1120659c6bb0SJed Brown    Input Parameter:
1121659c6bb0SJed Brown +  comm - Communicator for the new Mat
1122659c6bb0SJed Brown .  nr - number of nested row blocks
1123659c6bb0SJed Brown .  is_row - index sets for each nested row block, or PETSC_NULL to make contiguous
1124659c6bb0SJed Brown .  nc - number of nested column blocks
1125659c6bb0SJed Brown .  is_col - index sets for each nested column block, or PETSC_NULL to make contiguous
1126659c6bb0SJed Brown -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using PETSC_NULL
1127659c6bb0SJed Brown 
1128659c6bb0SJed Brown    Output Parameter:
1129659c6bb0SJed Brown .  B - new matrix
1130659c6bb0SJed Brown 
1131659c6bb0SJed Brown    Level: advanced
1132659c6bb0SJed Brown 
1133659c6bb0SJed Brown .seealso: MatCreate(), VecCreateNest(), DMGetMatrix(), MATNEST
1134659c6bb0SJed Brown @*/
11357087cfbeSBarry Smith PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1136d8588912SDave May {
1137d8588912SDave May   Mat            A;
1138d8588912SDave May   PetscErrorCode ierr;
1139d8588912SDave May 
1140d8588912SDave May   PetscFunctionBegin;
1141c8883902SJed Brown   *B = 0;
1142d8588912SDave May   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
1143c8883902SJed Brown   ierr = MatSetType(A,MATNEST);CHKERRQ(ierr);
1144c8883902SJed Brown   ierr = MatNestSetSubMats(A,nr,is_row,nc,is_col,a);CHKERRQ(ierr);
1145d8588912SDave May   *B = A;
1146d8588912SDave May   PetscFunctionReturn(0);
1147d8588912SDave May }
1148659c6bb0SJed Brown 
1149659c6bb0SJed Brown /*MC
1150659c6bb0SJed Brown   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.
1151659c6bb0SJed Brown 
1152659c6bb0SJed Brown   Level: intermediate
1153659c6bb0SJed Brown 
1154659c6bb0SJed Brown   Notes:
1155659c6bb0SJed Brown   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1156659c6bb0SJed Brown   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1157659c6bb0SJed Brown   It is usually used with DMComposite and DMGetMatrix()
1158659c6bb0SJed Brown 
1159659c6bb0SJed Brown .seealso: MatCreate(), MatType, MatCreateNest()
1160659c6bb0SJed Brown M*/
1161c8883902SJed Brown EXTERN_C_BEGIN
1162c8883902SJed Brown #undef __FUNCT__
1163c8883902SJed Brown #define __FUNCT__ "MatCreate_Nest"
1164c8883902SJed Brown PetscErrorCode MatCreate_Nest(Mat A)
1165c8883902SJed Brown {
1166c8883902SJed Brown   Mat_Nest       *s;
1167c8883902SJed Brown   PetscErrorCode ierr;
1168c8883902SJed Brown 
1169c8883902SJed Brown   PetscFunctionBegin;
1170c8883902SJed Brown   ierr = PetscNewLog(A,Mat_Nest,&s);CHKERRQ(ierr);
1171c8883902SJed Brown   A->data         = (void*)s;
1172c8883902SJed Brown   s->nr = s->nc   = -1;
1173c8883902SJed Brown   s->m            = PETSC_NULL;
1174c8883902SJed Brown 
1175c8883902SJed Brown   ierr = PetscMemzero(A->ops,sizeof(*A->ops));CHKERRQ(ierr);
1176c8883902SJed Brown   A->ops->mult                  = MatMult_Nest;
1177c8883902SJed Brown   A->ops->multadd               = 0;
1178c8883902SJed Brown   A->ops->multtranspose         = MatMultTranspose_Nest;
1179c8883902SJed Brown   A->ops->multtransposeadd      = 0;
1180c8883902SJed Brown   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1181c8883902SJed Brown   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1182c8883902SJed Brown   A->ops->zeroentries           = MatZeroEntries_Nest;
1183c8883902SJed Brown   A->ops->duplicate             = MatDuplicate_Nest;
1184c8883902SJed Brown   A->ops->getsubmatrix          = MatGetSubMatrix_Nest;
1185c8883902SJed Brown   A->ops->destroy               = MatDestroy_Nest;
1186c8883902SJed Brown   A->ops->view                  = MatView_Nest;
1187c8883902SJed Brown   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1188c8883902SJed Brown   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1189c8883902SJed Brown   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
11907874fa86SDave May   A->ops->getdiagonal           = MatGetDiagonal_Nest2; /* VECNEST version activated by calling MatNestSetVecType(A,VECNEST) */
11917874fa86SDave May   A->ops->diagonalscale         = MatDiagonalScale_Nest2; /* VECNEST version activated by calling MatNestSetVecType(A,VECNEST) */
1192c8883902SJed Brown 
1193c8883902SJed Brown   A->spptr        = 0;
1194c8883902SJed Brown   A->same_nonzero = PETSC_FALSE;
1195c8883902SJed Brown   A->assembled    = PETSC_FALSE;
1196c8883902SJed Brown 
1197c8883902SJed Brown   /* expose Nest api's */
1198c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C",   "MatNestGetSubMat_Nest",   MatNestGetSubMat_Nest);CHKERRQ(ierr);
1199c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C",  "MatNestGetSubMats_Nest",  MatNestGetSubMats_Nest);CHKERRQ(ierr);
1200c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C",     "MatNestGetSize_Nest",     MatNestGetSize_Nest);CHKERRQ(ierr);
1201c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetVecType_C",  "MatNestSetVecType_Nest",  MatNestSetVecType_Nest);CHKERRQ(ierr);
1202c8883902SJed Brown   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestSetSubMats_C",  "MatNestSetSubMats_Nest",  MatNestSetSubMats_Nest);CHKERRQ(ierr);
1203c8883902SJed Brown 
1204c8883902SJed Brown   ierr = PetscObjectChangeTypeName((PetscObject)A,MATNEST);CHKERRQ(ierr);
1205c8883902SJed Brown   PetscFunctionReturn(0);
1206c8883902SJed Brown }
120786e80854SHong Zhang EXTERN_C_END
1208