Actual source code: matnest.c

petsc-3.11.4 2019-09-28
Report Typos and Errors

  2:  #include <../src/mat/impls/nest/matnestimpl.h>
  3:  #include <../src/mat/impls/aij/seq/aij.h>
  4:  #include <petscsf.h>

  6: static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]);
  7: static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left);
  8: PETSC_INTERN PetscErrorCode MatConvert_Nest_IS(Mat,MatType,MatReuse,Mat*);

 10: /* private functions */
 11: static PetscErrorCode MatNestGetSizes_Private(Mat A,PetscInt *m,PetscInt *n,PetscInt *M,PetscInt *N)
 12: {
 13:   Mat_Nest       *bA = (Mat_Nest*)A->data;
 14:   PetscInt       i,j;

 18:   *m = *n = *M = *N = 0;
 19:   for (i=0; i<bA->nr; i++) {  /* rows */
 20:     PetscInt sm,sM;
 21:     ISGetLocalSize(bA->isglobal.row[i],&sm);
 22:     ISGetSize(bA->isglobal.row[i],&sM);
 23:     *m  += sm;
 24:     *M  += sM;
 25:   }
 26:   for (j=0; j<bA->nc; j++) {  /* cols */
 27:     PetscInt sn,sN;
 28:     ISGetLocalSize(bA->isglobal.col[j],&sn);
 29:     ISGetSize(bA->isglobal.col[j],&sN);
 30:     *n  += sn;
 31:     *N  += sN;
 32:   }
 33:   return(0);
 34: }

 36: /* operations */
 37: static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y)
 38: {
 39:   Mat_Nest       *bA = (Mat_Nest*)A->data;
 40:   Vec            *bx = bA->right,*by = bA->left;
 41:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;

 45:   for (i=0; i<nr; i++) {VecGetSubVector(y,bA->isglobal.row[i],&by[i]);}
 46:   for (i=0; i<nc; i++) {VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);}
 47:   for (i=0; i<nr; i++) {
 48:     VecZeroEntries(by[i]);
 49:     for (j=0; j<nc; j++) {
 50:       if (!bA->m[i][j]) continue;
 51:       /* y[i] <- y[i] + A[i][j] * x[j] */
 52:       MatMultAdd(bA->m[i][j],bx[j],by[i],by[i]);
 53:     }
 54:   }
 55:   for (i=0; i<nr; i++) {VecRestoreSubVector(y,bA->isglobal.row[i],&by[i]);}
 56:   for (i=0; i<nc; i++) {VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);}
 57:   return(0);
 58: }

 60: static PetscErrorCode MatMultAdd_Nest(Mat A,Vec x,Vec y,Vec z)
 61: {
 62:   Mat_Nest       *bA = (Mat_Nest*)A->data;
 63:   Vec            *bx = bA->right,*bz = bA->left;
 64:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;

 68:   for (i=0; i<nr; i++) {VecGetSubVector(z,bA->isglobal.row[i],&bz[i]);}
 69:   for (i=0; i<nc; i++) {VecGetSubVector(x,bA->isglobal.col[i],&bx[i]);}
 70:   for (i=0; i<nr; i++) {
 71:     if (y != z) {
 72:       Vec by;
 73:       VecGetSubVector(y,bA->isglobal.row[i],&by);
 74:       VecCopy(by,bz[i]);
 75:       VecRestoreSubVector(y,bA->isglobal.row[i],&by);
 76:     }
 77:     for (j=0; j<nc; j++) {
 78:       if (!bA->m[i][j]) continue;
 79:       /* y[i] <- y[i] + A[i][j] * x[j] */
 80:       MatMultAdd(bA->m[i][j],bx[j],bz[i],bz[i]);
 81:     }
 82:   }
 83:   for (i=0; i<nr; i++) {VecRestoreSubVector(z,bA->isglobal.row[i],&bz[i]);}
 84:   for (i=0; i<nc; i++) {VecRestoreSubVector(x,bA->isglobal.col[i],&bx[i]);}
 85:   return(0);
 86: }

 88: static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y)
 89: {
 90:   Mat_Nest       *bA = (Mat_Nest*)A->data;
 91:   Vec            *bx = bA->left,*by = bA->right;
 92:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;

 96:   for (i=0; i<nr; i++) {VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);}
 97:   for (i=0; i<nc; i++) {VecGetSubVector(y,bA->isglobal.col[i],&by[i]);}
 98:   for (j=0; j<nc; j++) {
 99:     VecZeroEntries(by[j]);
100:     for (i=0; i<nr; i++) {
101:       if (!bA->m[i][j]) continue;
102:       /* y[j] <- y[j] + (A[i][j])^T * x[i] */
103:       MatMultTransposeAdd(bA->m[i][j],bx[i],by[j],by[j]);
104:     }
105:   }
106:   for (i=0; i<nr; i++) {VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);}
107:   for (i=0; i<nc; i++) {VecRestoreSubVector(y,bA->isglobal.col[i],&by[i]);}
108:   return(0);
109: }

111: static PetscErrorCode MatMultTransposeAdd_Nest(Mat A,Vec x,Vec y,Vec z)
112: {
113:   Mat_Nest       *bA = (Mat_Nest*)A->data;
114:   Vec            *bx = bA->left,*bz = bA->right;
115:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;

119:   for (i=0; i<nr; i++) {VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);}
120:   for (i=0; i<nc; i++) {VecGetSubVector(z,bA->isglobal.col[i],&bz[i]);}
121:   for (j=0; j<nc; j++) {
122:     if (y != z) {
123:       Vec by;
124:       VecGetSubVector(y,bA->isglobal.col[j],&by);
125:       VecCopy(by,bz[j]);
126:       VecRestoreSubVector(y,bA->isglobal.col[j],&by);
127:     }
128:     for (i=0; i<nr; i++) {
129:       if (!bA->m[i][j]) continue;
130:       /* z[j] <- y[j] + (A[i][j])^T * x[i] */
131:       MatMultTransposeAdd(bA->m[i][j],bx[i],bz[j],bz[j]);
132:     }
133:   }
134:   for (i=0; i<nr; i++) {VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);}
135:   for (i=0; i<nc; i++) {VecRestoreSubVector(z,bA->isglobal.col[i],&bz[i]);}
136:   return(0);
137: }

139: static PetscErrorCode MatTranspose_Nest(Mat A,MatReuse reuse,Mat *B)
140: {
141:   Mat_Nest       *bA = (Mat_Nest*)A->data, *bC;
142:   Mat            C;
143:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;

147:   if (reuse == MAT_INPLACE_MATRIX && nr != nc) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_SIZ,"Square nested matrix only for in-place");

149:   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) {
150:     Mat *subs;
151:     IS  *is_row,*is_col;

153:     PetscCalloc1(nr * nc,&subs);
154:     PetscMalloc2(nr,&is_row,nc,&is_col);
155:     MatNestGetISs(A,is_row,is_col);
156:     if (reuse == MAT_INPLACE_MATRIX) {
157:       for (i=0; i<nr; i++) {
158:         for (j=0; j<nc; j++) {
159:           subs[i + nr * j] = bA->m[i][j];
160:         }
161:       }
162:     }

164:     MatCreateNest(PetscObjectComm((PetscObject)A),nc,is_col,nr,is_row,subs,&C);
165:     PetscFree(subs);
166:     PetscFree2(is_row,is_col);
167:   } else {
168:     C = *B;
169:   }

171:   bC = (Mat_Nest*)C->data;
172:   for (i=0; i<nr; i++) {
173:     for (j=0; j<nc; j++) {
174:       if (bA->m[i][j]) {
175:         MatTranspose(bA->m[i][j], reuse, &(bC->m[j][i]));
176:       } else {
177:         bC->m[j][i] = NULL;
178:       }
179:     }
180:   }

182:   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) {
183:     *B = C;
184:   } else {
185:     MatHeaderMerge(A, &C);
186:   }
187:   return(0);
188: }

190: static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list)
191: {
193:   IS             *lst = *list;
194:   PetscInt       i;

197:   if (!lst) return(0);
198:   for (i=0; i<n; i++) if (lst[i]) {ISDestroy(&lst[i]);}
199:   PetscFree(lst);
200:   *list = NULL;
201:   return(0);
202: }

