Actual source code: matnest.c

petsc-3.13.1 2020-05-02
Report Typos and Errors
  1:  #include <../src/mat/impls/nest/matnestimpl.h>
  2:  #include <../src/mat/impls/aij/seq/aij.h>
  3:  #include <petscsf.h>

  5: static PetscErrorCode MatSetUp_NestIS_Private(Mat,PetscInt,const IS[],PetscInt,const IS[]);
  6: static PetscErrorCode MatCreateVecs_Nest(Mat,Vec*,Vec*);
  7: static PetscErrorCode MatReset_Nest(Mat);

  9: PETSC_INTERN PetscErrorCode MatConvert_Nest_IS(Mat,MatType,MatReuse,Mat*);

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

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

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

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

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

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

 89: typedef struct {
 90:   Mat          *workC;    /* array of Mat with specific containers depending on the underlying MatMatMult implementation */
 91:   PetscScalar  *tarray;   /* buffer for storing all temporary products A[i][j] B[j] */
 92:   PetscInt     *dm,*dn,k; /* displacements and number of submatrices */
 93: } Nest_Dense;

 95: PETSC_INTERN PetscErrorCode MatMatMultNumeric_Nest_Dense(Mat A,Mat B,Mat C)
 96: {
 97:   Mat_Nest          *bA = (Mat_Nest*)A->data;
 98:   PetscContainer    container;
 99:   Nest_Dense        *contents;
100:   Mat               viewB,viewC,seq,productB,workC;
101:   const PetscScalar *barray;
102:   PetscScalar       *carray;
103:   PetscInt          i,j,M,N,nr = bA->nr,nc = bA->nc,ldb,ldc;
104:   PetscErrorCode    ierr;

107:   PetscObjectQuery((PetscObject)C,"workC",(PetscObject*)&container);
108:   if (!container) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exist");
109:   PetscContainerGetPointer(container,(void**)&contents);
110:   MatDenseGetLDA(B,&ldb);
111:   MatDenseGetLDA(C,&ldc);
112:   MatGetSize(B,NULL,&N);
113:   MatZeroEntries(C);
114:   MatDenseGetArrayRead(B,&barray);
115:   MatDenseGetArray(C,&carray);
116:   for (i=0; i<nr; i++) {
117:     ISGetSize(bA->isglobal.row[i],&M);
118:     MatCreateDense(PetscObjectComm((PetscObject)A),contents->dm[i+1]-contents->dm[i],PETSC_DECIDE,M,N,carray+contents->dm[i],&viewC);
119:     MatDenseGetLocalMatrix(viewC,&seq);
120:     MatSeqDenseSetLDA(seq,ldc);
121:     for (j=0; j<nc; j++) {
122:       if (!bA->m[i][j]) continue;
123:       ISGetSize(bA->isglobal.col[j],&M);
124:       MatCreateDense(PetscObjectComm((PetscObject)A),contents->dn[j+1]-contents->dn[j],PETSC_DECIDE,M,N,(PetscScalar*)(barray+contents->dn[j]),&viewB);
125:       MatDenseGetLocalMatrix(viewB,&seq);
126:       MatSeqDenseSetLDA(seq,ldb);

128:       /* MatMatMultNumeric(bA->m[i][j],viewB,contents->workC[i*nc + j]); */
129:       workC             = contents->workC[i*nc + j];
130:       productB          = workC->product->B;
131:       workC->product->B = viewB; /* use newly created dense matrix viewB */
132:       (workC->ops->productnumeric)(workC);
133:       MatDestroy(&viewB);
134:       workC->product->B = productB; /* resume original B */

136:       /* C[i] <- workC + C[i] */
137:       MatAXPY(viewC,1.0,contents->workC[i*nc + j],SAME_NONZERO_PATTERN);
138:     }
139:     MatDestroy(&viewC);
140:   }
141:   MatDenseRestoreArray(C,&carray);
142:   MatDenseRestoreArrayRead(B,&barray);

144:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
145:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
146:   return(0);
147: }

149: PetscErrorCode MatNest_DenseDestroy(void *ctx)
150: {
151:   Nest_Dense     *contents = (Nest_Dense*)ctx;
152:   PetscInt       i;

156:   PetscFree(contents->tarray);
157:   for (i=0; i<contents->k; i++) {
158:     MatDestroy(contents->workC + i);
159:   }
160:   PetscFree3(contents->dm,contents->dn,contents->workC);
161:   PetscFree(contents);
162:   return(0);
163: }

165: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_Nest_Dense(Mat A,Mat B,PetscReal fill,Mat C)
166: {
167:   Mat_Nest          *bA = (Mat_Nest*)A->data;
168:   Mat               viewB,viewSeq,workC;
169:   const PetscScalar *barray;
170:   PetscInt          i,j,M,N,m,nr = bA->nr,nc = bA->nc,maxm = 0,ldb;
171:   PetscContainer    container;
172:   Nest_Dense        *contents=NULL;
173:   PetscErrorCode    ierr;

176:   if (!C->assembled) {
177:     MatGetSize(B,NULL,&N);
178:     MatGetLocalSize(A,&m,NULL);
179:     MatGetSize(A,&M,NULL);

181:     MatSetSizes(C,m,PETSC_DECIDE,M,N);
182:     MatSetType(C,MATDENSE);
183:     MatSeqDenseSetPreallocation(C,NULL);
184:     MatMPIDenseSetPreallocation(C,NULL);
185:   }

187:   PetscNew(&contents);
188:   PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);
189:   PetscContainerSetPointer(container,contents);
190:   PetscContainerSetUserDestroy(container,MatNest_DenseDestroy);
191:   PetscObjectCompose((PetscObject)C,"workC",(PetscObject)container);
192:   PetscContainerDestroy(&container);
193:   PetscCalloc3(nr+1,&contents->dm,nc+1,&contents->dn,nr*nc,&contents->workC);
194:   contents->k = nr*nc;
195:   for (i=0; i<nr; i++) {
196:     ISGetLocalSize(bA->isglobal.row[i],contents->dm + i+1);
197:     maxm = PetscMax(maxm,contents->dm[i+1]);
198:     contents->dm[i+1] += contents->dm[i];
199:   }
200:   for (i=0; i<nc; i++) {
201:     ISGetLocalSize(bA->isglobal.col[i],contents->dn + i+1);
202:     contents->dn[i+1] += contents->dn[i];
203:   }
204:   PetscMalloc1(maxm*N,&contents->tarray);
205:   MatDenseGetLDA(B,&ldb);
206:   MatGetSize(B,NULL,&N);
207:   MatDenseGetArrayRead(B,&barray);
208:   /* loops are permuted compared to MatMatMultNumeric so that viewB is created only once per column of A */
209:   for (j=0; j<nc; j++) {
210:     ISGetSize(bA->isglobal.col[j],&M);
211:     MatCreateDense(PetscObjectComm((PetscObject)A),contents->dn[j+1]-contents->dn[j],PETSC_DECIDE,M,N,(PetscScalar*)(barray+contents->dn[j]),&viewB);
212:     MatDenseGetLocalMatrix(viewB,&viewSeq);
213:     MatSeqDenseSetLDA(viewSeq,ldb);
214:     for (i=0; i<nr; i++) {
215:       if (!bA->m[i][j]) continue;
216:       /* MatMatMultSymbolic may attach a specific container (depending on MatType of bA->m[i][j]) to workC[i][j] */

218:       MatProductCreate(bA->m[i][j],viewB,NULL,&contents->workC[i*nc + j]);
219:       workC = contents->workC[i*nc + j];
220:       MatProductSetType(workC,MATPRODUCT_AB);
221:       MatProductSetAlgorithm(workC,"default");
222:       MatProductSetFill(workC,fill);
223:       MatProductSetFromOptions(workC);
224:       MatProductSymbolic(workC);

226:       MatDenseGetLocalMatrix(workC,&viewSeq);
227:       /* free the memory allocated in MatMatMultSymbolic, since tarray will be shared by all Mat */
228:       MatSeqDenseSetPreallocation(viewSeq,contents->tarray);
229:     }
230:     MatDestroy(&viewB);
231:   }
232:   MatDenseRestoreArrayRead(B,&barray);

234:   C->ops->matmultnumeric = MatMatMultNumeric_Nest_Dense;
235:   return(0);
236: }

238: /* --------------------------------------------------------- */
239: static PetscErrorCode MatProductSetFromOptions_Nest_Dense_AB(Mat C)
240: {
242:   C->ops->matmultsymbolic = MatMatMultSymbolic_Nest_Dense;
243:   C->ops->productsymbolic = MatProductSymbolic_AB;
244:   return(0);
245: }

247: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_Nest_Dense(Mat C)
248: {
250:   Mat_Product    *product = C->product;

253:   MatSetType(C,MATDENSE);
254:   if (product->type == MATPRODUCT_AB) {
255:     MatProductSetFromOptions_Nest_Dense_AB(C);
256:   } else SETERRQ1(PetscObjectComm((PetscObject)C),PETSC_ERR_SUP,"MatProduct type %s is not supported for Nest and Dense matrices",MatProductTypes[product->type]);
257:   return(0);
258: }
259: /* --------------------------------------------------------- */

261: static PetscErrorCode MatMultTranspose_Nest(Mat A,Vec x,Vec y)
262: {
263:   Mat_Nest       *bA = (Mat_Nest*)A->data;
264:   Vec            *bx = bA->left,*by = bA->right;
265:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;

269:   for (i=0; i<nr; i++) {VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);}
270:   for (i=0; i<nc; i++) {VecGetSubVector(y,bA->isglobal.col[i],&by[i]);}
271:   for (j=0; j<nc; j++) {
272:     VecZeroEntries(by[j]);
273:     for (i=0; i<nr; i++) {
274:       if (!bA->m[i][j]) continue;
275:       /* y[j] <- y[j] + (A[i][j])^T * x[i] */
276:       MatMultTransposeAdd(bA->m[i][j],bx[i],by[j],by[j]);
277:     }
278:   }
279:   for (i=0; i<nr; i++) {VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);}
280:   for (i=0; i<nc; i++) {VecRestoreSubVector(y,bA->isglobal.col[i],&by[i]);}
281:   return(0);
282: }

284: static PetscErrorCode MatMultTransposeAdd_Nest(Mat A,Vec x,Vec y,Vec z)
285: {
286:   Mat_Nest       *bA = (Mat_Nest*)A->data;
287:   Vec            *bx = bA->left,*bz = bA->right;
288:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;

292:   for (i=0; i<nr; i++) {VecGetSubVector(x,bA->isglobal.row[i],&bx[i]);}
293:   for (i=0; i<nc; i++) {VecGetSubVector(z,bA->isglobal.col[i],&bz[i]);}
294:   for (j=0; j<nc; j++) {
295:     if (y != z) {
296:       Vec by;
297:       VecGetSubVector(y,bA->isglobal.col[j],&by);
298:       VecCopy(by,bz[j]);
299:       VecRestoreSubVector(y,bA->isglobal.col[j],&by);
300:     }
301:     for (i=0; i<nr; i++) {
302:       if (!bA->m[i][j]) continue;
303:       /* z[j] <- y[j] + (A[i][j])^T * x[i] */
304:       MatMultTransposeAdd(bA->m[i][j],bx[i],bz[j],bz[j]);
305:     }
306:   }
307:   for (i=0; i<nr; i++) {VecRestoreSubVector(x,bA->isglobal.row[i],&bx[i]);}
308:   for (i=0; i<nc; i++) {VecRestoreSubVector(z,bA->isglobal.col[i],&bz[i]);}
309:   return(0);
310: }

