Actual source code: matnest.c

petsc-dev 2014-04-15
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 MatGetVecs_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:   MPI_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:   }
399:   return(0);
400: }

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

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

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

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

438: static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
439: {
440:   Mat_Nest       *bA = (Mat_Nest*)A->data;
441:   PetscInt       i;

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

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

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

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

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

514: static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a)
515: {
516:   Mat_Nest       *bA = (Mat_Nest*)A->data;
517:   PetscInt       i;

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

530: static PetscErrorCode MatGetVecs_Nest(Mat A,Vec *right,Vec *left)
531: {
532:   Mat_Nest       *bA = (Mat_Nest*)A->data;
533:   Vec            *L,*R;
534:   MPI_Comm       comm;
535:   PetscInt       i,j;

539:   PetscObjectGetComm((PetscObject)A,&comm);
540:   if (right) {
541:     /* allocate R */
542:     PetscMalloc(sizeof(Vec) * bA->nc, &R);
543:     /* Create the right vectors */
544:     for (j=0; j<bA->nc; j++) {
545:       for (i=0; i<bA->nr; i++) {
546:         if (bA->m[i][j]) {
547:           MatGetVecs(bA->m[i][j],&R[j],NULL);
548:           break;
549:         }
550:       }
551:       if (i==bA->nr) {
552:         /* have an empty column */
553:         SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
554:       }
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:     PetscMalloc(sizeof(Vec) * 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:           MatGetVecs(bA->m[i][j],NULL,&L[i]);
572:           break;
573:         }
574:       }
575:       if (j==bA->nc) {
576:         /* have an empty row */
577:         SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
578:       }
579:     }

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

586:     PetscFree(L);
587:   }
588:   return(0);
589: }

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

601:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
602:   if (isascii) {

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

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

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

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

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

642: static PetscErrorCode MatZeroEntries_Nest(Mat A)
643: {
644:   Mat_Nest       *bA = (Mat_Nest*)A->data;
645:   PetscInt       i,j;

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

660: static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
661: {
662:   Mat_Nest       *bA = (Mat_Nest*)A->data;
663:   Mat            *b;
664:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;

668:   PetscMalloc1(nr*nc,&b);
669:   for (i=0; i<nr; i++) {
670:     for (j=0; j<nc; j++) {
671:       if (bA->m[i][j]) {
672:         MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);
673:       } else {
674:         b[i*nc+j] = NULL;
675:       }
676:     }
677:   }
678:   MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);
679:   /* Give the new MatNest exclusive ownership */
680:   for (i=0; i<nr*nc; i++) {
681:     MatDestroy(&b[i]);
682:   }
683:   PetscFree(b);

685:   MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
686:   MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
687:   return(0);
688: }

690: /* nest api */
693: PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
694: {
695:   Mat_Nest *bA = (Mat_Nest*)A->data;

698:   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
699:   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
700:   *mat = bA->m[idxm][jdxm];
701:   return(0);
702: }

706: /*@
707:  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.

709:  Not collective

711:  Input Parameters:
712: +   A  - nest matrix
713: .   idxm - index of the matrix within the nest matrix
714: -   jdxm - index of the matrix within the nest matrix

716:  Output Parameter:
717: .   sub - matrix at index idxm,jdxm within the nest matrix

719:  Level: developer

721: .seealso: MatNestGetSize(), MatNestGetSubMats()
722: @*/
723: PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
724: {

728:   PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));
729:   return(0);
730: }

734: PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat)
735: {
736:   Mat_Nest       *bA = (Mat_Nest*)A->data;
737:   PetscInt       m,n,M,N,mi,ni,Mi,Ni;

741:   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
742:   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
743:   MatGetLocalSize(mat,&m,&n);
744:   MatGetSize(mat,&M,&N);
745:   ISGetLocalSize(bA->isglobal.row[idxm],&mi);
746:   ISGetSize(bA->isglobal.row[idxm],&Mi);
747:   ISGetLocalSize(bA->isglobal.col[jdxm],&ni);
748:   ISGetSize(bA->isglobal.col[jdxm],&Ni);
749:   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);
750:   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);

752:   PetscObjectReference((PetscObject)mat);
753:   MatDestroy(&bA->m[idxm][jdxm]);
754:   bA->m[idxm][jdxm] = mat;
755:   return(0);
756: }