204: static PetscErrorCode MatDestroy_Nest(Mat A)
205: {
206:   Mat_Nest       *vs = (Mat_Nest*)A->data;
207:   PetscInt       i,j;

211:   /* release the matrices and the place holders */
212:   MatNestDestroyISList(vs->nr,&vs->isglobal.row);
213:   MatNestDestroyISList(vs->nc,&vs->isglobal.col);
214:   MatNestDestroyISList(vs->nr,&vs->islocal.row);
215:   MatNestDestroyISList(vs->nc,&vs->islocal.col);

217:   PetscFree(vs->row_len);
218:   PetscFree(vs->col_len);

220:   PetscFree2(vs->left,vs->right);

222:   /* release the matrices and the place holders */
223:   if (vs->m) {
224:     for (i=0; i<vs->nr; i++) {
225:       for (j=0; j<vs->nc; j++) {
226:         MatDestroy(&vs->m[i][j]);
227:       }
228:       PetscFree(vs->m[i]);
229:     }
230:     PetscFree(vs->m);
231:   }
232:   PetscFree(A->data);

234:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",0);
235:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",0);
236:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",0);
237:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",0);
238:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",0);
239:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",0);
240:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",0);
241:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",0);
242:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_mpiaij_C",0);
243:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_seqaij_C",0);
244:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",0);
245:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_is_C",0);
246:   return(0);
247: }

249: static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
250: {
251:   Mat_Nest       *vs = (Mat_Nest*)A->data;
252:   PetscInt       i,j;

256:   for (i=0; i<vs->nr; i++) {
257:     for (j=0; j<vs->nc; j++) {
258:       if (vs->m[i][j]) {
259:         MatAssemblyBegin(vs->m[i][j],type);
260:         if (!vs->splitassembly) {
261:           /* Note: split assembly will fail if the same block appears more than once (even indirectly through a nested
262:            * sub-block). This could be fixed by adding a flag to Mat so that there was a way to check if a Mat was
263:            * already performing an assembly, but the result would by more complicated and appears to offer less
264:            * potential for diagnostics and correctness checking. Split assembly should be fixed once there is an
265:            * interface for libraries to make asynchronous progress in "user-defined non-blocking collectives".
266:            */
267:           MatAssemblyEnd(vs->m[i][j],type);
268:         }
269:       }
270:     }
271:   }
272:   return(0);
273: }

275: static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
276: {
277:   Mat_Nest       *vs = (Mat_Nest*)A->data;
278:   PetscInt       i,j;

282:   for (i=0; i<vs->nr; i++) {
283:     for (j=0; j<vs->nc; j++) {
284:       if (vs->m[i][j]) {
285:         if (vs->splitassembly) {
286:           MatAssemblyEnd(vs->m[i][j],type);
287:         }
288:       }
289:     }
290:   }
291:   return(0);
292: }

294: static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
295: {
297:   Mat_Nest       *vs = (Mat_Nest*)A->data;
298:   PetscInt       j;
299:   Mat            sub;

302:   sub = (row < vs->nc) ? vs->m[row][row] : (Mat)NULL; /* Prefer to find on the diagonal */
303:   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
304:   if (sub) {MatSetUp(sub);}       /* Ensure that the sizes are available */
305:   *B = sub;
306:   return(0);
307: }

309: static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
310: {
312:   Mat_Nest       *vs = (Mat_Nest*)A->data;
313:   PetscInt       i;
314:   Mat            sub;

317:   sub = (col < vs->nr) ? vs->m[col][col] : (Mat)NULL; /* Prefer to find on the diagonal */
318:   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
319:   if (sub) {MatSetUp(sub);}       /* Ensure that the sizes are available */
320:   *B = sub;
321:   return(0);
322: }

324: static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
325: {
327:   PetscInt       i;
328:   PetscBool      flg;

334:   *found = -1;
335:   for (i=0; i<n; i++) {
336:     if (!list[i]) continue;
337:     ISEqual(list[i],is,&flg);
338:     if (flg) {
339:       *found = i;
340:       return(0);
341:     }
342:   }
343:   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Could not find index set");
344:   return(0);
345: }

347: /* Get a block row as a new MatNest */
348: static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B)
349: {
350:   Mat_Nest       *vs = (Mat_Nest*)A->data;
351:   char           keyname[256];

355:   *B   = NULL;
356:   PetscSNPrintf(keyname,sizeof(keyname),"NestRow_%D",row);
357:   PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);
358:   if (*B) return(0);

360:   MatCreateNest(PetscObjectComm((PetscObject)A),1,NULL,vs->nc,vs->isglobal.col,vs->m[row],B);

362:   (*B)->assembled = A->assembled;

364:   PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);
365:   PetscObjectDereference((PetscObject)*B); /* Leave the only remaining reference in the composition */
366:   return(0);
367: }

369: static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
370: {
371:   Mat_Nest       *vs = (Mat_Nest*)A->data;
373:   PetscInt       row,col;
374:   PetscBool      same,isFullCol,isFullColGlobal;

377:   /* Check if full column space. This is a hack */
378:   isFullCol = PETSC_FALSE;
379:   PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&same);
380:   if (same) {
381:     PetscInt n,first,step,i,an,am,afirst,astep;
382:     ISStrideGetInfo(iscol,&first,&step);
383:     ISGetLocalSize(iscol,&n);
384:     isFullCol = PETSC_TRUE;
385:     for (i=0,an=A->cmap->rstart; i<vs->nc; i++) {
386:       ISStrideGetInfo(is->col[i],&afirst,&astep);
387:       ISGetLocalSize(is->col[i],&am);
388:       if (afirst != an || astep != step) isFullCol = PETSC_FALSE;
389:       an += am;
390:     }
391:     if (an != A->cmap->rstart+n) isFullCol = PETSC_FALSE;
392:   }
393:   MPIU_Allreduce(&isFullCol,&isFullColGlobal,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)iscol));

395:   if (isFullColGlobal && vs->nc > 1) {
396:     PetscInt row;
397:     MatNestFindIS(A,vs->nr,is->row,isrow,&row);
398:     MatNestGetRow(A,row,B);
399:   } else {
400:     MatNestFindIS(A,vs->nr,is->row,isrow,&row);
401:     MatNestFindIS(A,vs->nc,is->col,iscol,&col);
402:     if (!vs->m[row][col]) {
403:       PetscInt lr,lc;

405:       MatCreate(PetscObjectComm((PetscObject)A),&vs->m[row][col]);
406:       ISGetLocalSize(vs->isglobal.row[row],&lr);
407:       ISGetLocalSize(vs->isglobal.col[col],&lc);
408:       MatSetSizes(vs->m[row][col],lr,lc,PETSC_DECIDE,PETSC_DECIDE);
409:       MatSetUp(vs->m[row][col]);
410:       MatAssemblyBegin(vs->m[row][col],MAT_FINAL_ASSEMBLY);
411:       MatAssemblyEnd(vs->m[row][col],MAT_FINAL_ASSEMBLY);
412:     }
413:     *B = vs->m[row][col];
414:   }
415:   return(0);
416: }

418: static PetscErrorCode MatCreateSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
419: {
421:   Mat_Nest       *vs = (Mat_Nest*)A->data;
422:   Mat            sub;

425:   MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);
426:   switch (reuse) {
427:   case MAT_INITIAL_MATRIX:
428:     if (sub) { PetscObjectReference((PetscObject)sub); }
429:     *B = sub;
430:     break;
431:   case MAT_REUSE_MATRIX:
432:     if (sub != *B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
433:     break;
434:   case MAT_IGNORE_MATRIX:       /* Nothing to do */
435:     break;
436:   case MAT_INPLACE_MATRIX:       /* Nothing to do */
437:     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX is not supported yet");
438:     break;
439:   }
440:   return(0);
441: }

443: PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
444: {
446:   Mat_Nest       *vs = (Mat_Nest*)A->data;
447:   Mat            sub;

450:   MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);
451:   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
452:   if (sub) {PetscObjectReference((PetscObject)sub);}
453:   *B = sub;
454:   return(0);
455: }

457: static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
458: {
460:   Mat_Nest       *vs = (Mat_Nest*)A->data;
461:   Mat            sub;

464:   MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);
465:   if (*B != sub) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
466:   if (sub) {
467:     if (((PetscObject)sub)->refct <= 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
468:     MatDestroy(B);
469:   }
470:   return(0);
471: }

473: static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
474: {
475:   Mat_Nest       *bA = (Mat_Nest*)A->data;
476:   PetscInt       i;

480:   for (i=0; i<bA->nr; i++) {
481:     Vec bv;
482:     VecGetSubVector(v,bA->isglobal.row[i],&bv);
483:     if (bA->m[i][i]) {
484:       MatGetDiagonal(bA->m[i][i],bv);
485:     } else {
486:       VecSet(bv,0.0);
487:     }
488:     VecRestoreSubVector(v,bA->isglobal.row[i],&bv);
489:   }
490:   return(0);
491: }

