Actual source code: matnest.c

petsc-3.6.4 2016-04-12
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:   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,0.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 MatCreateVecs_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:     PetscMalloc1(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:           MatCreateVecs(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:     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) {
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 MatCopy_Nest(Mat A,Mat B,MatStructure str)
661: {
662:   Mat_Nest       *bA = (Mat_Nest*)A->data,*bB = (Mat_Nest*)B->data;
663:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;

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

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

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

705:   MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
706:   MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
707:   return(0);
708: }

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

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

726: /*@
727:  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.

729:  Not collective

731:  Input Parameters:
732: +   A  - nest matrix
733: .   idxm - index of the matrix within the nest matrix
734: -   jdxm - index of the matrix within the nest matrix

736:  Output Parameter:
737: .   sub - matrix at index idxm,jdxm within the nest matrix

739:  Level: developer

741: .seealso: MatNestGetSize(), MatNestGetSubMats()
742: @*/
743: PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
744: {

748:   PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));
749:   return(0);
750: }

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

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

772:   PetscObjectReference((PetscObject)mat);
773:   MatDestroy(&bA->m[idxm][jdxm]);
774:   bA->m[idxm][jdxm] = mat;
775:   return(0);
776: }

780: /*@
781:  MatNestSetSubMat - Set a single submatrix in the nest matrix.

783:  Logically collective on the submatrix communicator

785:  Input Parameters:
786: +   A  - nest matrix
787: .   idxm - index of the matrix within the nest matrix
788: .   jdxm - index of the matrix within the nest matrix
789: -   sub - matrix at index idxm,jdxm within the nest matrix

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

794:  This increments the reference count of the submatrix.

796:  Level: developer

798: .seealso: MatNestSetSubMats(), MatNestGetSubMat()
799: @*/
800: PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
801: {

805:   PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));
806:   return(0);
807: }

811: PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
812: {
813:   Mat_Nest *bA = (Mat_Nest*)A->data;

816:   if (M)   *M   = bA->nr;
817:   if (N)   *N   = bA->nc;
818:   if (mat) *mat = bA->m;
819:   return(0);
820: }

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

827:  Not collective

829:  Input Parameters:
830: .   A  - nest matrix

832:  Output Parameter:
833: +   M - number of rows in the nest matrix
834: .   N - number of cols in the nest matrix
835: -   mat - 2d array of matrices

837:  Notes:

839:  The user should not free the array mat.

841:  Level: developer

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

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

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

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

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

871:  Not collective

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

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

880:  Notes:

882:  Level: developer

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

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

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

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

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

913:  Not collective

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

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

922:  Level: advanced

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

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

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

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

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

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

957:  Not collective

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

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

966:  Level: advanced

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

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

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

