Actual source code: matnest.c

petsc-master 2016-07-30
Report Typos and Errors
  2: #include <../src/mat/impls/nest/matnestimpl.h> /*I   "petscmat.h"   I*/
  3: #include <petscsf.h>

  5: static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]);
  6: static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left);

  8: /* 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 */
 39: static PetscErrorCode MatMult_Nest(Mat A,Vec x,Vec y)
 40: {
 41:   Mat_Nest       *bA = (Mat_Nest*)A->data;
 42:   Vec            *bx = bA->right,*by = bA->left;
 43:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;

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

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

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

 94: static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y)
 95: {
 96:   Mat_Nest       *bA = (Mat_Nest*)A->data;
 97:   Vec            *bx = bA->left,*by = bA->right;
 98:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;

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

119: static PetscErrorCode MatMultTransposeAdd_Nest(Mat A,Vec x,Vec y,Vec z)
120: {
121:   Mat_Nest       *bA = (Mat_Nest*)A->data;
122:   Vec            *bx = bA->left,*bz = bA->right;
123:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;

127:   for (i=0; i<nr; i++) {VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);}
128:   for (i=0; i<nc; i++) {VecGetSubVector(z,bA->isglobal.col[i],&bz[i]);}
129:   for (j=0; j<nc; j++) {
130:     if (y != z) {
131:       Vec by;
132:       VecGetSubVector(y,bA->isglobal.col[j],&by);
133:       VecCopy(by,bz[j]);
134:       VecRestoreSubVector(y,bA->isglobal.col[j],&by);
135:     }
136:     for (i=0; i<nr; i++) {
137:       if (!bA->m[i][j]) continue;
138:       /* z[j] <- y[j] + (A[i][j])^T * x[i] */
139:       MatMultTransposeAdd(bA->m[i][j],bx[i],bz[j],bz[j]);
140:     }
141:   }
142:   for (i=0; i<nr; i++) {VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);}
143:   for (i=0; i<nc; i++) {VecRestoreSubVector(z,bA->isglobal.col[i],&bz[i]);}
144:   return(0);
145: }

149: static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list)
150: {
152:   IS             *lst = *list;
153:   PetscInt       i;

156:   if (!lst) return(0);
157:   for (i=0; i<n; i++) if (lst[i]) {ISDestroy(&lst[i]);}
158:   PetscFree(lst);
159:   *list = NULL;
160:   return(0);
161: }

165: static PetscErrorCode MatDestroy_Nest(Mat A)
166: {
167:   Mat_Nest       *vs = (Mat_Nest*)A->data;
168:   PetscInt       i,j;

172:   /* release the matrices and the place holders */
173:   MatNestDestroyISList(vs->nr,&vs->isglobal.row);
174:   MatNestDestroyISList(vs->nc,&vs->isglobal.col);
175:   MatNestDestroyISList(vs->nr,&vs->islocal.row);
176:   MatNestDestroyISList(vs->nc,&vs->islocal.col);

178:   PetscFree(vs->row_len);
179:   PetscFree(vs->col_len);

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

183:   /* release the matrices and the place holders */
184:   if (vs->m) {
185:     for (i=0; i<vs->nr; i++) {
186:       for (j=0; j<vs->nc; j++) {
187:         MatDestroy(&vs->m[i][j]);
188:       }
189:       PetscFree(vs->m[i]);
190:     }
191:     PetscFree(vs->m);
192:   }
193:   PetscFree(A->data);

195:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",0);
196:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",0);
197:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",0);
198:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",0);
199:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",0);
200:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",0);
201:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",0);
202:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",0);
203:   return(0);
204: }

208: static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
209: {
210:   Mat_Nest       *vs = (Mat_Nest*)A->data;
211:   PetscInt       i,j;

215:   for (i=0; i<vs->nr; i++) {
216:     for (j=0; j<vs->nc; j++) {
217:       if (vs->m[i][j]) {
218:         MatAssemblyBegin(vs->m[i][j],type);
219:         if (!vs->splitassembly) {
220:           /* Note: split assembly will fail if the same block appears more than once (even indirectly through a nested
221:            * sub-block). This could be fixed by adding a flag to Mat so that there was a way to check if a Mat was
222:            * already performing an assembly, but the result would by more complicated and appears to offer less
223:            * potential for diagnostics and correctness checking. Split assembly should be fixed once there is an
224:            * interface for libraries to make asynchronous progress in "user-defined non-blocking collectives".
225:            */
226:           MatAssemblyEnd(vs->m[i][j],type);
227:         }
228:       }
229:     }
230:   }
231:   return(0);
232: }

236: static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
237: {
238:   Mat_Nest       *vs = (Mat_Nest*)A->data;
239:   PetscInt       i,j;

243:   for (i=0; i<vs->nr; i++) {
244:     for (j=0; j<vs->nc; j++) {
245:       if (vs->m[i][j]) {
246:         if (vs->splitassembly) {
247:           MatAssemblyEnd(vs->m[i][j],type);
248:         }
249:       }
250:     }
251:   }
252:   return(0);
253: }

257: static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
258: {
260:   Mat_Nest       *vs = (Mat_Nest*)A->data;
261:   PetscInt       j;
262:   Mat            sub;

265:   sub = (row < vs->nc) ? vs->m[row][row] : (Mat)NULL; /* Prefer to find on the diagonal */
266:   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
267:   if (sub) {MatSetUp(sub);}       /* Ensure that the sizes are available */
268:   *B = sub;
269:   return(0);
270: }