312: static PetscErrorCode MatTranspose_Nest(Mat A,MatReuse reuse,Mat *B)
313: {
314:   Mat_Nest       *bA = (Mat_Nest*)A->data, *bC;
315:   Mat            C;
316:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;

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

322:   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_INPLACE_MATRIX) {
323:     Mat *subs;
324:     IS  *is_row,*is_col;

326:     PetscCalloc1(nr * nc,&subs);
327:     PetscMalloc2(nr,&is_row,nc,&is_col);
328:     MatNestGetISs(A,is_row,is_col);
329:     if (reuse == MAT_INPLACE_MATRIX) {
330:       for (i=0; i<nr; i++) {
331:         for (j=0; j<nc; j++) {
332:           subs[i + nr * j] = bA->m[i][j];
333:         }
334:       }
335:     }

337:     MatCreateNest(PetscObjectComm((PetscObject)A),nc,is_col,nr,is_row,subs,&C);
338:     PetscFree(subs);
339:     PetscFree2(is_row,is_col);
340:   } else {
341:     C = *B;
342:   }

344:   bC = (Mat_Nest*)C->data;
345:   for (i=0; i<nr; i++) {
346:     for (j=0; j<nc; j++) {
347:       if (bA->m[i][j]) {
348:         MatTranspose(bA->m[i][j], reuse, &(bC->m[j][i]));
349:       } else {
350:         bC->m[j][i] = NULL;
351:       }
352:     }
353:   }

355:   if (reuse == MAT_INITIAL_MATRIX || reuse == MAT_REUSE_MATRIX) {
356:     *B = C;
357:   } else {
358:     MatHeaderMerge(A, &C);
359:   }
360:   return(0);
361: }

363: static PetscErrorCode MatNestDestroyISList(PetscInt n,IS **list)
364: {
366:   IS             *lst = *list;
367:   PetscInt       i;

370:   if (!lst) return(0);
371:   for (i=0; i<n; i++) if (lst[i]) {ISDestroy(&lst[i]);}
372:   PetscFree(lst);
373:   *list = NULL;
374:   return(0);
375: }

377: static PetscErrorCode MatReset_Nest(Mat A)
378: {
379:   Mat_Nest       *vs = (Mat_Nest*)A->data;
380:   PetscInt       i,j;

384:   /* release the matrices and the place holders */
385:   MatNestDestroyISList(vs->nr,&vs->isglobal.row);
386:   MatNestDestroyISList(vs->nc,&vs->isglobal.col);
387:   MatNestDestroyISList(vs->nr,&vs->islocal.row);
388:   MatNestDestroyISList(vs->nc,&vs->islocal.col);

390:   PetscFree(vs->row_len);
391:   PetscFree(vs->col_len);
392:   PetscFree(vs->nnzstate);

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

396:   /* release the matrices and the place holders */
397:   if (vs->m) {
398:     for (i=0; i<vs->nr; i++) {
399:       for (j=0; j<vs->nc; j++) {
400:         MatDestroy(&vs->m[i][j]);
401:       }
402:       PetscFree(vs->m[i]);
403:     }
404:     PetscFree(vs->m);
405:   }

407:   /* restore defaults */
408:   vs->nr = 0;
409:   vs->nc = 0;
410:   vs->splitassembly = PETSC_FALSE;
411:   return(0);
412: }

414: static PetscErrorCode MatDestroy_Nest(Mat A)
415: {

418:   MatReset_Nest(A);
419:   PetscFree(A->data);

421:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",0);
422:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",0);
423:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",0);
424:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",0);
425:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",0);
426:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",0);
427:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",0);
428:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",0);
429:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_mpiaij_C",0);
430:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_seqaij_C",0);
431:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",0);
432:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_is_C",0);
433:   PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_nest_seqdense_C",NULL);
434:   PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_nest_mpidense_C",NULL);
435:   PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_nest_dense_C",NULL);
436:   return(0);
437: }

439: static PetscErrorCode MatMissingDiagonal_Nest(Mat mat,PetscBool *missing,PetscInt *dd)
440: {
441:   Mat_Nest       *vs = (Mat_Nest*)mat->data;
442:   PetscInt       i;

446:   if (dd) *dd = 0;
447:   if (!vs->nr) {
448:     *missing = PETSC_TRUE;
449:     return(0);
450:   }
451:   *missing = PETSC_FALSE;
452:   for (i = 0; i < vs->nr && !(*missing); i++) {
453:     *missing = PETSC_TRUE;
454:     if (vs->m[i][i]) {
455:       MatMissingDiagonal(vs->m[i][i],missing,NULL);
456:       if (*missing && dd) SETERRQ(PetscObjectComm((PetscObject)mat),PETSC_ERR_SUP,"First missing entry not yet implemented");
457:     }
458:   }
459:   return(0);
460: }

462: static PetscErrorCode MatAssemblyBegin_Nest(Mat A,MatAssemblyType type)
463: {
464:   Mat_Nest       *vs = (Mat_Nest*)A->data;
465:   PetscInt       i,j;
467:   PetscBool      nnzstate = PETSC_FALSE;

470:   for (i=0; i<vs->nr; i++) {
471:     for (j=0; j<vs->nc; j++) {
472:       PetscObjectState subnnzstate = 0;
473:       if (vs->m[i][j]) {
474:         MatAssemblyBegin(vs->m[i][j],type);
475:         if (!vs->splitassembly) {
476:           /* Note: split assembly will fail if the same block appears more than once (even indirectly through a nested
477:            * sub-block). This could be fixed by adding a flag to Mat so that there was a way to check if a Mat was
478:            * already performing an assembly, but the result would by more complicated and appears to offer less
479:            * potential for diagnostics and correctness checking. Split assembly should be fixed once there is an
480:            * interface for libraries to make asynchronous progress in "user-defined non-blocking collectives".
481:            */
482:           MatAssemblyEnd(vs->m[i][j],type);
483:           MatGetNonzeroState(vs->m[i][j],&subnnzstate);
484:         }
485:       }
486:       nnzstate = (PetscBool)(nnzstate || vs->nnzstate[i*vs->nc+j] != subnnzstate);
487:       vs->nnzstate[i*vs->nc+j] = subnnzstate;
488:     }
489:   }
490:   if (nnzstate) A->nonzerostate++;
491:   return(0);
492: }

494: static PetscErrorCode MatAssemblyEnd_Nest(Mat A, MatAssemblyType type)
495: {
496:   Mat_Nest       *vs = (Mat_Nest*)A->data;
497:   PetscInt       i,j;

501:   for (i=0; i<vs->nr; i++) {
502:     for (j=0; j<vs->nc; j++) {
503:       if (vs->m[i][j]) {
504:         if (vs->splitassembly) {
505:           MatAssemblyEnd(vs->m[i][j],type);
506:         }
507:       }
508:     }
509:   }
510:   return(0);
511: }

513: static PetscErrorCode MatNestFindNonzeroSubMatRow(Mat A,PetscInt row,Mat *B)
514: {
516:   Mat_Nest       *vs = (Mat_Nest*)A->data;
517:   PetscInt       j;
518:   Mat            sub;

521:   sub = (row < vs->nc) ? vs->m[row][row] : (Mat)NULL; /* Prefer to find on the diagonal */
522:   for (j=0; !sub && j<vs->nc; j++) sub = vs->m[row][j];
523:   if (sub) {MatSetUp(sub);}       /* Ensure that the sizes are available */
524:   *B = sub;
525:   return(0);
526: }

528: static PetscErrorCode MatNestFindNonzeroSubMatCol(Mat A,PetscInt col,Mat *B)
529: {
531:   Mat_Nest       *vs = (Mat_Nest*)A->data;
532:   PetscInt       i;
533:   Mat            sub;

536:   sub = (col < vs->nr) ? vs->m[col][col] : (Mat)NULL; /* Prefer to find on the diagonal */
537:   for (i=0; !sub && i<vs->nr; i++) sub = vs->m[i][col];
538:   if (sub) {MatSetUp(sub);}       /* Ensure that the sizes are available */
539:   *B = sub;
540:   return(0);
541: }

543: static PetscErrorCode MatNestFindIS(Mat A,PetscInt n,const IS list[],IS is,PetscInt *found)
544: {
546:   PetscInt       i;
547:   PetscBool      flg;

553:   *found = -1;
554:   for (i=0; i<n; i++) {
555:     if (!list[i]) continue;
556:     ISEqualUnsorted(list[i],is,&flg);
557:     if (flg) {
558:       *found = i;
559:       return(0);
560:     }
561:   }
562:   SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_INCOMP,"Could not find index set");
563:   return(0);
564: }

566: /* Get a block row as a new MatNest */
567: static PetscErrorCode MatNestGetRow(Mat A,PetscInt row,Mat *B)
568: {
569:   Mat_Nest       *vs = (Mat_Nest*)A->data;
570:   char           keyname[256];

574:   *B   = NULL;
575:   PetscSNPrintf(keyname,sizeof(keyname),"NestRow_%D",row);
576:   PetscObjectQuery((PetscObject)A,keyname,(PetscObject*)B);
577:   if (*B) return(0);

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

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

583:   PetscObjectCompose((PetscObject)A,keyname,(PetscObject)*B);
584:   PetscObjectDereference((PetscObject)*B); /* Leave the only remaining reference in the composition */
585:   return(0);
586: }