760: /*@
761:  MatNestSetSubMat - Set a single submatrix in the nest matrix.

763:  Logically collective on the submatrix communicator

765:  Input Parameters:
766: +   A  - nest matrix
767: .   idxm - index of the matrix within the nest matrix
768: .   jdxm - index of the matrix within the nest matrix
769: -   sub - matrix at index idxm,jdxm within the nest matrix

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

774:  This increments the reference count of the submatrix.

776:  Level: developer

778: .seealso: MatNestSetSubMats(), MatNestGetSubMat()
779: @*/
780: PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
781: {

785:   PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));
786:   return(0);
787: }

791: PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
792: {
793:   Mat_Nest *bA = (Mat_Nest*)A->data;

796:   if (M)   *M   = bA->nr;
797:   if (N)   *N   = bA->nc;
798:   if (mat) *mat = bA->m;
799:   return(0);
800: }

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

807:  Not collective

809:  Input Parameters:
810: .   A  - nest matrix

812:  Output Parameter:
813: +   M - number of rows in the nest matrix
814: .   N - number of cols in the nest matrix
815: -   mat - 2d array of matrices

817:  Notes:

819:  The user should not free the array mat.

821:  Level: developer

823: .seealso: MatNestGetSize(), MatNestGetSubMat()
824: @*/
825: PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
826: {

830:   PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));
831:   return(0);
832: }

836: PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
837: {
838:   Mat_Nest *bA = (Mat_Nest*)A->data;

841:   if (M) *M = bA->nr;
842:   if (N) *N = bA->nc;
843:   return(0);
844: }

848: /*@
849:  MatNestGetSize - Returns the size of the nest matrix.

851:  Not collective

853:  Input Parameters:
854: .   A  - nest matrix

856:  Output Parameter:
857: +   M - number of rows in the nested mat
858: -   N - number of cols in the nested mat

860:  Notes:

862:  Level: developer

864: .seealso: MatNestGetSubMat(), MatNestGetSubMats()
865: @*/
866: PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
867: {

871:   PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));
872:   return(0);
873: }

877: static PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[])
878: {
879:   Mat_Nest *vs = (Mat_Nest*)A->data;
880:   PetscInt i;

883:   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i];
884:   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i];
885:   return(0);
886: }

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

893:  Not collective

895:  Input Parameters:
896: .   A  - nest matrix

898:  Output Parameter:
899: +   rows - array of row index sets
900: -   cols - array of column index sets

902:  Level: advanced

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

907: .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs()
908: @*/
909: PetscErrorCode  MatNestGetISs(Mat A,IS rows[],IS cols[])
910: {

915:   PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));
916:   return(0);
917: }

921: static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[])
922: {
923:   Mat_Nest *vs = (Mat_Nest*)A->data;
924:   PetscInt i;

927:   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i];
928:   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i];
929:   return(0);
930: }

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

937:  Not collective

939:  Input Parameters:
940: .   A  - nest matrix

942:  Output Parameter:
943: +   rows - array of row index sets (or NULL to ignore)
944: -   cols - array of column index sets (or NULL to ignore)

946:  Level: advanced

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

951: .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs()
952: @*/
953: PetscErrorCode  MatNestGetLocalISs(Mat A,IS rows[],IS cols[])
954: {

959:   PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));
960:   return(0);
961: }

965: PetscErrorCode  MatNestSetVecType_Nest(Mat A,VecType vtype)
966: {
968:   PetscBool      flg;

971:   PetscStrcmp(vtype,VECNEST,&flg);
972:   /* In reality, this only distinguishes VECNEST and "other" */
973:   if (flg) A->ops->getvecs = MatGetVecs_Nest;
974:   else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0;
975:   return(0);
976: }

980: /*@C
981:  MatNestSetVecType - Sets the type of Vec returned by MatGetVecs()

983:  Not collective

985:  Input Parameters:
986: +  A  - nest matrix
987: -  vtype - type to use for creating vectors

989:  Notes:

991:  Level: developer

993: .seealso: MatGetVecs()
994: @*/
995: PetscErrorCode  MatNestSetVecType(Mat A,VecType vtype)
996: {

1000:   PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));
1001:   return(0);
1002: }

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

1013:   s->nr = nr;
1014:   s->nc = nc;