274: static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
275: {
277:   Mat_Nest       *vs = (Mat_Nest*)A->data;
278:   PetscInt       i;
279:   Mat            sub;

282:   sub = (col < vs->nr) ? vs->m[col][col] : (Mat)NULL; /* Prefer to find on the diagonal */
283:   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
284:   if (sub) {MatSetUp(sub);}       /* Ensure that the sizes are available */
285:   *B = sub;
286:   return(0);
287: }

291: static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
292: {
294:   PetscInt       i;
295:   PetscBool      flg;

301:   *found = -1;
302:   for (i=0; i<n; i++) {
303:     if (!list[i]) continue;
304:     ISEqual(list[i],is,&flg);
305:     if (flg) {
306:       *found = i;
307:       return(0);
308:     }
309:   }
310:   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Could not find index set");
311:   return(0);
312: }

316: /* Get a block row as a new MatNest */
317: static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B)
318: {
319:   Mat_Nest       *vs = (Mat_Nest*)A->data;
320:   char           keyname[256];

324:   *B   = NULL;
325:   PetscSNPrintf(keyname,sizeof(keyname),"NestRow_%D",row);
326:   PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);
327:   if (*B) return(0);

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

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

333:   PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);
334:   PetscObjectDereference((PetscObject)*B); /* Leave the only remaining reference in the composition */
335:   return(0);
336: }

340: static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
341: {
342:   Mat_Nest       *vs = (Mat_Nest*)A->data;
344:   PetscInt       row,col;
345:   PetscBool      same,isFullCol,isFullColGlobal;

348:   /* Check if full column space. This is a hack */
349:   isFullCol = PETSC_FALSE;
350:   PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&same);
351:   if (same) {
352:     PetscInt n,first,step,i,an,am,afirst,astep;
353:     ISStrideGetInfo(iscol,&first,&step);
354:     ISGetLocalSize(iscol,&n);
355:     isFullCol = PETSC_TRUE;
356:     for (i=0,an=A->cmap->rstart; i<vs->nc; i++) {
357:       ISStrideGetInfo(is->col[i],&afirst,&astep);
358:       ISGetLocalSize(is->col[i],&am);
359:       if (afirst != an || astep != step) isFullCol = PETSC_FALSE;
360:       an += am;
361:     }
362:     if (an != A->cmap->rstart+n) isFullCol = PETSC_FALSE;
363:   }
364:   MPIU_Allreduce(&isFullCol,&isFullColGlobal,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)iscol));

366:   if (isFullColGlobal) {
367:     PetscInt row;
368:     MatNestFindIS(A,vs->nr,is->row,isrow,&row);
369:     MatNestGetRow(A,row,B);
370:   } else {
371:     MatNestFindIS(A,vs->nr,is->row,isrow,&row);
372:     MatNestFindIS(A,vs->nc,is->col,iscol,&col);
373:     *B   = vs->m[row][col];
374:   }
375:   return(0);
376: }

380: static PetscErrorCode MatGetSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
381: {
383:   Mat_Nest       *vs = (Mat_Nest*)A->data;
384:   Mat            sub;

387:   MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);
388:   switch (reuse) {
389:   case MAT_INITIAL_MATRIX:
390:     if (sub) { PetscObjectReference((PetscObject)sub); }
391:     *B = sub;
392:     break;
393:   case MAT_REUSE_MATRIX:
394:     if (sub != *B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
395:     break;
396:   case MAT_IGNORE_MATRIX:       /* Nothing to do */
397:     break;
398:   case MAT_INPLACE_MATRIX:       /* Nothing to do */
399:     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX is not supported yet");
400:     break;
401:   }
402:   return(0);
403: }

407: PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
408: {
410:   Mat_Nest       *vs = (Mat_Nest*)A->data;
411:   Mat            sub;

414:   MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);
415:   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
416:   if (sub) {PetscObjectReference((PetscObject)sub);}
417:   *B = sub;
418:   return(0);
419: }

423: static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
424: {
426:   Mat_Nest       *vs = (Mat_Nest*)A->data;
427:   Mat            sub;

430:   MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);
431:   if (*B != sub) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
432:   if (sub) {
433:     if (((PetscObject)sub)->refct <= 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
434:     MatDestroy(B);
435:   }
436:   return(0);
437: }

441: static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
442: {
443:   Mat_Nest       *bA = (Mat_Nest*)A->data;
444:   PetscInt       i;

448:   for (i=0; i<bA->nr; i++) {
449:     Vec bv;
450:     VecGetSubVector(v,bA->isglobal.row[i],&bv);
451:     if (bA->m[i][i]) {
452:       MatGetDiagonal(bA->m[i][i],bv);
453:     } else {
454:       VecSet(bv,0.0);
455:     }
456:     VecRestoreSubVector(v,bA->isglobal.row[i],&bv);
457:   }
458:   return(0);
459: }

463: static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
464: {
465:   Mat_Nest       *bA = (Mat_Nest*)A->data;
466:   Vec            bl,*br;
467:   PetscInt       i,j;

471:   PetscCalloc1(bA->nc,&br);
472:   if (r) {
473:     for (j=0; j<bA->nc; j++) {VecGetSubVector(r,bA->isglobal.col[j],&br[j]);}
474:   }
475:   bl = NULL;
476:   for (i=0; i<bA->nr; i++) {
477:     if (l) {
478:       VecGetSubVector(l,bA->isglobal.row[i],&bl);
479:     }
480:     for (j=0; j<bA->nc; j++) {
481:       if (bA->m[i][j]) {
482:         MatDiagonalScale(bA->m[i][j],bl,br[j]);
483:       }
484:     }
485:     if (l) {
486:       VecRestoreSubVector(l,bA->isglobal.row[i],&bl);
487:     }
488:   }
489:   if (r) {
490:     for (j=0; j<bA->nc; j++) {VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);}
491:   }
492:   PetscFree(br);
493:   return(0);
494: }

