xref: /petsc/src/mat/impls/nest/matnest.c (revision 841e96a3f21bd364d8ebd590337294c43808fb01)
1d8588912SDave May #define PETSCMAT_DLL
2d8588912SDave May 
3d8588912SDave May #include "matnestimpl.h"
4d8588912SDave May 
5d8588912SDave May /* private functions */
6d8588912SDave May #undef __FUNCT__
7d8588912SDave May #define __FUNCT__ "MatNestGetSize_Recursive"
8d8588912SDave May static PetscErrorCode MatNestGetSize_Recursive(Mat A,PetscBool globalsize,PetscInt *M,PetscInt *N)
9d8588912SDave May {
10d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
11d8588912SDave May   PetscInt       sizeM,sizeN,i,j,index;
12d8588912SDave May   PetscErrorCode ierr;
13d8588912SDave May 
14d8588912SDave May   PetscFunctionBegin;
15d8588912SDave May   /* rows */
16d8588912SDave May   /* now descend recursively */
17d8588912SDave May   for (i=0; i<bA->nr; i++) {
18d8588912SDave May     for (j=0; j<bA->nc; j++) {
19d8588912SDave May       if (bA->m[i][j]) { index = j; break; }
20d8588912SDave May     }
21d8588912SDave May     if (bA->m[i][index]) {
22d8588912SDave May       if (globalsize) { ierr = MatGetSize(bA->m[i][index],&sizeM,&sizeN);CHKERRQ(ierr); }
23d8588912SDave May       else {            ierr = MatGetLocalSize(bA->m[i][index],&sizeM,&sizeN);CHKERRQ(ierr); }
24d8588912SDave May       *M = *M + sizeM;
25d8588912SDave May     }
26d8588912SDave May   }
27d8588912SDave May 
28d8588912SDave May   /* cols */
29d8588912SDave May   /* now descend recursively */
30d8588912SDave May   for (j=0; j<bA->nc; j++) {
31d8588912SDave May     for (i=0; i<bA->nr; i++) {
32d8588912SDave May       if (bA->m[i][j]) { index = i; break; }
33d8588912SDave May     }
34d8588912SDave May     if (bA->m[index][j]) {
35d8588912SDave May       if (globalsize) { ierr = MatGetSize(bA->m[index][j],&sizeM,&sizeN);CHKERRQ(ierr); }
36d8588912SDave May       else {            ierr = MatGetLocalSize(bA->m[index][j],&sizeM,&sizeN);CHKERRQ(ierr); }
37d8588912SDave May       *N = *N + sizeN;
38d8588912SDave May     }
39d8588912SDave May   }
40d8588912SDave May   PetscFunctionReturn(0);
41d8588912SDave May }
42d8588912SDave May 
43d8588912SDave May #undef __FUNCT__
44d8588912SDave May #define __FUNCT__ "PETSc_MatNest_CheckMatVecCompatibility2"
45d8588912SDave May static PetscErrorCode PETSc_MatNest_CheckMatVecCompatibility2(Mat A,Vec x,Vec y)
46d8588912SDave May {
47d8588912SDave May   PetscBool      isANest, isxNest, isyNest;
48d8588912SDave May   PetscErrorCode ierr;
49d8588912SDave May 
50d8588912SDave May   PetscFunctionBegin;
51d8588912SDave May   isANest = isxNest = isyNest = PETSC_FALSE;
52d8588912SDave May   ierr = PetscTypeCompare( (PetscObject)A, MATNEST, &isANest );CHKERRQ(ierr);
53d8588912SDave May   ierr = PetscTypeCompare( (PetscObject)x, VECNEST, &isxNest );CHKERRQ(ierr);
54d8588912SDave May   ierr = PetscTypeCompare( (PetscObject)y, VECNEST, &isyNest );CHKERRQ(ierr);
55d8588912SDave May 
56d8588912SDave May   if (isANest && isxNest && isyNest) {
57d8588912SDave May     PetscFunctionReturn(0);
58d8588912SDave May   } else {
59d8588912SDave May     SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_SUP, "Operation only valid for Mat(Nest)-Vec(Nest) combinations");
60d8588912SDave May     PetscFunctionReturn(0);
61d8588912SDave May   }
62d8588912SDave May   PetscFunctionReturn(0);
63d8588912SDave May }
64d8588912SDave May 
65d8588912SDave May /*
66d8588912SDave May  This function is called every time we insert a sub matrix the Nest.
67d8588912SDave May  It ensures that every Nest along row r, and col c has the same dimensions
68d8588912SDave May  as the submat being inserted.
69d8588912SDave May */
70d8588912SDave May #undef __FUNCT__
71d8588912SDave May #define __FUNCT__ "PETSc_MatNest_CheckConsistency"
72d8588912SDave May static PetscErrorCode PETSc_MatNest_CheckConsistency(Mat A,Mat submat,PetscInt r,PetscInt c)
73d8588912SDave May {
74d8588912SDave May   Mat_Nest       *b = (Mat_Nest*)A->data;
75d8588912SDave May   PetscInt       i,j;
76d8588912SDave May   PetscInt       nr,nc;
77d8588912SDave May   PetscInt       sM,sN,M,N;
78d8588912SDave May   Mat            mat;
79d8588912SDave May   PetscErrorCode ierr;
80d8588912SDave May 
81d8588912SDave May   PetscFunctionBegin;
82d8588912SDave May   nr = b->nr;
83d8588912SDave May   nc = b->nc;
84d8588912SDave May   ierr = MatGetSize(submat,&sM,&sN);CHKERRQ(ierr);
85d8588912SDave May   for (i=0; i<nr; i++) {
86d8588912SDave May     mat = b->m[i][c];
87d8588912SDave May     if (mat) {
88d8588912SDave May       ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
89d8588912SDave May       /* Check columns match */
90d8588912SDave May       if (sN != N) {
91d8588912SDave May         SETERRQ3(((PetscObject)A)->comm,PETSC_ERR_SUP,"Inserted incompatible submatrix into Nest at (%D,%D). Submatrix must have %D rows",r,c,N );
92d8588912SDave May       }
93d8588912SDave May     }
94d8588912SDave May   }
95d8588912SDave May 
96d8588912SDave May   for (j=0; j<nc; j++) {
97d8588912SDave May     mat = b->m[r][j];
98d8588912SDave May     if (mat) {
99d8588912SDave May       ierr = MatGetSize(mat,&M,&N);CHKERRQ(ierr);
100d8588912SDave May       /* Check rows match */
101d8588912SDave May       if (sM != M) {
102d8588912SDave May         SETERRQ3(((PetscObject)A)->comm,PETSC_ERR_SUP,"Inserted incompatible submatrix into Nest at (%D,%D). Submatrix must have %D cols",r,c,M );
103d8588912SDave May       }
104d8588912SDave May     }
105d8588912SDave May   }
106d8588912SDave May   PetscFunctionReturn(0);
107d8588912SDave May }
108d8588912SDave May 
109d8588912SDave May /*
110d8588912SDave May  Checks if the row/col's contain a complete set of IS's.
111d8588912SDave May  Once they are filled, the offsets are computed. This saves having to update
112d8588912SDave May  the index set entries each time we insert something new.
113d8588912SDave May  */
114d8588912SDave May #undef __FUNCT__
115d8588912SDave May #define __FUNCT__ "PETSc_MatNest_UpdateStructure"
116d8588912SDave May static PetscErrorCode PETSc_MatNest_UpdateStructure(Mat A,PetscInt ridx,PetscInt cidx)
117d8588912SDave May {
118d8588912SDave May   Mat_Nest       *b = (Mat_Nest*)A->data;
119d8588912SDave May   PetscInt       m,n,i,j;
120d8588912SDave May   PetscBool      fullrow,fullcol;
121d8588912SDave May   PetscErrorCode ierr;
122d8588912SDave May 
123d8588912SDave May   PetscFunctionBegin;
124d8588912SDave May   ierr = MatGetLocalSize(b->m[ridx][cidx],&m,&n);CHKERRQ(ierr);
125d8588912SDave May   b->row_len[ridx] = m;
126d8588912SDave May   b->col_len[cidx] = n;
127d8588912SDave May 
128d8588912SDave May   fullrow = PETSC_TRUE;
129d8588912SDave May   for (i=0; i<b->nr; i++) {
130d8588912SDave May     if (b->row_len[i]==-1) { fullrow = PETSC_FALSE; break; }
131d8588912SDave May   }
132d8588912SDave May   if ( (fullrow) && (!b->is_row[0]) ){
133d8588912SDave May     PetscInt cnt = 0;
134d8588912SDave May     for (i=0; i<b->nr; i++) {
135d8588912SDave May       ierr = ISCreateStride(PETSC_COMM_SELF,b->row_len[i],cnt,1,&b->is_row[i]);CHKERRQ(ierr);
136d8588912SDave May       cnt = cnt + b->row_len[i];
137d8588912SDave May     }
138d8588912SDave May   }
139d8588912SDave May 
140d8588912SDave May   fullcol = PETSC_TRUE;
141d8588912SDave May   for (j=0; j<b->nc; j++) {
142d8588912SDave May     if (b->col_len[j]==-1) { fullcol = PETSC_FALSE; break; }
143d8588912SDave May   }
144d8588912SDave May   if( (fullcol) && (!b->is_col[0]) ){
145d8588912SDave May     PetscInt cnt = 0;
146d8588912SDave May     for (j=0; j<b->nc; j++) {
147d8588912SDave May       ierr = ISCreateStride(PETSC_COMM_SELF,b->col_len[j],cnt,1,&b->is_col[j]);CHKERRQ(ierr);
148d8588912SDave May       cnt = cnt + b->col_len[j];
149d8588912SDave May     }
150d8588912SDave May   }
151d8588912SDave May   PetscFunctionReturn(0);
152d8588912SDave May }
153d8588912SDave May 
154d8588912SDave May /* operations */
155d8588912SDave May #undef __FUNCT__
156d8588912SDave May #define __FUNCT__ "MatMult_Nest"
157d8588912SDave May PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y)
158d8588912SDave May {
159d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
160d8588912SDave May   Vec            *bx;
161d8588912SDave May   Vec            *by;
162d8588912SDave May   PetscInt       i,j;
163d8588912SDave May   PetscErrorCode ierr;
164d8588912SDave May 
165d8588912SDave May   PetscFunctionBegin;
166d8588912SDave May   ierr = PETSc_MatNest_CheckMatVecCompatibility2(A,x,y);CHKERRQ(ierr);
167d8588912SDave May   ierr = VecNestGetSubVecs(x,PETSC_NULL,&bx);CHKERRQ(ierr);
168d8588912SDave May   ierr = VecNestGetSubVecs(y,PETSC_NULL,&by);CHKERRQ(ierr);
169d8588912SDave May 
170d8588912SDave May   for (i=0; i<bA->nr; i++) {
171d8588912SDave May     ierr = VecZeroEntries(by[i]);CHKERRQ(ierr);
172d8588912SDave May     for (j=0; j<bA->nc; j++) {
173d8588912SDave May       if (!bA->m[i][j]) {
174d8588912SDave May         continue;
175d8588912SDave May       }
176d8588912SDave May       /* y[i] <- y[i] + A[i][j] * x[j] */
177d8588912SDave May       ierr = MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);CHKERRQ(ierr);
178d8588912SDave May     }
179d8588912SDave May   }
180d8588912SDave May   PetscFunctionReturn(0);
181d8588912SDave May }
182d8588912SDave May 
183d8588912SDave May #undef __FUNCT__
184d8588912SDave May #define __FUNCT__ "MatMultTranspose_Nest"
185d8588912SDave May PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y)
186d8588912SDave May {
187d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
188d8588912SDave May   Vec            *bx;
189d8588912SDave May   Vec            *by;
190d8588912SDave May   PetscInt       i,j;
191d8588912SDave May   PetscErrorCode ierr;
192d8588912SDave May 
193d8588912SDave May   PetscFunctionBegin;
194d8588912SDave May   ierr = PETSc_MatNest_CheckMatVecCompatibility2(A,x,y);CHKERRQ(ierr);
195d8588912SDave May   ierr = VecNestGetSubVecs(x,PETSC_NULL,&bx);CHKERRQ(ierr);
196d8588912SDave May   ierr = VecNestGetSubVecs(y,PETSC_NULL,&by);CHKERRQ(ierr);
197d8588912SDave May   if (A->symmetric) {
198d8588912SDave May     ierr = MatMult_Nest(A,x,y);CHKERRQ(ierr);
199d8588912SDave May     PetscFunctionReturn(0);
200d8588912SDave May   }
201d8588912SDave May   for (i=0; i<bA->nr; i++) {
202d8588912SDave May     ierr = VecZeroEntries(by[i]);CHKERRQ(ierr);
203d8588912SDave May     for (j=0; j<bA->nc; j++) {
204d8588912SDave May       if (!bA->m[j][i]) {
205d8588912SDave May         continue;
206d8588912SDave May       }
207d8588912SDave May       /* y[i] <- y[i] + A^T[i][j] * x[j], so we swap i,j in mat[][] */
208d8588912SDave May       ierr = MatMultTransposeAdd(bA->m[j][i],bx[j],by[i],by[i]);CHKERRQ(ierr);
209d8588912SDave May     }
210d8588912SDave May   }
211d8588912SDave May   PetscFunctionReturn(0);
212d8588912SDave May }
213d8588912SDave May 
214d8588912SDave May /* Returns the sum of the global size of all the consituent vectors in the nest */
215d8588912SDave May #undef __FUNCT__
216d8588912SDave May #define __FUNCT__ "MatGetSize_Nest"
217d8588912SDave May PetscErrorCode MatGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
218d8588912SDave May {
219d8588912SDave May   PetscFunctionBegin;
220d8588912SDave May   if (M) { *M = A->rmap->N; }
221d8588912SDave May   if (N) { *N = A->cmap->N; }
222d8588912SDave May   PetscFunctionReturn(0);
223d8588912SDave May }
224d8588912SDave May 
225d8588912SDave May /* Returns the sum of the local size of all the consituent vectors in the nest */
226d8588912SDave May #undef __FUNCT__
227d8588912SDave May #define __FUNCT__ "MatGetLocalSize_Nest"
228d8588912SDave May PetscErrorCode MatGetLocalSize_Nest(Mat A,PetscInt *m,PetscInt *n)
229d8588912SDave May {
230d8588912SDave May   PetscFunctionBegin;
231d8588912SDave May   if (m) { *m = A->rmap->n; }
232d8588912SDave May   if (n) { *n = A->cmap->n; }
233d8588912SDave May   PetscFunctionReturn(0);
234d8588912SDave May }
235d8588912SDave May 
236d8588912SDave May #undef __FUNCT__
237d8588912SDave May #define __FUNCT__ "MatDestroy_Nest"
238d8588912SDave May PetscErrorCode MatDestroy_Nest(Mat A)
239d8588912SDave May {
240d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
241d8588912SDave May   PetscInt       i,j;
242d8588912SDave May   PetscErrorCode ierr;
243d8588912SDave May 
244d8588912SDave May   PetscFunctionBegin;
245d8588912SDave May   /* release the matrices and the place holders */
246d8588912SDave May   if (vs->is_row) {
247d8588912SDave May     for (i=0; i<vs->nr; i++) {
248d8588912SDave May       if(vs->is_row[i]) ierr = ISDestroy(vs->is_row[i]);CHKERRQ(ierr);
249d8588912SDave May     }
250d8588912SDave May     ierr = PetscFree(vs->is_row);CHKERRQ(ierr);
251d8588912SDave May   }
252d8588912SDave May   if (vs->is_col) {
253d8588912SDave May     for (j=0; j<vs->nc; j++) {
254d8588912SDave May       if(vs->is_col[j]) ierr = ISDestroy(vs->is_col[j]);CHKERRQ(ierr);
255d8588912SDave May     }
256d8588912SDave May     ierr = PetscFree(vs->is_col);CHKERRQ(ierr);
257d8588912SDave May   }
258d8588912SDave May 
259d8588912SDave May   ierr = PetscFree(vs->row_len);CHKERRQ(ierr);
260d8588912SDave May   ierr = PetscFree(vs->col_len);CHKERRQ(ierr);
261d8588912SDave May 
262d8588912SDave May   /* release the matrices and the place holders */
263d8588912SDave May   if (vs->m) {
264d8588912SDave May     for (i=0; i<vs->nr; i++) {
265d8588912SDave May       for (j=0; j<vs->nc; j++) {
266d8588912SDave May 
267d8588912SDave May         if (vs->m[i][j]) {
268d8588912SDave May           ierr = MatDestroy(vs->m[i][j]);CHKERRQ(ierr);
269d8588912SDave May           vs->m[i][j] = PETSC_NULL;
270d8588912SDave May         }
271d8588912SDave May       }
272d8588912SDave May       ierr = PetscFree( vs->m[i] );CHKERRQ(ierr);
273d8588912SDave May       vs->m[i] = PETSC_NULL;
274d8588912SDave May     }
275d8588912SDave May     ierr = PetscFree(vs->m);CHKERRQ(ierr);
276d8588912SDave May     vs->m = PETSC_NULL;
277d8588912SDave May   }
278d8588912SDave May   ierr = PetscFree(vs);CHKERRQ(ierr);
279d8588912SDave May 
280d8588912SDave May   PetscFunctionReturn(0);
281d8588912SDave May }
282d8588912SDave May 
283d8588912SDave May #undef __FUNCT__
284d8588912SDave May #define __FUNCT__ "MatAssemblyBegin_Nest"
285d8588912SDave May PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
286d8588912SDave May {
287d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
288d8588912SDave May   PetscInt       i,j;
289d8588912SDave May   PetscErrorCode ierr;
290d8588912SDave May 
291d8588912SDave May   PetscFunctionBegin;
292d8588912SDave May   for (i=0; i<vs->nr; i++) {
293d8588912SDave May     for (j=0; j<vs->nc; j++) {
294d8588912SDave May       if (vs->m[i][j]) { ierr = MatAssemblyBegin(vs->m[i][j],type);CHKERRQ(ierr); }
295d8588912SDave May     }
296d8588912SDave May   }
297d8588912SDave May   PetscFunctionReturn(0);
298d8588912SDave May }
299d8588912SDave May 
300d8588912SDave May #undef __FUNCT__
301d8588912SDave May #define __FUNCT__ "MatAssemblyEnd_Nest"
302d8588912SDave May PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
303d8588912SDave May {
304d8588912SDave May   Mat_Nest       *vs = (Mat_Nest*)A->data;
305d8588912SDave May   PetscInt       i,j;
306d8588912SDave May   PetscErrorCode ierr;
307d8588912SDave May 
308d8588912SDave May   PetscFunctionBegin;
309d8588912SDave May   for (i=0; i<vs->nr; i++) {
310d8588912SDave May     for (j=0; j<vs->nc; j++) {
311d8588912SDave May       if (vs->m[i][j]) { ierr = MatAssemblyEnd(vs->m[i][j],type);CHKERRQ(ierr); }
312d8588912SDave May     }
313d8588912SDave May   }
314d8588912SDave May   PetscFunctionReturn(0);
315d8588912SDave May }
316d8588912SDave May 
317d8588912SDave May #undef __FUNCT__
318d8588912SDave May #define __FUNCT__ "MatGetLocalSubMatrix_Nest"
319d8588912SDave May PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS irow,IS icol,Mat *sub)
320d8588912SDave May {
321d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
322d8588912SDave May   PetscInt       i,j;
323d8588912SDave May   PetscBool      found_row,found_col;
324d8588912SDave May   PetscInt       row=-1,col=-1;
325d8588912SDave May   PetscErrorCode ierr;
326d8588912SDave May 
327d8588912SDave May   PetscFunctionBegin;
328d8588912SDave May   found_row = PETSC_FALSE;
329d8588912SDave May   for (i=0; i<bA->nr; i++) {
330d8588912SDave May     ierr = ISEqual(irow,bA->is_row[i],&found_row);CHKERRQ(ierr);
331d8588912SDave May     if(found_row){ row = i; break; }
332d8588912SDave May   }
333d8588912SDave May   found_col = PETSC_FALSE;
334d8588912SDave May   for (j=0; j<bA->nc; j++) {
335d8588912SDave May     ierr = ISEqual(icol,bA->is_col[j],&found_col);CHKERRQ(ierr);
336d8588912SDave May     if(found_col){ col = j; break; }
337d8588912SDave May   }
338d8588912SDave May   /* check valid i,j */
339d8588912SDave May   if ((row<0)||(col<0)) {
340d8588912SDave May     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Mat(Nest) must contain at least one matrix within each row and column");
341d8588912SDave May   }
342d8588912SDave May   if ((row>=bA->nr)||(col>=bA->nc)) {
343d8588912SDave May     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_USER,"Mat(Nest) row and column is too large");
344d8588912SDave May   }
345d8588912SDave May 
346d8588912SDave May   *sub = bA->m[row][col];
347d8588912SDave May   if (bA->m[row][col]) {
348d8588912SDave May     ierr = PetscObjectReference( (PetscObject)bA->m[row][col] );CHKERRQ(ierr);
349d8588912SDave May   }
350d8588912SDave May 
351d8588912SDave May   PetscFunctionReturn(0);
352d8588912SDave May }
353d8588912SDave May 
354d8588912SDave May #undef __FUNCT__
355d8588912SDave May #define __FUNCT__ "MatRestoreLocalSubMatrix_Nest"
356d8588912SDave May PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS irow,IS icol,Mat *sub)
357d8588912SDave May {
358d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
359d8588912SDave May   PetscInt       i,j;
360d8588912SDave May   PetscBool      found_row,found_col;
361d8588912SDave May   PetscInt       row=-1,col=-1;
362d8588912SDave May   PetscErrorCode ierr;
363d8588912SDave May 
364d8588912SDave May   PetscFunctionBegin;
365d8588912SDave May   found_row = PETSC_FALSE;
366d8588912SDave May   for (i=0; i<bA->nr; i++) {
367d8588912SDave May     ierr = ISEqual(irow,bA->is_row[i],&found_row);CHKERRQ(ierr);
368d8588912SDave May     if (found_row){ row = i; break; }
369d8588912SDave May   }
370d8588912SDave May   found_col = PETSC_FALSE;
371d8588912SDave May   for (j=0; j<bA->nc; j++) {
372d8588912SDave May     ierr = ISEqual(icol,bA->is_col[j],&found_col);CHKERRQ(ierr);
373d8588912SDave May     if (found_col){ col = j; break; }
374d8588912SDave May   }
375d8588912SDave May   /* check valid i,j */
376d8588912SDave May   if ((row<0)||(col<0)) {
377d8588912SDave May     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_SUP,"Mat(Nest) must contain at least one matrix within each row and column");
378d8588912SDave May   }
379d8588912SDave May   if ((row>=bA->nr)||(col>=bA->nc)) {
380d8588912SDave May     SETERRQ(((PetscObject)A)->comm,PETSC_ERR_USER,"Mat(Nest) row and column is too large");
381d8588912SDave May   }
382d8588912SDave May 
383d8588912SDave May   if (*sub) {
384d8588912SDave May     ierr = MatDestroy(*sub);CHKERRQ(ierr);
385d8588912SDave May   }
386d8588912SDave May 
387d8588912SDave May   PetscFunctionReturn(0);
388d8588912SDave May }
389d8588912SDave May 
390d8588912SDave May #undef __FUNCT__
391d8588912SDave May #define __FUNCT__ "MatGetVecs_Nest"
392d8588912SDave May PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left)
393d8588912SDave May {
394d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
395d8588912SDave May   Vec            *L,*R;
396d8588912SDave May   MPI_Comm       comm;
397d8588912SDave May   PetscInt       i,j;
398d8588912SDave May   PetscErrorCode ierr;
399d8588912SDave May 
400d8588912SDave May   PetscFunctionBegin;
401d8588912SDave May   comm = ((PetscObject)A)->comm;
402d8588912SDave May   if (right) {
403d8588912SDave May     /* allocate R */
404d8588912SDave May     ierr = PetscMalloc( sizeof(Vec) * bA->nc, &R );CHKERRQ(ierr);
405d8588912SDave May     /* Create the right vectors */
406d8588912SDave May     for (j=0; j<bA->nc; j++) {
407d8588912SDave May       for (i=0; i<bA->nr; i++) {
408d8588912SDave May         if (bA->m[i][j]) {
409d8588912SDave May           ierr = MatGetVecs(bA->m[i][j],&R[j],PETSC_NULL);CHKERRQ(ierr);
410d8588912SDave May           break;
411d8588912SDave May         }
412d8588912SDave May       }
413d8588912SDave May       if (i==bA->nr) {
414d8588912SDave May         /* have an empty column */
415d8588912SDave May         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
416d8588912SDave May       }
417d8588912SDave May     }
418d8588912SDave May     ierr = VecCreateNest(comm,bA->nc,bA->is_col,R,right);CHKERRQ(ierr);
419d8588912SDave May     /* hand back control to the nest vector */
420d8588912SDave May     for (j=0; j<bA->nc; j++) {
421d8588912SDave May       ierr = VecDestroy(R[j]);CHKERRQ(ierr);
422d8588912SDave May     }
423d8588912SDave May     ierr = PetscFree(R);CHKERRQ(ierr);
424d8588912SDave May   }
425d8588912SDave May 
426d8588912SDave May   if (left) {
427d8588912SDave May     /* allocate L */
428d8588912SDave May     ierr = PetscMalloc( sizeof(Vec) * bA->nr, &L );CHKERRQ(ierr);
429d8588912SDave May     /* Create the left vectors */
430d8588912SDave May     for (i=0; i<bA->nr; i++) {
431d8588912SDave May       for (j=0; j<bA->nc; j++) {
432d8588912SDave May         if (bA->m[i][j]) {
433d8588912SDave May           ierr = MatGetVecs(bA->m[i][j],PETSC_NULL,&L[i]);CHKERRQ(ierr);
434d8588912SDave May           break;
435d8588912SDave May         }
436d8588912SDave May       }
437d8588912SDave May       if (j==bA->nc) {
438d8588912SDave May         /* have an empty row */
439d8588912SDave May         SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
440d8588912SDave May       }
441d8588912SDave May     }
442d8588912SDave May 
443d8588912SDave May     ierr = VecCreateNest(comm,bA->nr,bA->is_row,L,left);CHKERRQ(ierr);
444d8588912SDave May     for (i=0; i<bA->nr; i++) {
445d8588912SDave May       ierr = VecDestroy(L[i]);CHKERRQ(ierr);
446d8588912SDave May     }
447d8588912SDave May 
448d8588912SDave May     ierr = PetscFree(L);CHKERRQ(ierr);
449d8588912SDave May   }
450d8588912SDave May   PetscFunctionReturn(0);
451d8588912SDave May }
452d8588912SDave May 
453d8588912SDave May #undef __FUNCT__
454d8588912SDave May #define __FUNCT__ "MatView_Nest"
455d8588912SDave May PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
456d8588912SDave May {
457d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
458d8588912SDave May   PetscBool      isascii;
459d8588912SDave May   PetscInt       i,j;
460d8588912SDave May   PetscErrorCode ierr;
461d8588912SDave May 
462d8588912SDave May   PetscFunctionBegin;
463d8588912SDave May   ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);CHKERRQ(ierr);
464d8588912SDave May   if (isascii) {
465d8588912SDave May 
466d8588912SDave May     PetscViewerASCIIPrintf(viewer,"Matrix object: \n" );
467d8588912SDave May     PetscViewerASCIIPushTab( viewer );    /* push0 */
468d8588912SDave May     PetscViewerASCIIPrintf( viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc);
469d8588912SDave May 
470d8588912SDave May     PetscViewerASCIIPrintf(viewer,"MatNest structure: \n" );
471d8588912SDave May     for (i=0; i<bA->nr; i++) {
472d8588912SDave May       for (j=0; j<bA->nc; j++) {
473d8588912SDave May         const MatType type;
474d8588912SDave May         const char *name;
475d8588912SDave May         PetscInt NR,NC;
476d8588912SDave May         PetscBool isNest = PETSC_FALSE;
477d8588912SDave May 
478d8588912SDave May         if (!bA->m[i][j]) {
479d8588912SDave May           PetscViewerASCIIPrintf( viewer, "(%D,%D) : PETSC_NULL \n",i,j);
480d8588912SDave May           continue;
481d8588912SDave May         }
482d8588912SDave May         ierr = MatGetSize(bA->m[i][j],&NR,&NC);CHKERRQ(ierr);
483d8588912SDave May         ierr = MatGetType( bA->m[i][j], &type );CHKERRQ(ierr);
484d8588912SDave May         name = ((PetscObject)bA->m[i][j])->prefix;
485d8588912SDave May         ierr = PetscTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);CHKERRQ(ierr);
486d8588912SDave May 
487d8588912SDave May         PetscViewerASCIIPrintf( viewer, "(%D,%D) : name=\"%s\", type=%s, rows=%D, cols=%D \n",i,j,name,type,NR,NC);
488d8588912SDave May 
489d8588912SDave May         if (isNest) {
490d8588912SDave May           PetscViewerASCIIPushTab( viewer );  /* push1 */
491d8588912SDave May           ierr = MatView( bA->m[i][j], viewer );CHKERRQ(ierr);
492d8588912SDave May           PetscViewerASCIIPopTab(viewer);    /* pop1 */
493d8588912SDave May         }
494d8588912SDave May       }
495d8588912SDave May     }
496d8588912SDave May     PetscViewerASCIIPopTab(viewer);    /* pop0 */
497d8588912SDave May   }
498d8588912SDave May   PetscFunctionReturn(0);
499d8588912SDave May }
500d8588912SDave May 
501d8588912SDave May #undef __FUNCT__
502d8588912SDave May #define __FUNCT__ "MatZeroEntries_Nest"
503d8588912SDave May PetscErrorCode MatZeroEntries_Nest(Mat A)
504d8588912SDave May {
505d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
506d8588912SDave May   PetscInt       i,j;
507d8588912SDave May   PetscErrorCode ierr;
508d8588912SDave May 
509d8588912SDave May   PetscFunctionBegin;
510d8588912SDave May   for (i=0; i<bA->nr; i++) {
511d8588912SDave May     for (j=0; j<bA->nc; j++) {
512d8588912SDave May       if (!bA->m[i][j]) continue;
513d8588912SDave May       ierr = MatZeroEntries(bA->m[i][j]);CHKERRQ(ierr);
514d8588912SDave May     }
515d8588912SDave May   }
516d8588912SDave May   PetscFunctionReturn(0);
517d8588912SDave May }
518d8588912SDave May 
519d8588912SDave May #undef __FUNCT__
520d8588912SDave May #define __FUNCT__ "MatDuplicate_Nest"
521d8588912SDave May PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
522d8588912SDave May {
523d8588912SDave May   Mat_Nest       *bA = (Mat_Nest*)A->data;
524*841e96a3SJed Brown   Mat            *b;
525*841e96a3SJed Brown   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
526d8588912SDave May   PetscErrorCode ierr;
527d8588912SDave May 
528d8588912SDave May   PetscFunctionBegin;
529*841e96a3SJed Brown   ierr = PetscMalloc(nr*nc*sizeof(Mat),&b);CHKERRQ(ierr);
530*841e96a3SJed Brown   for (i=0; i<nr; i++) {
531*841e96a3SJed Brown     for (j=0; j<nc; j++) {
532*841e96a3SJed Brown       if (bA->m[i][j]) {
533*841e96a3SJed Brown         ierr = MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);CHKERRQ(ierr);
534*841e96a3SJed Brown       } else {
535*841e96a3SJed Brown         b[i*nc+j] = PETSC_NULL;
536d8588912SDave May       }
537d8588912SDave May     }
538d8588912SDave May   }
539*841e96a3SJed Brown   ierr = MatCreateNest(((PetscObject)A)->comm,nr,bA->is_row,nc,bA->is_col,b,B);CHKERRQ(ierr);
540*841e96a3SJed Brown   /* Give the new MatNest exclusive ownership */
541*841e96a3SJed Brown   for (i=0; i<nr*nc; i++) {
542*841e96a3SJed Brown     if (b[i]) {ierr = MatDestroy(b[i]);CHKERRQ(ierr);}
543d8588912SDave May   }
544d8588912SDave May   ierr = PetscFree(b);CHKERRQ(ierr);
545d8588912SDave May 
546*841e96a3SJed Brown   ierr = MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
547*841e96a3SJed Brown   ierr = MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
548d8588912SDave May   PetscFunctionReturn(0);
549d8588912SDave May }
550d8588912SDave May 
551d8588912SDave May /* nest api */
552d8588912SDave May EXTERN_C_BEGIN
553d8588912SDave May #undef __FUNCT__
554d8588912SDave May #define __FUNCT__ "MatNestGetSubMat_Nest"
555d8588912SDave May PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
556d8588912SDave May {
557d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
558d8588912SDave May   PetscFunctionBegin;
559d8588912SDave May   if (idxm >= bA->nr) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
560d8588912SDave May   if (jdxm >= bA->nc) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
561d8588912SDave May   *mat = bA->m[idxm][jdxm];
562d8588912SDave May   PetscFunctionReturn(0);
563d8588912SDave May }
564d8588912SDave May EXTERN_C_END
565d8588912SDave May 
566d8588912SDave May #undef __FUNCT__
567d8588912SDave May #define __FUNCT__ "MatNestGetSubMat"
568d8588912SDave May /*@C
569d8588912SDave May  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.
570d8588912SDave May 
571d8588912SDave May  Not collective
572d8588912SDave May 
573d8588912SDave May  Input Parameters:
574d8588912SDave May  .  A  - nest matrix
575d8588912SDave May  .  idxm - index of the matrix within the nest matrix
576d8588912SDave May  .  jdxm - index of the matrix within the nest matrix
577d8588912SDave May 
578d8588912SDave May  Output Parameter:
579d8588912SDave May  .  sub - matrix at index idxm,jdxm within the nest matrix
580d8588912SDave May 
581d8588912SDave May  Notes:
582d8588912SDave May 
583d8588912SDave May  Level: developer
584d8588912SDave May 
585d8588912SDave May  .seealso: MatNestGetSize(), MatNestGetSubMats()
586d8588912SDave May @*/
587d8588912SDave May PetscErrorCode PETSCMAT_DLLEXPORT MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
588d8588912SDave May {
589d8588912SDave May   PetscErrorCode ierr,(*f)(Mat,PetscInt,PetscInt,Mat*);
590d8588912SDave May 
591d8588912SDave May   PetscFunctionBegin;
592d8588912SDave May   ierr = PetscObjectQueryFunction((PetscObject)A,"MatNestGetSubMat_C",(void (**)(void))&f);CHKERRQ(ierr);
593d8588912SDave May   if (f) {
594d8588912SDave May     ierr = (*f)(A,idxm,jdxm,sub);CHKERRQ(ierr);
595d8588912SDave May   }
596d8588912SDave May   PetscFunctionReturn(0);
597d8588912SDave May }
598d8588912SDave May 
599d8588912SDave May EXTERN_C_BEGIN
600d8588912SDave May #undef __FUNCT__
601d8588912SDave May #define __FUNCT__ "MatNestGetSubMats_Nest"
602d8588912SDave May PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
603d8588912SDave May {
604d8588912SDave May   Mat_Nest *bA = (Mat_Nest*)A->data;
605d8588912SDave May   PetscFunctionBegin;
606d8588912SDave May   if (M)   { *M   = bA->nr; }
607d8588912SDave May   if (N)   { *N   = bA->nc; }
608d8588912SDave May   if (mat) { *mat = bA->m;  }
609d8588912SDave May   PetscFunctionReturn(0);
610d8588912SDave May }
611d8588912SDave May EXTERN_C_END
612d8588912SDave May 
613d8588912SDave May #undef __FUNCT__
614d8588912SDave May #define __FUNCT__ "MatNestGetSubMats"
615d8588912SDave May /*@C
616d8588912SDave May  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.
617d8588912SDave May 
618d8588912SDave May  Not collective
619d8588912SDave May 
620d8588912SDave May  Input Parameters:
621d8588912SDave May  .  A  - nest matri
622d8588912SDave May 
623d8588912SDave May  Output Parameter:
624d8588912SDave May  .  M - number of rows in the nest matrix
625d8588912SDave May  .  N - number of cols in the nest matrix
626d8588912SDave May  .  mat - 2d array of matrices
627d8588912SDave May 
628d8588912SDave May  Notes:
629d8588912SDave May 
630d8588912SDave May  The user should not free the array mat.
631d8588912SDave May 
632d8588912SDave May  Level: developer
633d8588912SDave May 
634d8588912SDave May  .seealso: MatNestGetSize(), MatNestGetSubMat()
635d8588912SDave May @*/
636d8588912SDave May PetscErrorCode PETSCMAT_DLLEXPORT MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
637d8588912SDave May {
638d8588912SDave May   PetscErrorCode ierr,(*f)(Mat,PetscInt*,PetscInt*,Mat***);
639d8588912SDave May 
640d8588912SDave May   PetscFunctionBegin;
641d8588912SDave May   ierr = PetscObjectQueryFunction((PetscObject)A,"MatNestGetSubMats_C",(void (**)(void))&f);CHKERRQ(ierr);
642d8588912SDave May   if (f) {
643d8588912SDave May     ierr = (*f)(A,M,N,mat);CHKERRQ(ierr);
644d8588912SDave May   }
645d8588912SDave May   PetscFunctionReturn(0);
646d8588912SDave May }
647d8588912SDave May 
648d8588912SDave May EXTERN_C_BEGIN
649d8588912SDave May #undef __FUNCT__
650d8588912SDave May #define __FUNCT__ "MatNestGetSize_Nest"
651d8588912SDave May PetscErrorCode PETSCMAT_DLLEXPORT MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
652d8588912SDave May {
653d8588912SDave May   Mat_Nest  *bA = (Mat_Nest*)A->data;
654d8588912SDave May 
655d8588912SDave May   PetscFunctionBegin;
656d8588912SDave May   if (M) { *M  = bA->nr; }
657d8588912SDave May   if (N) { *N  = bA->nc; }
658d8588912SDave May   PetscFunctionReturn(0);
659d8588912SDave May }
660d8588912SDave May EXTERN_C_END
661d8588912SDave May 
662d8588912SDave May #undef __FUNCT__
663d8588912SDave May #define __FUNCT__ "MatNestGetSize"
664d8588912SDave May /*@C
665d8588912SDave May  MatNestGetSize - Returns the size of the nest matrix.
666d8588912SDave May 
667d8588912SDave May  Not collective
668d8588912SDave May 
669d8588912SDave May  Input Parameters:
670d8588912SDave May  .  A  - nest matrix
671d8588912SDave May 
672d8588912SDave May  Output Parameter:
673d8588912SDave May  .  M - number of rows in the nested mat
674d8588912SDave May  .  N - number of cols in the nested mat
675d8588912SDave May 
676d8588912SDave May  Notes:
677d8588912SDave May 
678d8588912SDave May  Level: developer
679d8588912SDave May 
680d8588912SDave May  .seealso: MatNestGetSubMat(), MatNestGetSubMats()
681d8588912SDave May @*/
682d8588912SDave May PetscErrorCode PETSCMAT_DLLEXPORT MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
683d8588912SDave May {
684d8588912SDave May   PetscErrorCode ierr,(*f)(Mat,PetscInt*,PetscInt*);
685d8588912SDave May 
686d8588912SDave May   PetscFunctionBegin;
687d8588912SDave May   ierr = PetscObjectQueryFunction((PetscObject)A,"MatNestGetSize_C",(void (**)(void))&f);CHKERRQ(ierr);
688d8588912SDave May   if (f) {
689d8588912SDave May     ierr = (*f)(A,M,N);CHKERRQ(ierr);
690d8588912SDave May   }
691d8588912SDave May   PetscFunctionReturn(0);
692d8588912SDave May }
693d8588912SDave May 
694d8588912SDave May /* constructor */
695d8588912SDave May #undef __FUNCT__
696d8588912SDave May #define __FUNCT__ "MatNestSetOps_Private"
697d8588912SDave May static PetscErrorCode MatNestSetOps_Private(struct _MatOps* ops)
698d8588912SDave May {
699d8588912SDave May   PetscFunctionBegin;
700d8588912SDave May   /* 0*/
701d8588912SDave May   ops->setvalues  = 0;
702d8588912SDave May   ops->getrow     = 0;
703d8588912SDave May   ops->restorerow = 0;
704d8588912SDave May   ops->mult       = MatMult_Nest;
705d8588912SDave May   ops->multadd    = 0;
706d8588912SDave May   /* 5*/
707d8588912SDave May   ops->multtranspose    = MatMultTranspose_Nest;
708d8588912SDave May   ops->multtransposeadd = 0;
709d8588912SDave May   ops->solve            = 0;
710d8588912SDave May   ops->solveadd         = 0;
711d8588912SDave May   ops->solvetranspose   = 0;
712d8588912SDave May   /*10*/
713d8588912SDave May   ops->solvetransposeadd = 0;
714d8588912SDave May   ops->lufactor          = 0;
715d8588912SDave May   ops->choleskyfactor    = 0;
716d8588912SDave May   ops->sor               = 0;
717d8588912SDave May   ops->transpose         = 0;
718d8588912SDave May   /*15*/
719d8588912SDave May   ops->getinfo       = 0;
720d8588912SDave May   ops->equal         = 0;
721d8588912SDave May   ops->getdiagonal   = 0;
722d8588912SDave May   ops->diagonalscale = 0;
723d8588912SDave May   ops->norm          = 0;
724d8588912SDave May   /*20*/
725d8588912SDave May   ops->assemblybegin = MatAssemblyBegin_Nest;
726d8588912SDave May   ops->assemblyend   = MatAssemblyEnd_Nest;
727d8588912SDave May   ops->setoption     = 0;
728d8588912SDave May   ops->zeroentries   = MatZeroEntries_Nest;
729d8588912SDave May   /*24*/
730d8588912SDave May   ops->zerorows               = 0;
731d8588912SDave May   ops->lufactorsymbolic       = 0;
732d8588912SDave May   ops->lufactornumeric        = 0;
733d8588912SDave May   ops->choleskyfactorsymbolic = 0;
734d8588912SDave May   ops->choleskyfactornumeric  = 0;
735d8588912SDave May   /*29*/
736d8588912SDave May   ops->setuppreallocation = 0;
737d8588912SDave May   ops->ilufactorsymbolic  = 0;
738d8588912SDave May   ops->iccfactorsymbolic  = 0;
739d8588912SDave May   ops->getarray           = 0;
740d8588912SDave May   ops->restorearray       = 0;
741d8588912SDave May   /*34*/
742d8588912SDave May   ops->duplicate     = MatDuplicate_Nest;
743d8588912SDave May   ops->forwardsolve  = 0;
744d8588912SDave May   ops->backwardsolve = 0;
745d8588912SDave May   ops->ilufactor     = 0;
746d8588912SDave May   ops->iccfactor     = 0;
747d8588912SDave May   /*39*/
748d8588912SDave May   ops->axpy            = 0;
749d8588912SDave May   ops->getsubmatrices  = 0;
750d8588912SDave May   ops->increaseoverlap = 0;
751d8588912SDave May   ops->getvalues       = 0;
752d8588912SDave May   ops->copy            = 0;
753d8588912SDave May   /*44*/
754d8588912SDave May   ops->getrowmax   = 0;
755d8588912SDave May   ops->scale       = 0;
756d8588912SDave May   ops->shift       = 0;
757d8588912SDave May   ops->diagonalset = 0;
758d8588912SDave May   ops->zerorowscolumns  = 0;
759d8588912SDave May   /*49*/
760d8588912SDave May   ops->setblocksize    = 0;
761d8588912SDave May   ops->getrowij        = 0;
762d8588912SDave May   ops->restorerowij    = 0;
763d8588912SDave May   ops->getcolumnij     = 0;
764d8588912SDave May   ops->restorecolumnij = 0;
765d8588912SDave May   /*54*/
766d8588912SDave May   ops->fdcoloringcreate = 0;
767d8588912SDave May   ops->coloringpatch    = 0;
768d8588912SDave May   ops->setunfactored    = 0;
769d8588912SDave May   ops->permute          = 0;
770d8588912SDave May   ops->setvaluesblocked = 0;
771d8588912SDave May   /*59*/
772d8588912SDave May   ops->getsubmatrix  = 0;
773d8588912SDave May   ops->destroy       = MatDestroy_Nest;
774d8588912SDave May   ops->view          = MatView_Nest;
775d8588912SDave May   ops->convertfrom   = 0;
776d8588912SDave May   ops->usescaledform = 0;
777d8588912SDave May   /*64*/
778d8588912SDave May   ops->scalesystem             = 0;
779d8588912SDave May   ops->unscalesystem           = 0;
780d8588912SDave May   ops->setlocaltoglobalmapping = 0;
781d8588912SDave May   ops->setvalueslocal          = 0;
782d8588912SDave May   ops->zerorowslocal           = 0;
783d8588912SDave May   /*69*/
784d8588912SDave May   ops->getrowmaxabs    = 0;
785d8588912SDave May   ops->getrowminabs    = 0;
786d8588912SDave May   ops->convert         = 0;/*MatConvert_Nest;*/
787d8588912SDave May   ops->setcoloring     = 0;
788d8588912SDave May   ops->setvaluesadic   = 0;
789d8588912SDave May   /* 74 */
790d8588912SDave May   ops->setvaluesadifor = 0;
791d8588912SDave May   ops->fdcoloringapply              = 0;
792d8588912SDave May   ops->setfromoptions               = 0;
793d8588912SDave May   ops->multconstrained              = 0;
794d8588912SDave May   ops->multtransposeconstrained     = 0;
795d8588912SDave May   /*79*/
796d8588912SDave May   ops->permutesparsify = 0;
797d8588912SDave May   ops->mults           = 0;
798d8588912SDave May   ops->solves          = 0;
799d8588912SDave May   ops->getinertia      = 0;
800d8588912SDave May   ops->load            = 0;
801d8588912SDave May   /*84*/
802d8588912SDave May   ops->issymmetric             = 0;
803d8588912SDave May   ops->ishermitian             = 0;
804d8588912SDave May   ops->isstructurallysymmetric = 0;
805d8588912SDave May   ops->dummy4                  = 0;
806d8588912SDave May   ops->getvecs                 = MatGetVecs_Nest;
807d8588912SDave May   /*89*/
808d8588912SDave May   ops->matmult         = 0;/*MatMatMult_Nest;*/
809d8588912SDave May   ops->matmultsymbolic = 0;
810d8588912SDave May   ops->matmultnumeric  = 0;
811d8588912SDave May   ops->ptap            = 0;
812d8588912SDave May   ops->ptapsymbolic    = 0;
813d8588912SDave May   /*94*/
814d8588912SDave May   ops->ptapnumeric              = 0;
815d8588912SDave May   ops->matmulttranspose         = 0;
816d8588912SDave May   ops->matmulttransposesymbolic = 0;
817d8588912SDave May   ops->matmulttransposenumeric  = 0;
818d8588912SDave May   ops->ptapsymbolic_seqaij      = 0;
819d8588912SDave May   /*99*/
820d8588912SDave May   ops->ptapnumeric_seqaij  = 0;
821d8588912SDave May   ops->ptapsymbolic_mpiaij = 0;
822d8588912SDave May   ops->ptapnumeric_mpiaij  = 0;
823d8588912SDave May   ops->conjugate           = 0;
824d8588912SDave May   ops->setsizes            = 0;
825d8588912SDave May   /*104*/
826d8588912SDave May   ops->setvaluesrow              = 0;
827d8588912SDave May   ops->realpart                  = 0;
828d8588912SDave May   ops->imaginarypart             = 0;
829d8588912SDave May   ops->getrowuppertriangular     = 0;
830d8588912SDave May   ops->restorerowuppertriangular = 0;
831d8588912SDave May   /*109*/
832d8588912SDave May   ops->matsolve           = 0;
833d8588912SDave May   ops->getredundantmatrix = 0;
834d8588912SDave May   ops->getrowmin          = 0;
835d8588912SDave May   ops->getcolumnvector    = 0;
836d8588912SDave May   ops->missingdiagonal    = 0;
837d8588912SDave May   /* 114 */
838d8588912SDave May   ops->getseqnonzerostructure = 0;
839d8588912SDave May   ops->create                 = 0;
840d8588912SDave May   ops->getghosts              = 0;
841d8588912SDave May   ops->getlocalsubmatrix      = MatGetLocalSubMatrix_Nest;
842d8588912SDave May   ops->restorelocalsubmatrix  = MatRestoreLocalSubMatrix_Nest;
843d8588912SDave May   /* 119 */
844d8588912SDave May   ops->multdiagonalblock         = 0;
845d8588912SDave May   ops->hermitiantranspose        = 0;
846d8588912SDave May   ops->multhermitiantranspose    = 0;
847d8588912SDave May   ops->multhermitiantransposeadd = 0;
848d8588912SDave May   ops->getmultiprocblock         = 0;
849d8588912SDave May   /* 124 */
850d8588912SDave May   ops->dummy1                 = 0;
851d8588912SDave May   ops->dummy2                 = 0;
852d8588912SDave May   ops->dummy3                 = 0;
853d8588912SDave May   ops->dummy4                 = 0;
854d8588912SDave May   ops->getsubmatricesparallel = 0;
855d8588912SDave May 
856d8588912SDave May   PetscFunctionReturn(0);
857d8588912SDave May }
858d8588912SDave May 
859d8588912SDave May #undef __FUNCT__
860d8588912SDave May #define __FUNCT__ "MatSetUp_Nest_Private"
861*841e96a3SJed Brown static PetscErrorCode MatSetUp_Nest_Private(Mat A,PetscInt nr,PetscInt nc,const Mat *sub)
862d8588912SDave May {
863d8588912SDave May   Mat_Nest       *ctx = (Mat_Nest*)A->data;
864d8588912SDave May   PetscInt       i,j;
865d8588912SDave May   PetscErrorCode ierr;
866d8588912SDave May 
867d8588912SDave May   PetscFunctionBegin;
868d8588912SDave May   if(ctx->setup_called) PetscFunctionReturn(0);
869d8588912SDave May 
870d8588912SDave May   ctx->nr = nr;
871d8588912SDave May   ctx->nc = nc;
872d8588912SDave May   if (ctx->nr < 0) {
873d8588912SDave May     SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Cannot create MAT_NEST with < 0 row Nests." );
874d8588912SDave May   }
875d8588912SDave May   if (ctx->nc < 0) {
876d8588912SDave May     SETERRQ( ((PetscObject)A)->comm, PETSC_ERR_ARG_WRONG, "Cannot create MAT_NEST with < 0 col Nests." );
877d8588912SDave May   }
878d8588912SDave May 
879d8588912SDave May   /* Create space */
880d8588912SDave May   ierr = PetscMalloc(sizeof(Mat*)*ctx->nr,&ctx->m);CHKERRQ(ierr);
881d8588912SDave May   for (i=0; i<ctx->nr; i++) {
882d8588912SDave May     ierr = PetscMalloc(sizeof(Mat)*ctx->nc,&ctx->m[i]);CHKERRQ(ierr);
883d8588912SDave May   }
884d8588912SDave May   for (i=0; i<ctx->nr; i++) {
885d8588912SDave May     for (j=0; j<ctx->nc; j++) {
886*841e96a3SJed Brown       ctx->m[i][j] = sub[i*nc+j];
887*841e96a3SJed Brown       if (sub[i*nc+j]) {
888*841e96a3SJed Brown         ierr = PetscObjectReference((PetscObject)sub[i*nc+j]);CHKERRQ(ierr);
889d8588912SDave May       }
890d8588912SDave May     }
891d8588912SDave May   }
892d8588912SDave May 
893d8588912SDave May   ierr = PetscMalloc(sizeof(IS)*ctx->nr,&ctx->is_row);CHKERRQ(ierr);
894d8588912SDave May   ierr = PetscMalloc(sizeof(IS)*ctx->nc,&ctx->is_col);CHKERRQ(ierr);
895d8588912SDave May 
896d8588912SDave May   ierr = PetscMalloc(sizeof(PetscInt)*ctx->nr,&ctx->row_len);CHKERRQ(ierr);
897d8588912SDave May   ierr = PetscMalloc(sizeof(PetscInt)*ctx->nc,&ctx->col_len);CHKERRQ(ierr);
898d8588912SDave May   for (i=0;i<ctx->nr;i++) ctx->row_len[i]=-1;
899d8588912SDave May   for (j=0;j<ctx->nc;j++) ctx->col_len[j]=-1;
900d8588912SDave May 
901d8588912SDave May   ctx->setup_called = PETSC_TRUE;
902d8588912SDave May 
903d8588912SDave May   PetscFunctionReturn(0);
904d8588912SDave May }
905d8588912SDave May 
906d8588912SDave May 
907d8588912SDave May /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
908d8588912SDave May /*
909d8588912SDave May   nprocessors = NP
910d8588912SDave May   Nest x^T = ( (g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1) )
911d8588912SDave May        proc 0: => (g_0,h_0,)
912d8588912SDave May        proc 1: => (g_1,h_1,)
913d8588912SDave May        ...
914d8588912SDave May        proc nprocs-1: => (g_NP-1,h_NP-1,)
915d8588912SDave May 
916d8588912SDave May             proc 0:                      proc 1:                    proc nprocs-1:
917d8588912SDave May     is[0] = ( 0,1,2,...,nlocal(g_0)-1 )  ( 0,1,...,nlocal(g_1)-1 )  ( 0,1,...,nlocal(g_NP-1) )
918d8588912SDave May 
919d8588912SDave May             proc 0:
920d8588912SDave May     is[1] = ( nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1 )
921d8588912SDave May             proc 1:
922d8588912SDave May     is[1] = ( nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1 )
923d8588912SDave May 
924d8588912SDave May             proc NP-1:
925d8588912SDave May     is[1] = ( nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1 )
926d8588912SDave May */
927d8588912SDave May #undef __FUNCT__
928d8588912SDave May #define __FUNCT__ "MatSetUp_NestIS_Private"
929*841e96a3SJed Brown static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
930d8588912SDave May {
931d8588912SDave May   Mat_Nest       *ctx = (Mat_Nest*)A->data;
932d8588912SDave May   PetscInt       i,j,offset,n,start,index;
933d8588912SDave May   PetscErrorCode ierr;
934d8588912SDave May 
935d8588912SDave May   PetscFunctionBegin;
936d8588912SDave May   if (is_row) { /* valid IS is passed in */
937d8588912SDave May     /* refs on is[] are incremeneted */
938d8588912SDave May     for (i=0; i<ctx->nr; i++) {
939d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_row[i]);CHKERRQ(ierr);
940d8588912SDave May       ctx->is_row[i] = is_row[i];
941d8588912SDave May     }
942d8588912SDave May   } else {
943d8588912SDave May     for (j=0; j<ctx->nc; j++) {
944d8588912SDave May       if (ctx->m[0][j]) { index = j; break; }
945d8588912SDave May     }
946d8588912SDave May     ierr = MatGetOwnershipRange(ctx->m[0][index],&start,PETSC_NULL);CHKERRQ(ierr);
947d8588912SDave May     offset = start;
948d8588912SDave May     for (i=0; i<ctx->nr; i++) {
949d8588912SDave May 
950d8588912SDave May      ierr = MatGetLocalSize(ctx->m[i][index],&n,PETSC_NULL);CHKERRQ(ierr);
951d8588912SDave May      ierr = ISCreateStride(((PetscObject)ctx->m[i][index])->comm,n,offset,1,&ctx->is_row[i]);CHKERRQ(ierr);
952d8588912SDave May      offset = offset + n;
953d8588912SDave May     }
954d8588912SDave May   }
955d8588912SDave May 
956d8588912SDave May   if (is_col) { /* valid IS is passed in */
957d8588912SDave May     /* refs on is[] are incremeneted */
958d8588912SDave May     for (j=0; j<ctx->nc; j++) {
959d8588912SDave May       ierr = PetscObjectReference((PetscObject)is_col[j]);CHKERRQ(ierr);
960d8588912SDave May       ctx->is_col[j] = is_col[j];
961d8588912SDave May     }
962d8588912SDave May   } else {
963d8588912SDave May     for (i=0; i<ctx->nr; i++) {
964d8588912SDave May       if (ctx->m[i][0]) { index = i; break; }
965d8588912SDave May     }
966d8588912SDave May     ierr = MatGetOwnershipRange(ctx->m[index][0],&start,PETSC_NULL);CHKERRQ(ierr);
967d8588912SDave May     offset = start;
968d8588912SDave May     for (j=0; j<ctx->nc; j++) {
969d8588912SDave May       for (i=0; i<ctx->nr; i++) {
970d8588912SDave May         if (ctx->m[i][j]) { index = i; break; }
971d8588912SDave May       }
972d8588912SDave May      ierr = MatGetLocalSize(ctx->m[index][j],PETSC_NULL,&n);CHKERRQ(ierr);
973d8588912SDave May      ierr = ISCreateStride(((PetscObject)ctx->m[index][j])->comm,n,offset,1,&ctx->is_col[j]);CHKERRQ(ierr);
974d8588912SDave May      offset = offset + n;
975d8588912SDave May     }
976d8588912SDave May   }
977d8588912SDave May   PetscFunctionReturn(0);
978d8588912SDave May }
979d8588912SDave May 
980d8588912SDave May #undef __FUNCT__
981d8588912SDave May #define __FUNCT__ "MatCreateNest"
982*841e96a3SJed Brown PetscErrorCode PETSCMAT_DLLEXPORT MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
983d8588912SDave May {
984d8588912SDave May   Mat            A;
985d8588912SDave May   Mat_Nest       *s;
986d8588912SDave May   PetscInt       m,n,M,N;
987d8588912SDave May   PetscErrorCode ierr;
988d8588912SDave May 
989d8588912SDave May   PetscFunctionBegin;
990*841e96a3SJed Brown   if (nr < 0) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
991*841e96a3SJed Brown   if (nr) PetscValidPointer(is_row,3);CHKERRQ(ierr);
992*841e96a3SJed Brown   if (nc < 0) SETERRQ(comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
993*841e96a3SJed Brown   if (nc) PetscValidPointer(is_col,5);CHKERRQ(ierr);
994*841e96a3SJed Brown   if (nr*nc) PetscValidPointer(a,6);
995*841e96a3SJed Brown   PetscValidPointer(B,7);
996*841e96a3SJed Brown 
997d8588912SDave May   ierr = MatCreate(comm,&A);CHKERRQ(ierr);
998d8588912SDave May 
999d8588912SDave May   ierr = PetscMalloc( sizeof(Mat_Nest), &s );CHKERRQ(ierr);
1000d8588912SDave May   ierr = PetscMemzero( s, sizeof(Mat_Nest) );CHKERRQ(ierr);
1001d8588912SDave May   A->data         = (void*)s;
1002d8588912SDave May   s->setup_called = PETSC_FALSE;
1003d8588912SDave May   s->nr = s->nc   = -1;
1004d8588912SDave May   s->m            = PETSC_NULL;
1005d8588912SDave May 
1006d8588912SDave May   /* define the operations */
1007d8588912SDave May   ierr = MatNestSetOps_Private(A->ops);CHKERRQ(ierr);
1008d8588912SDave May   A->spptr            = 0;
1009d8588912SDave May   A->same_nonzero     = PETSC_FALSE;
1010d8588912SDave May   A->assembled        = PETSC_FALSE;
1011d8588912SDave May 
1012d8588912SDave May   ierr = PetscObjectChangeTypeName((PetscObject) A, MATNEST );CHKERRQ(ierr);
1013d8588912SDave May   ierr = MatSetUp_Nest_Private(A,nr,nc,a);CHKERRQ(ierr);
1014d8588912SDave May 
1015d8588912SDave May   m = n = 0;
1016d8588912SDave May   M = N = 0;
1017d8588912SDave May   ierr = MatNestGetSize_Recursive(A,PETSC_TRUE,&M,&N);CHKERRQ(ierr);
1018d8588912SDave May   ierr = MatNestGetSize_Recursive(A,PETSC_FALSE,&m,&n);CHKERRQ(ierr);
1019d8588912SDave May 
1020d8588912SDave May   ierr = PetscLayoutSetSize(A->rmap,M);CHKERRQ(ierr);
1021d8588912SDave May   ierr = PetscLayoutSetLocalSize(A->rmap,m);CHKERRQ(ierr);
1022d8588912SDave May   ierr = PetscLayoutSetBlockSize(A->rmap,1);CHKERRQ(ierr);
1023d8588912SDave May   ierr = PetscLayoutSetSize(A->cmap,N);CHKERRQ(ierr);
1024d8588912SDave May   ierr = PetscLayoutSetLocalSize(A->cmap,n);CHKERRQ(ierr);
1025d8588912SDave May   ierr = PetscLayoutSetBlockSize(A->cmap,1);CHKERRQ(ierr);
1026d8588912SDave May 
1027d8588912SDave May   ierr = PetscLayoutSetUp(A->rmap);CHKERRQ(ierr);
1028d8588912SDave May   ierr = PetscLayoutSetUp(A->cmap);CHKERRQ(ierr);
1029d8588912SDave May 
1030*841e96a3SJed Brown   ierr = MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);CHKERRQ(ierr);
1031d8588912SDave May 
1032d8588912SDave May   /* expose Nest api's */
1033d8588912SDave May   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMat_C","MatNestGetSubMat_Nest",MatNestGetSubMat_Nest);CHKERRQ(ierr);
1034d8588912SDave May   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSubMats_C","MatNestGetSubMats_Nest",MatNestGetSubMats_Nest);CHKERRQ(ierr);
1035d8588912SDave May   ierr = PetscObjectComposeFunctionDynamic((PetscObject)A,"MatNestGetSize_C","MatNestGetSize_Nest",MatNestGetSize_Nest);CHKERRQ(ierr);
1036d8588912SDave May 
1037d8588912SDave May   *B = A;
1038d8588912SDave May   PetscFunctionReturn(0);
1039d8588912SDave May }
1040d8588912SDave May 
1041