1016:   /* Create space for submatrices */
1017:   PetscMalloc(sizeof(Mat*)*nr,&s->m);
1018:   for (i=0; i<nr; i++) {
1019:     PetscMalloc(sizeof(Mat)*nc,&s->m[i]);
1020:   }
1021:   for (i=0; i<nr; i++) {
1022:     for (j=0; j<nc; j++) {
1023:       s->m[i][j] = a[i*nc+j];
1024:       if (a[i*nc+j]) {
1025:         PetscObjectReference((PetscObject)a[i*nc+j]);
1026:       }
1027:     }
1028:   }

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

1032:   PetscMalloc(sizeof(PetscInt)*nr,&s->row_len);
1033:   PetscMalloc(sizeof(PetscInt)*nc,&s->col_len);
1034:   for (i=0; i<nr; i++) s->row_len[i]=-1;
1035:   for (j=0; j<nc; j++) s->col_len[j]=-1;

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

1039:   PetscLayoutSetSize(A->rmap,M);
1040:   PetscLayoutSetLocalSize(A->rmap,m);
1041:   PetscLayoutSetSize(A->cmap,N);
1042:   PetscLayoutSetLocalSize(A->cmap,n);

1044:   PetscLayoutSetUp(A->rmap);
1045:   PetscLayoutSetUp(A->cmap);

1047:   PetscCalloc2(nr,&s->left,nc,&s->right);
1048:   return(0);
1049: }

1053: /*@
1054:    MatNestSetSubMats - Sets the nested submatrices

1056:    Collective on Mat

1058:    Input Parameter:
1059: +  N - nested matrix
1060: .  nr - number of nested row blocks
1061: .  is_row - index sets for each nested row block, or NULL to make contiguous
1062: .  nc - number of nested column blocks
1063: .  is_col - index sets for each nested column block, or NULL to make contiguous
1064: -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL

1066:    Level: advanced

1068: .seealso: MatCreateNest(), MATNEST
1069: @*/
1070: PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1071: {
1073:   PetscInt       i;

1077:   if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
1078:   if (nr && is_row) {
1081:   }
1082:   if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
1083:   if (nc && is_col) {
1086:   }
1088:   PetscUseMethod(A,"MatNestSetSubMats_C",(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]),(A,nr,is_row,nc,is_col,a));
1089:   return(0);
1090: }

1094: static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog,ISLocalToGlobalMapping *ltogb)
1095: {
1097:   PetscBool      flg;
1098:   PetscInt       i,j,m,mi,*ix;

1101:   for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) {
1102:     if (islocal[i]) {
1103:       ISGetSize(islocal[i],&mi);
1104:       flg  = PETSC_TRUE;      /* We found a non-trivial entry */
1105:     } else {
1106:       ISGetSize(isglobal[i],&mi);
1107:     }
1108:     m += mi;
1109:   }
1110:   if (flg) {
1111:     PetscMalloc1(m,&ix);
1112:     for (i=0,n=0; i<n; i++) {
1113:       ISLocalToGlobalMapping smap = NULL;
1114:       VecScatter             scat;
1115:       IS                     isreq;
1116:       Vec                    lvec,gvec;
1117:       union {char padding[sizeof(PetscScalar)]; PetscInt integer;} *x;
1118:       Mat sub;

1120:       if (sizeof(*x) != sizeof(PetscScalar)) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"No support when scalars smaller than integers");
1121:       if (colflg) {
1122:         MatNestFindNonzeroSubMatRow(A,i,&sub);
1123:       } else {
1124:         MatNestFindNonzeroSubMatCol(A,i,&sub);
1125:       }
1126:       if (sub) {MatGetLocalToGlobalMapping(sub,&smap,NULL);}
1127:       if (islocal[i]) {
1128:         ISGetSize(islocal[i],&mi);
1129:       } else {
1130:         ISGetSize(isglobal[i],&mi);
1131:       }
1132:       for (j=0; j<mi; j++) ix[m+j] = j;
1133:       if (smap) {ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);}
1134:       /*
1135:         Now we need to extract the monolithic global indices that correspond to the given split global indices.
1136:         In many/most cases, we only want MatGetLocalSubMatrix() to work, in which case we only need to know the size of the local spaces.
1137:         The approach here is ugly because it uses VecScatter to move indices.
1138:        */
1139:       VecCreateSeq(PETSC_COMM_SELF,mi,&lvec);
1140:       VecCreateMPI(((PetscObject)isglobal[i])->comm,mi,PETSC_DECIDE,&gvec);
1141:       ISCreateGeneral(((PetscObject)isglobal[i])->comm,mi,ix+m,PETSC_COPY_VALUES,&isreq);
1142:       VecScatterCreate(gvec,isreq,lvec,NULL,&scat);
1143:       VecGetArray(gvec,(PetscScalar**)&x);
1144:       for (j=0; j<mi; j++) x[j].integer = ix[m+j];
1145:       VecRestoreArray(gvec,(PetscScalar**)&x);
1146:       VecScatterBegin(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);
1147:       VecScatterEnd(scat,gvec,lvec,INSERT_VALUES,SCATTER_FORWARD);
1148:       VecGetArray(lvec,(PetscScalar**)&x);
1149:       for (j=0; j<mi; j++) ix[m+j] = x[j].integer;
1150:       VecRestoreArray(lvec,(PetscScalar**)&x);
1151:       VecDestroy(&lvec);
1152:       VecDestroy(&gvec);
1153:       ISDestroy(&isreq);
1154:       VecScatterDestroy(&scat);
1155:       m   += mi;
1156:     }
1157:     ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),m,ix,PETSC_OWN_POINTER,ltog);
1158:     *ltogb = NULL;
1159:   } else {
1160:     *ltog  = NULL;
1161:     *ltogb = NULL;
1162:   }
1163:   return(0);
1164: }