493: static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
494: {
495:   Mat_Nest       *bA = (Mat_Nest*)A->data;
496:   Vec            bl,*br;
497:   PetscInt       i,j;

501:   PetscCalloc1(bA->nc,&br);
502:   if (r) {
503:     for (j=0; j<bA->nc; j++) {VecGetSubVector(r,bA->isglobal.col[j],&br[j]);}
504:   }
505:   bl = NULL;
506:   for (i=0; i<bA->nr; i++) {
507:     if (l) {
508:       VecGetSubVector(l,bA->isglobal.row[i],&bl);
509:     }
510:     for (j=0; j<bA->nc; j++) {
511:       if (bA->m[i][j]) {
512:         MatDiagonalScale(bA->m[i][j],bl,br[j]);
513:       }
514:     }
515:     if (l) {
516:       VecRestoreSubVector(l,bA->isglobal.row[i],&bl);
517:     }
518:   }
519:   if (r) {
520:     for (j=0; j<bA->nc; j++) {VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);}
521:   }
522:   PetscFree(br);
523:   return(0);
524: }

526: static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a)
527: {
528:   Mat_Nest       *bA = (Mat_Nest*)A->data;
529:   PetscInt       i,j;

533:   for (i=0; i<bA->nr; i++) {
534:     for (j=0; j<bA->nc; j++) {
535:       if (bA->m[i][j]) {
536:         MatScale(bA->m[i][j],a);
537:       }
538:     }
539:   }
540:   return(0);
541: }

543: static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a)
544: {
545:   Mat_Nest       *bA = (Mat_Nest*)A->data;
546:   PetscInt       i;

550:   for (i=0; i<bA->nr; i++) {
551:     if (!bA->m[i][i]) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support for shifting an empty diagonal block, insert a matrix in block (%D,%D)",i,i);
552:     MatShift(bA->m[i][i],a);
553:   }
554:   return(0);
555: }

557: static PetscErrorCode MatDiagonalSet_Nest(Mat A,Vec D,InsertMode is)
558: {
559:   Mat_Nest       *bA = (Mat_Nest*)A->data;
560:   PetscInt       i;

564:   for (i=0; i<bA->nr; i++) {
565:     Vec bv;
566:     VecGetSubVector(D,bA->isglobal.row[i],&bv);
567:     if (bA->m[i][i]) {
568:       MatDiagonalSet(bA->m[i][i],bv,is);
569:     }
570:     VecRestoreSubVector(D,bA->isglobal.row[i],&bv);
571:   }
572:   return(0);
573: }

575: static PetscErrorCode MatSetRandom_Nest(Mat A,PetscRandom rctx)
576: {
577:   Mat_Nest       *bA = (Mat_Nest*)A->data;
578:   PetscInt       i,j;

582:   for (i=0; i<bA->nr; i++) {
583:     for (j=0; j<bA->nc; j++) {
584:       if (bA->m[i][j]) {
585:         MatSetRandom(bA->m[i][j],rctx);
586:       }
587:     }
588:   }
589:   return(0);
590: }

592: static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left)
593: {
594:   Mat_Nest       *bA = (Mat_Nest*)A->data;
595:   Vec            *L,*R;
596:   MPI_Comm       comm;
597:   PetscInt       i,j;

601:   PetscObjectGetComm((PetscObject)A,&comm);
602:   if (right) {
603:     /* allocate R */
604:     PetscMalloc1(bA->nc, &R);
605:     /* Create the right vectors */
606:     for (j=0; j<bA->nc; j++) {
607:       for (i=0; i<bA->nr; i++) {
608:         if (bA->m[i][j]) {
609:           MatCreateVecs(bA->m[i][j],&R[j],NULL);
610:           break;
611:         }
612:       }
613:       if (i==bA->nr) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
614:     }
615:     VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);
616:     /* hand back control to the nest vector */
617:     for (j=0; j<bA->nc; j++) {
618:       VecDestroy(&R[j]);
619:     }
620:     PetscFree(R);
621:   }

623:   if (left) {
624:     /* allocate L */
625:     PetscMalloc1(bA->nr, &L);
626:     /* Create the left vectors */
627:     for (i=0; i<bA->nr; i++) {
628:       for (j=0; j<bA->nc; j++) {
629:         if (bA->m[i][j]) {
630:           MatCreateVecs(bA->m[i][j],NULL,&L[i]);
631:           break;
632:         }
633:       }
634:       if (j==bA->nc) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
635:     }

637:     VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);
638:     for (i=0; i<bA->nr; i++) {
639:       VecDestroy(&L[i]);
640:     }

642:     PetscFree(L);
643:   }
644:   return(0);
645: }

647: static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
648: {
649:   Mat_Nest       *bA = (Mat_Nest*)A->data;
650:   PetscBool      isascii;
651:   PetscInt       i,j;

655:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
656:   if (isascii) {

658:     PetscViewerASCIIPrintf(viewer,"Matrix object: \n");
659:     PetscViewerASCIIPushTab(viewer);
660:     PetscViewerASCIIPrintf(viewer, "type=nest, rows=%D, cols=%D \n",bA->nr,bA->nc);

662:     PetscViewerASCIIPrintf(viewer,"MatNest structure: \n");
663:     for (i=0; i<bA->nr; i++) {
664:       for (j=0; j<bA->nc; j++) {
665:         MatType   type;
666:         char      name[256] = "",prefix[256] = "";
667:         PetscInt  NR,NC;
668:         PetscBool isNest = PETSC_FALSE;

670:         if (!bA->m[i][j]) {
671:           PetscViewerASCIIPrintf(viewer, "(%D,%D) : NULL \n",i,j);
672:           continue;
673:         }
674:         MatGetSize(bA->m[i][j],&NR,&NC);
675:         MatGetType(bA->m[i][j], &type);
676:         if (((PetscObject)bA->m[i][j])->name) {PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);}
677:         if (((PetscObject)bA->m[i][j])->prefix) {PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);}
678:         PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);

680:         PetscViewerASCIIPrintf(viewer,"(%D,%D) : %s%stype=%s, rows=%D, cols=%D \n",i,j,name,prefix,type,NR,NC);

682:         if (isNest) {
683:           PetscViewerASCIIPushTab(viewer);  /* push1 */
684:           MatView(bA->m[i][j],viewer);
685:           PetscViewerASCIIPopTab(viewer);    /* pop1 */
686:         }
687:       }
688:     }
689:     PetscViewerASCIIPopTab(viewer);    /* pop0 */
690:   }
691:   return(0);
692: }

694: static PetscErrorCode MatZeroEntries_Nest(Mat A)
695: {
696:   Mat_Nest       *bA = (Mat_Nest*)A->data;
697:   PetscInt       i,j;

701:   for (i=0; i<bA->nr; i++) {
702:     for (j=0; j<bA->nc; j++) {
703:       if (!bA->m[i][j]) continue;
704:       MatZeroEntries(bA->m[i][j]);
705:     }
706:   }
707:   return(0);
708: }