588: static PetscErrorCode MatNestFindSubMat(Mat A,struct MatNestISPair *is,IS isrow,IS iscol,Mat *B)
589: {
590:   Mat_Nest       *vs = (Mat_Nest*)A->data;
592:   PetscInt       row,col;
593:   PetscBool      same,isFullCol,isFullColGlobal;

596:   /* Check if full column space. This is a hack */
597:   isFullCol = PETSC_FALSE;
598:   PetscObjectTypeCompare((PetscObject)iscol,ISSTRIDE,&same);
599:   if (same) {
600:     PetscInt n,first,step,i,an,am,afirst,astep;
601:     ISStrideGetInfo(iscol,&first,&step);
602:     ISGetLocalSize(iscol,&n);
603:     isFullCol = PETSC_TRUE;
604:     for (i=0,an=A->cmap->rstart; i<vs->nc; i++) {
605:       PetscObjectTypeCompare((PetscObject)is->col[i],ISSTRIDE,&same);
606:       ISGetLocalSize(is->col[i],&am);
607:       if (same) {
608:         ISStrideGetInfo(is->col[i],&afirst,&astep);
609:         if (afirst != an || astep != step) isFullCol = PETSC_FALSE;
610:       } else isFullCol = PETSC_FALSE;
611:       an += am;
612:     }
613:     if (an != A->cmap->rstart+n) isFullCol = PETSC_FALSE;
614:   }
615:   MPIU_Allreduce(&isFullCol,&isFullColGlobal,1,MPIU_BOOL,MPI_LAND,PetscObjectComm((PetscObject)iscol));

617:   if (isFullColGlobal && vs->nc > 1) {
618:     PetscInt row;
619:     MatNestFindIS(A,vs->nr,is->row,isrow,&row);
620:     MatNestGetRow(A,row,B);
621:   } else {
622:     MatNestFindIS(A,vs->nr,is->row,isrow,&row);
623:     MatNestFindIS(A,vs->nc,is->col,iscol,&col);
624:     if (!vs->m[row][col]) {
625:       PetscInt lr,lc;

627:       MatCreate(PetscObjectComm((PetscObject)A),&vs->m[row][col]);
628:       ISGetLocalSize(vs->isglobal.row[row],&lr);
629:       ISGetLocalSize(vs->isglobal.col[col],&lc);
630:       MatSetSizes(vs->m[row][col],lr,lc,PETSC_DECIDE,PETSC_DECIDE);
631:       MatSetType(vs->m[row][col],MATAIJ);
632:       MatSeqAIJSetPreallocation(vs->m[row][col],0,NULL);
633:       MatMPIAIJSetPreallocation(vs->m[row][col],0,NULL,0,NULL);
634:       MatSetUp(vs->m[row][col]);
635:       MatAssemblyBegin(vs->m[row][col],MAT_FINAL_ASSEMBLY);
636:       MatAssemblyEnd(vs->m[row][col],MAT_FINAL_ASSEMBLY);
637:     }
638:     *B = vs->m[row][col];
639:   }
640:   return(0);
641: }

643: /*
644:    TODO: This does not actually returns a submatrix we can modify
645: */
646: static PetscErrorCode MatCreateSubMatrix_Nest(Mat A,IS isrow,IS iscol,MatReuse reuse,Mat *B)
647: {
649:   Mat_Nest       *vs = (Mat_Nest*)A->data;
650:   Mat            sub;

653:   MatNestFindSubMat(A,&vs->isglobal,isrow,iscol,&sub);
654:   switch (reuse) {
655:   case MAT_INITIAL_MATRIX:
656:     if (sub) { PetscObjectReference((PetscObject)sub); }
657:     *B = sub;
658:     break;
659:   case MAT_REUSE_MATRIX:
660:     if (sub != *B) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Submatrix was not used before in this call");
661:     break;
662:   case MAT_IGNORE_MATRIX:       /* Nothing to do */
663:     break;
664:   case MAT_INPLACE_MATRIX:       /* Nothing to do */
665:     SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MAT_INPLACE_MATRIX is not supported yet");
666:     break;
667:   }
668:   return(0);
669: }

671: PetscErrorCode MatGetLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
672: {
674:   Mat_Nest       *vs = (Mat_Nest*)A->data;
675:   Mat            sub;

678:   MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);
679:   /* We allow the submatrix to be NULL, perhaps it would be better for the user to return an empty matrix instead */
680:   if (sub) {PetscObjectReference((PetscObject)sub);}
681:   *B = sub;
682:   return(0);
683: }

685: static PetscErrorCode MatRestoreLocalSubMatrix_Nest(Mat A,IS isrow,IS iscol,Mat *B)
686: {
688:   Mat_Nest       *vs = (Mat_Nest*)A->data;
689:   Mat            sub;

692:   MatNestFindSubMat(A,&vs->islocal,isrow,iscol,&sub);
693:   if (*B != sub) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has not been gotten");
694:   if (sub) {
695:     if (((PetscObject)sub)->refct <= 1) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Local submatrix has had reference count decremented too many times");
696:     MatDestroy(B);
697:   }
698:   return(0);
699: }

701: static PetscErrorCode MatGetDiagonal_Nest(Mat A,Vec v)
702: {
703:   Mat_Nest       *bA = (Mat_Nest*)A->data;
704:   PetscInt       i;

708:   for (i=0; i<bA->nr; i++) {
709:     Vec bv;
710:     VecGetSubVector(v,bA->isglobal.row[i],&bv);
711:     if (bA->m[i][i]) {
712:       MatGetDiagonal(bA->m[i][i],bv);
713:     } else {
714:       VecSet(bv,0.0);
715:     }
716:     VecRestoreSubVector(v,bA->isglobal.row[i],&bv);
717:   }
718:   return(0);
719: }

721: static PetscErrorCode MatDiagonalScale_Nest(Mat A,Vec l,Vec r)
722: {
723:   Mat_Nest       *bA = (Mat_Nest*)A->data;
724:   Vec            bl,*br;
725:   PetscInt       i,j;

729:   PetscCalloc1(bA->nc,&br);
730:   if (r) {
731:     for (j=0; j<bA->nc; j++) {VecGetSubVector(r,bA->isglobal.col[j],&br[j]);}
732:   }
733:   bl = NULL;
734:   for (i=0; i<bA->nr; i++) {
735:     if (l) {
736:       VecGetSubVector(l,bA->isglobal.row[i],&bl);
737:     }
738:     for (j=0; j<bA->nc; j++) {
739:       if (bA->m[i][j]) {
740:         MatDiagonalScale(bA->m[i][j],bl,br[j]);
741:       }
742:     }
743:     if (l) {
744:       VecRestoreSubVector(l,bA->isglobal.row[i],&bl);
745:     }
746:   }
747:   if (r) {
748:     for (j=0; j<bA->nc; j++) {VecRestoreSubVector(r,bA->isglobal.col[j],&br[j]);}
749:   }
750:   PetscFree(br);
751:   return(0);
752: }

754: static PetscErrorCode MatScale_Nest(Mat A,PetscScalar a)
755: {
756:   Mat_Nest       *bA = (Mat_Nest*)A->data;
757:   PetscInt       i,j;

761:   for (i=0; i<bA->nr; i++) {
762:     for (j=0; j<bA->nc; j++) {
763:       if (bA->m[i][j]) {
764:         MatScale(bA->m[i][j],a);
765:       }
766:     }
767:   }
768:   return(0);
769: }

771: static PetscErrorCode MatShift_Nest(Mat A,PetscScalar a)
772: {
773:   Mat_Nest       *bA = (Mat_Nest*)A->data;
774:   PetscInt       i;
776:   PetscBool      nnzstate = PETSC_FALSE;

779:   for (i=0; i<bA->nr; i++) {
780:     PetscObjectState subnnzstate = 0;
781:     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);
782:     MatShift(bA->m[i][i],a);
783:     MatGetNonzeroState(bA->m[i][i],&subnnzstate);
784:     nnzstate = (PetscBool)(nnzstate || bA->nnzstate[i*bA->nc+i] != subnnzstate);
785:     bA->nnzstate[i*bA->nc+i] = subnnzstate;
786:   }
787:   if (nnzstate) A->nonzerostate++;
788:   return(0);
789: }

791: static PetscErrorCode MatDiagonalSet_Nest(Mat A,Vec D,InsertMode is)
792: {
793:   Mat_Nest       *bA = (Mat_Nest*)A->data;
794:   PetscInt       i;
796:   PetscBool      nnzstate = PETSC_FALSE;

799:   for (i=0; i<bA->nr; i++) {
800:     PetscObjectState subnnzstate = 0;
801:     Vec              bv;
802:     VecGetSubVector(D,bA->isglobal.row[i],&bv);
803:     if (bA->m[i][i]) {
804:       MatDiagonalSet(bA->m[i][i],bv,is);
805:       MatGetNonzeroState(bA->m[i][i],&subnnzstate);
806:     }
807:     VecRestoreSubVector(D,bA->isglobal.row[i],&bv);
808:     nnzstate = (PetscBool)(nnzstate || bA->nnzstate[i*bA->nc+i] != subnnzstate);
809:     bA->nnzstate[i*bA->nc+i] = subnnzstate;
810:   }
811:   if (nnzstate) A->nonzerostate++;
812:   return(0);
813: }

815: static PetscErrorCode MatSetRandom_Nest(Mat A,PetscRandom rctx)
816: {
817:   Mat_Nest       *bA = (Mat_Nest*)A->data;
818:   PetscInt       i,j;

822:   for (i=0; i<bA->nr; i++) {
823:     for (j=0; j<bA->nc; j++) {
824:       if (bA->m[i][j]) {
825:         MatSetRandom(bA->m[i][j],rctx);
826:       }
827:     }
828:   }
829:   return(0);
830: }

832: static PetscErrorCode MatCreateVecs_Nest(Mat A,Vec *right,Vec *left)
833: {
834:   Mat_Nest       *bA = (Mat_Nest*)A->data;
835:   Vec            *L,*R;
836:   MPI_Comm       comm;
837:   PetscInt       i,j;

841:   PetscObjectGetComm((PetscObject)A,&comm);
842:   if (right) {
843:     /* allocate R */
844:     PetscMalloc1(bA->nc, &R);
845:     /* Create the right vectors */
846:     for (j=0; j<bA->nc; j++) {
847:       for (i=0; i<bA->nr; i++) {
848:         if (bA->m[i][j]) {
849:           MatCreateVecs(bA->m[i][j],&R[j],NULL);
850:           break;
851:         }
852:       }
853:       if (i==bA->nr) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null column.");
854:     }
855:     VecCreateNest(comm,bA->nc,bA->isglobal.col,R,right);
856:     /* hand back control to the nest vector */
857:     for (j=0; j<bA->nc; j++) {
858:       VecDestroy(&R[j]);
859:     }
860:     PetscFree(R);
861:   }

863:   if (left) {
864:     /* allocate L */
865:     PetscMalloc1(bA->nr, &L);
866:     /* Create the left vectors */
867:     for (i=0; i<bA->nr; i++) {
868:       for (j=0; j<bA->nc; j++) {
869:         if (bA->m[i][j]) {
870:           MatCreateVecs(bA->m[i][j],NULL,&L[i]);
871:           break;
872:         }
873:       }
874:       if (j==bA->nc) SETERRQ(PetscObjectComm((PetscObject)A), PETSC_ERR_ARG_WRONG, "Mat(Nest) contains a null row.");
875:     }

877:     VecCreateNest(comm,bA->nr,bA->isglobal.row,L,left);
878:     for (i=0; i<bA->nr; i++) {
879:       VecDestroy(&L[i]);
880:     }

882:     PetscFree(L);
883:   }
884:   return(0);
885: }

