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