498: static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a)
499: {
500:   Mat_Nest       *bA = (Mat_Nest*)A->data;
501:   PetscInt       i,j;

505:   for (i=0; i<bA->nr; i++) {
506:     for (j=0; j<bA->nc; j++) {
507:       if (bA->m[i][j]) {
508:         MatScale(bA->m[i][j],a);
509:       }
510:     }
511:   }
512:   return(0);
513: }

517: static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a)
518: {
519:   Mat_Nest       *bA = (Mat_Nest*)A->data;
520:   PetscInt       i;

524:   for (i=0; i<bA->nr; i++) {
525:     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);
526:     MatShift(bA->m[i][i],a);
527:   }
528:   return(0);
529: }

533: static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left)
534: {
535:   Mat_Nest       *bA = (Mat_Nest*)A->data;
536:   Vec            *L,*R;
537:   MPI_Comm       comm;
538:   PetscInt       i,j;

542:   PetscObjectGetComm((PetscObject)A,&comm);
543:   if (right) {
544:     /* allocate R */
545:     PetscMalloc1(bA->nc, &R);
546:     /* Create the right vectors */
547:     for (j=0; j<bA->nc; j++) {
548:       for (i=0; i<bA->nr; i++) {
549:         if (bA->m[i][j]) {
550:           MatCreateVecs(bA->m[i][j],&R[j],NULL);
551:           break;
552:         }
553:       }
554:       if (i==bA->nr) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
555:     }
556:     VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);
557:     /* hand back control to the nest vector */
558:     for (j=0; j<bA->nc; j++) {
559:       VecDestroy(&R[j]);
560:     }
561:     PetscFree(R);
562:   }

564:   if (left) {
565:     /* allocate L */
566:     PetscMalloc1(bA->nr, &L);
567:     /* Create the left vectors */
568:     for (i=0; i<bA->nr; i++) {
569:       for (j=0; j<bA->nc; j++) {
570:         if (bA->m[i][j]) {
571:           MatCreateVecs(bA->m[i][j],NULL,&L[i]);
572:           break;
573:         }
574:       }
575:       if (j==bA->nc) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
576:     }

578:     VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);
579:     for (i=0; i<bA->nr; i++) {
580:       VecDestroy(&L[i]);
581:     }

583:     PetscFree(L);
584:   }
585:   return(0);
586: }

590: static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
591: {
592:   Mat_Nest       *bA = (Mat_Nest*)A->data;
593:   PetscBool      isascii;
594:   PetscInt       i,j;

598:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
599:   if (isascii) {

601:     PetscViewerASCIIPrintf(viewer,"Matrix object: \n");
602:     PetscViewerASCIIPushTab(viewer);    /* push0 */
603:     PetscViewerASCIIPrintf(viewer, "type=nest, rows=%d, cols=%d \n",bA->nr,bA->nc);

605:     PetscViewerASCIIPrintf(viewer,"MatNest structure: \n");
606:     for (i=0; i<bA->nr; i++) {
607:       for (j=0; j<bA->nc; j++) {
608:         MatType   type;
609:         char      name[256] = "",prefix[256] = "";
610:         PetscInt  NR,NC;
611:         PetscBool isNest = PETSC_FALSE;

613:         if (!bA->m[i][j]) {
614:           PetscViewerASCIIPrintf(viewer, "(%D,%D) : NULL \n",i,j);
615:           continue;
616:         }
617:         MatGetSize(bA->m[i][j],&NR,&NC);
618:         MatGetType(bA->m[i][j], &type);
619:         if (((PetscObject)bA->m[i][j])->name) {PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);}
620:         if (((PetscObject)bA->m[i][j])->prefix) {PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);}
621:         PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);

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

625:         if (isNest) {
626:           PetscViewerASCIIPushTab(viewer);  /* push1 */
627:           MatView(bA->m[i][j],viewer);
628:           PetscViewerASCIIPopTab(viewer);    /* pop1 */
629:         }
630:       }
631:     }
632:     PetscViewerASCIIPopTab(viewer);    /* pop0 */
633:   }
634:   return(0);
635: }

639: static PetscErrorCode MatZeroEntries_Nest(Mat A)
640: {
641:   Mat_Nest       *bA = (Mat_Nest*)A->data;
642:   PetscInt       i,j;

646:   for (i=0; i<bA->nr; i++) {
647:     for (j=0; j<bA->nc; j++) {
648:       if (!bA->m[i][j]) continue;
649:       MatZeroEntries(bA->m[i][j]);
650:     }
651:   }
652:   return(0);
653: }

657: static PetscErrorCode MatCopy_Nest(Mat A,Mat B,MatStructure str)
658: {
659:   Mat_Nest       *bA = (Mat_Nest*)A->data,*bB = (Mat_Nest*)B->data;
660:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;

664:   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);
665:   for (i=0; i<nr; i++) {
666:     for (j=0; j<nc; j++) {
667:       if (bA->m[i][j] && bB->m[i][j]) {
668:         MatCopy(bA->m[i][j],bB->m[i][j],str);
669:       } 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);
670:     }
671:   }
672:   return(0);
673: }