1167: /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1168: /*
1169:   nprocessors = NP
1170:   Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1))
1171:        proc 0: => (g_0,h_0,)
1172:        proc 1: => (g_1,h_1,)
1173:        ...
1174:        proc nprocs-1: => (g_NP-1,h_NP-1,)

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

1179:             proc 0:
1180:     is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1)
1181:             proc 1:
1182:     is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1)

1184:             proc NP-1:
1185:     is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1)
1186: */
1189: static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1190: {
1191:   Mat_Nest       *vs = (Mat_Nest*)A->data;
1192:   PetscInt       i,j,offset,n,nsum,bs;
1194:   Mat            sub = NULL;

1197:   PetscMalloc(sizeof(IS)*nr,&vs->isglobal.row);
1198:   PetscMalloc(sizeof(IS)*nc,&vs->isglobal.col);
1199:   if (is_row) { /* valid IS is passed in */
1200:     /* refs on is[] are incremeneted */
1201:     for (i=0; i<vs->nr; i++) {
1202:       PetscObjectReference((PetscObject)is_row[i]);

1204:       vs->isglobal.row[i] = is_row[i];
1205:     }
1206:   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
1207:     nsum = 0;
1208:     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
1209:       MatNestFindNonzeroSubMatRow(A,i,&sub);
1210:       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
1211:       MatGetLocalSize(sub,&n,NULL);
1212:       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
1213:       nsum += n;
1214:     }
1215:     MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));
1216:     offset -= nsum;
1217:     for (i=0; i<vs->nr; i++) {
1218:       MatNestFindNonzeroSubMatRow(A,i,&sub);
1219:       MatGetLocalSize(sub,&n,NULL);
1220:       MatGetBlockSize(sub,&bs);
1221:       ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);
1222:       ISSetBlockSize(vs->isglobal.row[i],bs);
1223:       offset += n;
1224:     }
1225:   }

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

1232:       vs->isglobal.col[j] = is_col[j];
1233:     }
1234:   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
1235:     offset = A->cmap->rstart;
1236:     nsum   = 0;
1237:     for (j=0; j<vs->nc; j++) {
1238:       MatNestFindNonzeroSubMatCol(A,j,&sub);
1239:       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
1240:       MatGetLocalSize(sub,NULL,&n);
1241:       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
1242:       nsum += n;
1243:     }
1244:     MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));
1245:     offset -= nsum;
1246:     for (j=0; j<vs->nc; j++) {
1247:       MatNestFindNonzeroSubMatCol(A,j,&sub);
1248:       MatGetLocalSize(sub,NULL,&n);
1249:       MatGetBlockSize(sub,&bs);
1250:       ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);
1251:       ISSetBlockSize(vs->isglobal.col[j],bs);
1252:       offset += n;
1253:     }
1254:   }