710: static PetscErrorCode MatCopy_Nest(Mat A,Mat B,MatStructure str)
711: {
712:   Mat_Nest       *bA = (Mat_Nest*)A->data,*bB = (Mat_Nest*)B->data;
713:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;

717:   if (nr != bB->nr || nc != bB->nc) SETERRQ4(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Cannot copy a Mat_Nest of block size (%D,%D) to a Mat_Nest of block size (%D,%D)",bB->nr,bB->nc,nr,nc);
718:   for (i=0; i<nr; i++) {
719:     for (j=0; j<nc; j++) {
720:       if (bA->m[i][j] && bB->m[i][j]) {
721:         MatCopy(bA->m[i][j],bB->m[i][j],str);
722:       } else if (bA->m[i][j] || bB->m[i][j]) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Matrix block does not exist at %D,%D",i,j);
723:     }
724:   }
725:   PetscObjectStateIncrease((PetscObject)B);
726:   return(0);
727: }

729: static PetscErrorCode MatAXPY_Nest(Mat Y,PetscScalar a,Mat X,MatStructure str)
730: {
731:   Mat_Nest       *bY = (Mat_Nest*)Y->data,*bX = (Mat_Nest*)X->data;
732:   PetscInt       i,j,nr = bY->nr,nc = bY->nc;

736:   if (nr != bX->nr || nc != bX->nc) SETERRQ4(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_INCOMP,"Cannot AXPY a MatNest of block size (%D,%D) with a MatNest of block size (%D,%D)",bX->nr,bX->nc,nr,nc);
737:   for (i=0; i<nr; i++) {
738:     for (j=0; j<nc; j++) {
739:       if (bY->m[i][j] && bX->m[i][j]) {
740:         MatAXPY(bY->m[i][j],a,bX->m[i][j],str);
741:       } else if (bX->m[i][j]) SETERRQ2(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_INCOMP,"Matrix block does not exist at %D,%D",i,j);
742:     }
743:   }
744:   return(0);
745: }

747: static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
748: {
749:   Mat_Nest       *bA = (Mat_Nest*)A->data;
750:   Mat            *b;
751:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;

755:   PetscMalloc1(nr*nc,&b);
756:   for (i=0; i<nr; i++) {
757:     for (j=0; j<nc; j++) {
758:       if (bA->m[i][j]) {
759:         MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);
760:       } else {
761:         b[i*nc+j] = NULL;
762:       }
763:     }
764:   }
765:   MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);
766:   /* Give the new MatNest exclusive ownership */
767:   for (i=0; i<nr*nc; i++) {
768:     MatDestroy(&b[i]);
769:   }
770:   PetscFree(b);

772:   MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
773:   MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
774:   return(0);
775: }

777: /* nest api */
778: PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
779: {
780:   Mat_Nest *bA = (Mat_Nest*)A->data;

783:   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
784:   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
785:   *mat = bA->m[idxm][jdxm];
786:   return(0);
787: }

789: /*@
790:  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.

792:  Not collective

794:  Input Parameters:
795: +   A  - nest matrix
796: .   idxm - index of the matrix within the nest matrix
797: -   jdxm - index of the matrix within the nest matrix

799:  Output Parameter:
800: .   sub - matrix at index idxm,jdxm within the nest matrix

802:  Level: developer

804: .seealso: MatNestGetSize(), MatNestGetSubMats()
805: @*/
806: PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
807: {

811:   PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));
812:   return(0);
813: }

815: PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat)
816: {
817:   Mat_Nest       *bA = (Mat_Nest*)A->data;
818:   PetscInt       m,n,M,N,mi,ni,Mi,Ni;

822:   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
823:   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
824:   MatGetLocalSize(mat,&m,&n);
825:   MatGetSize(mat,&M,&N);
826:   ISGetLocalSize(bA->isglobal.row[idxm],&mi);
827:   ISGetSize(bA->isglobal.row[idxm],&Mi);
828:   ISGetLocalSize(bA->isglobal.col[jdxm],&ni);
829:   ISGetSize(bA->isglobal.col[jdxm],&Ni);
830:   if (M != Mi || N != Ni) SETERRQ4(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"Submatrix dimension (%D,%D) incompatible with nest block (%D,%D)",M,N,Mi,Ni);
831:   if (m != mi || n != ni) SETERRQ4(PetscObjectComm((PetscObject)mat),PETSC_ERR_ARG_INCOMP,"Submatrix local dimension (%D,%D) incompatible with nest block (%D,%D)",m,n,mi,ni);

833:   PetscObjectReference((PetscObject)mat);
834:   MatDestroy(&bA->m[idxm][jdxm]);
835:   bA->m[idxm][jdxm] = mat;
836:   return(0);
837: }

839: /*@
840:  MatNestSetSubMat - Set a single submatrix in the nest matrix.

842:  Logically collective on the submatrix communicator

844:  Input Parameters:
845: +   A  - nest matrix
846: .   idxm - index of the matrix within the nest matrix
847: .   jdxm - index of the matrix within the nest matrix
848: -   sub - matrix at index idxm,jdxm within the nest matrix

850:  Notes:
851:  The new submatrix must have the same size and communicator as that block of the nest.

853:  This increments the reference count of the submatrix.

855:  Level: developer

857: .seealso: MatNestSetSubMats(), MatNestGetSubMats()
858: @*/
859: PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
860: {

864:   PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));
865:   return(0);
866: }

868: PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
869: {
870:   Mat_Nest *bA = (Mat_Nest*)A->data;

873:   if (M)   *M   = bA->nr;
874:   if (N)   *N   = bA->nc;
875:   if (mat) *mat = bA->m;
876:   return(0);
877: }

879: /*@C
880:  MatNestGetSubMats - Returns the entire two dimensional array of matrices defining a nest matrix.

882:  Not collective

884:  Input Parameters:
885: .   A  - nest matrix

887:  Output Parameter:
888: +   M - number of rows in the nest matrix
889: .   N - number of cols in the nest matrix
890: -   mat - 2d array of matrices

892:  Notes:

894:  The user should not free the array mat.

896:  In Fortran, this routine has a calling sequence
897: $   call MatNestGetSubMats(A, M, N, mat, ierr)
898:  where the space allocated for the optional argument mat is assumed large enough (if provided).

900:  Level: developer

902: .seealso: MatNestGetSize(), MatNestGetSubMat()
903: @*/
904: PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
905: {

909:   PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));
910:   return(0);
911: }

913: PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
914: {
915:   Mat_Nest *bA = (Mat_Nest*)A->data;

918:   if (M) *M = bA->nr;
919:   if (N) *N = bA->nc;
920:   return(0);
921: }

923: /*@
924:  MatNestGetSize - Returns the size of the nest matrix.

926:  Not collective

928:  Input Parameters:
929: .   A  - nest matrix

931:  Output Parameter:
932: +   M - number of rows in the nested mat
933: -   N - number of cols in the nested mat

935:  Notes:

937:  Level: developer

939: .seealso: MatNestGetSubMat(), MatNestGetSubMats()
940: @*/
941: PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
942: {

946:   PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));
947:   return(0);
948: }

950: static PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[])
951: {
952:   Mat_Nest *vs = (Mat_Nest*)A->data;
953:   PetscInt i;

956:   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i];
957:   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i];
958:   return(0);
959: }

961: /*@C
962:  MatNestGetISs - Returns the index sets partitioning the row and column spaces

964:  Not collective

966:  Input Parameters:
967: .   A  - nest matrix

969:  Output Parameter:
970: +   rows - array of row index sets
971: -   cols - array of column index sets

973:  Level: advanced

975:  Notes:
976:  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.

978: .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs()
979: @*/
980: PetscErrorCode  MatNestGetISs(Mat A,IS rows[],IS cols[])
981: {

986:   PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));
987:   return(0);
988: }

990: static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[])
991: {
992:   Mat_Nest *vs = (Mat_Nest*)A->data;
993:   PetscInt i;

996:   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i];
997:   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i];
998:   return(0);
999: }

1001: /*@C
1002:  MatNestGetLocalISs - Returns the index sets partitioning the row and column spaces

1004:  Not collective

1006:  Input Parameters:
1007: .   A  - nest matrix

1009:  Output Parameter:
1010: +   rows - array of row index sets (or NULL to ignore)
1011: -   cols - array of column index sets (or NULL to ignore)

1013:  Level: advanced

1015:  Notes:
1016:  The user must have allocated arrays of the correct size. The reference count is not increased on the returned ISs.

1018: .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs()
1019: @*/
1020: PetscErrorCode  MatNestGetLocalISs(Mat A,IS rows[],IS cols[])
1021: {

1026:   PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));
1027:   return(0);
1028: }

1030: PetscErrorCode  MatNestSetVecType_Nest(Mat A,VecType vtype)
1031: {
1033:   PetscBool      flg;

1036:   PetscStrcmp(vtype,VECNEST,&flg);
1037:   /* In reality, this only distinguishes VECNEST and "other" */
1038:   if (flg) A->ops->getvecs = MatCreateVecs_Nest;
1039:   else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0;
1040:   return(0);
1041: }

1043: /*@C
1044:  MatNestSetVecType - Sets the type of Vec returned by MatCreateVecs()

1046:  Not collective

1048:  Input Parameters:
1049: +  A  - nest matrix
1050: -  vtype - type to use for creating vectors

1052:  Notes:

1054:  Level: developer

1056: .seealso: MatCreateVecs()
1057: @*/
1058: PetscErrorCode  MatNestSetVecType(Mat A,VecType vtype)
1059: {

1063:   PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));
1064:   return(0);
1065: }