677: static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
678: {
679:   Mat_Nest       *bA = (Mat_Nest*)A->data;
680:   Mat            *b;
681:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;

685:   PetscMalloc1(nr*nc,&b);
686:   for (i=0; i<nr; i++) {
687:     for (j=0; j<nc; j++) {
688:       if (bA->m[i][j]) {
689:         MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);
690:       } else {
691:         b[i*nc+j] = NULL;
692:       }
693:     }
694:   }
695:   MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);
696:   /* Give the new MatNest exclusive ownership */
697:   for (i=0; i<nr*nc; i++) {
698:     MatDestroy(&b[i]);
699:   }
700:   PetscFree(b);

702:   MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
703:   MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
704:   return(0);
705: }

707: /* nest api */
710: PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
711: {
712:   Mat_Nest *bA = (Mat_Nest*)A->data;

715:   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
716:   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
717:   *mat = bA->m[idxm][jdxm];
718:   return(0);
719: }

723: /*@
724:  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.

726:  Not collective

728:  Input Parameters:
729: +   A  - nest matrix
730: .   idxm - index of the matrix within the nest matrix
731: -   jdxm - index of the matrix within the nest matrix

733:  Output Parameter:
734: .   sub - matrix at index idxm,jdxm within the nest matrix

736:  Level: developer

738: .seealso: MatNestGetSize(), MatNestGetSubMats()
739: @*/
740: PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
741: {

745:   PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));
746:   return(0);
747: }

751: PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat)
752: {
753:   Mat_Nest       *bA = (Mat_Nest*)A->data;
754:   PetscInt       m,n,M,N,mi,ni,Mi,Ni;

758:   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
759:   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
760:   MatGetLocalSize(mat,&m,&n);
761:   MatGetSize(mat,&M,&N);
762:   ISGetLocalSize(bA->isglobal.row[idxm],&mi);
763:   ISGetSize(bA->isglobal.row[idxm],&Mi);
764:   ISGetLocalSize(bA->isglobal.col[jdxm],&ni);
765:   ISGetSize(bA->isglobal.col[jdxm],&Ni);
766:   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);
767:   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);

769:   PetscObjectReference((PetscObject)mat);
770:   MatDestroy(&bA->m[idxm][jdxm]);
771:   bA->m[idxm][jdxm] = mat;
772:   return(0);
773: }

777: /*@
778:  MatNestSetSubMat - Set a single submatrix in the nest matrix.

780:  Logically collective on the submatrix communicator

782:  Input Parameters:
783: +   A  - nest matrix
784: .   idxm - index of the matrix within the nest matrix
785: .   jdxm - index of the matrix within the nest matrix
786: -   sub - matrix at index idxm,jdxm within the nest matrix

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

791:  This increments the reference count of the submatrix.

793:  Level: developer

795: .seealso: MatNestSetSubMats(), MatNestGetSubMat()
796: @*/
797: PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
798: {

802:   PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));
803:   return(0);
804: }

808: PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
809: {
810:   Mat_Nest *bA = (Mat_Nest*)A->data;

813:   if (M)   *M   = bA->nr;
814:   if (N)   *N   = bA->nc;
815:   if (mat) *mat = bA->m;
816:   return(0);
817: }

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

824:  Not collective

826:  Input Parameters:
827: .   A  - nest matrix

829:  Output Parameter:
830: +   M - number of rows in the nest matrix
831: .   N - number of cols in the nest matrix
832: -   mat - 2d array of matrices

834:  Notes:

836:  The user should not free the array mat.

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

842:  Level: developer

844: .seealso: MatNestGetSize(), MatNestGetSubMat()
845: @*/
846: PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
847: {

851:   PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));
852:   return(0);
853: }

857: PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
858: {
859:   Mat_Nest *bA = (Mat_Nest*)A->data;

862:   if (M) *M = bA->nr;
863:   if (N) *N = bA->nc;
864:   return(0);
865: }

869: /*@
870:  MatNestGetSize - Returns the size of the nest matrix.

872:  Not collective

874:  Input Parameters:
875: .   A  - nest matrix

877:  Output Parameter:
878: +   M - number of rows in the nested mat
879: -   N - number of cols in the nested mat

881:  Notes:

883:  Level: developer

885: .seealso: MatNestGetSubMat(), MatNestGetSubMats()
886: @*/
887: PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
888: {

892:   PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));
893:   return(0);
894: }

898: static PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[])
899: {
900:   Mat_Nest *vs = (Mat_Nest*)A->data;
901:   PetscInt i;

904:   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i];
905:   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i];
906:   return(0);
907: }

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

914:  Not collective

916:  Input Parameters:
917: .   A  - nest matrix

919:  Output Parameter:
920: +   rows - array of row index sets
921: -   cols - array of column index sets

923:  Level: advanced

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

928: .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs()
929: @*/
930: PetscErrorCode  MatNestGetISs(Mat A,IS rows[],IS cols[])
931: {

936:   PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));
937:   return(0);
938: }

942: static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[])
943: {
944:   Mat_Nest *vs = (Mat_Nest*)A->data;
945:   PetscInt i;

948:   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i];
949:   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i];
950:   return(0);
951: }

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

958:  Not collective

960:  Input Parameters:
961: .   A  - nest matrix

963:  Output Parameter:
964: +   rows - array of row index sets (or NULL to ignore)
965: -   cols - array of column index sets (or NULL to ignore)

967:  Level: advanced

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

972: .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs()
973: @*/
974: PetscErrorCode  MatNestGetLocalISs(Mat A,IS rows[],IS cols[])
975: {

980:   PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));
981:   return(0);
982: }