1256:   /* Set up the local ISs */
1257:   PetscMalloc1(vs->nr,&vs->islocal.row);
1258:   PetscMalloc1(vs->nc,&vs->islocal.col);
1259:   for (i=0,offset=0; i<vs->nr; i++) {
1260:     IS                     isloc;
1261:     ISLocalToGlobalMapping rmap = NULL;
1262:     PetscInt               nlocal,bs;
1263:     MatNestFindNonzeroSubMatRow(A,i,&sub);
1264:     if (sub) {MatGetLocalToGlobalMapping(sub,&rmap,NULL);}
1265:     if (rmap) {
1266:       MatGetBlockSize(sub,&bs);
1267:       ISLocalToGlobalMappingGetSize(rmap,&nlocal);
1268:       ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);
1269:       ISSetBlockSize(isloc,bs);
1270:     } else {
1271:       nlocal = 0;
1272:       isloc  = NULL;
1273:     }
1274:     vs->islocal.row[i] = isloc;
1275:     offset            += nlocal;
1276:   }
1277:   for (i=0,offset=0; i<vs->nc; i++) {
1278:     IS                     isloc;
1279:     ISLocalToGlobalMapping cmap = NULL;
1280:     PetscInt               nlocal,bs;
1281:     MatNestFindNonzeroSubMatCol(A,i,&sub);
1282:     if (sub) {MatGetLocalToGlobalMapping(sub,NULL,&cmap);}
1283:     if (cmap) {
1284:       MatGetBlockSize(sub,&bs);
1285:       ISLocalToGlobalMappingGetSize(cmap,&nlocal);
1286:       ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);
1287:       ISSetBlockSize(isloc,bs);
1288:     } else {
1289:       nlocal = 0;
1290:       isloc  = NULL;
1291:     }
1292:     vs->islocal.col[i] = isloc;
1293:     offset            += nlocal;
1294:   }

1296:   /* Set up the aggregate ISLocalToGlobalMapping */
1297:   {
1298:     ISLocalToGlobalMapping rmap,rmapb,cmap,cmapb;
1299:     MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap,&rmapb);
1300:     MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap,&cmapb);
1301:     if (rmap && cmap) {MatSetLocalToGlobalMapping(A,rmap,cmap);}
1302:     if (rmapb && cmapb) {MatSetLocalToGlobalMappingBlock(A,rmapb,cmapb);}
1303:     ISLocalToGlobalMappingDestroy(&rmap);
1304:     ISLocalToGlobalMappingDestroy(&rmapb);
1305:     ISLocalToGlobalMappingDestroy(&cmap);
1306:     ISLocalToGlobalMappingDestroy(&cmapb);
1307:   }

1309: #if defined(PETSC_USE_DEBUG)
1310:   for (i=0; i<vs->nr; i++) {
1311:     for (j=0; j<vs->nc; j++) {
1312:       PetscInt m,n,M,N,mi,ni,Mi,Ni;
1313:       Mat      B = vs->m[i][j];
1314:       if (!B) continue;
1315:       MatGetSize(B,&M,&N);
1316:       MatGetLocalSize(B,&m,&n);
1317:       ISGetSize(vs->isglobal.row[i],&Mi);
1318:       ISGetSize(vs->isglobal.col[j],&Ni);
1319:       ISGetLocalSize(vs->isglobal.row[i],&mi);
1320:       ISGetLocalSize(vs->isglobal.col[j],&ni);
1321:       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);
1322:       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);
1323:     }
1324:   }
1325: #endif

1327:   /* Set A->assembled if all non-null blocks are currently assembled */
1328:   for (i=0; i<vs->nr; i++) {
1329:     for (j=0; j<vs->nc; j++) {
1330:       if (vs->m[i][j] && !vs->m[i][j]->assembled) return(0);
1331:     }
1332:   }
1333:   A->assembled = PETSC_TRUE;
1334:   return(0);
1335: }

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

1342:    Collective on Mat

1344:    Input Parameter:
1345: +  comm - Communicator for the new Mat
1346: .  nr - number of nested row blocks
1347: .  is_row - index sets for each nested row block, or NULL to make contiguous
1348: .  nc - number of nested column blocks
1349: .  is_col - index sets for each nested column block, or NULL to make contiguous
1350: -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL

1352:    Output Parameter:
1353: .  B - new matrix

1355:    Level: advanced

1357: .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST
1358: @*/
1359: PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1360: {
1361:   Mat            A;

1365:   *B   = 0;
1366:   MatCreate(comm,&A);
1367:   MatSetType(A,MATNEST);
1368:   MatSetUp(A);
1369:   MatNestSetSubMats(A,nr,is_row,nc,is_col,a);
1370:   *B   = A;
1371:   return(0);
1372: }