887: static PetscErrorCode MatView_Nest(Mat A,PetscViewer viewer)
888: {
889:   Mat_Nest       *bA = (Mat_Nest*)A->data;
890:   PetscBool      isascii,viewSub = PETSC_FALSE;
891:   PetscInt       i,j;

895:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
896:   if (isascii) {

898:     PetscOptionsGetBool(((PetscObject)A)->options,((PetscObject)A)->prefix,"-mat_view_nest_sub",&viewSub,NULL);
899:     PetscViewerASCIIPrintf(viewer,"Matrix object: \n");
900:     PetscViewerASCIIPushTab(viewer);
901:     PetscViewerASCIIPrintf(viewer, "type=nest, rows=%D, cols=%D \n",bA->nr,bA->nc);

903:     PetscViewerASCIIPrintf(viewer,"MatNest structure: \n");
904:     for (i=0; i<bA->nr; i++) {
905:       for (j=0; j<bA->nc; j++) {
906:         MatType   type;
907:         char      name[256] = "",prefix[256] = "";
908:         PetscInt  NR,NC;
909:         PetscBool isNest = PETSC_FALSE;

911:         if (!bA->m[i][j]) {
912:           PetscViewerASCIIPrintf(viewer, "(%D,%D) : NULL \n",i,j);
913:           continue;
914:         }
915:         MatGetSize(bA->m[i][j],&NR,&NC);
916:         MatGetType(bA->m[i][j], &type);
917:         if (((PetscObject)bA->m[i][j])->name) {PetscSNPrintf(name,sizeof(name),"name=\"%s\", ",((PetscObject)bA->m[i][j])->name);}
918:         if (((PetscObject)bA->m[i][j])->prefix) {PetscSNPrintf(prefix,sizeof(prefix),"prefix=\"%s\", ",((PetscObject)bA->m[i][j])->prefix);}
919:         PetscObjectTypeCompare((PetscObject)bA->m[i][j],MATNEST,&isNest);

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

923:         if (isNest || viewSub) {
924:           PetscViewerASCIIPushTab(viewer);  /* push1 */
925:           MatView(bA->m[i][j],viewer);
926:           PetscViewerASCIIPopTab(viewer);    /* pop1 */
927:         }
928:       }
929:     }
930:     PetscViewerASCIIPopTab(viewer);    /* pop0 */
931:   }
932:   return(0);
933: }

935: static PetscErrorCode MatZeroEntries_Nest(Mat A)
936: {
937:   Mat_Nest       *bA = (Mat_Nest*)A->data;
938:   PetscInt       i,j;

942:   for (i=0; i<bA->nr; i++) {
943:     for (j=0; j<bA->nc; j++) {
944:       if (!bA->m[i][j]) continue;
945:       MatZeroEntries(bA->m[i][j]);
946:     }
947:   }
948:   return(0);
949: }

951: static PetscErrorCode MatCopy_Nest(Mat A,Mat B,MatStructure str)
952: {
953:   Mat_Nest       *bA = (Mat_Nest*)A->data,*bB = (Mat_Nest*)B->data;
954:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
956:   PetscBool      nnzstate = PETSC_FALSE;

959:   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);
960:   for (i=0; i<nr; i++) {
961:     for (j=0; j<nc; j++) {
962:       PetscObjectState subnnzstate = 0;
963:       if (bA->m[i][j] && bB->m[i][j]) {
964:         MatCopy(bA->m[i][j],bB->m[i][j],str);
965:       } 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);
966:       MatGetNonzeroState(bB->m[i][j],&subnnzstate);
967:       nnzstate = (PetscBool)(nnzstate || bB->nnzstate[i*nc+j] != subnnzstate);
968:       bB->nnzstate[i*nc+j] = subnnzstate;
969:     }
970:   }
971:   if (nnzstate) B->nonzerostate++;
972:   return(0);
973: }

975: static PetscErrorCode MatAXPY_Nest(Mat Y,PetscScalar a,Mat X,MatStructure str)
976: {
977:   Mat_Nest       *bY = (Mat_Nest*)Y->data,*bX = (Mat_Nest*)X->data;
978:   PetscInt       i,j,nr = bY->nr,nc = bY->nc;
980:   PetscBool      nnzstate = PETSC_FALSE;

983:   if (nr != bX->nr || nc != bX->nc) SETERRQ4(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_INCOMP,"Cannot AXPY a MatNest of block size (%D,%D) with a MatNest of block size (%D,%D)",bX->nr,bX->nc,nr,nc);
984:   for (i=0; i<nr; i++) {
985:     for (j=0; j<nc; j++) {
986:       PetscObjectState subnnzstate = 0;
987:       if (bY->m[i][j] && bX->m[i][j]) {
988:         MatAXPY(bY->m[i][j],a,bX->m[i][j],str);
989:       } else if (bX->m[i][j]) {
990:         Mat M;

992:         if (str != DIFFERENT_NONZERO_PATTERN) SETERRQ2(PetscObjectComm((PetscObject)Y),PETSC_ERR_ARG_INCOMP,"Matrix block does not exist at %D,%D. Use DIFFERENT_NONZERO_PATTERN",i,j);
993:         MatDuplicate(bX->m[i][j],MAT_COPY_VALUES,&M);
994:         MatNestSetSubMat(Y,i,j,M);
995:         MatDestroy(&M);
996:       }
997:       if (bY->m[i][j]) { MatGetNonzeroState(bY->m[i][j],&subnnzstate); }
998:       nnzstate = (PetscBool)(nnzstate || bY->nnzstate[i*nc+j] != subnnzstate);
999:       bY->nnzstate[i*nc+j] = subnnzstate;
1000:     }
1001:   }
1002:   if (nnzstate) Y->nonzerostate++;
1003:   return(0);
1004: }

1006: static PetscErrorCode MatDuplicate_Nest(Mat A,MatDuplicateOption op,Mat *B)
1007: {
1008:   Mat_Nest       *bA = (Mat_Nest*)A->data;
1009:   Mat            *b;
1010:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;

1014:   PetscMalloc1(nr*nc,&b);
1015:   for (i=0; i<nr; i++) {
1016:     for (j=0; j<nc; j++) {
1017:       if (bA->m[i][j]) {
1018:         MatDuplicate(bA->m[i][j],op,&b[i*nc+j]);
1019:       } else {
1020:         b[i*nc+j] = NULL;
1021:       }
1022:     }
1023:   }
1024:   MatCreateNest(PetscObjectComm((PetscObject)A),nr,bA->isglobal.row,nc,bA->isglobal.col,b,B);
1025:   /* Give the new MatNest exclusive ownership */
1026:   for (i=0; i<nr*nc; i++) {
1027:     MatDestroy(&b[i]);
1028:   }
1029:   PetscFree(b);

1031:   MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY);
1032:   MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY);
1033:   return(0);
1034: }

1036: /* nest api */
1037: PetscErrorCode MatNestGetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat *mat)
1038: {
1039:   Mat_Nest *bA = (Mat_Nest*)A->data;

1042:   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
1043:   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
1044:   *mat = bA->m[idxm][jdxm];
1045:   return(0);
1046: }

1048: /*@
1049:  MatNestGetSubMat - Returns a single, sub-matrix from a nest matrix.

1051:  Not collective

1053:  Input Parameters:
1054: +   A  - nest matrix
1055: .   idxm - index of the matrix within the nest matrix
1056: -   jdxm - index of the matrix within the nest matrix

1058:  Output Parameter:
1059: .   sub - matrix at index idxm,jdxm within the nest matrix

1061:  Level: developer

1063: .seealso: MatNestGetSize(), MatNestGetSubMats(), MatCreateNest(), MATNEST, MatNestSetSubMat(),
1064:           MatNestGetLocalISs(), MatNestGetISs()
1065: @*/
1066: PetscErrorCode  MatNestGetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat *sub)
1067: {

1071:   PetscUseMethod(A,"MatNestGetSubMat_C",(Mat,PetscInt,PetscInt,Mat*),(A,idxm,jdxm,sub));
1072:   return(0);
1073: }

1075: PetscErrorCode MatNestSetSubMat_Nest(Mat A,PetscInt idxm,PetscInt jdxm,Mat mat)
1076: {
1077:   Mat_Nest       *bA = (Mat_Nest*)A->data;
1078:   PetscInt       m,n,M,N,mi,ni,Mi,Ni;

1082:   if (idxm >= bA->nr) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Row too large: row %D max %D",idxm,bA->nr-1);
1083:   if (jdxm >= bA->nc) SETERRQ2(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Col too large: row %D max %D",jdxm,bA->nc-1);
1084:   MatGetLocalSize(mat,&m,&n);
1085:   MatGetSize(mat,&M,&N);
1086:   ISGetLocalSize(bA->isglobal.row[idxm],&mi);
1087:   ISGetSize(bA->isglobal.row[idxm],&Mi);
1088:   ISGetLocalSize(bA->isglobal.col[jdxm],&ni);
1089:   ISGetSize(bA->isglobal.col[jdxm],&Ni);
1090:   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);
1091:   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);

1093:   /* do not increase object state */
1094:   if (mat == bA->m[idxm][jdxm]) return(0);

1096:   PetscObjectReference((PetscObject)mat);
1097:   MatDestroy(&bA->m[idxm][jdxm]);
1098:   bA->m[idxm][jdxm] = mat;
1099:   PetscObjectStateIncrease((PetscObject)A);
1100:   MatGetNonzeroState(mat,&bA->nnzstate[idxm*bA->nc+jdxm]);
1101:   A->nonzerostate++;
1102:   return(0);
1103: }

1105: /*@
1106:  MatNestSetSubMat - Set a single submatrix in the nest matrix.

1108:  Logically collective on the submatrix communicator

1110:  Input Parameters:
1111: +   A  - nest matrix
1112: .   idxm - index of the matrix within the nest matrix
1113: .   jdxm - index of the matrix within the nest matrix
1114: -   sub - matrix at index idxm,jdxm within the nest matrix

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

1119:  This increments the reference count of the submatrix.

1121:  Level: developer

1123: .seealso: MatNestSetSubMats(), MatNestGetSubMats(), MatNestGetLocalISs(), MATNEST, MatCreateNest(),
1124:           MatNestGetSubMat(), MatNestGetISs(), MatNestGetSize()
1125: @*/
1126: PetscErrorCode  MatNestSetSubMat(Mat A,PetscInt idxm,PetscInt jdxm,Mat sub)
1127: {

1131:   PetscUseMethod(A,"MatNestSetSubMat_C",(Mat,PetscInt,PetscInt,Mat),(A,idxm,jdxm,sub));
1132:   return(0);
1133: }