1067: PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1068: {
1069:   Mat_Nest       *s = (Mat_Nest*)A->data;
1070:   PetscInt       i,j,m,n,M,N;

1074:   s->nr = nr;
1075:   s->nc = nc;

1077:   /* Create space for submatrices */
1078:   PetscMalloc1(nr,&s->m);
1079:   for (i=0; i<nr; i++) {
1080:     PetscMalloc1(nc,&s->m[i]);
1081:   }
1082:   for (i=0; i<nr; i++) {
1083:     for (j=0; j<nc; j++) {
1084:       s->m[i][j] = a[i*nc+j];
1085:       if (a[i*nc+j]) {
1086:         PetscObjectReference((PetscObject)a[i*nc+j]);
1087:       }
1088:     }
1089:   }

1091:   MatSetUp_NestIS_Private(A,nr,is_row,nc,is_col);

1093:   PetscMalloc1(nr,&s->row_len);
1094:   PetscMalloc1(nc,&s->col_len);
1095:   for (i=0; i<nr; i++) s->row_len[i]=-1;
1096:   for (j=0; j<nc; j++) s->col_len[j]=-1;

1098:   MatNestGetSizes_Private(A,&m,&n,&M,&N);

1100:   PetscLayoutSetSize(A->rmap,M);
1101:   PetscLayoutSetLocalSize(A->rmap,m);
1102:   PetscLayoutSetSize(A->cmap,N);
1103:   PetscLayoutSetLocalSize(A->cmap,n);

1105:   PetscLayoutSetUp(A->rmap);
1106:   PetscLayoutSetUp(A->cmap);

1108:   PetscCalloc2(nr,&s->left,nc,&s->right);
1109:   return(0);
1110: }

1112: /*@
1113:    MatNestSetSubMats - Sets the nested submatrices

1115:    Collective on Mat

1117:    Input Parameter:
1118: +  N - nested matrix
1119: .  nr - number of nested row blocks
1120: .  is_row - index sets for each nested row block, or NULL to make contiguous
1121: .  nc - number of nested column blocks
1122: .  is_col - index sets for each nested column block, or NULL to make contiguous
1123: -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL

1125:    Level: advanced

1127: .seealso: MatCreateNest(), MATNEST
1128: @*/
1129: PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1130: {
1132:   PetscInt       i,nr_nc;

1136:   if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
1137:   if (nr && is_row) {
1140:   }
1141:   if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
1142:   if (nc && is_col) {
1145:   }
1146:   nr_nc=nr*nc;
1148:   PetscUseMethod(A,"MatNestSetSubMats_C",(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]),(A,nr,is_row,nc,is_col,a));
1149:   return(0);
1150: }

1152: static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog)
1153: {
1155:   PetscBool      flg;
1156:   PetscInt       i,j,m,mi,*ix;

1159:   for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) {
1160:     if (islocal[i]) {
1161:       ISGetSize(islocal[i],&mi);
1162:       flg  = PETSC_TRUE;      /* We found a non-trivial entry */
1163:     } else {
1164:       ISGetSize(isglobal[i],&mi);
1165:     }
1166:     m += mi;
1167:   }
1168:   if (flg) {
1169:     PetscMalloc1(m,&ix);
1170:     for (i=0,m=0; i<n; i++) {
1171:       ISLocalToGlobalMapping smap = NULL;
1172:       Mat                    sub = NULL;
1173:       PetscSF                sf;
1174:       PetscLayout            map;
1175:       PetscInt               *ix2;

1177:       if (!colflg) {
1178:         MatNestFindNonzeroSubMatRow(A,i,&sub);
1179:       } else {
1180:         MatNestFindNonzeroSubMatCol(A,i,&sub);
1181:       }
1182:       if (sub) {
1183:         if (!colflg) {
1184:           MatGetLocalToGlobalMapping(sub,&smap,NULL);
1185:         } else {
1186:           MatGetLocalToGlobalMapping(sub,NULL,&smap);
1187:         }
1188:       }
1189:       if (islocal[i]) {
1190:         ISGetSize(islocal[i],&mi);
1191:       } else {
1192:         ISGetSize(isglobal[i],&mi);
1193:       }
1194:       for (j=0; j<mi; j++) ix[m+j] = j;
1195:       if (smap) {ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);}

1197:       /*
1198:         Now we need to extract the monolithic global indices that correspond to the given split global indices.
1199:         In many/most cases, we only want MatGetLocalSubMatrix() to work, in which case we only need to know the size of the local spaces.
1200:        */
1201:       PetscMalloc1(mi,&ix2);
1202:       PetscSFCreate(((PetscObject)isglobal[i])->comm,&sf);
1203:       PetscLayoutCreate(((PetscObject)isglobal[i])->comm,&map);
1204:       PetscLayoutSetLocalSize(map,mi);
1205:       PetscLayoutSetUp(map);
1206:       PetscSFSetGraphLayout(sf,map,mi,NULL,PETSC_USE_POINTER,ix+m);
1207:       PetscLayoutDestroy(&map);
1208:       for (j=0; j<mi; j++) ix2[j] = ix[m+j];
1209:       PetscSFBcastBegin(sf,MPIU_INT,ix2,ix + m);
1210:       PetscSFBcastEnd(sf,MPIU_INT,ix2,ix + m);
1211:       PetscSFDestroy(&sf);
1212:       PetscFree(ix2);
1213:       m   += mi;
1214:     }
1215:     ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),1,m,ix,PETSC_OWN_POINTER,ltog);
1216:   } else {
1217:     *ltog = NULL;
1218:   }
1219:   return(0);
1220: }


1223: /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1224: /*
1225:   nprocessors = NP
1226:   Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1))
1227:        proc 0: => (g_0,h_0,)
1228:        proc 1: => (g_1,h_1,)
1229:        ...
1230:        proc nprocs-1: => (g_NP-1,h_NP-1,)

1232:             proc 0:                      proc 1:                    proc nprocs-1:
1233:     is[0] = (0,1,2,...,nlocal(g_0)-1)  (0,1,...,nlocal(g_1)-1)  (0,1,...,nlocal(g_NP-1))

1235:             proc 0:
1236:     is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1)
1237:             proc 1:
1238:     is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1)

1240:             proc NP-1:
1241:     is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1)
1242: */
1243: static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1244: {
1245:   Mat_Nest       *vs = (Mat_Nest*)A->data;
1246:   PetscInt       i,j,offset,n,nsum,bs;
1248:   Mat            sub = NULL;

1251:   PetscMalloc1(nr,&vs->isglobal.row);
1252:   PetscMalloc1(nc,&vs->isglobal.col);
1253:   if (is_row) { /* valid IS is passed in */
1254:     /* refs on is[] are incremeneted */
1255:     for (i=0; i<vs->nr; i++) {
1256:       PetscObjectReference((PetscObject)is_row[i]);

1258:       vs->isglobal.row[i] = is_row[i];
1259:     }
1260:   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
1261:     nsum = 0;
1262:     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
1263:       MatNestFindNonzeroSubMatRow(A,i,&sub);
1264:       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
1265:       MatGetLocalSize(sub,&n,NULL);
1266:       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
1267:       nsum += n;
1268:     }
1269:     MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));
1270:     offset -= nsum;
1271:     for (i=0; i<vs->nr; i++) {
1272:       MatNestFindNonzeroSubMatRow(A,i,&sub);
1273:       MatGetLocalSize(sub,&n,NULL);
1274:       MatGetBlockSize(sub,&bs);
1275:       ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);
1276:       ISSetBlockSize(vs->isglobal.row[i],bs);
1277:       offset += n;
1278:     }
1279:   }

1281:   if (is_col) { /* valid IS is passed in */
1282:     /* refs on is[] are incremeneted */
1283:     for (j=0; j<vs->nc; j++) {
1284:       PetscObjectReference((PetscObject)is_col[j]);

1286:       vs->isglobal.col[j] = is_col[j];
1287:     }
1288:   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
1289:     offset = A->cmap->rstart;
1290:     nsum   = 0;
1291:     for (j=0; j<vs->nc; j++) {
1292:       MatNestFindNonzeroSubMatCol(A,j,&sub);
1293:       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
1294:       MatGetLocalSize(sub,NULL,&n);
1295:       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
1296:       nsum += n;
1297:     }
1298:     MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));
1299:     offset -= nsum;
1300:     for (j=0; j<vs->nc; j++) {
1301:       MatNestFindNonzeroSubMatCol(A,j,&sub);
1302:       MatGetLocalSize(sub,NULL,&n);
1303:       MatGetBlockSize(sub,&bs);
1304:       ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);
1305:       ISSetBlockSize(vs->isglobal.col[j],bs);
1306:       offset += n;
1307:     }
1308:   }