985: PetscErrorCode  MatNestSetVecType_Nest(Mat A,VecType vtype)
986: {
988:   PetscBool      flg;

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

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

1003:  Not collective

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

1009:  Notes:

1011:  Level: developer

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

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

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

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

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

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

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

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

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

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

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

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

1076:    Collective on Mat

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

1086:    Level: advanced

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

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

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

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

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


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

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

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

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

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

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

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

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

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

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

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

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

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

1357:    Collective on Mat

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

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

1370:    Level: advanced

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

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

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

1399:   MatGetSize(A,&M,&N);
1400:   MatGetLocalSize(A,&m,&n);
1401:   switch (reuse) {
1402:   case MAT_INITIAL_MATRIX:
1403:     MatCreate(PetscObjectComm((PetscObject)A),&C);
1404:     MatSetType(C,newtype);
1405:     MatSetSizes(C,m,n,M,N);
1406:     *newmat = C;
1407:     break;
1408:   case MAT_REUSE_MATRIX:
1409:     C = *newmat;
1410:     break;
1411:   default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse");
1412:   }

1414:   /* Preallocation */
1415:   PetscMalloc1(2*m,&dnnz);
1416:   onnz = dnnz + m;
1417:   for (k=0; k<m; k++) {
1418:     dnnz[k] = 0;
1419:     onnz[k] = 0;
1420:   }
1421:   for (j=0; j<nest->nc; ++j) {
1422:     IS             bNis;
1423:     PetscInt       bN;
1424:     const PetscInt *bNindices;
1425:     /* Using global column indices and ISAllGather() is not scalable. */
1426:     ISAllGather(nest->isglobal.col[j], &bNis);
1427:     ISGetSize(bNis, &bN);
1428:     ISGetIndices(bNis,&bNindices);
1429:     for (i=0; i<nest->nr; ++i) {
1430:       PetscSF        bmsf;
1431:       PetscSFNode    *bmedges;
1432:       Mat            B;
1433:       PetscInt       bm, *bmdnnz, br;
1434:       const PetscInt *bmindices;
1435:       B = nest->m[i][j];
1436:       if (!B) continue;
1437:       ISGetLocalSize(nest->isglobal.row[i],&bm);
1438:       ISGetIndices(nest->isglobal.row[i],&bmindices);
1439:       PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);
1440:       PetscMalloc1(bm,&bmedges);
1441:       PetscMalloc1(2*bm,&bmdnnz);
1442:       for (k = 0; k < 2*bm; ++k) bmdnnz[k] = 0;
1443:       /*
1444:        Locate the owners for all of the locally-owned global row indices for this row block.
1445:        These determine the roots of PetscSF used to communicate preallocation data to row owners.
1446:        The roots correspond to the dnnz and onnz entries; thus, there are two roots per row.
1447:        */
1448:       MatGetOwnershipRange(B,&rstart,NULL);
1449:       for (br = 0; br < bm; ++br) {
1450:         PetscInt       row = bmindices[br], rowowner = 0, brncols, col, colowner = 0;
1451:         const PetscInt *brcols;
1452:         PetscInt       rowrel = 0; /* row's relative index on its owner rank */
1453:         PetscInt       rowownerm; /* local row size on row's owning rank. */
1454:         PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);

1456:         rowownerm = A->rmap->range[rowowner+1]-A->rmap->range[rowowner];

1458:         bmedges[br].rank = rowowner; bmedges[br].index = rowrel;           /* edge from bmdnnz to dnnz */
1459:         bmedges[br].rank = rowowner; bmedges[br].index = rowrel+rowownerm; /* edge from bmonnz to onnz */
1460:         /* Now actually compute the data -- bmdnnz and bmonnz by looking at the global columns in the br row of this block. */
1461:         /* Note that this is not a pessimistic bound only because we assume the index sets embedding the blocks do not overlap. */
1462:         MatGetRow(B,br+rstart,&brncols,&brcols,NULL);
1463:         for (k=0; k<brncols; k++) {
1464:           col  = bNindices[brcols[k]];
1465:           PetscLayoutFindOwnerIndex(A->cmap,col,&colowner,NULL);
1466:           if (colowner == rowowner) bmdnnz[br]++;
1467:           else onnz[br]++;
1468:         }
1469:         MatRestoreRow(B,br+rstart,&brncols,&brcols,NULL);
1470:       }
1471:       ISRestoreIndices(nest->isglobal.row[i],&bmindices);
1472:       /* bsf will have to take care of disposing of bedges. */
1473:       PetscSFSetGraph(bmsf,m,bm,NULL,PETSC_COPY_VALUES,bmedges,PETSC_OWN_POINTER);
1474:       PetscSFReduceBegin(bmsf,MPIU_INT,bmdnnz,dnnz,MPIU_SUM);
1475:       PetscSFReduceEnd(bmsf,MPIU_INT,bmdnnz,dnnz,MPIU_SUM);
1476:       PetscFree(bmdnnz);
1477:       PetscSFDestroy(&bmsf);
1478:     }
1479:     ISRestoreIndices(bNis,&bNindices);
1480:     ISDestroy(&bNis);
1481:   }
1482:   /* dnnz is not correct */
1483:   MatSeqAIJSetPreallocation(C,0,dnnz);
1484:   MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);
1485:   PetscFree(dnnz);

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

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

1533:   Level: intermediate

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

1540: .seealso: MatCreate(), MatType, MatCreateNest()
1541: M*/
1544: PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A)
1545: {
1546:   Mat_Nest       *s;

1550:   PetscNewLog(A,&s);
1551:   A->data = (void*)s;

1553:   s->nr            = -1;
1554:   s->nc            = -1;
1555:   s->m             = NULL;
1556:   s->splitassembly = PETSC_FALSE;

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

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

1580:   A->spptr        = 0;
1581:   A->assembled    = PETSC_FALSE;

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

1594:   PetscObjectChangeTypeName((PetscObject)A,MATNEST);
1595:   return(0);
1596: }