1135: PetscErrorCode MatNestGetSubMats_Nest(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
1136: {
1137:   Mat_Nest *bA = (Mat_Nest*)A->data;

1140:   if (M)   *M   = bA->nr;
1141:   if (N)   *N   = bA->nc;
1142:   if (mat) *mat = bA->m;
1143:   return(0);
1144: }

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

1149:  Not collective

1151:  Input Parameters:
1152: .   A  - nest matrix

1154:  Output Parameter:
1155: +   M - number of rows in the nest matrix
1156: .   N - number of cols in the nest matrix
1157: -   mat - 2d array of matrices

1159:  Notes:

1161:  The user should not free the array mat.

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

1167:  Level: developer

1169: .seealso: MatNestGetSize(), MatNestGetSubMat(), MatNestGetLocalISs(), MATNEST, MatCreateNest(),
1170:           MatNestSetSubMats(), MatNestGetISs(), MatNestSetSubMat()
1171: @*/
1172: PetscErrorCode  MatNestGetSubMats(Mat A,PetscInt *M,PetscInt *N,Mat ***mat)
1173: {

1177:   PetscUseMethod(A,"MatNestGetSubMats_C",(Mat,PetscInt*,PetscInt*,Mat***),(A,M,N,mat));
1178:   return(0);
1179: }

1181: PetscErrorCode  MatNestGetSize_Nest(Mat A,PetscInt *M,PetscInt *N)
1182: {
1183:   Mat_Nest *bA = (Mat_Nest*)A->data;

1186:   if (M) *M = bA->nr;
1187:   if (N) *N = bA->nc;
1188:   return(0);
1189: }

1191: /*@
1192:  MatNestGetSize - Returns the size of the nest matrix.

1194:  Not collective

1196:  Input Parameters:
1197: .   A  - nest matrix

1199:  Output Parameter:
1200: +   M - number of rows in the nested mat
1201: -   N - number of cols in the nested mat

1203:  Notes:

1205:  Level: developer

1207: .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MATNEST, MatCreateNest(), MatNestGetLocalISs(),
1208:           MatNestGetISs()
1209: @*/
1210: PetscErrorCode  MatNestGetSize(Mat A,PetscInt *M,PetscInt *N)
1211: {

1215:   PetscUseMethod(A,"MatNestGetSize_C",(Mat,PetscInt*,PetscInt*),(A,M,N));
1216:   return(0);
1217: }

1219: static PetscErrorCode MatNestGetISs_Nest(Mat A,IS rows[],IS cols[])
1220: {
1221:   Mat_Nest *vs = (Mat_Nest*)A->data;
1222:   PetscInt i;

1225:   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->isglobal.row[i];
1226:   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->isglobal.col[i];
1227:   return(0);
1228: }

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

1233:  Not collective

1235:  Input Parameters:
1236: .   A  - nest matrix

1238:  Output Parameter:
1239: +   rows - array of row index sets
1240: -   cols - array of column index sets

1242:  Level: advanced

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

1247: .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetLocalISs(), MATNEST,
1248:           MatCreateNest(), MatNestGetSubMats(), MatNestSetSubMats()
1249: @*/
1250: PetscErrorCode  MatNestGetISs(Mat A,IS rows[],IS cols[])
1251: {

1256:   PetscUseMethod(A,"MatNestGetISs_C",(Mat,IS[],IS[]),(A,rows,cols));
1257:   return(0);
1258: }

1260: static PetscErrorCode MatNestGetLocalISs_Nest(Mat A,IS rows[],IS cols[])
1261: {
1262:   Mat_Nest *vs = (Mat_Nest*)A->data;
1263:   PetscInt i;

1266:   if (rows) for (i=0; i<vs->nr; i++) rows[i] = vs->islocal.row[i];
1267:   if (cols) for (i=0; i<vs->nc; i++) cols[i] = vs->islocal.col[i];
1268:   return(0);
1269: }

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

1274:  Not collective

1276:  Input Parameters:
1277: .   A  - nest matrix

1279:  Output Parameter:
1280: +   rows - array of row index sets (or NULL to ignore)
1281: -   cols - array of column index sets (or NULL to ignore)

1283:  Level: advanced

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

1288: .seealso: MatNestGetSubMat(), MatNestGetSubMats(), MatNestGetSize(), MatNestGetISs(), MatCreateNest(),
1289:           MATNEST, MatNestSetSubMats(), MatNestSetSubMat()
1290: @*/
1291: PetscErrorCode  MatNestGetLocalISs(Mat A,IS rows[],IS cols[])
1292: {

1297:   PetscUseMethod(A,"MatNestGetLocalISs_C",(Mat,IS[],IS[]),(A,rows,cols));
1298:   return(0);
1299: }

1301: PetscErrorCode  MatNestSetVecType_Nest(Mat A,VecType vtype)
1302: {
1304:   PetscBool      flg;

1307:   PetscStrcmp(vtype,VECNEST,&flg);
1308:   /* In reality, this only distinguishes VECNEST and "other" */
1309:   if (flg) A->ops->getvecs = MatCreateVecs_Nest;
1310:   else A->ops->getvecs = (PetscErrorCode (*)(Mat,Vec*,Vec*)) 0;
1311:   return(0);
1312: }

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

1317:  Not collective

1319:  Input Parameters:
1320: +  A  - nest matrix
1321: -  vtype - type to use for creating vectors

1323:  Notes:

1325:  Level: developer

1327: .seealso: MatCreateVecs(), MATNEST, MatCreateNest()
1328: @*/
1329: PetscErrorCode  MatNestSetVecType(Mat A,VecType vtype)
1330: {

1334:   PetscTryMethod(A,"MatNestSetVecType_C",(Mat,VecType),(A,vtype));
1335:   return(0);
1336: }

1338: PetscErrorCode MatNestSetSubMats_Nest(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1339: {
1340:   Mat_Nest       *s = (Mat_Nest*)A->data;
1341:   PetscInt       i,j,m,n,M,N;
1343:   PetscBool      cong;

1346:   MatReset_Nest(A);

1348:   s->nr = nr;
1349:   s->nc = nc;

1351:   /* Create space for submatrices */
1352:   PetscMalloc1(nr,&s->m);
1353:   for (i=0; i<nr; i++) {
1354:     PetscMalloc1(nc,&s->m[i]);
1355:   }
1356:   for (i=0; i<nr; i++) {
1357:     for (j=0; j<nc; j++) {
1358:       s->m[i][j] = a[i*nc+j];
1359:       if (a[i*nc+j]) {
1360:         PetscObjectReference((PetscObject)a[i*nc+j]);
1361:       }
1362:     }
1363:   }

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

1367:   PetscMalloc1(nr,&s->row_len);
1368:   PetscMalloc1(nc,&s->col_len);
1369:   for (i=0; i<nr; i++) s->row_len[i]=-1;
1370:   for (j=0; j<nc; j++) s->col_len[j]=-1;

1372:   PetscCalloc1(nr*nc,&s->nnzstate);
1373:   for (i=0; i<nr; i++) {
1374:     for (j=0; j<nc; j++) {
1375:       if (s->m[i][j]) {
1376:         MatGetNonzeroState(s->m[i][j],&s->nnzstate[i*nc+j]);
1377:       }
1378:     }
1379:   }

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

1383:   PetscLayoutSetSize(A->rmap,M);
1384:   PetscLayoutSetLocalSize(A->rmap,m);
1385:   PetscLayoutSetSize(A->cmap,N);
1386:   PetscLayoutSetLocalSize(A->cmap,n);

1388:   PetscLayoutSetUp(A->rmap);
1389:   PetscLayoutSetUp(A->cmap);

1391:   /* disable operations that are not supported for non-square matrices,
1392:      or matrices for which is_row != is_col  */
1393:   MatHasCongruentLayouts(A,&cong);
1394:   if (cong && nr != nc) cong = PETSC_FALSE;
1395:   if (cong) {
1396:     for (i = 0; cong && i < nr; i++) {
1397:       ISEqualUnsorted(s->isglobal.row[i],s->isglobal.col[i],&cong);
1398:     }
1399:   }
1400:   if (!cong) {
1401:     A->ops->missingdiagonal = NULL;
1402:     A->ops->getdiagonal     = NULL;
1403:     A->ops->shift           = NULL;
1404:     A->ops->diagonalset     = NULL;
1405:   }

1407:   PetscCalloc2(nr,&s->left,nc,&s->right);
1408:   PetscObjectStateIncrease((PetscObject)A);
1409:   A->nonzerostate++;
1410:   return(0);
1411: }

1413: /*@
1414:    MatNestSetSubMats - Sets the nested submatrices

1416:    Collective on Mat

1418:    Input Parameter:
1419: +  A - nested matrix
1420: .  nr - number of nested row blocks
1421: .  is_row - index sets for each nested row block, or NULL to make contiguous
1422: .  nc - number of nested column blocks
1423: .  is_col - index sets for each nested column block, or NULL to make contiguous
1424: -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL

1426:    Notes: this always resets any submatrix information previously set

1428:    Level: advanced

1430: .seealso: MatCreateNest(), MATNEST, MatNestSetSubMat(), MatNestGetSubMat(), MatNestGetSubMats()
1431: @*/
1432: PetscErrorCode MatNestSetSubMats(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[])
1433: {
1435:   PetscInt       i;

1439:   if (nr < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of rows cannot be negative");
1440:   if (nr && is_row) {
1443:   }
1444:   if (nc < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_OUTOFRANGE,"Number of columns cannot be negative");
1445:   if (nc && is_col) {
1448:   }
1450:   PetscUseMethod(A,"MatNestSetSubMats_C",(Mat,PetscInt,const IS[],PetscInt,const IS[],const Mat[]),(A,nr,is_row,nc,is_col,a));
1451:   return(0);
1452: }

1454: static PetscErrorCode MatNestCreateAggregateL2G_Private(Mat A,PetscInt n,const IS islocal[],const IS isglobal[],PetscBool colflg,ISLocalToGlobalMapping *ltog)
1455: {
1457:   PetscBool      flg;
1458:   PetscInt       i,j,m,mi,*ix;

1461:   *ltog = NULL;
1462:   for (i=0,m=0,flg=PETSC_FALSE; i<n; i++) {
1463:     if (islocal[i]) {
1464:       ISGetLocalSize(islocal[i],&mi);
1465:       flg  = PETSC_TRUE;      /* We found a non-trivial entry */
1466:     } else {
1467:       ISGetLocalSize(isglobal[i],&mi);
1468:     }
1469:     m += mi;
1470:   }
1471:   if (!flg) return(0);

1473:   PetscMalloc1(m,&ix);
1474:   for (i=0,m=0; i<n; i++) {
1475:     ISLocalToGlobalMapping smap = NULL;
1476:     Mat                    sub = NULL;
1477:     PetscSF                sf;
1478:     PetscLayout            map;
1479:     const PetscInt         *ix2;

1481:     if (!colflg) {
1482:       MatNestFindNonzeroSubMatRow(A,i,&sub);
1483:     } else {
1484:       MatNestFindNonzeroSubMatCol(A,i,&sub);
1485:     }
1486:     if (sub) {
1487:       if (!colflg) {
1488:         MatGetLocalToGlobalMapping(sub,&smap,NULL);
1489:       } else {
1490:         MatGetLocalToGlobalMapping(sub,NULL,&smap);
1491:       }
1492:     }
1493:     /*
1494:        Now we need to extract the monolithic global indices that correspond to the given split global indices.
1495:        In many/most cases, we only want MatGetLocalSubMatrix() to work, in which case we only need to know the size of the local spaces.
1496:     */
1497:     ISGetIndices(isglobal[i],&ix2);
1498:     if (islocal[i]) {
1499:       PetscInt *ilocal,*iremote;
1500:       PetscInt mil,nleaves;

1502:       ISGetLocalSize(islocal[i],&mi);
1503:       if (!smap) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"Missing local to global map");
1504:       for (j=0; j<mi; j++) ix[m+j] = j;
1505:       ISLocalToGlobalMappingApply(smap,mi,ix+m,ix+m);

1507:       /* PetscSFSetGraphLayout does not like negative indices */
1508:       PetscMalloc2(mi,&ilocal,mi,&iremote);
1509:       for (j=0, nleaves = 0; j<mi; j++) {
1510:         if (ix[m+j] < 0) continue;
1511:         ilocal[nleaves]  = j;
1512:         iremote[nleaves] = ix[m+j];
1513:         nleaves++;
1514:       }
1515:       ISGetLocalSize(isglobal[i],&mil);
1516:       PetscSFCreate(PetscObjectComm((PetscObject)A),&sf);
1517:       PetscLayoutCreate(PetscObjectComm((PetscObject)A),&map);
1518:       PetscLayoutSetLocalSize(map,mil);
1519:       PetscLayoutSetUp(map);
1520:       PetscSFSetGraphLayout(sf,map,nleaves,ilocal,PETSC_USE_POINTER,iremote);
1521:       PetscLayoutDestroy(&map);
1522:       PetscSFBcastBegin(sf,MPIU_INT,ix2,ix + m);
1523:       PetscSFBcastEnd(sf,MPIU_INT,ix2,ix + m);
1524:       PetscSFDestroy(&sf);
1525:       PetscFree2(ilocal,iremote);
1526:     } else {
1527:       ISGetLocalSize(isglobal[i],&mi);
1528:       for (j=0; j<mi; j++) ix[m+j] = ix2[i];
1529:     }
1530:     ISRestoreIndices(isglobal[i],&ix2);
1531:     m   += mi;
1532:   }
1533:   ISLocalToGlobalMappingCreate(PetscObjectComm((PetscObject)A),1,m,ix,PETSC_OWN_POINTER,ltog);
1534:   return(0);
1535: }


1538: /* If an IS was provided, there is nothing Nest needs to do, otherwise Nest will build a strided IS */
1539: /*
1540:   nprocessors = NP
1541:   Nest x^T = ((g_0,g_1,...g_nprocs-1), (h_0,h_1,...h_NP-1))
1542:        proc 0: => (g_0,h_0,)
1543:        proc 1: => (g_1,h_1,)
1544:        ...
1545:        proc nprocs-1: => (g_NP-1,h_NP-1,)

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

1550:             proc 0:
1551:     is[1] = (nlocal(g_0),nlocal(g_0)+1,...,nlocal(g_0)+nlocal(h_0)-1)
1552:             proc 1:
1553:     is[1] = (nlocal(g_1),nlocal(g_1)+1,...,nlocal(g_1)+nlocal(h_1)-1)

1555:             proc NP-1:
1556:     is[1] = (nlocal(g_NP-1),nlocal(g_NP-1)+1,...,nlocal(g_NP-1)+nlocal(h_NP-1)-1)
1557: */
1558: static PetscErrorCode MatSetUp_NestIS_Private(Mat A,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[])
1559: {
1560:   Mat_Nest       *vs = (Mat_Nest*)A->data;
1561:   PetscInt       i,j,offset,n,nsum,bs;
1563:   Mat            sub = NULL;

1566:   PetscMalloc1(nr,&vs->isglobal.row);
1567:   PetscMalloc1(nc,&vs->isglobal.col);
1568:   if (is_row) { /* valid IS is passed in */
1569:     /* refs on is[] are incremeneted */
1570:     for (i=0; i<vs->nr; i++) {
1571:       PetscObjectReference((PetscObject)is_row[i]);

1573:       vs->isglobal.row[i] = is_row[i];
1574:     }
1575:   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each row */
1576:     nsum = 0;
1577:     for (i=0; i<vs->nr; i++) {  /* Add up the local sizes to compute the aggregate offset */
1578:       MatNestFindNonzeroSubMatRow(A,i,&sub);
1579:       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in row %D",i);
1580:       MatGetLocalSize(sub,&n,NULL);
1581:       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
1582:       nsum += n;
1583:     }
1584:     MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));
1585:     offset -= nsum;
1586:     for (i=0; i<vs->nr; i++) {
1587:       MatNestFindNonzeroSubMatRow(A,i,&sub);
1588:       MatGetLocalSize(sub,&n,NULL);
1589:       MatGetBlockSize(sub,&bs);
1590:       ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.row[i]);
1591:       ISSetBlockSize(vs->isglobal.row[i],bs);
1592:       offset += n;
1593:     }
1594:   }

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