986: PetscErrorCode  MatNestSetVecType_Nest(Mat A,VecType vtype)
987: {
989:   PetscBool      flg;

992:   PetscStrcmp(vtype,VECNEST,&flg);
993:   /* In reality, this only distinguishes VECNEST and "other" */
994:   if (flg) A->ops->getvecs = MatCreateVecs_Nest;
995:   else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0;
996:   return(0);
997: }

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

1004:  Not collective

1006:  Input Parameters:
1007: +  A  - nest matrix
1008: -  vtype - type to use for creating vectors

1010:  Notes:

1012:  Level: developer

1014: .seealso: MatCreateVecs()
1015: @*/
1016: PetscErrorCode  MatNestSetVecType(Mat A,VecType vtype)
1017: {

1021:   PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));
1022:   return(0);
1023: }

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

1034:   s->nr = nr;
1035:   s->nc = nc;

1037:   /* Create space for submatrices */
1038:   PetscMalloc1(nr,&s->m);
1039:   for (i=0; i<nr; i++) {
1040:     PetscMalloc1(nc,&s->m[i]);
1041:   }
1042:   for (i=0; i<nr; i++) {
1043:     for (j=0; j<nc; j++) {
1044:       s->m[i][j] = a[i*nc+j];
1045:       if (a[i*nc+j]) {
1046:         PetscObjectReference((PetscObject)a[i*nc+j]);
1047:       }
1048:     }
1049:   }

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

1053:   PetscMalloc1(nr,&s->row_len);
1054:   PetscMalloc1(nc,&s->col_len);
1055:   for (i=0; i<nr; i++) s->row_len[i]=-1;
1056:   for (j=0; j<nc; j++) s->col_len[j]=-1;

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

1060:   PetscLayoutSetSize(A->rmap,M);
1061:   PetscLayoutSetLocalSize(A->rmap,m);
1062:   PetscLayoutSetSize(A->cmap,N);
1063:   PetscLayoutSetLocalSize(A->cmap,n);

1065:   PetscLayoutSetUp(A->rmap);
1066:   PetscLayoutSetUp(A->cmap);

1068:   PetscCalloc2(nr,&s->left,nc,&s->right);
1069:   return(0);
1070: }

1074: /*@
1075:    MatNestSetSubMats - Sets the nested submatrices

1077:    Collective on Mat

1079:    Input Parameter:
1080: +  N - nested matrix
1081: .  nr - number of nested row blocks
1082: .  is_row - index sets for each nested row block, or NULL to make contiguous
1083: .  nc - number of nested column blocks
1084: .  is_col - index sets for each nested column block, or NULL to make contiguous
1085: -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL

1087:    Level: advanced

1089: .seealso: MatCreateNest(), MATNEST
1090: @*/
1091: PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1092: {
1094:   PetscInt       i;

1098:   if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
1099:   if (nr && is_row) {
1102:   }
1103:   if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
1104:   if (nc && is_col) {
1107:   }
1109:   PetscUseMethod(A,"MatNestSetSubMats_C",(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]),(A,nr,is_row,nc,is_col,a));
1110:   return(0);
1111: }

1115: static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog)
1116: {
1118:   PetscBool      flg;
1119:   PetscInt       i,j,m,mi,*ix;

1122:   for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) {
1123:     if (islocal[i]) {
1124:       ISGetSize(islocal[i],&mi);
1125:       flg  = PETSC_TRUE;      /* We found a non-trivial entry */
1126:     } else {
1127:       ISGetSize(isglobal[i],&mi);
1128:     }
1129:     m += mi;
1130:   }
1131:   if (flg) {
1132:     PetscMalloc1(m,&ix);
1133:     for (i=0,n=0; i<n; i++) {
1134:       ISLocalToGlobalMapping smap = NULL;
1135:       VecScatter             scat;
1136:       IS                     isreq;
1137:       Vec                    lvec,gvec;
1138:       union {char padding[sizeof(PetscScalar)]; PetscInt integer;} *x;
1139:       Mat sub;

1141:       if (sizeof(*x) != sizeof(PetscScalar)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support when scalars smaller than integers");
1142:       if (colflg) {
1143:         MatNestFindNonzeroSubMatRow(A,i,&sub);
1144:       } else {
1145:         MatNestFindNonzeroSubMatCol(A,i,&sub);
1146:       }
1147:       if (sub) {MatGetLocalToGlobalMapping(sub,&smap,NULL);}
1148:       if (islocal[i]) {
1149:         ISGetSize(islocal[i],&mi);
1150:       } else {
1151:         ISGetSize(isglobal[i],&mi);
1152:       }
1153:       for (j=0; j<mi; j++) ix[m+j] = j;
1154:       if (smap) {ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);}
1155:       /*
1156:         Now we need to extract the monolithic global indices that correspond to the given split global indices.
1157:         In many/most cases, we only want MatGetLocalSubMatrix() to work, in which case we only need to know the size of the local spaces.
1158:         The approach here is ugly because it uses VecScatter to move indices.
1159:        */
1160:       VecCreateSeq(PETSC_COMM_SELF,mi,&lvec);
1161:       VecCreateMPI(((PetscObject)isglobal[i])->comm,mi,PETSC_DECIDE,&gvec);
1162:       ISCreateGeneral(((PetscObject)isglobal[i])->comm,mi,ix+m,PETSC_COPY_VALUES,&isreq);
1163:       VecScatterCreate(gvec,isreq,lvec,NULL,&scat);
1164:       VecGetArray(gvec,(PetscScalar**)&x);
1165:       for (j=0; j<mi; j++) x[j].integer = ix[m+j];
1166:       VecRestoreArray(gvec,(PetscScalar**)&x);
1167:       VecScatterBegin(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);
1168:       VecScatterEnd(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);
1169:       VecGetArray(lvec,(PetscScalar**)&x);
1170:       for (j=0; j<mi; j++) ix[m+j] = x[j].integer;
1171:       VecRestoreArray(lvec,(PetscScalar**)&x);
1172:       VecDestroy(&lvec);
1173:       VecDestroy(&gvec);
1174:       ISDestroy(&isreq);
1175:       VecScatterDestroy(&scat);
1176:       m   += mi;
1177:     }
1178:     ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),1,m,ix,PETSC_OWN_POINTER,ltog);
1179:   } else {
1180:     *ltog  = NULL;
1181:   }
1182:   return(0);
1183: }