1376: PETSC_EXTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1377: {
1379:   Mat_Nest       *nest = (Mat_Nest*)A->data;
1380:   PetscInt       m,n,M,N,i,j,k,*dnnz,*onnz;
1381:   Mat            C;

1384:   MatGetSize(A,&M,&N);
1385:   MatGetLocalSize(A,&m,&n);
1386:   switch (reuse) {
1387:   case MAT_INITIAL_MATRIX:
1388:     MatCreate(PetscObjectComm((PetscObject)A),&C);
1389:     MatSetType(C,newtype);
1390:     MatSetSizes(C,m,n,M,N);
1391:     *newmat = C;
1392:     break;
1393:   case MAT_REUSE_MATRIX:
1394:     C = *newmat;
1395:     break;
1396:   default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse");
1397:   }

1399:   /* Preallocation */
1400:   PetscMalloc1(2*m,&dnnz);
1401:   onnz = dnnz + m;
1402:   for (k=0; k<m; k++) {
1403:     dnnz[k] = 0;
1404:     onnz[k] = 0;
1405:   }
1406:   for (j=0; j<nest->nc; ++j) {
1407:     IS             bNis;
1408:     PetscInt       bN;
1409:     const PetscInt *bNindices;
1410:     /* Using global column indices and ISAllGather() is not scalable. */
1411:     ISAllGather(nest->isglobal.col[j], &bNis);
1412:     ISGetSize(bNis, &bN);
1413:     ISGetIndices(bNis,&bNindices);
1414:     for (i=0; i<nest->nr; ++i) {
1415:       PetscSF        bmsf;
1416:       PetscSFNode    *bmedges;
1417:       Mat            B;
1418:       PetscInt       bm, *bmdnnz, br;
1419:       const PetscInt *bmindices;
1420:       B = nest->m[i][j];
1421:       if (!B) continue;
1422:       ISGetLocalSize(nest->isglobal.row[i],&bm);
1423:       ISGetIndices(nest->isglobal.row[i],&bmindices);
1424:       PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);
1425:       PetscMalloc1(2*bm,&bmedges);
1426:       PetscMalloc1(2*bm,&bmdnnz);
1427:       for (k = 0; k < 2*bm; ++k) bmdnnz[k] = 0;
1428:       /*
1429:        Locate the owners for all of the locally-owned global row indices for this row block.
1430:        These determine the roots of PetscSF used to communicate preallocation data to row owners.
1431:        The roots correspond to the dnnz and onnz entries; thus, there are two roots per row.
1432:        */
1433:       for (br = 0; br < bm; ++br) {
1434:         PetscInt       row = bmindices[br], rowowner = 0, brncols, col, colowner = 0;
1435:         const PetscInt *brcols;
1436:         PetscInt       rowrel = 0; /* row's relative index on its owner rank */
1437:         PetscInt       rowownerm; /* local row size on row's owning rank. */

1439:         PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);
1440:         rowownerm = A->rmap->range[rowowner+1]-A->rmap->range[rowowner];

1442:         bmedges[br].rank = rowowner; bmedges[br].index = rowrel;           /* edge from bmdnnz to dnnz */
1443:         bmedges[br].rank = rowowner; bmedges[br].index = rowrel+rowownerm; /* edge from bmonnz to onnz */
1444:         /* Now actually compute the data -- bmdnnz and bmonnz by looking at the global columns in the br row of this block. */
1445:         /* Note that this is not a pessimistic bound only because we assume the index sets embedding the blocks do not overlap. */
1446:         MatGetRow(B,br,&brncols,&brcols,NULL);
1447:         for (k=0; k<brncols; k++) {
1448:           col  = bNindices[brcols[k]];
1449:           PetscLayoutFindOwnerIndex(A->cmap,col,&colowner,NULL);
1450:           if (colowner == rowowner) bmdnnz[br]++;
1451:           else onnz[br]++;
1452:         }
1453:         MatRestoreRow(B,br,&brncols,&brcols,NULL);
1454:       }
1455:       ISRestoreIndices(nest->isglobal.row[i],&bmindices);
1456:       /* bsf will have to take care of disposing of bedges. */
1457:       PetscSFSetGraph(bmsf,m,2*bm,NULL,PETSC_COPY_VALUES,bmedges,PETSC_OWN_POINTER);
1458:       PetscSFReduceBegin(bmsf,MPIU_INT,bmdnnz,dnnz,MPIU_SUM);
1459:       PetscSFReduceEnd(bmsf,MPIU_INT,bmdnnz,dnnz,MPIU_SUM);
1460:       PetscFree(bmdnnz);
1461:       PetscSFDestroy(&bmsf);
1462:     }
1463:     ISRestoreIndices(bNis,&bNindices);
1464:     ISDestroy(&bNis);
1465:   }
1466:   MatSeqAIJSetPreallocation(C,0,dnnz);
1467:   MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);
1468:   PetscFree(dnnz);