1601:       vs->isglobal.col[j] = is_col[j];
1602:     }
1603:   } else {                      /* Create the ISs by inspecting sizes of a submatrix in each column */
1604:     offset = A->cmap->rstart;
1605:     nsum   = 0;
1606:     for (j=0; j<vs->nc; j++) {
1607:       MatNestFindNonzeroSubMatCol(A,j,&sub);
1608:       if (!sub) SETERRQ1(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONG,"No nonzero submatrix in column %D",i);
1609:       MatGetLocalSize(sub,NULL,&n);
1610:       if (n < 0) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_ARG_WRONGSTATE,"Sizes have not yet been set for submatrix");
1611:       nsum += n;
1612:     }
1613:     MPI_Scan(&nsum,&offset,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)A));
1614:     offset -= nsum;
1615:     for (j=0; j<vs->nc; j++) {
1616:       MatNestFindNonzeroSubMatCol(A,j,&sub);
1617:       MatGetLocalSize(sub,NULL,&n);
1618:       MatGetBlockSize(sub,&bs);
1619:       ISCreateStride(PetscObjectComm((PetscObject)sub),n,offset,1,&vs->isglobal.col[j]);
1620:       ISSetBlockSize(vs->isglobal.col[j],bs);
1621:       offset += n;
1622:     }
1623:   }

1625:   /* Set up the local ISs */
1626:   PetscMalloc1(vs->nr,&vs->islocal.row);
1627:   PetscMalloc1(vs->nc,&vs->islocal.col);
1628:   for (i=0,offset=0; i<vs->nr; i++) {
1629:     IS                     isloc;
1630:     ISLocalToGlobalMapping rmap = NULL;
1631:     PetscInt               nlocal,bs;
1632:     MatNestFindNonzeroSubMatRow(A,i,&sub);
1633:     if (sub) {MatGetLocalToGlobalMapping(sub,&rmap,NULL);}
1634:     if (rmap) {
1635:       MatGetBlockSize(sub,&bs);
1636:       ISLocalToGlobalMappingGetSize(rmap,&nlocal);
1637:       ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);
1638:       ISSetBlockSize(isloc,bs);
1639:     } else {
1640:       nlocal = 0;
1641:       isloc  = NULL;
1642:     }
1643:     vs->islocal.row[i] = isloc;
1644:     offset            += nlocal;
1645:   }
1646:   for (i=0,offset=0; i<vs->nc; i++) {
1647:     IS                     isloc;
1648:     ISLocalToGlobalMapping cmap = NULL;
1649:     PetscInt               nlocal,bs;
1650:     MatNestFindNonzeroSubMatCol(A,i,&sub);
1651:     if (sub) {MatGetLocalToGlobalMapping(sub,NULL,&cmap);}
1652:     if (cmap) {
1653:       MatGetBlockSize(sub,&bs);
1654:       ISLocalToGlobalMappingGetSize(cmap,&nlocal);
1655:       ISCreateStride(PETSC_COMM_SELF,nlocal,offset,1,&isloc);
1656:       ISSetBlockSize(isloc,bs);
1657:     } else {
1658:       nlocal = 0;
1659:       isloc  = NULL;
1660:     }
1661:     vs->islocal.col[i] = isloc;
1662:     offset            += nlocal;
1663:   }

1665:   /* Set up the aggregate ISLocalToGlobalMapping */
1666:   {
1667:     ISLocalToGlobalMapping rmap,cmap;
1668:     MatNestCreateAggregateL2G_Private(A,vs->nr,vs->islocal.row,vs->isglobal.row,PETSC_FALSE,&rmap);
1669:     MatNestCreateAggregateL2G_Private(A,vs->nc,vs->islocal.col,vs->isglobal.col,PETSC_TRUE,&cmap);
1670:     if (rmap && cmap) {MatSetLocalToGlobalMapping(A,rmap,cmap);}
1671:     ISLocalToGlobalMappingDestroy(&rmap);
1672:     ISLocalToGlobalMappingDestroy(&cmap);
1673:   }

1675: #if defined(PETSC_USE_DEBUG)
1676:   for (i=0; i<vs->nr; i++) {
1677:     for (j=0; j<vs->nc; j++) {
1678:       PetscInt m,n,M,N,mi,ni,Mi,Ni;
1679:       Mat      B = vs->m[i][j];
1680:       if (!B) continue;
1681:       MatGetSize(B,&M,&N);
1682:       MatGetLocalSize(B,&m,&n);
1683:       ISGetSize(vs->isglobal.row[i],&Mi);
1684:       ISGetSize(vs->isglobal.col[j],&Ni);
1685:       ISGetLocalSize(vs->isglobal.row[i],&mi);
1686:       ISGetLocalSize(vs->isglobal.col[j],&ni);
1687:       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);
1688:       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);
1689:     }
1690:   }
1691: #endif

1693:   /* Set A->assembled if all non-null blocks are currently assembled */
1694:   for (i=0; i<vs->nr; i++) {
1695:     for (j=0; j<vs->nc; j++) {
1696:       if (vs->m[i][j] && !vs->m[i][j]->assembled) return(0);
1697:     }
1698:   }
1699:   A->assembled = PETSC_TRUE;
1700:   return(0);
1701: }

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

1706:    Collective on Mat

1708:    Input Parameter:
1709: +  comm - Communicator for the new Mat
1710: .  nr - number of nested row blocks
1711: .  is_row - index sets for each nested row block, or NULL to make contiguous
1712: .  nc - number of nested column blocks
1713: .  is_col - index sets for each nested column block, or NULL to make contiguous
1714: -  a - row-aligned array of nr*nc submatrices, empty submatrices can be passed using NULL

1716:    Output Parameter:
1717: .  B - new matrix

1719:    Level: advanced

1721: .seealso: MatCreate(), VecCreateNest(), DMCreateMatrix(), MATNEST, MatNestSetSubMat(),
1722:           MatNestGetSubMat(), MatNestGetLocalISs(), MatNestGetSize(),
1723:           MatNestGetISs(), MatNestSetSubMats(), MatNestGetSubMats()
1724: @*/
1725: PetscErrorCode MatCreateNest(MPI_Comm comm,PetscInt nr,const IS is_row[],PetscInt nc,const IS is_col[],const Mat a[],Mat *B)
1726: {
1727:   Mat            A;

1731:   *B   = 0;
1732:   MatCreate(comm,&A);
1733:   MatSetType(A,MATNEST);
1734:   A->preallocated = PETSC_TRUE;
1735:   MatNestSetSubMats(A,nr,is_row,nc,is_col,a);
1736:   *B   = A;
1737:   return(0);
1738: }