1186: /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1187: /*
1188:   nprocessors = NP
1189:   Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1))
1190:        proc 0: => (g_0,h_0,)
1191:        proc 1: => (g_1,h_1,)
1192:        ...
1193:        proc nprocs-1: => (g_NP-1,h_NP-1,)

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

1198:             proc 0:
1199:     is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1)
1200:             proc 1:
1201:     is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1)

1203:             proc NP-1:
1204:     is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1)
1205: */
1208: static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1209: {
1210:   Mat_Nest       *vs = (Mat_Nest*)A->data;
1211:   PetscInt       i,j,offset,n,nsum,bs;
1213:   Mat            sub = NULL;

1216:   PetscMalloc1(nr,&vs->isglobal.row);
1217:   PetscMalloc1(nc,&vs->isglobal.col);
1218:   if (is_row) { /* valid IS is passed in */
1219:     /* refs on is[] are incremeneted */
1220:     for (i=0; i<vs->nr; i++) {
1221:       PetscObjectReference((PetscObject)is_row[i]);

1223:       vs->isglobal.row[i] = is_row[i];
1224:     }
1225:   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
1226:     nsum = 0;
1227:     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
1228:       MatNestFindNonzeroSubMatRow(A,i,&sub);
1229:       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
1230:       MatGetLocalSize(sub,&n,NULL);
1231:       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
1232:       nsum += n;
1233:     }
1234:     MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));
1235:     offset -= nsum;
1236:     for (i=0; i<vs->nr; i++) {
1237:       MatNestFindNonzeroSubMatRow(A,i,&sub);
1238:       MatGetLocalSize(sub,&n,NULL);
1239:       MatGetBlockSize(sub,&bs);
1240:       ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);
1241:       ISSetBlockSize(vs->isglobal.row[i],bs);
1242:       offset += n;
1243:     }
1244:   }

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

1251:       vs->isglobal.col[j] = is_col[j];
1252:     }
1253:   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
1254:     offset = A->cmap->rstart;
1255:     nsum   = 0;
1256:     for (j=0; j<vs->nc; j++) {
1257:       MatNestFindNonzeroSubMatCol(A,j,&sub);
1258:       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
1259:       MatGetLocalSize(sub,NULL,&n);
1260:       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
1261:       nsum += n;
1262:     }
1263:     MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));
1264:     offset -= nsum;
1265:     for (j=0; j<vs->nc; j++) {
1266:       MatNestFindNonzeroSubMatCol(A,j,&sub);
1267:       MatGetLocalSize(sub,NULL,&n);
1268:       MatGetBlockSize(sub,&bs);
1269:       ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);
1270:       ISSetBlockSize(vs->isglobal.col[j],bs);
1271:       offset += n;
1272:     }
1273:   }

1275:   /* Set up the local ISs */
1276:   PetscMalloc1(vs->nr,&vs->islocal.row);
1277:   PetscMalloc1(vs->nc,&vs->islocal.col);
1278:   for (i=0,offset=0; i<vs->nr; i++) {
1279:     IS                     isloc;
1280:     ISLocalToGlobalMapping rmap = NULL;
1281:     PetscInt               nlocal,bs;
1282:     MatNestFindNonzeroSubMatRow(A,i,&sub);
1283:     if (sub) {MatGetLocalToGlobalMapping(sub,&rmap,NULL);}
1284:     if (rmap) {
1285:       MatGetBlockSize(sub,&bs);
1286:       ISLocalToGlobalMappingGetSize(rmap,&nlocal);
1287:       ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);
1288:       ISSetBlockSize(isloc,bs);
1289:     } else {
1290:       nlocal = 0;
1291:       isloc  = NULL;
1292:     }
1293:     vs->islocal.row[i] = isloc;
1294:     offset            += nlocal;
1295:   }
1296:   for (i=0,offset=0; i<vs->nc; i++) {
1297:     IS                     isloc;
1298:     ISLocalToGlobalMapping cmap = NULL;
1299:     PetscInt               nlocal,bs;
1300:     MatNestFindNonzeroSubMatCol(A,i,&sub);
1301:     if (sub) {MatGetLocalToGlobalMapping(sub,NULL,&cmap);}
1302:     if (cmap) {
1303:       MatGetBlockSize(sub,&bs);
1304:       ISLocalToGlobalMappingGetSize(cmap,&nlocal);
1305:       ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);
1306:       ISSetBlockSize(isloc,bs);
1307:     } else {
1308:       nlocal = 0;
1309:       isloc  = NULL;
1310:     }
1311:     vs->islocal.col[i] = isloc;
1312:     offset            += nlocal;
1313:   }