1470:   /* Fill by row */
1471:   for (j=0; j<nest->nc; ++j) {
1472:     /* Using global column indices and ISAllGather() is not scalable. */
1473:     IS             bNis;
1474:     PetscInt       bN;
1475:     const PetscInt *bNindices;
1476:     ISAllGather(nest->isglobal.col[j], &bNis);
1477:     ISGetSize(bNis,&bN);
1478:     ISGetIndices(bNis,&bNindices);
1479:     for (i=0; i<nest->nr; ++i) {
1480:       Mat            B;
1481:       PetscInt       bm, br;
1482:       const PetscInt *bmindices;
1483:       B = nest->m[i][j];
1484:       if (!B) continue;
1485:       ISGetLocalSize(nest->isglobal.row[i],&bm);
1486:       ISGetIndices(nest->isglobal.row[i],&bmindices);
1487:       for (br = 0; br < bm; ++br) {
1488:         PetscInt          row = bmindices[br], brncols,  *cols;
1489:         const PetscInt    *brcols;
1490:         const PetscScalar *brcoldata;
1491:         MatGetRow(B,br,&brncols,&brcols,&brcoldata);
1492:         PetscMalloc1(brncols,&cols);
1493:         for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]];
1494:         /*
1495:          Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match.
1496:          Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES.
1497:          */
1498:         MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);
1499:         MatRestoreRow(B,br,&brncols,&brcols,&brcoldata);
1500:         PetscFree(cols);
1501:       }
1502:       ISRestoreIndices(nest->isglobal.row[i],&bmindices);
1503:     }
1504:     ISRestoreIndices(bNis,&bNindices);
1505:     ISDestroy(&bNis);
1506:   }
1507:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1508:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1509:   return(0);
1510: }

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

1515:   Level: intermediate

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

1522: .seealso: MatCreate(), MatType, MatCreateNest()
1523: M*/
1526: PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A)
1527: {
1528:   Mat_Nest       *s;

1532:   PetscNewLog(A,&s);
1533:   A->data = (void*)s;

1535:   s->nr            = -1;
1536:   s->nc            = -1;
1537:   s->m             = NULL;
1538:   s->splitassembly = PETSC_FALSE;

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

1542:   A->ops->mult                  = MatMult_Nest;
1543:   A->ops->multadd               = MatMultAdd_Nest;
1544:   A->ops->multtranspose         = MatMultTranspose_Nest;
1545:   A->ops->multtransposeadd      = MatMultTransposeAdd_Nest;
1546:   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
1547:   A->ops->assemblyend           = MatAssemblyEnd_Nest;
1548:   A->ops->zeroentries           = MatZeroEntries_Nest;
1549:   A->ops->duplicate             = MatDuplicate_Nest;
1550:   A->ops->getsubmatrix          = MatGetSubMatrix_Nest;
1551:   A->ops->destroy               = MatDestroy_Nest;
1552:   A->ops->view                  = MatView_Nest;
1553:   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
1554:   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
1555:   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
1556:   A->ops->getdiagonal           = MatGetDiagonal_Nest;
1557:   A->ops->diagonalscale         = MatDiagonalScale_Nest;
1558:   A->ops->scale                 = MatScale_Nest;
1559:   A->ops->shift                 = MatShift_Nest;

1561:   A->spptr        = 0;
1562:   A->assembled    = PETSC_FALSE;

1564:   /* expose Nest api's */
1565:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",   MatNestGetSubMat_Nest);
1566:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",   MatNestSetSubMat_Nest);
1567:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",  MatNestGetSubMats_Nest);
1568:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",     MatNestGetSize_Nest);
1569:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",      MatNestGetISs_Nest);
1570:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C", MatNestGetLocalISs_Nest);
1571:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",  MatNestSetVecType_Nest);
1572:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",  MatNestSetSubMats_Nest);

1574:   PetscObjectChangeTypeName((PetscObject)A,MATNEST);
1575:   return(0);
1576: }