1740: static PetscErrorCode MatConvert_Nest_SeqAIJ_fast(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1741: {
1742:   Mat_Nest       *nest = (Mat_Nest*)A->data;
1743:   Mat            *trans;
1744:   PetscScalar    **avv;
1745:   PetscScalar    *vv;
1746:   PetscInt       **aii,**ajj;
1747:   PetscInt       *ii,*jj,*ci;
1748:   PetscInt       nr,nc,nnz,i,j;
1749:   PetscBool      done;

1753:   MatGetSize(A,&nr,&nc);
1754:   if (reuse == MAT_REUSE_MATRIX) {
1755:     PetscInt rnr;

1757:     MatGetRowIJ(*newmat,0,PETSC_FALSE,PETSC_FALSE,&rnr,(const PetscInt**)&ii,(const PetscInt**)&jj,&done);
1758:     if (!done) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_PLIB,"MatGetRowIJ");
1759:     if (rnr != nr) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Cannot reuse matrix, wrong number of rows");
1760:     MatSeqAIJGetArray(*newmat,&vv);
1761:   }
1762:   /* extract CSR for nested SeqAIJ matrices */
1763:   nnz  = 0;
1764:   PetscCalloc4(nest->nr*nest->nc,&aii,nest->nr*nest->nc,&ajj,nest->nr*nest->nc,&avv,nest->nr*nest->nc,&trans);
1765:   for (i=0; i<nest->nr; ++i) {
1766:     for (j=0; j<nest->nc; ++j) {
1767:       Mat B = nest->m[i][j];
1768:       if (B) {
1769:         PetscScalar *naa;
1770:         PetscInt    *nii,*njj,nnr;
1771:         PetscBool   istrans;

1773:         PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&istrans);
1774:         if (istrans) {
1775:           Mat Bt;

1777:           MatTransposeGetMat(B,&Bt);
1778:           MatTranspose(Bt,MAT_INITIAL_MATRIX,&trans[i*nest->nc+j]);
1779:           B    = trans[i*nest->nc+j];
1780:         }
1781:         MatGetRowIJ(B,0,PETSC_FALSE,PETSC_FALSE,&nnr,(const PetscInt**)&nii,(const PetscInt**)&njj,&done);
1782:         if (!done) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"MatGetRowIJ");
1783:         MatSeqAIJGetArray(B,&naa);
1784:         nnz += nii[nnr];

1786:         aii[i*nest->nc+j] = nii;
1787:         ajj[i*nest->nc+j] = njj;
1788:         avv[i*nest->nc+j] = naa;
1789:       }
1790:     }
1791:   }
1792:   if (reuse != MAT_REUSE_MATRIX) {
1793:     PetscMalloc1(nr+1,&ii);
1794:     PetscMalloc1(nnz,&jj);
1795:     PetscMalloc1(nnz,&vv);
1796:   } else {
1797:     if (nnz != ii[nr]) SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_USER,"Cannot reuse matrix, wrong number of nonzeros");
1798:   }

1800:   /* new row pointer */
1801:   PetscArrayzero(ii,nr+1);
1802:   for (i=0; i<nest->nr; ++i) {
1803:     PetscInt       ncr,rst;

1805:     ISStrideGetInfo(nest->isglobal.row[i],&rst,NULL);
1806:     ISGetLocalSize(nest->isglobal.row[i],&ncr);
1807:     for (j=0; j<nest->nc; ++j) {
1808:       if (aii[i*nest->nc+j]) {
1809:         PetscInt    *nii = aii[i*nest->nc+j];
1810:         PetscInt    ir;

1812:         for (ir=rst; ir<ncr+rst; ++ir) {
1813:           ii[ir+1] += nii[1]-nii[0];
1814:           nii++;
1815:         }
1816:       }
1817:     }
1818:   }
1819:   for (i=0; i<nr; i++) ii[i+1] += ii[i];

1821:   /* construct CSR for the new matrix */
1822:   PetscCalloc1(nr,&ci);
1823:   for (i=0; i<nest->nr; ++i) {
1824:     PetscInt       ncr,rst;

1826:     ISStrideGetInfo(nest->isglobal.row[i],&rst,NULL);
1827:     ISGetLocalSize(nest->isglobal.row[i],&ncr);
1828:     for (j=0; j<nest->nc; ++j) {
1829:       if (aii[i*nest->nc+j]) {
1830:         PetscScalar *nvv = avv[i*nest->nc+j];
1831:         PetscInt    *nii = aii[i*nest->nc+j];
1832:         PetscInt    *njj = ajj[i*nest->nc+j];
1833:         PetscInt    ir,cst;

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

1839:           for (ij=0;ij<rsize;ij++) {
1840:             jj[ist+ij] = *njj+cst;
1841:             vv[ist+ij] = *nvv;
1842:             njj++;
1843:             nvv++;
1844:           }
1845:           ci[ir] += rsize;
1846:           nii++;
1847:         }
1848:       }
1849:     }
1850:   }
1851:   PetscFree(ci);

1853:   /* restore info */
1854:   for (i=0; i<nest->nr; ++i) {
1855:     for (j=0; j<nest->nc; ++j) {
1856:       Mat B = nest->m[i][j];
1857:       if (B) {
1858:         PetscInt nnr = 0, k = i*nest->nc+j;

1860:         B    = (trans[k] ? trans[k] : B);
1861:         MatRestoreRowIJ(B,0,PETSC_FALSE,PETSC_FALSE,&nnr,(const PetscInt**)&aii[k],(const PetscInt**)&ajj[k],&done);
1862:         if (!done) SETERRQ(PetscObjectComm((PetscObject)B),PETSC_ERR_PLIB,"MatRestoreRowIJ");
1863:         MatSeqAIJRestoreArray(B,&avv[k]);
1864:         MatDestroy(&trans[k]);
1865:       }
1866:     }
1867:   }
1868:   PetscFree4(aii,ajj,avv,trans);

1870:   /* finalize newmat */
1871:   if (reuse == MAT_INITIAL_MATRIX) {
1872:     MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),nr,nc,ii,jj,vv,newmat);
1873:   } else if (reuse == MAT_INPLACE_MATRIX) {
1874:     Mat B;

1876:     MatCreateSeqAIJWithArrays(PetscObjectComm((PetscObject)A),nr,nc,ii,jj,vv,&B);
1877:     MatHeaderReplace(A,&B);
1878:   }
1879:   MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);
1880:   MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);
1881:   {
1882:     Mat_SeqAIJ *a = (Mat_SeqAIJ*)((*newmat)->data);
1883:     a->free_a     = PETSC_TRUE;
1884:     a->free_ij    = PETSC_TRUE;
1885:   }
1886:   return(0);
1887: }

1889: PETSC_INTERN PetscErrorCode MatConvert_Nest_AIJ(Mat A,MatType newtype,MatReuse reuse,Mat *newmat)
1890: {
1892:   Mat_Nest       *nest = (Mat_Nest*)A->data;
1893:   PetscInt       m,n,M,N,i,j,k,*dnnz,*onnz,rstart;
1894:   PetscInt       cstart,cend;
1895:   PetscMPIInt    size;
1896:   Mat            C;

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

1904:     PetscStrcmp(newtype,MATAIJ,&fast);
1905:     if (!fast) {
1906:       PetscStrcmp(newtype,MATSEQAIJ,&fast);
1907:     }
1908:     for (i=0; i<nest->nr && fast; ++i) {
1909:       for (j=0; j<nest->nc && fast; ++j) {
1910:         Mat B = nest->m[i][j];
1911:         if (B) {
1912:           PetscObjectTypeCompare((PetscObject)B,MATSEQAIJ,&fast);
1913:           if (!fast) {
1914:             PetscBool istrans;

1916:             PetscObjectTypeCompare((PetscObject)B,MATTRANSPOSEMAT,&istrans);
1917:             if (istrans) {
1918:               Mat Bt;

1920:               MatTransposeGetMat(B,&Bt);
1921:               PetscObjectTypeCompare((PetscObject)Bt,MATSEQAIJ,&fast);
1922:             }
1923:           }
1924:         }
1925:       }
1926:     }
1927:     for (i=0, nf=0; i<nest->nr && fast; ++i) {
1928:       PetscObjectTypeCompare((PetscObject)nest->isglobal.row[i],ISSTRIDE,&fast);
1929:       if (fast) {
1930:         PetscInt f,s;

1932:         ISStrideGetInfo(nest->isglobal.row[i],&f,&s);
1933:         if (f != nf || s != 1) { fast = PETSC_FALSE; }
1934:         else {
1935:           ISGetSize(nest->isglobal.row[i],&f);
1936:           nf  += f;
1937:         }
1938:       }
1939:     }
1940:     for (i=0, nf=0; i<nest->nc && fast; ++i) {
1941:       PetscObjectTypeCompare((PetscObject)nest->isglobal.col[i],ISSTRIDE,&fast);
1942:       if (fast) {
1943:         PetscInt f,s;

1945:         ISStrideGetInfo(nest->isglobal.col[i],&f,&s);
1946:         if (f != nf || s != 1) { fast = PETSC_FALSE; }
1947:         else {
1948:           ISGetSize(nest->isglobal.col[i],&f);
1949:           nf  += f;
1950:         }
1951:       }
1952:     }
1953:     if (fast) {
1954:       MatConvert_Nest_SeqAIJ_fast(A,newtype,reuse,newmat);
1955:       return(0);
1956:     }
1957:   }
1958:   MatGetSize(A,&M,&N);
1959:   MatGetLocalSize(A,&m,&n);
1960:   MatGetOwnershipRangeColumn(A,&cstart,&cend);
1961:   switch (reuse) {
1962:   case MAT_INITIAL_MATRIX:
1963:     MatCreate(PetscObjectComm((PetscObject)A),&C);
1964:     MatSetType(C,newtype);
1965:     MatSetSizes(C,m,n,M,N);
1966:     *newmat = C;
1967:     break;
1968:   case MAT_REUSE_MATRIX:
1969:     C = *newmat;
1970:     break;
1971:   default: SETERRQ(PetscObjectComm((PetscObject)A),PETSC_ERR_SUP,"MatReuse");
1972:   }
1973:   PetscMalloc1(2*m,&dnnz);
1974:   onnz = dnnz + m;
1975:   for (k=0; k<m; k++) {
1976:     dnnz[k] = 0;
1977:     onnz[k] = 0;
1978:   }
1979:   for (j=0; j<nest->nc; ++j) {
1980:     IS             bNis;
1981:     PetscInt       bN;
1982:     const PetscInt *bNindices;
1983:     /* Using global column indices and ISAllGather() is not scalable. */
1984:     ISAllGather(nest->isglobal.col[j], &bNis);
1985:     ISGetSize(bNis, &bN);
1986:     ISGetIndices(bNis,&bNindices);
1987:     for (i=0; i<nest->nr; ++i) {
1988:       PetscSF        bmsf;
1989:       PetscSFNode    *iremote;
1990:       Mat            B;
1991:       PetscInt       bm, *sub_dnnz,*sub_onnz, br;
1992:       const PetscInt *bmindices;
1993:       B = nest->m[i][j];
1994:       if (!B) continue;
1995:       ISGetLocalSize(nest->isglobal.row[i],&bm);
1996:       ISGetIndices(nest->isglobal.row[i],&bmindices);
1997:       PetscSFCreate(PetscObjectComm((PetscObject)A), &bmsf);
1998:       PetscMalloc1(bm,&iremote);
1999:       PetscMalloc1(bm,&sub_dnnz);
2000:       PetscMalloc1(bm,&sub_onnz);
2001:       for (k = 0; k < bm; ++k){
2002:             sub_dnnz[k] = 0;
2003:             sub_onnz[k] = 0;
2004:       }
2005:       /*
2006:        Locate the owners for all of the locally-owned global row indices for this row block.
2007:        These determine the roots of PetscSF used to communicate preallocation data to row owners.
2008:        The roots correspond to the dnnz and onnz entries; thus, there are two roots per row.
2009:        */
2010:       MatGetOwnershipRange(B,&rstart,NULL);
2011:       for (br = 0; br < bm; ++br) {
2012:         PetscInt       row = bmindices[br], brncols, col;
2013:         const PetscInt *brcols;
2014:         PetscInt       rowrel = 0; /* row's relative index on its owner rank */
2015:         PetscMPIInt    rowowner = 0;
2016:         PetscLayoutFindOwnerIndex(A->rmap,row,&rowowner,&rowrel);
2017:         /* how many roots  */
2018:         iremote[br].rank = rowowner; iremote[br].index = rowrel;           /* edge from bmdnnz to dnnz */
2019:         /* get nonzero pattern */
2020:         MatGetRow(B,br+rstart,&brncols,&brcols,NULL);
2021:         for (k=0; k<brncols; k++) {
2022:           col  = bNindices[brcols[k]];
2023:           if (col>=A->cmap->range[rowowner] && col<A->cmap->range[rowowner+1]) {
2024:             sub_dnnz[br]++;
2025:           } else {
2026:             sub_onnz[br]++;
2027:           }
2028:         }
2029:         MatRestoreRow(B,br+rstart,&brncols,&brcols,NULL);
2030:       }
2031:       ISRestoreIndices(nest->isglobal.row[i],&bmindices);
2032:       /* bsf will have to take care of disposing of bedges. */
2033:       PetscSFSetGraph(bmsf,m,bm,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);
2034:       PetscSFReduceBegin(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);
2035:       PetscSFReduceEnd(bmsf,MPIU_INT,sub_dnnz,dnnz,MPI_SUM);
2036:       PetscSFReduceBegin(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);
2037:       PetscSFReduceEnd(bmsf,MPIU_INT,sub_onnz,onnz,MPI_SUM);
2038:       PetscFree(sub_dnnz);
2039:       PetscFree(sub_onnz);
2040:       PetscSFDestroy(&bmsf);
2041:     }
2042:     ISRestoreIndices(bNis,&bNindices);
2043:     ISDestroy(&bNis);
2044:   }
2045:   /* Resize preallocation if overestimated */
2046:   for (i=0;i<m;i++) {
2047:     dnnz[i] = PetscMin(dnnz[i],A->cmap->n);
2048:     onnz[i] = PetscMin(onnz[i],A->cmap->N - A->cmap->n);
2049:   }
2050:   MatSeqAIJSetPreallocation(C,0,dnnz);
2051:   MatMPIAIJSetPreallocation(C,0,dnnz,0,onnz);
2052:   PetscFree(dnnz);