1315:   /* Set up the aggregate ISLocalToGlobalMapping */
1316:   {
1317:     ISLocalToGlobalMapping rmap,cmap;
1318:     MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap);
1319:     MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap);
1320:     if (rmap && cmap) {MatSetLocalToGlobalMapping(A,rmap,cmap);}
1321:     ISLocalToGlobalMappingDestroy(&rmap);
1322:     ISLocalToGlobalMappingDestroy(&cmap);
1323:   }

1325: #if defined(PETSC_USE_DEBUG)
1326:   for (i=0; i<vs->nr; i++) {
1327:     for (j=0; j<vs->nc; j++) {
1328:       PetscInt m,n,M,N,mi,ni,Mi,Ni;
1329:       Mat      B = vs->m[i][j];
1330:       if (!B) continue;
1331:       MatGetSize(B,&M,&N);
1332:       MatGetLocalSize(B,&m,&n);
1333:       ISGetSize(vs->isglobal.row[i],&Mi);
1334:       ISGetSize(vs->isglobal.col[j],&Ni);
1335:       ISGetLocalSize(vs->isglobal.row[i],&mi);
1336:       ISGetLocalSize(vs->isglobal.col[j],&ni);
1337:       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);
1338:       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);
1339:     }
1340:   }
1341: #endif

1343:   /* Set A->assembled if all non-null blocks are currently assembled */
1344:   for (i=0; i<vs->nr; i++) {
1345:     for (j=0; j<vs->nc; j++) {
1346:       if (vs->m[i][j] && !vs->m[i][j]->assembled) return(0);
1347:     }
1348:   }
1349:   A->assembled = PETSC_TRUE;
1350:   return(0);
1351: }

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

1358:    Collective on Mat

1360:    Input Parameter:
1361: +  comm - Communicator for the new Mat
1362: .  nr - number of nested row blocks
1363: .  is_row - index sets for each nested row block, or NULL to make contiguous
1364: .  nc - number of nested column blocks
1365: .  is_col - index sets for each nested column block, or NULL to make contiguous
1366: -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL

1368:    Output Parameter:
1369: .  B - new matrix

1371:    Level: advanced

1373: .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST
1374: @*/
1375: PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1376: {
1377:   Mat            A;

1381:   *B   = 0;
1382:   MatCreate(comm,&A);
1383:   MatSetType(A,MATNEST);
1384:   A->preallocated = PETSC_TRUE;
1385:   MatNestSetSubMats(A,nr,is_row,nc,is_col,a);
1386:   *B   = A;
1387:   return(0);
1388: }

1392: PETSC_INTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1393: {
1395:   Mat_Nest       *nest = (Mat_Nest*)A->data;
1396:   PetscInt       m,n,M,N,i,j,k,*dnnz,*onnz,rstart;
1397:   PetscInt       cstart,cend;
1398:   Mat            C;

1401:   MatGetSize(A,&M,&N);
1402:   MatGetLocalSize(A,&m,&n);
1403:   MatGetOwnershipRangeColumn(A,&cstart,&cend);
1404:   switch (reuse) {
1405:   case MAT_INITIAL_MATRIX:
1406:     MatCreate(PetscObjectComm((PetscObject)A),&C);
1407:     MatSetType(C,newtype);
1408:     MatSetSizes(C,m,n,M,N);
1409:     *newmat = C;
1410:     break;
1411:   case MAT_REUSE_MATRIX:
1412:     C = *newmat;
1413:     break;
1414:   default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse");
1415:   }
1416:   PetscMalloc1(2*m,&dnnz);
1417:   onnz = dnnz + m;
1418:   for (k=0; k<m; k++) {
1419:     dnnz[k] = 0;
1420:     onnz[k] = 0;
1421:   }
1422:   for (j=0; j<nest->nc; ++j) {
1423:     IS             bNis;
1424:     PetscInt       bN;
1425:     const PetscInt *bNindices;
1426:     /* Using global column indices and ISAllGather() is not scalable. */
1427:     ISAllGather(nest->isglobal.col[j], &bNis);
1428:     ISGetSize(bNis, &bN);
1429:     ISGetIndices(bNis,&bNindices);
1430:     for (i=0; i<nest->nr; ++i) {
1431:       PetscSF        bmsf;
1432:       PetscSFNode    *iremote;
1433:       Mat            B;
1434:       PetscInt       bm, *sub_dnnz,*sub_onnz, br;
1435:       const PetscInt *bmindices;
1436:       B = nest->m[i][j];
1437:       if (!B) continue;
1438:       ISGetLocalSize(nest->isglobal.row[i],&bm);
1439:       ISGetIndices(nest->isglobal.row[i],&bmindices);
1440:       PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);
1441:       PetscMalloc1(bm,&iremote);
1442:       PetscMalloc1(bm,&sub_dnnz);
1443:       PetscMalloc1(bm,&sub_onnz);
1444:       for (k = 0; k < bm; ++k){
1445:             sub_dnnz[k] = 0;
1446:             sub_onnz[k] = 0;
1447:       }
1448:       /*
1449:        Locate the owners for all of the locally-owned global row indices for this row block.
1450:        These determine the roots of PetscSF used to communicate preallocation data to row owners.
1451:        The roots correspond to the dnnz and onnz entries; thus, there are two roots per row.
1452:        */
1453:       MatGetOwnershipRange(B,&rstart,NULL);
1454:       for (br = 0; br < bm; ++br) {
1455:         PetscInt       row = bmindices[br], rowowner = 0, brncols, col;
1456:         const PetscInt *brcols;
1457:         PetscInt       rowrel = 0; /* row's relative index on its owner rank */
1458:         PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);
1459:         /* how many roots  */
1460:         iremote[br].rank = rowowner; iremote[br].index = rowrel;           /* edge from bmdnnz to dnnz */
1461:         /* get nonzero pattern */
1462:         MatGetRow(B,br+rstart,&brncols,&brcols,NULL);
1463:         for (k=0; k<brncols; k++) {
1464:           col  = bNindices[brcols[k]];
1465:           if(col>=A->cmap->range[rowowner] && col<A->cmap->range[rowowner+1]){
1466:                 sub_dnnz[br]++;
1467:           }else{
1468:                 sub_onnz[br]++;
1469:           }
1470:         }
1471:         MatRestoreRow(B,br+rstart,&brncols,&brcols,NULL);
1472:       }
1473:       ISRestoreIndices(nest->isglobal.row[i],&bmindices);
1474:       /* bsf will have to take care of disposing of bedges. */
1475:       PetscSFSetGraph(bmsf,m,bm,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);
1476:       PetscSFReduceBegin(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);
1477:       PetscSFReduceEnd(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);
1478:       PetscSFReduceBegin(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);
1479:       PetscSFReduceEnd(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);
1480:       PetscFree(sub_dnnz);
1481:       PetscFree(sub_onnz);
1482:       PetscSFDestroy(&bmsf);
1483:     }
1484:     ISRestoreIndices(bNis,&bNindices);
1485:     ISDestroy(&bNis);
1486:   }
1487:   MatSeqAIJSetPreallocation(C,0,dnnz);
1488:   MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);
1489:   PetscFree(dnnz);