1310:   /* Set up the local ISs */
1311:   PetscMalloc1(vs->nr,&vs->islocal.row);
1312:   PetscMalloc1(vs->nc,&vs->islocal.col);
1313:   for (i=0,offset=0; i<vs->nr; i++) {
1314:     IS                     isloc;
1315:     ISLocalToGlobalMapping rmap = NULL;
1316:     PetscInt               nlocal,bs;
1317:     MatNestFindNonzeroSubMatRow(A,i,&sub);
1318:     if (sub) {MatGetLocalToGlobalMapping(sub,&rmap,NULL);}
1319:     if (rmap) {
1320:       MatGetBlockSize(sub,&bs);
1321:       ISLocalToGlobalMappingGetSize(rmap,&nlocal);
1322:       ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);
1323:       ISSetBlockSize(isloc,bs);
1324:     } else {
1325:       nlocal = 0;
1326:       isloc  = NULL;
1327:     }
1328:     vs->islocal.row[i] = isloc;
1329:     offset            += nlocal;
1330:   }
1331:   for (i=0,offset=0; i<vs->nc; i++) {
1332:     IS                     isloc;
1333:     ISLocalToGlobalMapping cmap = NULL;
1334:     PetscInt               nlocal,bs;
1335:     MatNestFindNonzeroSubMatCol(A,i,&sub);
1336:     if (sub) {MatGetLocalToGlobalMapping(sub,NULL,&cmap);}
1337:     if (cmap) {
1338:       MatGetBlockSize(sub,&bs);
1339:       ISLocalToGlobalMappingGetSize(cmap,&nlocal);
1340:       ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);
1341:       ISSetBlockSize(isloc,bs);
1342:     } else {
1343:       nlocal = 0;
1344:       isloc  = NULL;
1345:     }
1346:     vs->islocal.col[i] = isloc;
1347:     offset            += nlocal;
1348:   }

1350:   /* Set up the aggregate ISLocalToGlobalMapping */
1351:   {
1352:     ISLocalToGlobalMapping rmap,cmap;
1353:     MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap);
1354:     MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap);
1355:     if (rmap && cmap) {MatSetLocalToGlobalMapping(A,rmap,cmap);}
1356:     ISLocalToGlobalMappingDestroy(&rmap);
1357:     ISLocalToGlobalMappingDestroy(&cmap);
1358:   }

1360: #if defined(PETSC_USE_DEBUG)
1361:   for (i=0; i<vs->nr; i++) {
1362:     for (j=0; j<vs->nc; j++) {
1363:       PetscInt m,n,M,N,mi,ni,Mi,Ni;
1364:       Mat      B = vs->m[i][j];
1365:       if (!B) continue;
1366:       MatGetSize(B,&M,&N);
1367:       MatGetLocalSize(B,&m,&n);
1368:       ISGetSize(vs->isglobal.row[i],&Mi);
1369:       ISGetSize(vs->isglobal.col[j],&Ni);
1370:       ISGetLocalSize(vs->isglobal.row[i],&mi);
1371:       ISGetLocalSize(vs->isglobal.col[j],&ni);
1372:       if (M != Mi || N != Ni) SETERRQ6(PetscObjectComm((PetscObject)sub),PETSC_ERR_ARG_INCOMP,"Global sizes (%D,%D) of nested submatrix (%D,%D) do not agree with space defined by index sets (%D,%D)",M,N,i,j,Mi,Ni);
1373:       if (m != mi || n != ni) SETERRQ6(PetscObjectComm((PetscObject)sub),PETSC_ERR_ARG_INCOMP,"Local sizes (%D,%D) of nested submatrix (%D,%D) do not agree with space defined by index sets (%D,%D)",m,n,i,j,mi,ni);
1374:     }
1375:   }
1376: #endif

1378:   /* Set A->assembled if all non-null blocks are currently assembled */
1379:   for (i=0; i<vs->nr; i++) {
1380:     for (j=0; j<vs->nc; j++) {
1381:       if (vs->m[i][j] && !vs->m[i][j]->assembled) return(0);
1382:     }
1383:   }
1384:   A->assembled = PETSC_TRUE;
1385:   return(0);
1386: }

1388: /*@C
1389:    MatCreateNest - Creates a new matrix containing several nested submatrices, each stored separately

1391:    Collective on Mat

1393:    Input Parameter:
1394: +  comm - Communicator for the new Mat
1395: .  nr - number of nested row blocks
1396: .  is_row - index sets for each nested row block, or NULL to make contiguous
1397: .  nc - number of nested column blocks
1398: .  is_col - index sets for each nested column block, or NULL to make contiguous
1399: -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL

1401:    Output Parameter:
1402: .  B - new matrix

1404:    Level: advanced

1406: .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST
1407: @*/
1408: PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1409: {
1410:   Mat            A;

1414:   *B   = 0;
1415:   MatCreate(comm,&A);
1416:   MatSetType(A,MATNEST);
1417:   A->preallocated = PETSC_TRUE;
1418:   MatNestSetSubMats(A,nr,is_row,nc,is_col,a);
1419:   *B   = A;
1420:   return(0);
1421: }

1423: static PetscErrorCode MatConvert_Nest_SeqAIJ_fast(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1424: {
1425:   Mat_Nest       *nest = (Mat_Nest*)A->data;
1426:   Mat            *trans;
1427:   PetscScalar    **avv;
1428:   PetscScalar    *vv;
1429:   PetscInt       **aii,**ajj;
1430:   PetscInt       *ii,*jj,*ci;
1431:   PetscInt       nr,nc,nnz,i,j;
1432:   PetscBool      done;

1436:   MatGetSize(A,&nr,&nc);
1437:   if (reuse == MAT_REUSE_MATRIX) {
1438:     PetscInt rnr;

1440:     MatGetRowIJ(*newmat,0,PETSC_FALSE,PETSC_FALSE,&rnr,(const PetscInt**)&ii,(const PetscInt**)&jj,&done);
1441:     if (!done) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"MatGetRowIJ");
1442:     if (rnr != nr) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Cannot reuse matrix, wrong number of rows");
1443:     MatSeqAIJGetArray(*newmat,&vv);
1444:   }
1445:   /* extract CSR for nested SeqAIJ matrices */
1446:   nnz  = 0;
1447:   PetscCalloc4(nest->nr*nest->nc,&aii,nest->nr*nest->nc,&ajj,nest->nr*nest->nc,&avv,nest->nr*nest->nc,&trans);
1448:   for (i=0; i<nest->nr; ++i) {
1449:     for (j=0; j<nest->nc; ++j) {
1450:       Mat B = nest->m[i][j];
1451:       if (B) {
1452:         PetscScalar *naa;
1453:         PetscInt    *nii,*njj,nnr;
1454:         PetscBool   istrans;

1456:         PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&istrans);
1457:         if (istrans) {
1458:           Mat Bt;

1460:           MatTransposeGetMat(B,&Bt);
1461:           MatTranspose(Bt,MAT_INITIAL_MATRIX,&trans[i*nest->nc+j]);
1462:           B    = trans[i*nest->nc+j];
1463:         }
1464:         MatGetRowIJ(B,0,PETSC_FALSE,PETSC_FALSE,&nnr,(const PetscInt**)&nii,(const PetscInt**)&njj,&done);
1465:         if (!done) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"MatGetRowIJ");
1466:         MatSeqAIJGetArray(B,&naa);
1467:         nnz += nii[nnr];

1469:         aii[i*nest->nc+j] = nii;
1470:         ajj[i*nest->nc+j] = njj;
1471:         avv[i*nest->nc+j] = naa;
1472:       }
1473:     }
1474:   }
1475:   if (reuse != MAT_REUSE_MATRIX) {
1476:     PetscMalloc1(nr+1,&ii);
1477:     PetscMalloc1(nnz,&jj);
1478:     PetscMalloc1(nnz,&vv);
1479:   } else {
1480:     if (nnz != ii[nr]) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Cannot reuse matrix, wrong number of nonzeros");
1481:   }