2054:   /* Fill by row */
2055:   for (j=0; j<nest->nc; ++j) {
2056:     /* Using global column indices and ISAllGather() is not scalable. */
2057:     IS             bNis;
2058:     PetscInt       bN;
2059:     const PetscInt *bNindices;
2060:     ISAllGather(nest->isglobal.col[j], &bNis);
2061:     ISGetSize(bNis,&bN);
2062:     ISGetIndices(bNis,&bNindices);
2063:     for (i=0; i<nest->nr; ++i) {
2064:       Mat            B;
2065:       PetscInt       bm, br;
2066:       const PetscInt *bmindices;
2067:       B = nest->m[i][j];
2068:       if (!B) continue;
2069:       ISGetLocalSize(nest->isglobal.row[i],&bm);
2070:       ISGetIndices(nest->isglobal.row[i],&bmindices);
2071:       MatGetOwnershipRange(B,&rstart,NULL);
2072:       for (br = 0; br < bm; ++br) {
2073:         PetscInt          row = bmindices[br], brncols,  *cols;
2074:         const PetscInt    *brcols;
2075:         const PetscScalar *brcoldata;
2076:         MatGetRow(B,br+rstart,&brncols,&brcols,&brcoldata);
2077:         PetscMalloc1(brncols,&cols);
2078:         for (k=0; k<brncols; k++) cols[k] = bNindices[brcols[k]];
2079:         /*
2080:           Nest blocks are required to be nonoverlapping -- otherwise nest and monolithic index layouts wouldn't match.
2081:           Thus, we could use INSERT_VALUES, but I prefer ADD_VALUES.
2082:          */
2083:         MatSetValues(C,1,&row,brncols,cols,brcoldata,ADD_VALUES);
2084:         MatRestoreRow(B,br+rstart,&brncols,&brcols,&brcoldata);
2085:         PetscFree(cols);
2086:       }
2087:       ISRestoreIndices(nest->isglobal.row[i],&bmindices);
2088:     }
2089:     ISRestoreIndices(bNis,&bNindices);
2090:     ISDestroy(&bNis);
2091:   }
2092:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
2093:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
2094:   return(0);
2095: }

2097: PetscErrorCode MatHasOperation_Nest(Mat mat,MatOperation op,PetscBool *has)
2098: {
2099:   Mat_Nest       *bA = (Mat_Nest*)mat->data;
2100:   PetscInt       i,j,nr = bA->nr,nc = bA->nc;
2101:   PetscBool      flg;

2105:   *has = PETSC_FALSE;
2106:   if (op == MATOP_MULT_TRANSPOSE || op == MATOP_MAT_MULT) {
2107:     for (j=0; j<nc; j++) {
2108:       for (i=0; i<nr; i++) {
2109:         if (!bA->m[i][j]) continue;
2110:         MatHasOperation(bA->m[i][j],op,&flg);
2111:         if (!flg) return(0);
2112:       }
2113:     }
2114:   }
2115:   if (((void**)mat->ops)[op] || (op == MATOP_MAT_MULT && flg)) *has = PETSC_TRUE;
2116:   return(0);
2117: }

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

2122:   Level: intermediate

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

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

2133: .seealso: MatCreate(), MatType, MatCreateNest(), MatNestSetSubMat(), MatNestGetSubMat(),
2134:           VecCreateNest(), DMCreateMatrix(), DMCOMPOSITE, MatNestSetVecType(), MatNestGetLocalISs(),
2135:           MatNestGetISs(), MatNestSetSubMats(), MatNestGetSubMats()
2136: M*/
2137: PETSC_EXTERN PetscErrorCode MatCreate_Nest(Mat A)
2138: {
2139:   Mat_Nest       *s;

2143:   PetscNewLog(A,&s);
2144:   A->data = (void*)s;

2146:   s->nr            = -1;
2147:   s->nc            = -1;
2148:   s->m             = NULL;
2149:   s->splitassembly = PETSC_FALSE;

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

2153:   A->ops->mult                  = MatMult_Nest;
2154:   A->ops->multadd               = MatMultAdd_Nest;
2155:   A->ops->multtranspose         = MatMultTranspose_Nest;
2156:   A->ops->multtransposeadd      = MatMultTransposeAdd_Nest;
2157:   A->ops->transpose             = MatTranspose_Nest;
2158:   A->ops->assemblybegin         = MatAssemblyBegin_Nest;
2159:   A->ops->assemblyend           = MatAssemblyEnd_Nest;
2160:   A->ops->zeroentries           = MatZeroEntries_Nest;
2161:   A->ops->copy                  = MatCopy_Nest;
2162:   A->ops->axpy                  = MatAXPY_Nest;
2163:   A->ops->duplicate             = MatDuplicate_Nest;
2164:   A->ops->createsubmatrix       = MatCreateSubMatrix_Nest;
2165:   A->ops->destroy               = MatDestroy_Nest;
2166:   A->ops->view                  = MatView_Nest;
2167:   A->ops->getvecs               = 0; /* Use VECNEST by calling MatNestSetVecType(A,VECNEST) */
2168:   A->ops->getlocalsubmatrix     = MatGetLocalSubMatrix_Nest;
2169:   A->ops->restorelocalsubmatrix = MatRestoreLocalSubMatrix_Nest;
2170:   A->ops->getdiagonal           = MatGetDiagonal_Nest;
2171:   A->ops->diagonalscale         = MatDiagonalScale_Nest;
2172:   A->ops->scale                 = MatScale_Nest;
2173:   A->ops->shift                 = MatShift_Nest;
2174:   A->ops->diagonalset           = MatDiagonalSet_Nest;
2175:   A->ops->setrandom             = MatSetRandom_Nest;
2176:   A->ops->hasoperation          = MatHasOperation_Nest;
2177:   A->ops->missingdiagonal       = MatMissingDiagonal_Nest;

2179:   A->spptr        = 0;
2180:   A->assembled    = PETSC_FALSE;

2182:   /* expose Nest api's */
2183:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMat_C",        MatNestGetSubMat_Nest);
2184:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMat_C",        MatNestSetSubMat_Nest);
2185:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSubMats_C",       MatNestGetSubMats_Nest);
2186:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetSize_C",          MatNestGetSize_Nest);
2187:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetISs_C",           MatNestGetISs_Nest);
2188:   PetscObjectComposeFunction((PetscObject)A,"MatNestGetLocalISs_C",      MatNestGetLocalISs_Nest);
2189:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetVecType_C",       MatNestSetVecType_Nest);
2190:   PetscObjectComposeFunction((PetscObject)A,"MatNestSetSubMats_C",       MatNestSetSubMats_Nest);
2191:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_mpiaij_C",  MatConvert_Nest_AIJ);
2192:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_seqaij_C",  MatConvert_Nest_AIJ);
2193:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_aij_C",     MatConvert_Nest_AIJ);
2194:   PetscObjectComposeFunction((PetscObject)A,"MatConvert_nest_is_C",      MatConvert_Nest_IS);
2195:   PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_nest_seqdense_C",MatProductSetFromOptions_Nest_Dense);
2196:   PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_nest_mpidense_C",MatProductSetFromOptions_Nest_Dense);
2197:   PetscObjectComposeFunction((PetscObject)A,"MatProductSetFromOptions_nest_dense_C",MatProductSetFromOptions_Nest_Dense);

2199:   PetscObjectChangeTypeName((PetscObject)A,MATNEST);
2200:   return(0);
2201: }