1491:   /* Fill by row */
1492:   for (j=0; j<nest->nc; ++j) {
1493:     /* Using global column indices and ISAllGather() is not scalable. */
1494:     IS             bNis;
1495:     PetscInt       bN;
1496:     const PetscInt *bNindices;
1497:     ISAllGather(nest->isglobal.col[j], &bNis);
1498:     ISGetSize(bNis,&bN);
1499:     ISGetIndices(bNis,&bNindices);
1500:     for (i=0; i<nest->nr; ++i) {
1501:       Mat            B;
1502:       PetscInt       bm, br;
1503:       const PetscInt *bmindices;
1504:       B = nest->m[i][j];
1505:       if (!B) continue;
1506:       ISGetLocalSize(nest->isglobal.row[i],&bm);
1507:       ISGetIndices(nest->isglobal.row[i],&bmindices);
1508:       MatGetOwnershipRange(B,&rstart,NULL);
1509:       for (br = 0; br < bm; ++br) {
1510:         PetscInt          row = bmindices[br], brncols,  *cols;
1511:         const PetscInt    *brcols;
1512:         const PetscScalar *brcoldata;
1513:         MatGetRow(B,br+rstart,&brncols,&brcols,&brcoldata);
1514:         PetscMalloc1(brncols,&cols);
1515:         for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]];
1516:         /*
1517:           Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match.
1518:           Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES.
1519:          */
1520:         MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);
1521:         MatRestoreRow(B,br+rstart,&brncols,&brcols,&brcoldata);
1522:         PetscFree(cols);
1523:       }
1524:       ISRestoreIndices(nest->isglobal.row[i],&bmindices);
1525:     }
1526:     ISRestoreIndices(bNis,&bNindices);
1527:     ISDestroy(&bNis);
1528:   }
1529:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1530:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1531:   return(0);
1532: }

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

1537:   Level: intermediate

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

1544: .seealso: MatCreate(), MatType, MatCreateNest()
1545: M*/
1548: PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A)
1549: {
1550:   Mat_Nest       *s;

1554:   PetscNewLog(A,&s);
1555:   A->data = (void*)s;

1557:   s->nr            = -1;
1558:   s->nc            = -1;
1559:   s->m             = NULL;
1560:   s->splitassembly = PETSC_FALSE;

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

1564:   A->ops->mult                  = MatMult_Nest;
1565:   A->ops->multadd               = MatMultAdd_Nest;
1566:   A->ops->multtranspose         = MatMultTranspose_Nest;
1567:   A->ops->multtransposeadd      = MatMultTransposeAdd_Nest;
1568:   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1569:   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1570:   A->ops->zeroentries           = MatZeroEntries_Nest;
1571:   A->ops->copy                  = MatCopy_Nest;
1572:   A->ops->duplicate             = MatDuplicate_Nest;
1573:   A->ops->getsubmatrix          = MatGetSubMatrix_Nest;
1574:   A->ops->destroy               = MatDestroy_Nest;
1575:   A->ops->view                  = MatView_Nest;
1576:   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1577:   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1578:   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
1579:   A->ops->getdiagonal           = MatGetDiagonal_Nest;
1580:   A->ops->diagonalscale         = MatDiagonalScale_Nest;
1581:   A->ops->scale                 = MatScale_Nest;
1582:   A->ops->shift                 = MatShift_Nest;

1584:   A->spptr        = 0;
1585:   A->assembled    = PETSC_FALSE;

1587:   /* expose Nest api's */
1588:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",   MatNestGetSubMat_Nest);
1589:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",   MatNestSetSubMat_Nest);
1590:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",  MatNestGetSubMats_Nest);
1591:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",     MatNestGetSize_Nest);
1592:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",      MatNestGetISs_Nest);
1593:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C", MatNestGetLocalISs_Nest);
1594:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",  MatNestSetVecType_Nest);
1595:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",  MatNestSetSubMats_Nest);
1596:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",MatConvert_Nest_AIJ);

1598:   PetscObjectChangeTypeName((PetscObject)A,MATNEST);
1599:   return(0);
1600: }