1483:   /* new row pointer */
1484:   PetscMemzero(ii,(nr+1)*sizeof(PetscInt));
1485:   for (i=0; i<nest->nr; ++i) {
1486:     PetscInt       ncr,rst;

1488:     ISStrideGetInfo(nest->isglobal.row[i],&rst,NULL);
1489:     ISGetLocalSize(nest->isglobal.row[i],&ncr);
1490:     for (j=0; j<nest->nc; ++j) {
1491:       if (aii[i*nest->nc+j]) {
1492:         PetscInt    *nii = aii[i*nest->nc+j];
1493:         PetscInt    ir;

1495:         for (ir=rst; ir<ncr+rst; ++ir) {
1496:           ii[ir+1] += nii[1]-nii[0];
1497:           nii++;
1498:         }
1499:       }
1500:     }
1501:   }
1502:   for (i=0; i<nr; i++) ii[i+1] += ii[i];

1504:   /* construct CSR for the new matrix */
1505:   PetscCalloc1(nr,&ci);
1506:   for (i=0; i<nest->nr; ++i) {
1507:     PetscInt       ncr,rst;

1509:     ISStrideGetInfo(nest->isglobal.row[i],&rst,NULL);
1510:     ISGetLocalSize(nest->isglobal.row[i],&ncr);
1511:     for (j=0; j<nest->nc; ++j) {
1512:       if (aii[i*nest->nc+j]) {
1513:         PetscScalar *nvv = avv[i*nest->nc+j];
1514:         PetscInt    *nii = aii[i*nest->nc+j];
1515:         PetscInt    *njj = ajj[i*nest->nc+j];
1516:         PetscInt    ir,cst;

1518:         ISStrideGetInfo(nest->isglobal.col[j],&cst,NULL);
1519:         for (ir=rst; ir<ncr+rst; ++ir) {
1520:           PetscInt ij,rsize = nii[1]-nii[0],ist = ii[ir]+ci[ir];

1522:           for (ij=0;ij<rsize;ij++) {
1523:             jj[ist+ij] = *njj+cst;
1524:             vv[ist+ij] = *nvv;
1525:             njj++;
1526:             nvv++;
1527:           }
1528:           ci[ir] += rsize;
1529:           nii++;
1530:         }
1531:       }
1532:     }
1533:   }
1534:   PetscFree(ci);

1536:   /* restore info */
1537:   for (i=0; i<nest->nr; ++i) {
1538:     for (j=0; j<nest->nc; ++j) {
1539:       Mat B = nest->m[i][j];
1540:       if (B) {
1541:         PetscInt nnr = 0, k = i*nest->nc+j;

1543:         B    = (trans[k] ? trans[k] : B);
1544:         MatRestoreRowIJ(B,0,PETSC_FALSE,PETSC_FALSE,&nnr,(const PetscInt**)&aii[k],(const PetscInt**)&ajj[k],&done);
1545:         if (!done) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"MatRestoreRowIJ");
1546:         MatSeqAIJRestoreArray(B,&avv[k]);
1547:         MatDestroy(&trans[k]);
1548:       }
1549:     }
1550:   }
1551:   PetscFree4(aii,ajj,avv,trans);

1553:   /* finalize newmat */
1554:   if (reuse == MAT_INITIAL_MATRIX) {
1555:     MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),nr,nc,ii,jj,vv,newmat);
1556:   } else if (reuse == MAT_INPLACE_MATRIX) {
1557:     Mat B;

1559:     MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),nr,nc,ii,jj,vv,&B);
1560:     MatHeaderReplace(A,&B);
1561:   }
1562:   MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);
1563:   MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);
1564:   {
1565:     Mat_SeqAIJ *a = (Mat_SeqAIJ*)((*newmat)->data);
1566:     a->free_a     = PETSC_TRUE;
1567:     a->free_ij    = PETSC_TRUE;
1568:   }
1569:   return(0);
1570: }

1572: PETSC_INTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1573: {
1575:   Mat_Nest       *nest = (Mat_Nest*)A->data;
1576:   PetscInt       m,n,M,N,i,j,k,*dnnz,*onnz,rstart;
1577:   PetscInt       cstart,cend;
1578:   PetscMPIInt    size;
1579:   Mat            C;

1582:   MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);
1583:   if (size == 1) { /* look for a special case with SeqAIJ matrices and strided-1, contiguous, blocks */
1584:     PetscInt  nf;
1585:     PetscBool fast;

1587:     PetscStrcmp(newtype,MATAIJ,&fast);
1588:     if (!fast) {
1589:       PetscStrcmp(newtype,MATSEQAIJ,&fast);
1590:     }
1591:     for (i=0; i<nest->nr && fast; ++i) {
1592:       for (j=0; j<nest->nc && fast; ++j) {
1593:         Mat B = nest->m[i][j];
1594:         if (B) {
1595:           PetscObjectTypeCompare((PetscObject)B,MATSEQAIJ,&fast);
1596:           if (!fast) {
1597:             PetscBool istrans;

1599:             PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&istrans);
1600:             if (istrans) {
1601:               Mat Bt;

1603:               MatTransposeGetMat(B,&Bt);
1604:               PetscObjectTypeCompare((PetscObject)Bt,MATSEQAIJ,&fast);
1605:             }
1606:           }
1607:         }
1608:       }
1609:     }
1610:     for (i=0, nf=0; i<nest->nr && fast; ++i) {
1611:       PetscObjectTypeCompare((PetscObject)nest->isglobal.row[i],ISSTRIDE,&fast);
1612:       if (fast) {
1613:         PetscInt f,s;

1615:         ISStrideGetInfo(nest->isglobal.row[i],&f,&s);
1616:         if (f != nf || s != 1) { fast = PETSC_FALSE; }
1617:         else {
1618:           ISGetSize(nest->isglobal.row[i],&f);
1619:           nf  += f;
1620:         }
1621:       }
1622:     }
1623:     for (i=0, nf=0; i<nest->nc && fast; ++i) {
1624:       PetscObjectTypeCompare((PetscObject)nest->isglobal.col[i],ISSTRIDE,&fast);
1625:       if (fast) {
1626:         PetscInt f,s;

1628:         ISStrideGetInfo(nest->isglobal.col[i],&f,&s);
1629:         if (f != nf || s != 1) { fast = PETSC_FALSE; }
1630:         else {
1631:           ISGetSize(nest->isglobal.col[i],&f);
1632:           nf  += f;
1633:         }
1634:       }
1635:     }
1636:     if (fast) {
1637:       MatConvert_Nest_SeqAIJ_fast(A,newtype,reuse,newmat);
1638:       return(0);
1639:     }
1640:   }
1641:   MatGetSize(A,&M,&N);
1642:   MatGetLocalSize(A,&m,&n);
1643:   MatGetOwnershipRangeColumn(A,&cstart,&cend);
1644:   switch (reuse) {
1645:   case MAT_INITIAL_MATRIX:
1646:     MatCreate(PetscObjectComm((PetscObject)A),&C);
1647:     MatSetType(C,newtype);
1648:     MatSetSizes(C,m,n,M,N);
1649:     *newmat = C;
1650:     break;
1651:   case MAT_REUSE_MATRIX:
1652:     C = *newmat;
1653:     break;
1654:   default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse");
1655:   }
1656:   PetscMalloc1(2*m,&dnnz);
1657:   onnz = dnnz + m;
1658:   for (k=0; k<m; k++) {
1659:     dnnz[k] = 0;
1660:     onnz[k] = 0;
1661:   }
1662:   for (j=0; j<nest->nc; ++j) {
1663:     IS             bNis;
1664:     PetscInt       bN;
1665:     const PetscInt *bNindices;
1666:     /* Using global column indices and ISAllGather() is not scalable. */
1667:     ISAllGather(nest->isglobal.col[j], &bNis);
1668:     ISGetSize(bNis, &bN);
1669:     ISGetIndices(bNis,&bNindices);
1670:     for (i=0; i<nest->nr; ++i) {
1671:       PetscSF        bmsf;
1672:       PetscSFNode    *iremote;
1673:       Mat            B;
1674:       PetscInt       bm, *sub_dnnz,*sub_onnz, br;
1675:       const PetscInt *bmindices;
1676:       B = nest->m[i][j];
1677:       if (!B) continue;
1678:       ISGetLocalSize(nest->isglobal.row[i],&bm);
1679:       ISGetIndices(nest->isglobal.row[i],&bmindices);
1680:       PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);
1681:       PetscMalloc1(bm,&iremote);
1682:       PetscMalloc1(bm,&sub_dnnz);
1683:       PetscMalloc1(bm,&sub_onnz);
1684:       for (k = 0; k < bm; ++k){
1685:             sub_dnnz[k] = 0;
1686:             sub_onnz[k] = 0;
1687:       }
1688:       /*
1689:        Locate the owners for all of the locally-owned global row indices for this row block.
1690:        These determine the roots of PetscSF used to communicate preallocation data to row owners.
1691:        The roots correspond to the dnnz and onnz entries; thus, there are two roots per row.
1692:        */
1693:       MatGetOwnershipRange(B,&rstart,NULL);
1694:       for (br = 0; br < bm; ++br) {
1695:         PetscInt       row = bmindices[br], rowowner = 0, brncols, col;
1696:         const PetscInt *brcols;
1697:         PetscInt       rowrel = 0; /* row's relative index on its owner rank */
1698:         PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);
1699:         /* how many roots  */
1700:         iremote[br].rank = rowowner; iremote[br].index = rowrel;           /* edge from bmdnnz to dnnz */
1701:         /* get nonzero pattern */
1702:         MatGetRow(B,br+rstart,&brncols,&brcols,NULL);
1703:         for (k=0; k<brncols; k++) {
1704:           col  = bNindices[brcols[k]];
1705:           if (col>=A->cmap->range[rowowner] && col<A->cmap->range[rowowner+1]) {
1706:             sub_dnnz[br]++;
1707:           } else {
1708:             sub_onnz[br]++;
1709:           }
1710:         }
1711:         MatRestoreRow(B,br+rstart,&brncols,&brcols,NULL);
1712:       }
1713:       ISRestoreIndices(nest->isglobal.row[i],&bmindices);
1714:       /* bsf will have to take care of disposing of bedges. */
1715:       PetscSFSetGraph(bmsf,m,bm,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);
1716:       PetscSFReduceBegin(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);
1717:       PetscSFReduceEnd(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);
1718:       PetscSFReduceBegin(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);
1719:       PetscSFReduceEnd(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);
1720:       PetscFree(sub_dnnz);
1721:       PetscFree(sub_onnz);
1722:       PetscSFDestroy(&bmsf);
1723:     }
1724:     ISRestoreIndices(bNis,&bNindices);
1725:     ISDestroy(&bNis);
1726:   }
1727:   /* Resize preallocation if overestimated */
1728:   for (i=0;i<m;i++) {
1729:     dnnz[i] = PetscMin(dnnz[i],A->cmap->n);
1730:     onnz[i] = PetscMin(onnz[i],A->cmap->N - A->cmap->n);
1731:   }
1732:   MatSeqAIJSetPreallocation(C,0,dnnz);
1733:   MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);
1734:   PetscFree(dnnz);

1736:   /* Fill by row */
1737:   for (j=0; j<nest->nc; ++j) {
1738:     /* Using global column indices and ISAllGather() is not scalable. */
1739:     IS             bNis;
1740:     PetscInt       bN;
1741:     const PetscInt *bNindices;
1742:     ISAllGather(nest->isglobal.col[j], &bNis);
1743:     ISGetSize(bNis,&bN);
1744:     ISGetIndices(bNis,&bNindices);
1745:     for (i=0; i<nest->nr; ++i) {
1746:       Mat            B;
1747:       PetscInt       bm, br;
1748:       const PetscInt *bmindices;
1749:       B = nest->m[i][j];
1750:       if (!B) continue;
1751:       ISGetLocalSize(nest->isglobal.row[i],&bm);
1752:       ISGetIndices(nest->isglobal.row[i],&bmindices);
1753:       MatGetOwnershipRange(B,&rstart,NULL);
1754:       for (br = 0; br < bm; ++br) {
1755:         PetscInt          row = bmindices[br], brncols,  *cols;
1756:         const PetscInt    *brcols;
1757:         const PetscScalar *brcoldata;
1758:         MatGetRow(B,br+rstart,&brncols,&brcols,&brcoldata);
1759:         PetscMalloc1(brncols,&cols);
1760:         for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]];
1761:         /*
1762:           Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match.
1763:           Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES.
1764:          */
1765:         MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);
1766:         MatRestoreRow(B,br+rstart,&brncols,&brcols,&brcoldata);
1767:         PetscFree(cols);
1768:       }
1769:       ISRestoreIndices(nest->isglobal.row[i],&bmindices);
1770:     }
1771:     ISRestoreIndices(bNis,&bNindices);
1772:     ISDestroy(&bNis);
1773:   }
1774:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1775:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1776:   return(0);
1777: }

1779: PetscErrorCode  MatHasOperation_Nest(Mat mat,MatOperation op,PetscBool  *has)
1780: {
1782:   *has = PETSC_FALSE;
1783:   if (op == MATOP_MULT_TRANSPOSE) {
1784:     Mat_Nest       *bA = (Mat_Nest*)mat->data;
1785:     PetscInt       i,j,nr = bA->nr,nc = bA->nc;
1787:     PetscBool      flg;

1789:     for (j=0; j<nc; j++) {
1790:       for (i=0; i<nr; i++) {
1791:         if (!bA->m[i][j]) continue;
1792:         MatHasOperation(bA->m[i][j],MATOP_MULT_TRANSPOSE_ADD,&flg);
1793:         if (!flg) return(0);
1794:       }
1795:     }
1796:   }
1797:   if (((void**)mat->ops)[op]) *has =  PETSC_TRUE;
1798:   return(0);
1799: }

1801: /*MC
1802:   MATNEST - MATNEST = "nest" - Matrix type consisting of nested submatrices, each stored separately.

1804:   Level: intermediate

1806:   Notes:
1807:   This matrix type permits scalable use of PCFieldSplit and avoids the large memory costs of extracting submatrices.
1808:   It allows the use of symmetric and block formats for parts of multi-physics simulations.
1809:   It is usually used with DMComposite and DMCreateMatrix()

1811:   Each of the submatrices lives on the same MPI communicator as the original nest matrix (though they can have zero
1812:   rows/columns on some processes.) Thus this is not meant for cases where the submatrices live on far fewer processes
1813:   than the nest matrix.

1815: .seealso: MatCreate(), MatType, MatCreateNest()
1816: M*/
1817: PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A)
1818: {
1819:   Mat_Nest       *s;

1823:   PetscNewLog(A,&s);
1824:   A->data = (void*)s;

1826:   s->nr            = -1;
1827:   s->nc            = -1;
1828:   s->m             = NULL;
1829:   s->splitassembly = PETSC_FALSE;

1831:   PetscMemzero(A->ops,sizeof(*A->ops));

1833:   A->ops->mult                  = MatMult_Nest;
1834:   A->ops->multadd               = MatMultAdd_Nest;
1835:   A->ops->multtranspose         = MatMultTranspose_Nest;
1836:   A->ops->multtransposeadd      = MatMultTransposeAdd_Nest;
1837:   A->ops->transpose             = MatTranspose_Nest;
1838:   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1839:   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1840:   A->ops->zeroentries           = MatZeroEntries_Nest;
1841:   A->ops->copy                  = MatCopy_Nest;
1842:   A->ops->axpy                  = MatAXPY_Nest;
1843:   A->ops->duplicate             = MatDuplicate_Nest;
1844:   A->ops->createsubmatrix       = MatCreateSubMatrix_Nest;
1845:   A->ops->destroy               = MatDestroy_Nest;
1846:   A->ops->view                  = MatView_Nest;
1847:   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1848:   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1849:   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
1850:   A->ops->getdiagonal           = MatGetDiagonal_Nest;
1851:   A->ops->diagonalscale         = MatDiagonalScale_Nest;
1852:   A->ops->scale                 = MatScale_Nest;
1853:   A->ops->shift                 = MatShift_Nest;
1854:   A->ops->diagonalset           = MatDiagonalSet_Nest;
1855:   A->ops->setrandom             = MatSetRandom_Nest;
1856:   A->ops->hasoperation          = MatHasOperation_Nest;

1858:   A->spptr        = 0;
1859:   A->assembled    = PETSC_FALSE;

1861:   /* expose Nest api's */
1862:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",      MatNestGetSubMat_Nest);
1863:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",      MatNestSetSubMat_Nest);
1864:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",     MatNestGetSubMats_Nest);
1865:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",        MatNestGetSize_Nest);
1866:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",         MatNestGetISs_Nest);
1867:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",    MatNestGetLocalISs_Nest);
1868:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",     MatNestSetVecType_Nest);
1869:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",     MatNestSetSubMats_Nest);
1870:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_mpiaij_C",MatConvert_Nest_AIJ);
1871:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_seqaij_C",MatConvert_Nest_AIJ);
1872:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",   MatConvert_Nest_AIJ);
1873:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_is_C",    MatConvert_Nest_IS);

1875:   PetscObjectChangeTypeName((PetscObject)A,MATNEST);
1876:   return(0);
1877: }