Actual source code: mpimatmatmult.c

petsc-3.13.0 2020-03-29
Report Typos and Errors

  2: /*
  3:   Defines matrix-matrix product routines for pairs of MPIAIJ matrices
  4:           C = A * B
  5: */
  6:  #include <../src/mat/impls/aij/seq/aij.h>
  7:  #include <../src/mat/utils/freespace.h>
  8:  #include <../src/mat/impls/aij/mpi/mpiaij.h>
  9:  #include <petscbt.h>
 10:  #include <../src/mat/impls/dense/mpi/mpidense.h>
 11:  #include <petsc/private/vecimpl.h>
 12:  #include <petsc/private/vecscatterimpl.h>

 14: #if defined(PETSC_HAVE_HYPRE)
 15: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat);
 16: #endif

 18: PETSC_INTERN PetscErrorCode MatProductSymbolic_AB_MPIAIJ_MPIAIJ(Mat C)
 19: {
 20:   PetscErrorCode      ierr;
 21:   Mat_Product         *product = C->product;
 22:   Mat                 A=product->A,B=product->B;
 23:   MatProductAlgorithm alg=product->alg;
 24:   PetscReal           fill=product->fill;
 25:   PetscBool           flg;

 28:   /* scalable */
 29:   PetscStrcmp(alg,"scalable",&flg);
 30:   if (flg) {
 31:     MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);
 32:     goto next;
 33:   }

 35:   /* nonscalable */
 36:   PetscStrcmp(alg,"nonscalable",&flg);
 37:   if (flg) {
 38:     MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);
 39:     goto next;
 40:   }

 42:   /* seqmpi */
 43:   PetscStrcmp(alg,"seqmpi",&flg);
 44:   if (flg) {
 45:     MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(A,B,fill,C);
 46:     goto next;
 47:   }

 49: #if defined(PETSC_HAVE_HYPRE)
 50:   PetscStrcmp(alg,"hypre",&flg);
 51:   if (flg) {
 52:     MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);
 53:     return(0);
 54:   }
 55: #endif
 56:   SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_SUP,"Mat Product Algorithm is not supported");

 58: next:
 59:   {
 60:     Mat_MPIAIJ *c  = (Mat_MPIAIJ*)C->data;
 61:     Mat_APMPI  *ap = c->ap;
 62:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatFreeIntermediateDataStructures","Mat");
 63:     ap->freestruct = PETSC_FALSE;
 64:     PetscOptionsBool("-mat_freeintermediatedatastructures","Free intermediate data structures", "MatFreeIntermediateDataStructures",ap->freestruct,&ap->freestruct, NULL);
 65:     PetscOptionsEnd();
 66:   }
 67:   return(0);
 68: }

 70: PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
 71: {
 73:   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
 74:   Mat_APMPI      *ptap = a->ap;

 77:   PetscFree2(ptap->startsj_s,ptap->startsj_r);
 78:   PetscFree(ptap->bufa);
 79:   MatDestroy(&ptap->P_loc);
 80:   MatDestroy(&ptap->P_oth);
 81:   MatDestroy(&ptap->Pt);
 82:   PetscFree(ptap->api);
 83:   PetscFree(ptap->apj);
 84:   PetscFree(ptap->apa);
 85:   ptap->destroy(A);
 86:   PetscFree(ptap);
 87:   return(0);
 88: }

 90: PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C)
 91: {
 93:   Mat_MPIAIJ     *a  =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
 94:   Mat_SeqAIJ     *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
 95:   Mat_SeqAIJ     *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
 96:   PetscScalar    *cda=cd->a,*coa=co->a;
 97:   Mat_SeqAIJ     *p_loc,*p_oth;
 98:   PetscScalar    *apa,*ca;
 99:   PetscInt       cm   =C->rmap->n;
100:   Mat_APMPI      *ptap=c->ap;
101:   PetscInt       *api,*apj,*apJ,i,k;
102:   PetscInt       cstart=C->cmap->rstart;
103:   PetscInt       cdnz,conz,k0,k1;
104:   MPI_Comm       comm;
105:   PetscMPIInt    size;

108:   PetscObjectGetComm((PetscObject)A,&comm);
109:   MPI_Comm_size(comm,&size);

111:   if (!ptap->P_oth && size>1) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"AP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");

113:   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
114:   /*-----------------------------------------------------*/
115:   /* update numerical values of P_oth and P_loc */
116:   MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);
117:   MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);

119:   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
120:   /*----------------------------------------------------------*/
121:   /* get data from symbolic products */
122:   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
123:   p_oth = NULL;
124:   if (size >1) {
125:     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
126:   }

128:   /* get apa for storing dense row A[i,:]*P */
129:   apa = ptap->apa;

131:   api = ptap->api;
132:   apj = ptap->apj;
133:   for (i=0; i<cm; i++) {
134:     /* compute apa = A[i,:]*P */
135:     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);

137:     /* set values in C */
138:     apJ  = apj + api[i];
139:     cdnz = cd->i[i+1] - cd->i[i];
140:     conz = co->i[i+1] - co->i[i];

142:     /* 1st off-diagoanl part of C */
143:     ca = coa + co->i[i];
144:     k  = 0;
145:     for (k0=0; k0<conz; k0++) {
146:       if (apJ[k] >= cstart) break;
147:       ca[k0]      = apa[apJ[k]];
148:       apa[apJ[k++]] = 0.0;
149:     }

151:     /* diagonal part of C */
152:     ca = cda + cd->i[i];
153:     for (k1=0; k1<cdnz; k1++) {
154:       ca[k1]      = apa[apJ[k]];
155:       apa[apJ[k++]] = 0.0;
156:     }

158:     /* 2nd off-diagoanl part of C */
159:     ca = coa + co->i[i];
160:     for (; k0<conz; k0++) {
161:       ca[k0]      = apa[apJ[k]];
162:       apa[apJ[k++]] = 0.0;
163:     }
164:   }
165:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
166:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

168:   if (ptap->freestruct) {
169:     MatFreeIntermediateDataStructures(C);
170:   }
171:   return(0);
172: }

174: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat C)
175: {
176:   PetscErrorCode     ierr;
177:   MPI_Comm           comm;
178:   PetscMPIInt        size;
179:   Mat_APMPI          *ptap;
180:   PetscFreeSpaceList free_space=NULL,current_space=NULL;
181:   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data,*c;
182:   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
183:   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
184:   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
185:   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
186:   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
187:   PetscBT            lnkbt;
188:   PetscReal          afill;
189:   MatType            mtype;

192:   PetscObjectGetComm((PetscObject)A,&comm);
193:   MPI_Comm_size(comm,&size);

195:   /* create struct Mat_APMPI and attached it to C later */
196:   PetscNew(&ptap);

198:   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
199:   MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);

201:   /* get P_loc by taking all local rows of P */
202:   MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);

204:   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
205:   pi_loc = p_loc->i; pj_loc = p_loc->j;
206:   if (size > 1) {
207:     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
208:     pi_oth = p_oth->i; pj_oth = p_oth->j;
209:   } else {
210:     p_oth = NULL;
211:     pi_oth = NULL; pj_oth = NULL;
212:   }

214:   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
215:   /*-------------------------------------------------------------------*/
216:   PetscMalloc1(am+2,&api);
217:   ptap->api = api;
218:   api[0]    = 0;

220:   /* create and initialize a linked list */
221:   PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);

223:   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
224:   PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);
225:   current_space = free_space;

227:   MatPreallocateInitialize(comm,am,pn,dnz,onz);
228:   for (i=0; i<am; i++) {
229:     /* diagonal portion of A */
230:     nzi = adi[i+1] - adi[i];
231:     for (j=0; j<nzi; j++) {
232:       row  = *adj++;
233:       pnz  = pi_loc[row+1] - pi_loc[row];
234:       Jptr = pj_loc + pi_loc[row];
235:       /* add non-zero cols of P into the sorted linked list lnk */
236:       PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);
237:     }
238:     /* off-diagonal portion of A */
239:     nzi = aoi[i+1] - aoi[i];
240:     for (j=0; j<nzi; j++) {
241:       row  = *aoj++;
242:       pnz  = pi_oth[row+1] - pi_oth[row];
243:       Jptr = pj_oth + pi_oth[row];
244:       PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);
245:     }

247:     apnz     = lnk[0];
248:     api[i+1] = api[i] + apnz;

250:     /* if free space is not available, double the total space in the list */
251:     if (current_space->local_remaining<apnz) {
252:       PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);
253:       nspacedouble++;
254:     }

256:     /* Copy data into free space, then initialize lnk */
257:     PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);
258:     MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);

260:     current_space->array           += apnz;
261:     current_space->local_used      += apnz;
262:     current_space->local_remaining -= apnz;
263:   }

265:   /* Allocate space for apj, initialize apj, and */
266:   /* destroy list of free space and other temporary array(s) */
267:   PetscMalloc1(api[am]+1,&ptap->apj);
268:   apj  = ptap->apj;
269:   PetscFreeSpaceContiguous(&free_space,ptap->apj);
270:   PetscLLDestroy(lnk,lnkbt);

272:   /* malloc apa to store dense row A[i,:]*P */
273:   PetscCalloc1(pN,&ptap->apa);

275:   /* set and assemble symbolic parallel matrix C */
276:   /*---------------------------------------------*/
277:   MatSetSizes(C,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);
278:   MatSetBlockSizesFromMats(C,A,P);

280:   MatGetType(A,&mtype);
281:   MatSetType(C,mtype);
282:   MatMPIAIJSetPreallocation(C,0,dnz,0,onz);

284:   MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(C, apj, api);
285:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
286:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
287:   MatPreallocateFinalize(dnz,onz);

289:   ptap->destroy        = C->ops->destroy;
290:   ptap->duplicate      = C->ops->duplicate;
291:   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
292:   C->ops->productnumeric = MatProductNumeric_AB;
293:   C->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
294:   C->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;

296:   /* attach the supporting struct to C for reuse */
297:   c     = (Mat_MPIAIJ*)C->data;
298:   c->ap = ptap;

300:   /* set MatInfo */
301:   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
302:   if (afill < 1.0) afill = 1.0;
303:   C->info.mallocs           = nspacedouble;
304:   C->info.fill_ratio_given  = fill;
305:   C->info.fill_ratio_needed = afill;

307: #if defined(PETSC_USE_INFO)
308:   if (api[am]) {
309:     PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
310:     PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
311:   } else {
312:     PetscInfo(C,"Empty matrix product\n");
313:   }
314: #endif
315:   return(0);
316: }

318: /* ------------------------------------------------------- */
319: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense_AB(Mat C)
320: {
321:   Mat_Product *product = C->product;
322:   Mat         A = product->A,B=product->B;

325:   if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend)
326:     SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,B->rmap->rstart,B->rmap->rend);

328:   C->ops->matmultsymbolic = MatMatMultSymbolic_MPIAIJ_MPIDense;
329:   C->ops->productsymbolic = MatProductSymbolic_AB;
330:   C->ops->productnumeric  = MatProductNumeric_AB;
331:   return(0);
332: }
333: /* -------------------------------------------------------------------- */
334: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense_AtB(Mat C)
335: {
336:   Mat_Product *product = C->product;
337:   Mat         A = product->A,B=product->B;

340:   if (A->rmap->rstart != B->rmap->rstart || A->rmap->rend != B->rmap->rend)
341:     SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->rmap->rstart,A->rmap->rend,B->rmap->rstart,B->rmap->rend);

343:   C->ops->transposematmultsymbolic = MatTransposeMatMultSymbolic_MPIAIJ_MPIDense;
344:   C->ops->productsymbolic          = MatProductSymbolic_AtB;
345:   return(0);
346: }

348: /* --------------------------------------------------------------------- */
349: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense(Mat C)
350: {
352:   Mat_Product    *product = C->product;

355:   switch (product->type) {
356:   case MATPRODUCT_AB:
357:     MatProductSetFromOptions_MPIAIJ_MPIDense_AB(C);
358:     break;
359:   case MATPRODUCT_AtB:
360:     MatProductSetFromOptions_MPIAIJ_MPIDense_AtB(C);
361:     break;
362:   default:
363:     /* Use MatProduct_Basic() if there is no specific implementation */
364:     C->ops->productsymbolic = MatProductSymbolic_Basic;
365:   }
366:   return(0);
367: }
368: /* ------------------------------------------------------- */

370: typedef struct {
371:   Mat          workB,Bb,Cb,workB1,Bb1,Cb1;
372:   MPI_Request  *rwaits,*swaits;
373:   PetscInt     numBb;  /* num of Bb matrices */
374:   PetscInt     nsends,nrecvs;
375:   MPI_Datatype *stype,*rtype;
376: } MPIAIJ_MPIDense;

378: PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
379: {
380:   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*)ctx;
381:   PetscErrorCode  ierr;
382:   PetscInt        i;

385:   MatDestroy(&contents->workB);

387:   if (contents->numBb) {
388:     MatDestroy(&contents->Bb);
389:     MatDestroy(&contents->Cb);

391:     MatDestroy(&contents->workB1);
392:     MatDestroy(&contents->Bb1);
393:     MatDestroy(&contents->Cb1);
394:   }
395:   for (i=0; i<contents->nsends; i++) {
396:     MPI_Type_free(&contents->stype[i]);
397:   }
398:   for (i=0; i<contents->nrecvs; i++) {
399:     MPI_Type_free(&contents->rtype[i]);
400:   }
401:   PetscFree4(contents->stype,contents->rtype,contents->rwaits,contents->swaits);
402:   PetscFree(contents);
403:   return(0);
404: }

406: /*
407:     This is a "dummy function" that handles the case where matrix C was created as a dense matrix
408:   directly by the user and passed to MatMatMult() with the MAT_REUSE_MATRIX option

410:   It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C
411: */
412: PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C)
413: {
414:   PetscBool      flg;

418:   PetscObjectTypeCompare((PetscObject)A,MATNEST,&flg);
419:   if (flg) {
420:     PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
421:     MatMatMultSymbolic_Nest_Dense(A,B,PETSC_DEFAULT,&C);
422:     PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
423:     C->ops->matmultnumeric = MatMatMultNumeric_Nest_Dense;
424:   } else {
425:     MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,PETSC_DEFAULT,C);
426:   }
427:   (*C->ops->matmultnumeric)(A,B,C);
428:   return(0);
429: }

431: /*
432:   Create Bb, Cb, Bb1 and Cb1 matrices to be used by MatMatMultSymbolic_MPIAIJ_MPIDense().
433:   These matrices are used as wrappers for sub-columns of B and C, thus their own matrix operations are not used.
434:   Modified from MatCreateDense().
435: */
436: PETSC_STATIC_INLINE PetscErrorCode MatCreateSubMPIDense_private(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt rbs,PetscInt cbs,PetscScalar *data,Mat *A)
437: {

441:   MatCreate(comm,A);
442:   MatSetSizes(*A,m,n,M,N);
443:   MatSetBlockSizes(*A,rbs,cbs);
444:   MatSetType(*A,MATMPIDENSE);
445:   MatMPIDenseSetPreallocation(*A,data);
446:   (*A)->assembled = PETSC_TRUE;
447:   return(0);
448: }

450: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat C)
451: {
452:   PetscErrorCode  ierr;
453:   Mat_MPIAIJ      *aij=(Mat_MPIAIJ*)A->data;
454:   Mat_MPIDense    *b=(Mat_MPIDense*)B->data;
455:   Mat_SeqDense    *bseq=(Mat_SeqDense*)(b->A)->data;
456:   PetscInt        nz=aij->B->cmap->n,nsends,nrecvs,i,nrows_to,j,lda=bseq->lda;
457:   PetscContainer  container;
458:   MPIAIJ_MPIDense *contents;
459:   VecScatter      ctx=aij->Mvctx;
460:   PetscInt        Am=A->rmap->n,Bm=B->rmap->n,BN=B->cmap->N,Bbn,Bbn1,bs,nrows_from;
461:   MPI_Comm        comm;
462:   MPI_Datatype    type1,*stype,*rtype;
463:   const PetscInt  *sindices,*sstarts,*rstarts;
464:   PetscMPIInt     *disp;

467:   PetscObjectGetComm((PetscObject)A,&comm);
468:   if (!C->preallocated) {
469:     MatSetSizes(C,Am,B->cmap->n,A->rmap->N,BN);
470:     MatSetBlockSizesFromMats(C,A,B);
471:     MatSetType(C,MATMPIDENSE);
472:     MatMPIDenseSetPreallocation(C,NULL);
473:   }

475:   PetscNew(&contents);
476:   contents->numBb = 0;

478:   VecScatterGetRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,NULL,NULL);
479:   VecScatterGetRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,NULL,NULL);

481:   /* Create column block of B and C for memory scalability when BN is too large */
482:   /* Estimate Bbn, column size of Bb */
483:   if (nz) {
484:     Bbn1 = 2*Am*BN/nz;
485:   } else Bbn1 = BN;

487:   bs = PetscAbs(B->cmap->bs);
488:   Bbn1 = Bbn1/bs *bs; /* Bbn1 is a multiple of bs */
489:   if (Bbn1 > BN) Bbn1 = BN;
490:   MPI_Allreduce(&Bbn1,&Bbn,1,MPIU_INT,MPI_MAX,comm);

492:   /* Enable runtime option for Bbn */
493:   PetscOptionsBegin(comm,((PetscObject)C)->prefix,"MatMatMult","Mat");
494:   PetscOptionsInt("-matmatmult_Bbn","Number of columns in Bb","MatMatMult",Bbn,&Bbn,NULL);
495:   if (Bbn > BN) SETERRQ2(comm,PETSC_ERR_ARG_SIZ,"Bbn=%D cannot be larger than %D, column size of B",Bbn,BN);
496:   PetscOptionsEnd();

498:   if (Bbn < BN) {
499:     contents->numBb = BN/Bbn;
500:     Bbn1 = BN - contents->numBb*Bbn;
501:   }

503:   if (contents->numBb) {
504:     PetscScalar data[1]; /* fake array for Bb and Cb */
505:     PetscInfo3(C,"use Bb, BN=%D, Bbn=%D; numBb=%D\n",BN,Bbn,contents->numBb);
506:     MatCreateSubMPIDense_private(comm,B->rmap->n,PETSC_DECIDE,A->rmap->N,Bbn,B->rmap->bs,B->cmap->bs,data,&contents->Bb);
507:     MatCreateSubMPIDense_private(comm,Am,PETSC_DECIDE,A->rmap->N,Bbn,C->rmap->bs,C->cmap->bs,data,&contents->Cb);

509:     if (Bbn1) { /* Create Bb1 and Cb1 for the remaining columns */
510:       PetscInfo2(C,"use Bb1, BN=%D, Bbn1=%D\n",BN,Bbn1);
511:       MatCreateSubMPIDense_private(comm,B->rmap->n,PETSC_DECIDE,A->rmap->N,Bbn1,B->rmap->bs,B->cmap->bs,data,&contents->Bb1);
512:       MatCreateSubMPIDense_private(comm,Am,PETSC_DECIDE,A->rmap->N,Bbn1,C->rmap->bs,C->cmap->bs,data,&contents->Cb1);

514:       /* Create work matrix used to store off processor rows of B needed for local product */
515:       MatCreateSeqDense(PETSC_COMM_SELF,nz,Bbn1,NULL,&contents->workB1);
516:     }
517:   }

519:   /* Create work matrix used to store off processor rows of B needed for local product */
520:   MatCreateSeqDense(PETSC_COMM_SELF,nz,Bbn,NULL,&contents->workB);

522:   /* Use MPI derived data type to reduce memory required by the send/recv buffers */
523:   PetscMalloc4(nsends,&stype,nrecvs,&rtype,nrecvs,&contents->rwaits,nsends,&contents->swaits);
524:   contents->stype  = stype;
525:   contents->nsends = nsends;

527:   contents->rtype  = rtype;
528:   contents->nrecvs = nrecvs;

530:   PetscMalloc1(Bm+1,&disp);
531:   for (i=0; i<nsends; i++) {
532:     nrows_to = sstarts[i+1]-sstarts[i];
533:     for (j=0; j<nrows_to; j++){
534:       disp[j] = sindices[sstarts[i]+j]; /* rowB to be sent */
535:     }
536:     MPI_Type_create_indexed_block(nrows_to,1,(const PetscMPIInt *)disp,MPIU_SCALAR,&type1);

538:     MPI_Type_create_resized(type1,0,lda*sizeof(PetscScalar),&stype[i]);
539:     MPI_Type_commit(&stype[i]);
540:     MPI_Type_free(&type1);
541:   }

543:   for (i=0; i<nrecvs; i++) {
544:     /* received values from a process form a (nrows_from x Bbn) row block in workB (column-wise) */
545:     nrows_from = rstarts[i+1]-rstarts[i];
546:     disp[0] = 0;
547:     MPI_Type_create_indexed_block(1, nrows_from, (const PetscMPIInt *)disp, MPIU_SCALAR, &type1);
548:     MPI_Type_create_resized(type1, 0, nz*sizeof(PetscScalar), &rtype[i]);
549:     MPI_Type_commit(&rtype[i]);
550:     MPI_Type_free(&type1);
551:   }

553:   PetscFree(disp);
554:   VecScatterRestoreRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,NULL,NULL);
555:   VecScatterRestoreRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,NULL,NULL);

557:   PetscContainerCreate(comm,&container);
558:   PetscContainerSetPointer(container,contents);
559:   PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);
560:   PetscObjectCompose((PetscObject)C,"workB",(PetscObject)container);
561:   PetscContainerDestroy(&container);
562:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
563:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
564:   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
565:   return(0);
566: }

568: extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
569: /*
570:     Performs an efficient scatter on the rows of B needed by this process; this is
571:     a modification of the VecScatterBegin_() routines.

573:     Input: Bbidx = 0: B = Bb
574:                  = 1: B = Bb1, see MatMatMultSymbolic_MPIAIJ_MPIDense()
575: */
576: PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,PetscInt Bbidx,Mat C,Mat *outworkB)
577: {
578:   Mat_MPIAIJ        *aij = (Mat_MPIAIJ*)A->data;
579:   PetscErrorCode    ierr;
580:   const PetscScalar *b;
581:   PetscScalar       *rvalues;
582:   VecScatter        ctx = aij->Mvctx;
583:   const PetscInt    *sindices,*sstarts,*rstarts;
584:   const PetscMPIInt *sprocs,*rprocs;
585:   PetscInt          i,nsends,nrecvs,nrecvs2;
586:   MPI_Request       *swaits,*rwaits;
587:   MPI_Comm          comm;
588:   PetscMPIInt       tag=((PetscObject)ctx)->tag,ncols=B->cmap->N,nrows=aij->B->cmap->n,imdex,nsends_mpi,nrecvs_mpi;
589:   MPI_Status        status;
590:   MPIAIJ_MPIDense   *contents;
591:   PetscContainer    container;
592:   Mat               workB;
593:   MPI_Datatype      *stype,*rtype;

596:   PetscObjectGetComm((PetscObject)A,&comm);
597:   PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);
598:   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
599:   PetscContainerGetPointer(container,(void**)&contents);

601:   VecScatterGetRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,&sprocs,NULL/*bs*/);
602:   VecScatterGetRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,&rprocs,NULL/*bs*/);
603:   PetscMPIIntCast(nsends,&nsends_mpi);
604:   PetscMPIIntCast(nrecvs,&nrecvs_mpi);
605:   if (Bbidx == 0) {
606:     workB = *outworkB = contents->workB;
607:   } else {
608:     workB = *outworkB = contents->workB1;
609:   }
610:   if (nrows != workB->rmap->n) SETERRQ2(comm,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",workB->cmap->n,nrows);
611:   swaits  = contents->swaits;
612:   rwaits  = contents->rwaits;

614:   MatDenseGetArrayRead(B,&b);
615:   MatDenseGetArray(workB,&rvalues);

617:   /* Post recv, use MPI derived data type to save memory */
618:   rtype = contents->rtype;
619:   for (i=0; i<nrecvs; i++) {
620:     MPI_Irecv(rvalues+(rstarts[i]-rstarts[0]),ncols,rtype[i],rprocs[i],tag,comm,rwaits+i);
621:   }

623:   stype = contents->stype;
624:   for (i=0; i<nsends; i++) {
625:     MPI_Isend(b,ncols,stype[i],sprocs[i],tag,comm,swaits+i);
626:   }

628:   nrecvs2 = nrecvs;
629:   while (nrecvs2) {
630:     MPI_Waitany(nrecvs_mpi,rwaits,&imdex,&status);
631:     nrecvs2--;
632:   }
633:   if (nsends) {MPI_Waitall(nsends_mpi,swaits,MPI_STATUSES_IGNORE);}

635:   VecScatterRestoreRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,&sprocs,NULL);
636:   VecScatterRestoreRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,&rprocs,NULL);
637:   MatDenseRestoreArrayRead(B,&b);
638:   MatDenseRestoreArray(workB,&rvalues);
639:   MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);
640:   MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);
641:   return(0);
642: }

644: /*
645:   Compute Cb = A*Bb
646: */
647: PETSC_STATIC_INLINE PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense_private(Mat A,Mat Bb,PetscInt Bbidx,PetscInt start,Mat C,const PetscScalar *barray,PetscScalar *carray,Mat Cb)
648: {
649:   PetscErrorCode  ierr;
650:   PetscInt        start1;
651:   Mat             workB;
652:   Mat_MPIAIJ      *aij = (Mat_MPIAIJ*)A->data;
653:   Mat_MPIDense    *cbdense = (Mat_MPIDense*)Cb->data;

656:   /* Place barray to Bb */
657:   start1 = start*Bb->rmap->n;
658:   MatDensePlaceArray(Bb,barray+start1);

660:   /* get off processor parts of Bb needed to complete Cb=A*Bb */
661:   MatMPIDenseScatter(A,Bb,Bbidx,C,&workB);
662:   MatDenseResetArray(Bb);

664:   /* off-diagonal block of A times nonlocal rows of Bb */
665:   /* Place carray to Cb */
666:   start1 = start*Cb->rmap->n;
667:   MatDensePlaceArray(Cb,carray+start1);
668:   MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cbdense->A);
669:   MatDenseResetArray(Cb);
670:   return(0);
671: }

673: PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
674: {
675:   PetscErrorCode  ierr;
676:   Mat_MPIAIJ      *aij    = (Mat_MPIAIJ*)A->data;
677:   Mat_MPIDense    *bdense = (Mat_MPIDense*)B->data;
678:   Mat_MPIDense    *cdense = (Mat_MPIDense*)C->data;
679:   Mat             workB;
680:   MPIAIJ_MPIDense *contents;
681:   PetscContainer  container;
682:   MPI_Comm        comm;
683:   PetscInt        numBb;

686:   /* diagonal block of A times all local rows of B*/
687:   MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);

689:   PetscObjectGetComm((PetscObject)A,&comm);
690:   PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);
691:   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
692:   PetscContainerGetPointer(container,(void**)&contents);
693:   numBb = contents->numBb;

695:   if (!numBb) {
696:     /* get off processor parts of B needed to complete C=A*B */
697:     MatMPIDenseScatter(A,B,0,C,&workB);

699:     /* off-diagonal block of A times nonlocal rows of B */
700:     MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);

702:   } else {
703:     const PetscScalar *barray;
704:     PetscScalar       *carray;
705:     Mat               Bb=contents->Bb,Cb=contents->Cb;
706:     PetscInt          BbN=Bb->cmap->N,start,i;

708:     MatDenseGetArrayRead(B,&barray);
709:     MatDenseGetArray(C,&carray);
710:     for (i=0; i<numBb; i++) {
711:       start = i*BbN;
712:       MatMatMultNumeric_MPIAIJ_MPIDense_private(A,Bb,0,start,C,barray,carray,Cb);
713:     }

715:     if (contents->Bb1) {
716:       Bb = contents->Bb1; Cb = contents->Cb1;
717:       start = i*BbN;
718:       MatMatMultNumeric_MPIAIJ_MPIDense_private(A,Bb,1,start,C,barray,carray,Cb);
719:     }
720:     MatDenseRestoreArrayRead(B,&barray);
721:     MatDenseRestoreArray(C,&carray);
722:   }
723:   return(0);
724: }

726: PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
727: {
729:   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
730:   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
731:   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
732:   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
733:   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
734:   Mat_SeqAIJ     *p_loc,*p_oth;
735:   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
736:   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
737:   PetscInt       cm    = C->rmap->n,anz,pnz;
738:   Mat_APMPI      *ptap = c->ap;
739:   PetscScalar    *apa_sparse;
740:   PetscInt       *api,*apj,*apJ,i,j,k,row;
741:   PetscInt       cstart = C->cmap->rstart;
742:   PetscInt       cdnz,conz,k0,k1,nextp;
743:   MPI_Comm       comm;
744:   PetscMPIInt    size;

747:   PetscObjectGetComm((PetscObject)C,&comm);
748:   MPI_Comm_size(comm,&size);

750:   if (!ptap->P_oth && size>1) {
751:     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"AP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
752:   }
753:   apa_sparse = ptap->apa;

755:   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
756:   /*-----------------------------------------------------*/
757:   /* update numerical values of P_oth and P_loc */
758:   MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);
759:   MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);

761:   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
762:   /*----------------------------------------------------------*/
763:   /* get data from symbolic products */
764:   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
765:   pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a;
766:   if (size >1) {
767:     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
768:     pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a;
769:   } else {
770:     p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL;
771:   }

773:   api = ptap->api;
774:   apj = ptap->apj;
775:   for (i=0; i<cm; i++) {
776:     apJ = apj + api[i];

778:     /* diagonal portion of A */
779:     anz = adi[i+1] - adi[i];
780:     adj = ad->j + adi[i];
781:     ada = ad->a + adi[i];
782:     for (j=0; j<anz; j++) {
783:       row = adj[j];
784:       pnz = pi_loc[row+1] - pi_loc[row];
785:       pj  = pj_loc + pi_loc[row];
786:       pa  = pa_loc + pi_loc[row];
787:       /* perform sparse axpy */
788:       valtmp = ada[j];
789:       nextp  = 0;
790:       for (k=0; nextp<pnz; k++) {
791:         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
792:           apa_sparse[k] += valtmp*pa[nextp++];
793:         }
794:       }
795:       PetscLogFlops(2.0*pnz);
796:     }

798:     /* off-diagonal portion of A */
799:     anz = aoi[i+1] - aoi[i];
800:     aoj = ao->j + aoi[i];
801:     aoa = ao->a + aoi[i];
802:     for (j=0; j<anz; j++) {
803:       row = aoj[j];
804:       pnz = pi_oth[row+1] - pi_oth[row];
805:       pj  = pj_oth + pi_oth[row];
806:       pa  = pa_oth + pi_oth[row];
807:       /* perform sparse axpy */
808:       valtmp = aoa[j];
809:       nextp  = 0;
810:       for (k=0; nextp<pnz; k++) {
811:         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
812:           apa_sparse[k] += valtmp*pa[nextp++];
813:         }
814:       }
815:       PetscLogFlops(2.0*pnz);
816:     }

818:     /* set values in C */
819:     cdnz = cd->i[i+1] - cd->i[i];
820:     conz = co->i[i+1] - co->i[i];

822:     /* 1st off-diagoanl part of C */
823:     ca = coa + co->i[i];
824:     k  = 0;
825:     for (k0=0; k0<conz; k0++) {
826:       if (apJ[k] >= cstart) break;
827:       ca[k0]        = apa_sparse[k];
828:       apa_sparse[k] = 0.0;
829:       k++;
830:     }

832:     /* diagonal part of C */
833:     ca = cda + cd->i[i];
834:     for (k1=0; k1<cdnz; k1++) {
835:       ca[k1]        = apa_sparse[k];
836:       apa_sparse[k] = 0.0;
837:       k++;
838:     }

840:     /* 2nd off-diagoanl part of C */
841:     ca = coa + co->i[i];
842:     for (; k0<conz; k0++) {
843:       ca[k0]        = apa_sparse[k];
844:       apa_sparse[k] = 0.0;
845:       k++;
846:     }
847:   }
848:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
849:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

851:   if (ptap->freestruct) {
852:     MatFreeIntermediateDataStructures(C);
853:   }
854:   return(0);
855: }

857: /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
858: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat C)
859: {
860:   PetscErrorCode     ierr;
861:   MPI_Comm           comm;
862:   PetscMPIInt        size;
863:   Mat_APMPI          *ptap;
864:   PetscFreeSpaceList free_space = NULL,current_space=NULL;
865:   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
866:   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
867:   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
868:   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
869:   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
870:   PetscInt           am=A->rmap->n,pn=P->cmap->n,pm=P->rmap->n,lsize=pn+20;
871:   PetscReal          afill;
872:   MatType            mtype;

875:   PetscObjectGetComm((PetscObject)A,&comm);
876:   MPI_Comm_size(comm,&size);

878:   /* create struct Mat_APMPI and attached it to C later */
879:   PetscNew(&ptap);

881:   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
882:   MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);

884:   /* get P_loc by taking all local rows of P */
885:   MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);

887:   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
888:   pi_loc = p_loc->i; pj_loc = p_loc->j;
889:   if (size > 1) {
890:     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
891:     pi_oth = p_oth->i; pj_oth = p_oth->j;
892:   } else {
893:     p_oth  = NULL;
894:     pi_oth = NULL; pj_oth = NULL;
895:   }

897:   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
898:   /*-------------------------------------------------------------------*/
899:   PetscMalloc1(am+2,&api);
900:   ptap->api = api;
901:   api[0]    = 0;

903:   PetscLLCondensedCreate_Scalable(lsize,&lnk);

905:   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
906:   PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);
907:   current_space = free_space;
908:   MatPreallocateInitialize(comm,am,pn,dnz,onz);
909:   for (i=0; i<am; i++) {
910:     /* diagonal portion of A */
911:     nzi = adi[i+1] - adi[i];
912:     for (j=0; j<nzi; j++) {
913:       row  = *adj++;
914:       pnz  = pi_loc[row+1] - pi_loc[row];
915:       Jptr = pj_loc + pi_loc[row];
916:       /* Expand list if it is not long enough */
917:       if (pnz+apnz_max > lsize) {
918:         lsize = pnz+apnz_max;
919:         PetscLLCondensedExpand_Scalable(lsize, &lnk);
920:       }
921:       /* add non-zero cols of P into the sorted linked list lnk */
922:       PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);
923:       apnz     = *lnk; /* The first element in the list is the number of items in the list */
924:       api[i+1] = api[i] + apnz;
925:       if (apnz > apnz_max) apnz_max = apnz;
926:     }
927:     /* off-diagonal portion of A */
928:     nzi = aoi[i+1] - aoi[i];
929:     for (j=0; j<nzi; j++) {
930:       row  = *aoj++;
931:       pnz  = pi_oth[row+1] - pi_oth[row];
932:       Jptr = pj_oth + pi_oth[row];
933:       /* Expand list if it is not long enough */
934:       if (pnz+apnz_max > lsize) {
935:         lsize = pnz + apnz_max;
936:         PetscLLCondensedExpand_Scalable(lsize, &lnk);
937:       }
938:       /* add non-zero cols of P into the sorted linked list lnk */
939:       PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);
940:       apnz     = *lnk;  /* The first element in the list is the number of items in the list */
941:       api[i+1] = api[i] + apnz;
942:       if (apnz > apnz_max) apnz_max = apnz;
943:     }
944:     apnz     = *lnk;
945:     api[i+1] = api[i] + apnz;
946:     if (apnz > apnz_max) apnz_max = apnz;

948:     /* if free space is not available, double the total space in the list */
949:     if (current_space->local_remaining<apnz) {
950:       PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);
951:       nspacedouble++;
952:     }

954:     /* Copy data into free space, then initialize lnk */
955:     PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);
956:     MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);

958:     current_space->array           += apnz;
959:     current_space->local_used      += apnz;
960:     current_space->local_remaining -= apnz;
961:   }

963:   /* Allocate space for apj, initialize apj, and */
964:   /* destroy list of free space and other temporary array(s) */
965:   PetscMalloc1(api[am]+1,&ptap->apj);
966:   apj  = ptap->apj;
967:   PetscFreeSpaceContiguous(&free_space,ptap->apj);
968:   PetscLLCondensedDestroy_Scalable(lnk);

970:   /* create and assemble symbolic parallel matrix C */
971:   /*----------------------------------------------------*/
972:   MatSetSizes(C,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);
973:   MatSetBlockSizesFromMats(C,A,P);
974:   MatGetType(A,&mtype);
975:   MatSetType(C,mtype);
976:   MatMPIAIJSetPreallocation(C,0,dnz,0,onz);

978:   /* malloc apa for assembly C */
979:   PetscCalloc1(apnz_max,&ptap->apa);

981:   MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(C, apj, api);
982:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
983:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
984:   MatPreallocateFinalize(dnz,onz);

986:   ptap->destroy             = C->ops->destroy;
987:   ptap->duplicate           = C->ops->duplicate;
988:   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
989:   C->ops->productnumeric = MatProductNumeric_AB;
990:   C->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
991:   C->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;

993:   /* attach the supporting struct to C for reuse */
994:   c     = (Mat_MPIAIJ*)C->data;
995:   c->ap = ptap;

997:   /* set MatInfo */
998:   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
999:   if (afill < 1.0) afill = 1.0;
1000:   C->info.mallocs           = nspacedouble;
1001:   C->info.fill_ratio_given  = fill;
1002:   C->info.fill_ratio_needed = afill;

1004: #if defined(PETSC_USE_INFO)
1005:   if (api[am]) {
1006:     PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
1007:     PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
1008:   } else {
1009:     PetscInfo(C,"Empty matrix product\n");
1010:   }
1011: #endif
1012:   return(0);
1013: }

1015: /* This function is needed for the seqMPI matrix-matrix multiplication.  */
1016: /* Three input arrays are merged to one output array. The size of the    */
1017: /* output array is also output. Duplicate entries only show up once.     */
1018: static void Merge3SortedArrays(PetscInt  size1, PetscInt *in1,
1019:                                PetscInt  size2, PetscInt *in2,
1020:                                PetscInt  size3, PetscInt *in3,
1021:                                PetscInt *size4, PetscInt *out)
1022: {
1023:   int i = 0, j = 0, k = 0, l = 0;

1025:   /* Traverse all three arrays */
1026:   while (i<size1 && j<size2 && k<size3) {
1027:     if (in1[i] < in2[j] && in1[i] < in3[k]) {
1028:       out[l++] = in1[i++];
1029:     }
1030:     else if(in2[j] < in1[i] && in2[j] < in3[k]) {
1031:       out[l++] = in2[j++];
1032:     }
1033:     else if(in3[k] < in1[i] && in3[k] < in2[j]) {
1034:       out[l++] = in3[k++];
1035:     }
1036:     else if(in1[i] == in2[j] && in1[i] < in3[k]) {
1037:       out[l++] = in1[i];
1038:       i++, j++;
1039:     }
1040:     else if(in1[i] == in3[k] && in1[i] < in2[j]) {
1041:       out[l++] = in1[i];
1042:       i++, k++;
1043:     }
1044:     else if(in3[k] == in2[j] && in2[j] < in1[i])  {
1045:       out[l++] = in2[j];
1046:       k++, j++;
1047:     }
1048:     else if(in1[i] == in2[j] && in1[i] == in3[k]) {
1049:       out[l++] = in1[i];
1050:       i++, j++, k++;
1051:     }
1052:   }

1054:   /* Traverse two remaining arrays */
1055:   while (i<size1 && j<size2) {
1056:     if (in1[i] < in2[j]) {
1057:       out[l++] = in1[i++];
1058:     }
1059:     else if(in1[i] > in2[j]) {
1060:       out[l++] = in2[j++];
1061:     }
1062:     else {
1063:       out[l++] = in1[i];
1064:       i++, j++;
1065:     }
1066:   }

1068:   while (i<size1 && k<size3) {
1069:     if (in1[i] < in3[k]) {
1070:       out[l++] = in1[i++];
1071:     }
1072:     else if(in1[i] > in3[k]) {
1073:       out[l++] = in3[k++];
1074:     }
1075:     else {
1076:       out[l++] = in1[i];
1077:       i++, k++;
1078:     }
1079:   }

1081:   while (k<size3 && j<size2)  {
1082:     if (in3[k] < in2[j]) {
1083:       out[l++] = in3[k++];
1084:     }
1085:     else if(in3[k] > in2[j]) {
1086:       out[l++] = in2[j++];
1087:     }
1088:     else {
1089:       out[l++] = in3[k];
1090:       k++, j++;
1091:     }
1092:   }

1094:   /* Traverse one remaining array */
1095:   while (i<size1) out[l++] = in1[i++];
1096:   while (j<size2) out[l++] = in2[j++];
1097:   while (k<size3) out[l++] = in3[k++];

1099:   *size4 = l;
1100: }

1102: /* This matrix-matrix multiplication algorithm divides the multiplication into three multiplications and  */
1103: /* adds up the products. Two of these three multiplications are performed with existing (sequential)      */
1104: /* matrix-matrix multiplications.  */
1105: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat A, Mat P, PetscReal fill, Mat C)
1106: {
1107:   PetscErrorCode     ierr;
1108:   MPI_Comm           comm;
1109:   PetscMPIInt        size;
1110:   Mat_APMPI          *ptap;
1111:   PetscFreeSpaceList free_space_diag=NULL, current_space=NULL;
1112:   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data;
1113:   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc;
1114:   Mat_MPIAIJ         *p        =(Mat_MPIAIJ*)P->data;
1115:   Mat_MPIAIJ         *c;
1116:   Mat_SeqAIJ         *adpd_seq, *p_off, *aopoth_seq;
1117:   PetscInt           adponz, adpdnz;
1118:   PetscInt           *pi_loc,*dnz,*onz;
1119:   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,rstart=A->rmap->rstart;
1120:   PetscInt           *lnk,i, i1=0,pnz,row,*adpoi,*adpoj, *api, *adpoJ, *aopJ, *apJ,*Jptr, aopnz, nspacedouble=0,j,nzi,
1121:                      *apj,apnz, *adpdi, *adpdj, *adpdJ, *poff_i, *poff_j, *j_temp, *aopothi, *aopothj;
1122:   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n, p_colstart, p_colend;
1123:   PetscBT            lnkbt;
1124:   PetscReal          afill;
1125:   PetscMPIInt        rank;
1126:   Mat                adpd, aopoth;
1127:   MatType            mtype;
1128:   const char         *prefix;

1131:   PetscObjectGetComm((PetscObject)A,&comm);
1132:   MPI_Comm_size(comm,&size);
1133:   MPI_Comm_rank(comm, &rank);
1134:   MatGetOwnershipRangeColumn(P, &p_colstart, &p_colend);

1136:   /* create struct Mat_APMPI and attached it to C later */
1137:   PetscNew(&ptap);

1139:   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
1140:   MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);

1142:   /* get P_loc by taking all local rows of P */
1143:   MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);


1146:   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1147:   pi_loc = p_loc->i;

1149:   /* Allocate memory for the i arrays of the matrices A*P, A_diag*P_off and A_offd * P */
1150:   PetscMalloc1(am+2,&api);
1151:   PetscMalloc1(am+2,&adpoi);

1153:   adpoi[0]    = 0;
1154:   ptap->api = api;
1155:   api[0] = 0;

1157:   /* create and initialize a linked list, will be used for both A_diag * P_loc_off and A_offd * P_oth */
1158:   PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);
1159:   MatPreallocateInitialize(comm,am,pn,dnz,onz);

1161:   /* Symbolic calc of A_loc_diag * P_loc_diag */
1162:   MatGetOptionsPrefix(A,&prefix);
1163:   MatProductCreate(a->A,p->A,NULL,&adpd);
1164:   MatGetOptionsPrefix(A,&prefix);
1165:   MatSetOptionsPrefix(adpd,prefix);
1166:   MatAppendOptionsPrefix(adpd,"inner_diag_");

1168:   MatProductSetType(adpd,MATPRODUCT_AB);
1169:   MatProductSetAlgorithm(adpd,"sorted");
1170:   MatProductSetFill(adpd,fill);
1171:   MatProductSetFromOptions(adpd);
1172:   MatProductSymbolic(adpd);

1174:   adpd_seq = (Mat_SeqAIJ*)((adpd)->data);
1175:   adpdi = adpd_seq->i; adpdj = adpd_seq->j;
1176:   p_off = (Mat_SeqAIJ*)((p->B)->data);
1177:   poff_i = p_off->i; poff_j = p_off->j;

1179:   /* j_temp stores indices of a result row before they are added to the linked list */
1180:   PetscMalloc1(pN+2,&j_temp);


1183:   /* Symbolic calc of the A_diag * p_loc_off */
1184:   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
1185:   PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space_diag);
1186:   current_space = free_space_diag;

1188:   for (i=0; i<am; i++) {
1189:     /* A_diag * P_loc_off */
1190:     nzi = adi[i+1] - adi[i];
1191:     for (j=0; j<nzi; j++) {
1192:       row  = *adj++;
1193:       pnz  = poff_i[row+1] - poff_i[row];
1194:       Jptr = poff_j + poff_i[row];
1195:       for(i1 = 0; i1 < pnz; i1++) {
1196:         j_temp[i1] = p->garray[Jptr[i1]];
1197:       }
1198:       /* add non-zero cols of P into the sorted linked list lnk */
1199:       PetscLLCondensedAddSorted(pnz,j_temp,lnk,lnkbt);
1200:     }

1202:     adponz     = lnk[0];
1203:     adpoi[i+1] = adpoi[i] + adponz;

1205:     /* if free space is not available, double the total space in the list */
1206:     if (current_space->local_remaining<adponz) {
1207:       PetscFreeSpaceGet(PetscIntSumTruncate(adponz,current_space->total_array_size),&current_space);
1208:       nspacedouble++;
1209:     }

1211:     /* Copy data into free space, then initialize lnk */
1212:     PetscLLCondensedClean(pN,adponz,current_space->array,lnk,lnkbt);

1214:     current_space->array           += adponz;
1215:     current_space->local_used      += adponz;
1216:     current_space->local_remaining -= adponz;
1217:   }

1219:   /* Symbolic calc of A_off * P_oth */
1220:   MatSetOptionsPrefix(a->B,prefix);
1221:   MatAppendOptionsPrefix(a->B,"inner_offdiag_");
1222:   MatCreate(PETSC_COMM_SELF,&aopoth);
1223:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->B, ptap->P_oth, fill, aopoth);
1224:   aopoth_seq = (Mat_SeqAIJ*)((aopoth)->data);
1225:   aopothi = aopoth_seq->i; aopothj = aopoth_seq->j;

1227:   /* Allocate space for apj, adpj, aopj, ... */
1228:   /* destroy lists of free space and other temporary array(s) */

1230:   PetscMalloc1(aopothi[am] + adpoi[am] + adpdi[am]+2, &ptap->apj);
1231:   PetscMalloc1(adpoi[am]+2, &adpoj);

1233:   /* Copy from linked list to j-array */
1234:   PetscFreeSpaceContiguous(&free_space_diag,adpoj);
1235:   PetscLLDestroy(lnk,lnkbt);

1237:   adpoJ = adpoj;
1238:   adpdJ = adpdj;
1239:   aopJ = aopothj;
1240:   apj  = ptap->apj;
1241:   apJ = apj; /* still empty */

1243:   /* Merge j-arrays of A_off * P, A_diag * P_loc_off, and */
1244:   /* A_diag * P_loc_diag to get A*P */
1245:   for (i = 0; i < am; i++) {
1246:     aopnz  =  aopothi[i+1] -  aopothi[i];
1247:     adponz = adpoi[i+1] - adpoi[i];
1248:     adpdnz = adpdi[i+1] - adpdi[i];

1250:     /* Correct indices from A_diag*P_diag */
1251:     for(i1 = 0; i1 < adpdnz; i1++) {
1252:       adpdJ[i1] += p_colstart;
1253:     }
1254:     /* Merge j-arrays of A_diag * P_loc_off and A_diag * P_loc_diag and A_off * P_oth */
1255:     Merge3SortedArrays(adponz, adpoJ, adpdnz, adpdJ, aopnz, aopJ, &apnz, apJ);
1256:     MatPreallocateSet(i+rstart, apnz, apJ, dnz, onz);

1258:     aopJ += aopnz;
1259:     adpoJ += adponz;
1260:     adpdJ += adpdnz;
1261:     apJ += apnz;
1262:     api[i+1] = api[i] + apnz;
1263:   }

1265:   /* malloc apa to store dense row A[i,:]*P */
1266:   PetscCalloc1(pN+2,&ptap->apa);

1268:   /* create and assemble symbolic parallel matrix C */
1269:   MatSetSizes(C,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);
1270:   MatSetBlockSizesFromMats(C,A,P);
1271:   MatGetType(A,&mtype);
1272:   MatSetType(C,mtype);
1273:   MatMPIAIJSetPreallocation(C,0,dnz,0,onz);


1276:   MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(C, apj, api);
1277:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1278:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1279:   MatPreallocateFinalize(dnz,onz);


1282:   ptap->destroy        = C->ops->destroy;
1283:   ptap->duplicate      = C->ops->duplicate;
1284:   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1285:   C->ops->productnumeric = MatProductNumeric_AB;
1286:   C->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;

1288:   /* attach the supporting struct to C for reuse */
1289:   c       = (Mat_MPIAIJ*)C->data;
1290:   c->ap = ptap;

1292:   /* set MatInfo */
1293:   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
1294:   if (afill < 1.0) afill = 1.0;
1295:   C->info.mallocs           = nspacedouble;
1296:   C->info.fill_ratio_given  = fill;
1297:   C->info.fill_ratio_needed = afill;

1299: #if defined(PETSC_USE_INFO)
1300:   if (api[am]) {
1301:     PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
1302:     PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
1303:   } else {
1304:     PetscInfo(C,"Empty matrix product\n");
1305:   }
1306: #endif

1308:   MatDestroy(&aopoth);
1309:   MatDestroy(&adpd);
1310:   PetscFree(j_temp);
1311:   PetscFree(adpoj);
1312:   PetscFree(adpoi);
1313:   return(0);
1314: }

1316: /*-------------------------------------------------------------------------*/
1317: /* This routine only works when scall=MAT_REUSE_MATRIX! */
1318: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
1319: {
1321:   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
1322:   Mat_APMPI      *ptap= c->ap;
1323:   Mat            Pt;

1326:   if (!ptap->Pt) {
1327:     MPI_Comm comm;
1328:     PetscObjectGetComm((PetscObject)C,&comm);
1329:     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
1330:   }

1332:   Pt = ptap->Pt;
1333:   MatTranspose(P,MAT_REUSE_MATRIX,&Pt);
1334:   MatMatMultNumeric_MPIAIJ_MPIAIJ(Pt,A,C);

1336:   /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */
1337:   if (ptap->freestruct) {
1338:     MatFreeIntermediateDataStructures(C);
1339:   }
1340:   return(0);
1341: }

1343: /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1344: PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat C)
1345: {
1346:   PetscErrorCode      ierr;
1347:   Mat_APMPI           *ptap;
1348:   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*c;
1349:   MPI_Comm            comm;
1350:   PetscMPIInt         size,rank;
1351:   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1352:   PetscInt            pn=P->cmap->n,aN=A->cmap->N,an=A->cmap->n;
1353:   PetscInt            *lnk,i,k,nsend;
1354:   PetscBT             lnkbt;
1355:   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
1356:   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1357:   PetscInt            len,proc,*dnz,*onz,*owners,nzi;
1358:   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1359:   MPI_Request         *swaits,*rwaits;
1360:   MPI_Status          *sstatus,rstatus;
1361:   PetscLayout         rowmap;
1362:   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
1363:   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
1364:   PetscInt            *Jptr,*prmap=p->garray,con,j,Crmax;
1365:   Mat_SeqAIJ          *a_loc,*c_loc,*c_oth;
1366:   PetscTable          ta;
1367:   MatType             mtype;
1368:   const char          *prefix;

1371:   PetscObjectGetComm((PetscObject)A,&comm);
1372:   MPI_Comm_size(comm,&size);
1373:   MPI_Comm_rank(comm,&rank);

1375:   /* create symbolic parallel matrix C */
1376:   MatGetType(A,&mtype);
1377:   MatSetType(C,mtype);

1379:   C->ops->transposematmultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;

1381:   /* create struct Mat_APMPI and attached it to C later */
1382:   PetscNew(&ptap);
1383:   ptap->reuse = MAT_INITIAL_MATRIX;

1385:   /* (0) compute Rd = Pd^T, Ro = Po^T  */
1386:   /* --------------------------------- */
1387:   MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);
1388:   MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);

1390:   /* (1) compute symbolic A_loc */
1391:   /* ---------------------------*/
1392:   MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&ptap->A_loc);

1394:   /* (2-1) compute symbolic C_oth = Ro*A_loc  */
1395:   /* ------------------------------------ */
1396:   MatGetOptionsPrefix(A,&prefix);
1397:   MatSetOptionsPrefix(ptap->Ro,prefix);
1398:   MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");
1399:   MatCreate(PETSC_COMM_SELF,&ptap->C_oth);
1400:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->A_loc,fill,ptap->C_oth);

1402:   /* (3) send coj of C_oth to other processors  */
1403:   /* ------------------------------------------ */
1404:   /* determine row ownership */
1405:   PetscLayoutCreate(comm,&rowmap);
1406:   rowmap->n  = pn;
1407:   rowmap->bs = 1;
1408:   PetscLayoutSetUp(rowmap);
1409:   owners = rowmap->range;

1411:   /* determine the number of messages to send, their lengths */
1412:   PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);
1413:   PetscArrayzero(len_s,size);
1414:   PetscArrayzero(len_si,size);

1416:   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
1417:   coi   = c_oth->i; coj = c_oth->j;
1418:   con   = ptap->C_oth->rmap->n;
1419:   proc  = 0;
1420:   for (i=0; i<con; i++) {
1421:     while (prmap[i] >= owners[proc+1]) proc++;
1422:     len_si[proc]++;               /* num of rows in Co(=Pt*A) to be sent to [proc] */
1423:     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
1424:   }

1426:   len          = 0; /* max length of buf_si[], see (4) */
1427:   owners_co[0] = 0;
1428:   nsend        = 0;
1429:   for (proc=0; proc<size; proc++) {
1430:     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1431:     if (len_s[proc]) {
1432:       nsend++;
1433:       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
1434:       len         += len_si[proc];
1435:     }
1436:   }

1438:   /* determine the number and length of messages to receive for coi and coj  */
1439:   PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);
1440:   PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);

1442:   /* post the Irecv and Isend of coj */
1443:   PetscCommGetNewTag(comm,&tagj);
1444:   PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);
1445:   PetscMalloc1(nsend+1,&swaits);
1446:   for (proc=0, k=0; proc<size; proc++) {
1447:     if (!len_s[proc]) continue;
1448:     i    = owners_co[proc];
1449:     MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);
1450:     k++;
1451:   }

1453:   /* (2-2) compute symbolic C_loc = Rd*A_loc */
1454:   /* ---------------------------------------- */
1455:   MatSetOptionsPrefix(ptap->Rd,prefix);
1456:   MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");
1457:   MatCreate(PETSC_COMM_SELF,&ptap->C_loc);
1458:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->A_loc,fill,ptap->C_loc);
1459:   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;

1461:   /* receives coj are complete */
1462:   for (i=0; i<nrecv; i++) {
1463:     MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);
1464:   }
1465:   PetscFree(rwaits);
1466:   if (nsend) {MPI_Waitall(nsend,swaits,sstatus);}

1468:   /* add received column indices into ta to update Crmax */
1469:   a_loc = (Mat_SeqAIJ*)(ptap->A_loc)->data;

1471:   /* create and initialize a linked list */
1472:   PetscTableCreate(an,aN,&ta); /* for compute Crmax */
1473:   MatRowMergeMax_SeqAIJ(a_loc,ptap->A_loc->rmap->N,ta);

1475:   for (k=0; k<nrecv; k++) {/* k-th received message */
1476:     Jptr = buf_rj[k];
1477:     for (j=0; j<len_r[k]; j++) {
1478:       PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);
1479:     }
1480:   }
1481:   PetscTableGetCount(ta,&Crmax);
1482:   PetscTableDestroy(&ta);

1484:   /* (4) send and recv coi */
1485:   /*-----------------------*/
1486:   PetscCommGetNewTag(comm,&tagi);
1487:   PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);
1488:   PetscMalloc1(len+1,&buf_s);
1489:   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1490:   for (proc=0,k=0; proc<size; proc++) {
1491:     if (!len_s[proc]) continue;
1492:     /* form outgoing message for i-structure:
1493:          buf_si[0]:                 nrows to be sent
1494:                [1:nrows]:           row index (global)
1495:                [nrows+1:2*nrows+1]: i-structure index
1496:     */
1497:     /*-------------------------------------------*/
1498:     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
1499:     buf_si_i    = buf_si + nrows+1;
1500:     buf_si[0]   = nrows;
1501:     buf_si_i[0] = 0;
1502:     nrows       = 0;
1503:     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1504:       nzi = coi[i+1] - coi[i];
1505:       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1506:       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1507:       nrows++;
1508:     }
1509:     MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);
1510:     k++;
1511:     buf_si += len_si[proc];
1512:   }
1513:   for (i=0; i<nrecv; i++) {
1514:     MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);
1515:   }
1516:   PetscFree(rwaits);
1517:   if (nsend) {MPI_Waitall(nsend,swaits,sstatus);}

1519:   PetscFree4(len_s,len_si,sstatus,owners_co);
1520:   PetscFree(len_ri);
1521:   PetscFree(swaits);
1522:   PetscFree(buf_s);

1524:   /* (5) compute the local portion of C      */
1525:   /* ------------------------------------------ */
1526:   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of C */
1527:   PetscFreeSpaceGet(Crmax,&free_space);
1528:   current_space = free_space;

1530:   PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);
1531:   for (k=0; k<nrecv; k++) {
1532:     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1533:     nrows       = *buf_ri_k[k];
1534:     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1535:     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1536:   }

1538:   MatPreallocateInitialize(comm,pn,an,dnz,onz);
1539:   PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);
1540:   for (i=0; i<pn; i++) {
1541:     /* add C_loc into C */
1542:     nzi  = c_loc->i[i+1] - c_loc->i[i];
1543:     Jptr = c_loc->j + c_loc->i[i];
1544:     PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);

1546:     /* add received col data into lnk */
1547:     for (k=0; k<nrecv; k++) { /* k-th received message */
1548:       if (i == *nextrow[k]) { /* i-th row */
1549:         nzi  = *(nextci[k]+1) - *nextci[k];
1550:         Jptr = buf_rj[k] + *nextci[k];
1551:         PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);
1552:         nextrow[k]++; nextci[k]++;
1553:       }
1554:     }
1555:     nzi = lnk[0];

1557:     /* copy data into free space, then initialize lnk */
1558:     PetscLLCondensedClean(aN,nzi,current_space->array,lnk,lnkbt);
1559:     MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);
1560:   }
1561:   PetscFree3(buf_ri_k,nextrow,nextci);
1562:   PetscLLDestroy(lnk,lnkbt);
1563:   PetscFreeSpaceDestroy(free_space);

1565:   /* local sizes and preallocation */
1566:   MatSetSizes(C,pn,an,PETSC_DETERMINE,PETSC_DETERMINE);
1567:   if (P->cmap->bs > 0) {PetscLayoutSetBlockSize(C->rmap,P->cmap->bs);}
1568:   if (A->cmap->bs > 0) {PetscLayoutSetBlockSize(C->cmap,A->cmap->bs);}
1569:   MatMPIAIJSetPreallocation(C,0,dnz,0,onz);
1570:   MatPreallocateFinalize(dnz,onz);

1572:   /* members in merge */
1573:   PetscFree(id_r);
1574:   PetscFree(len_r);
1575:   PetscFree(buf_ri[0]);
1576:   PetscFree(buf_ri);
1577:   PetscFree(buf_rj[0]);
1578:   PetscFree(buf_rj);
1579:   PetscLayoutDestroy(&rowmap);

1581:   /* attach the supporting struct to C for reuse */
1582:   c = (Mat_MPIAIJ*)C->data;
1583:   c->ap         = ptap;
1584:   ptap->destroy = C->ops->destroy;

1586:   /* C is not ready for use - assembly will be done by MatPtAPNumeric() */
1587:   C->assembled        = PETSC_FALSE;
1588:   C->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
1589:   C->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;
1590:   return(0);
1591: }

1593: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
1594: {
1595:   PetscErrorCode    ierr;
1596:   Mat_MPIAIJ        *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1597:   Mat_SeqAIJ        *c_seq;
1598:   Mat_APMPI         *ptap = c->ap;
1599:   Mat               A_loc,C_loc,C_oth;
1600:   PetscInt          i,rstart,rend,cm,ncols,row;
1601:   const PetscInt    *cols;
1602:   const PetscScalar *vals;

1605:   if (!ptap->A_loc) {
1606:     MPI_Comm comm;
1607:     PetscObjectGetComm((PetscObject)C,&comm);
1608:     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
1609:   }

1611:   MatZeroEntries(C);

1613:   if (ptap->reuse == MAT_REUSE_MATRIX) {
1614:     /* These matrices are obtained in MatTransposeMatMultSymbolic() */
1615:     /* 1) get R = Pd^T, Ro = Po^T */
1616:     /*----------------------------*/
1617:     MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);
1618:     MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);

1620:     /* 2) compute numeric A_loc */
1621:     /*--------------------------*/
1622:     MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&ptap->A_loc);
1623:   }

1625:   /* 3) C_loc = Rd*A_loc, C_oth = Ro*A_loc */
1626:   A_loc = ptap->A_loc;
1627:   ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,A_loc,ptap->C_loc);
1628:   ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,A_loc,ptap->C_oth);
1629:   C_loc = ptap->C_loc;
1630:   C_oth = ptap->C_oth;

1632:   /* add C_loc and Co to to C */
1633:   MatGetOwnershipRange(C,&rstart,&rend);

1635:   /* C_loc -> C */
1636:   cm    = C_loc->rmap->N;
1637:   c_seq = (Mat_SeqAIJ*)C_loc->data;
1638:   cols = c_seq->j;
1639:   vals = c_seq->a;
1640:   for (i=0; i<cm; i++) {
1641:     ncols = c_seq->i[i+1] - c_seq->i[i];
1642:     row = rstart + i;
1643:     MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);
1644:     cols += ncols; vals += ncols;
1645:   }

1647:   /* Co -> C, off-processor part */
1648:   cm    = C_oth->rmap->N;
1649:   c_seq = (Mat_SeqAIJ*)C_oth->data;
1650:   cols  = c_seq->j;
1651:   vals  = c_seq->a;
1652:   for (i=0; i<cm; i++) {
1653:     ncols = c_seq->i[i+1] - c_seq->i[i];
1654:     row = p->garray[i];
1655:     MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);
1656:     cols += ncols; vals += ncols;
1657:   }
1658:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1659:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

1661:   ptap->reuse = MAT_REUSE_MATRIX;

1663:   /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */
1664:   if (ptap->freestruct) {
1665:     MatFreeIntermediateDataStructures(C);
1666:   }
1667:   return(0);
1668: }

1670: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1671: {
1672:   PetscErrorCode      ierr;
1673:   Mat_Merge_SeqsToMPI *merge;
1674:   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1675:   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1676:   Mat_APMPI           *ptap;
1677:   PetscInt            *adj;
1678:   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1679:   MatScalar           *ada,*ca,valtmp;
1680:   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1681:   MPI_Comm            comm;
1682:   PetscMPIInt         size,rank,taga,*len_s;
1683:   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1684:   PetscInt            **buf_ri,**buf_rj;
1685:   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1686:   MPI_Request         *s_waits,*r_waits;
1687:   MPI_Status          *status;
1688:   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1689:   PetscInt            *ai,*aj,*coi,*coj,*poJ,*pdJ;
1690:   Mat                 A_loc;
1691:   Mat_SeqAIJ          *a_loc;

1694:   PetscObjectGetComm((PetscObject)C,&comm);
1695:   MPI_Comm_size(comm,&size);
1696:   MPI_Comm_rank(comm,&rank);

1698:   ptap  = c->ap;
1699:   if (!ptap->A_loc) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
1700:   merge = ptap->merge;

1702:   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1703:   /*------------------------------------------*/
1704:   /* get data from symbolic products */
1705:   coi    = merge->coi; coj = merge->coj;
1706:   PetscCalloc1(coi[pon]+1,&coa);
1707:   bi     = merge->bi; bj = merge->bj;
1708:   owners = merge->rowmap->range;
1709:   PetscCalloc1(bi[cm]+1,&ba);

1711:   /* get A_loc by taking all local rows of A */
1712:   A_loc = ptap->A_loc;
1713:   MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);
1714:   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1715:   ai    = a_loc->i;
1716:   aj    = a_loc->j;

1718:   for (i=0; i<am; i++) {
1719:     anz = ai[i+1] - ai[i];
1720:     adj = aj + ai[i];
1721:     ada = a_loc->a + ai[i];

1723:     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1724:     /*-------------------------------------------------------------*/
1725:     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1726:     pnz = po->i[i+1] - po->i[i];
1727:     poJ = po->j + po->i[i];
1728:     pA  = po->a + po->i[i];
1729:     for (j=0; j<pnz; j++) {
1730:       row = poJ[j];
1731:       cj  = coj + coi[row];
1732:       ca  = coa + coi[row];
1733:       /* perform sparse axpy */
1734:       nexta  = 0;
1735:       valtmp = pA[j];
1736:       for (k=0; nexta<anz; k++) {
1737:         if (cj[k] == adj[nexta]) {
1738:           ca[k] += valtmp*ada[nexta];
1739:           nexta++;
1740:         }
1741:       }
1742:       PetscLogFlops(2.0*anz);
1743:     }

1745:     /* put the value into Cd (diagonal part) */
1746:     pnz = pd->i[i+1] - pd->i[i];
1747:     pdJ = pd->j + pd->i[i];
1748:     pA  = pd->a + pd->i[i];
1749:     for (j=0; j<pnz; j++) {
1750:       row = pdJ[j];
1751:       cj  = bj + bi[row];
1752:       ca  = ba + bi[row];
1753:       /* perform sparse axpy */
1754:       nexta  = 0;
1755:       valtmp = pA[j];
1756:       for (k=0; nexta<anz; k++) {
1757:         if (cj[k] == adj[nexta]) {
1758:           ca[k] += valtmp*ada[nexta];
1759:           nexta++;
1760:         }
1761:       }
1762:       PetscLogFlops(2.0*anz);
1763:     }
1764:   }

1766:   /* 3) send and recv matrix values coa */
1767:   /*------------------------------------*/
1768:   buf_ri = merge->buf_ri;
1769:   buf_rj = merge->buf_rj;
1770:   len_s  = merge->len_s;
1771:   PetscCommGetNewTag(comm,&taga);
1772:   PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);

1774:   PetscMalloc2(merge->nsend+1,&s_waits,size,&status);
1775:   for (proc=0,k=0; proc<size; proc++) {
1776:     if (!len_s[proc]) continue;
1777:     i    = merge->owners_co[proc];
1778:     MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);
1779:     k++;
1780:   }
1781:   if (merge->nrecv) {MPI_Waitall(merge->nrecv,r_waits,status);}
1782:   if (merge->nsend) {MPI_Waitall(merge->nsend,s_waits,status);}

1784:   PetscFree2(s_waits,status);
1785:   PetscFree(r_waits);
1786:   PetscFree(coa);

1788:   /* 4) insert local Cseq and received values into Cmpi */
1789:   /*----------------------------------------------------*/
1790:   PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);
1791:   for (k=0; k<merge->nrecv; k++) {
1792:     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1793:     nrows       = *(buf_ri_k[k]);
1794:     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1795:     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1796:   }

1798:   for (i=0; i<cm; i++) {
1799:     row  = owners[rank] + i; /* global row index of C_seq */
1800:     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1801:     ba_i = ba + bi[i];
1802:     bnz  = bi[i+1] - bi[i];
1803:     /* add received vals into ba */
1804:     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1805:       /* i-th row */
1806:       if (i == *nextrow[k]) {
1807:         cnz    = *(nextci[k]+1) - *nextci[k];
1808:         cj     = buf_rj[k] + *(nextci[k]);
1809:         ca     = abuf_r[k] + *(nextci[k]);
1810:         nextcj = 0;
1811:         for (j=0; nextcj<cnz; j++) {
1812:           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1813:             ba_i[j] += ca[nextcj++];
1814:           }
1815:         }
1816:         nextrow[k]++; nextci[k]++;
1817:         PetscLogFlops(2.0*cnz);
1818:       }
1819:     }
1820:     MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);
1821:   }
1822:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1823:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

1825:   PetscFree(ba);
1826:   PetscFree(abuf_r[0]);
1827:   PetscFree(abuf_r);
1828:   PetscFree3(buf_ri_k,nextrow,nextci);

1830:   if (ptap->freestruct) {
1831:     MatFreeIntermediateDataStructures(C);
1832:   }
1833:   return(0);
1834: }

1836: PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat C)
1837: {
1838:   PetscErrorCode      ierr;
1839:   Mat                 A_loc,POt,PDt;
1840:   Mat_APMPI           *ptap;
1841:   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1842:   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*a=(Mat_MPIAIJ*)A->data,*c;
1843:   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1844:   PetscInt            nnz;
1845:   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1846:   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1847:   MPI_Comm            comm;
1848:   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1849:   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1850:   PetscInt            len,proc,*dnz,*onz,*owners;
1851:   PetscInt            nzi,*bi,*bj;
1852:   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1853:   MPI_Request         *swaits,*rwaits;
1854:   MPI_Status          *sstatus,rstatus;
1855:   Mat_Merge_SeqsToMPI *merge;
1856:   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1857:   PetscReal           afill  =1.0,afill_tmp;
1858:   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax;
1859:   Mat_SeqAIJ          *a_loc,*pdt,*pot;
1860:   PetscTable          ta;
1861:   MatType             mtype;

1864:   PetscObjectGetComm((PetscObject)A,&comm);
1865:   /* check if matrix local sizes are compatible */
1866:   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != P (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);

1868:   MPI_Comm_size(comm,&size);
1869:   MPI_Comm_rank(comm,&rank);

1871:   /* create struct Mat_APMPI and attached it to C later */
1872:   PetscNew(&ptap);

1874:   /* get A_loc by taking all local rows of A */
1875:   MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);

1877:   ptap->A_loc = A_loc;
1878:   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1879:   ai          = a_loc->i;
1880:   aj          = a_loc->j;

1882:   /* determine symbolic Co=(p->B)^T*A - send to others */
1883:   /*----------------------------------------------------*/
1884:   MatTransposeSymbolic_SeqAIJ(p->A,&PDt);
1885:   pdt  = (Mat_SeqAIJ*)PDt->data;
1886:   pdti = pdt->i; pdtj = pdt->j;

1888:   MatTransposeSymbolic_SeqAIJ(p->B,&POt);
1889:   pot  = (Mat_SeqAIJ*)POt->data;
1890:   poti = pot->i; potj = pot->j;

1892:   /* then, compute symbolic Co = (p->B)^T*A */
1893:   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1894:                          >= (num of nonzero rows of C_seq) - pn */
1895:   PetscMalloc1(pon+1,&coi);
1896:   coi[0] = 0;

1898:   /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1899:   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am]));
1900:   PetscFreeSpaceGet(nnz,&free_space);
1901:   current_space = free_space;

1903:   /* create and initialize a linked list */
1904:   PetscTableCreate(A->cmap->n + a->B->cmap->N,aN,&ta);
1905:   MatRowMergeMax_SeqAIJ(a_loc,am,ta);
1906:   PetscTableGetCount(ta,&Armax);

1908:   PetscLLCondensedCreate_Scalable(Armax,&lnk);

1910:   for (i=0; i<pon; i++) {
1911:     pnz = poti[i+1] - poti[i];
1912:     ptJ = potj + poti[i];
1913:     for (j=0; j<pnz; j++) {
1914:       row  = ptJ[j]; /* row of A_loc == col of Pot */
1915:       anz  = ai[row+1] - ai[row];
1916:       Jptr = aj + ai[row];
1917:       /* add non-zero cols of AP into the sorted linked list lnk */
1918:       PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);
1919:     }
1920:     nnz = lnk[0];

1922:     /* If free space is not available, double the total space in the list */
1923:     if (current_space->local_remaining<nnz) {
1924:       PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);
1925:       nspacedouble++;
1926:     }

1928:     /* Copy data into free space, and zero out denserows */
1929:     PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);

1931:     current_space->array           += nnz;
1932:     current_space->local_used      += nnz;
1933:     current_space->local_remaining -= nnz;

1935:     coi[i+1] = coi[i] + nnz;
1936:   }

1938:   PetscMalloc1(coi[pon]+1,&coj);
1939:   PetscFreeSpaceContiguous(&free_space,coj);
1940:   PetscLLCondensedDestroy_Scalable(lnk); /* must destroy to get a new one for C */

1942:   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1943:   if (afill_tmp > afill) afill = afill_tmp;

1945:   /* send j-array (coj) of Co to other processors */
1946:   /*----------------------------------------------*/
1947:   /* determine row ownership */
1948:   PetscNew(&merge);
1949:   PetscLayoutCreate(comm,&merge->rowmap);

1951:   merge->rowmap->n  = pn;
1952:   merge->rowmap->bs = 1;

1954:   PetscLayoutSetUp(merge->rowmap);
1955:   owners = merge->rowmap->range;

1957:   /* determine the number of messages to send, their lengths */
1958:   PetscCalloc1(size,&len_si);
1959:   PetscCalloc1(size,&merge->len_s);

1961:   len_s        = merge->len_s;
1962:   merge->nsend = 0;

1964:   PetscMalloc1(size+2,&owners_co);

1966:   proc = 0;
1967:   for (i=0; i<pon; i++) {
1968:     while (prmap[i] >= owners[proc+1]) proc++;
1969:     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1970:     len_s[proc] += coi[i+1] - coi[i];
1971:   }

1973:   len          = 0; /* max length of buf_si[] */
1974:   owners_co[0] = 0;
1975:   for (proc=0; proc<size; proc++) {
1976:     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1977:     if (len_si[proc]) {
1978:       merge->nsend++;
1979:       len_si[proc] = 2*(len_si[proc] + 1);
1980:       len         += len_si[proc];
1981:     }
1982:   }

1984:   /* determine the number and length of messages to receive for coi and coj  */
1985:   PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);
1986:   PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);

1988:   /* post the Irecv and Isend of coj */
1989:   PetscCommGetNewTag(comm,&tagj);
1990:   PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);
1991:   PetscMalloc1(merge->nsend+1,&swaits);
1992:   for (proc=0, k=0; proc<size; proc++) {
1993:     if (!len_s[proc]) continue;
1994:     i    = owners_co[proc];
1995:     MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);
1996:     k++;
1997:   }

1999:   /* receives and sends of coj are complete */
2000:   PetscMalloc1(size,&sstatus);
2001:   for (i=0; i<merge->nrecv; i++) {
2002:     PetscMPIInt icompleted;
2003:     MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);
2004:   }
2005:   PetscFree(rwaits);
2006:   if (merge->nsend) {MPI_Waitall(merge->nsend,swaits,sstatus);}

2008:   /* add received column indices into table to update Armax */
2009:   /* Armax can be as large as aN if a P[row,:] is dense, see src/ksp/ksp/tutorials/ex56.c! */
2010:   for (k=0; k<merge->nrecv; k++) {/* k-th received message */
2011:     Jptr = buf_rj[k];
2012:     for (j=0; j<merge->len_r[k]; j++) {
2013:       PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);
2014:     }
2015:   }
2016:   PetscTableGetCount(ta,&Armax);
2017:   /* printf("Armax %d, an %d + Bn %d = %d, aN %d\n",Armax,A->cmap->n,a->B->cmap->N,A->cmap->n+a->B->cmap->N,aN); */

2019:   /* send and recv coi */
2020:   /*-------------------*/
2021:   PetscCommGetNewTag(comm,&tagi);
2022:   PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);
2023:   PetscMalloc1(len+1,&buf_s);
2024:   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
2025:   for (proc=0,k=0; proc<size; proc++) {
2026:     if (!len_s[proc]) continue;
2027:     /* form outgoing message for i-structure:
2028:          buf_si[0]:                 nrows to be sent
2029:                [1:nrows]:           row index (global)
2030:                [nrows+1:2*nrows+1]: i-structure index
2031:     */
2032:     /*-------------------------------------------*/
2033:     nrows       = len_si[proc]/2 - 1;
2034:     buf_si_i    = buf_si + nrows+1;
2035:     buf_si[0]   = nrows;
2036:     buf_si_i[0] = 0;
2037:     nrows       = 0;
2038:     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
2039:       nzi               = coi[i+1] - coi[i];
2040:       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
2041:       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
2042:       nrows++;
2043:     }
2044:     MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);
2045:     k++;
2046:     buf_si += len_si[proc];
2047:   }
2048:   i = merge->nrecv;
2049:   while (i--) {
2050:     PetscMPIInt icompleted;
2051:     MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);
2052:   }
2053:   PetscFree(rwaits);
2054:   if (merge->nsend) {MPI_Waitall(merge->nsend,swaits,sstatus);}
2055:   PetscFree(len_si);
2056:   PetscFree(len_ri);
2057:   PetscFree(swaits);
2058:   PetscFree(sstatus);
2059:   PetscFree(buf_s);

2061:   /* compute the local portion of C (mpi mat) */
2062:   /*------------------------------------------*/
2063:   /* allocate bi array and free space for accumulating nonzero column info */
2064:   PetscMalloc1(pn+1,&bi);
2065:   bi[0] = 0;

2067:   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
2068:   nnz           = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am])));
2069:   PetscFreeSpaceGet(nnz,&free_space);
2070:   current_space = free_space;

2072:   PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);
2073:   for (k=0; k<merge->nrecv; k++) {
2074:     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
2075:     nrows       = *buf_ri_k[k];
2076:     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
2077:     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
2078:   }

2080:   PetscLLCondensedCreate_Scalable(Armax,&lnk);
2081:   MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);
2082:   rmax = 0;
2083:   for (i=0; i<pn; i++) {
2084:     /* add pdt[i,:]*AP into lnk */
2085:     pnz = pdti[i+1] - pdti[i];
2086:     ptJ = pdtj + pdti[i];
2087:     for (j=0; j<pnz; j++) {
2088:       row  = ptJ[j];  /* row of AP == col of Pt */
2089:       anz  = ai[row+1] - ai[row];
2090:       Jptr = aj + ai[row];
2091:       /* add non-zero cols of AP into the sorted linked list lnk */
2092:       PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);
2093:     }

2095:     /* add received col data into lnk */
2096:     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
2097:       if (i == *nextrow[k]) { /* i-th row */
2098:         nzi  = *(nextci[k]+1) - *nextci[k];
2099:         Jptr = buf_rj[k] + *nextci[k];
2100:         PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);
2101:         nextrow[k]++; nextci[k]++;
2102:       }
2103:     }
2104:     nnz = lnk[0];

2106:     /* if free space is not available, make more free space */
2107:     if (current_space->local_remaining<nnz) {
2108:       PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);
2109:       nspacedouble++;
2110:     }
2111:     /* copy data into free space, then initialize lnk */
2112:     PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);
2113:     MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);

2115:     current_space->array           += nnz;
2116:     current_space->local_used      += nnz;
2117:     current_space->local_remaining -= nnz;

2119:     bi[i+1] = bi[i] + nnz;
2120:     if (nnz > rmax) rmax = nnz;
2121:   }
2122:   PetscFree3(buf_ri_k,nextrow,nextci);

2124:   PetscMalloc1(bi[pn]+1,&bj);
2125:   PetscFreeSpaceContiguous(&free_space,bj);
2126:   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
2127:   if (afill_tmp > afill) afill = afill_tmp;
2128:   PetscLLCondensedDestroy_Scalable(lnk);
2129:   PetscTableDestroy(&ta);

2131:   MatDestroy(&POt);
2132:   MatDestroy(&PDt);

2134:   /* create symbolic parallel matrix C - why cannot be assembled in Numeric part   */
2135:   /*-------------------------------------------------------------------------------*/
2136:   MatSetSizes(C,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);
2137:   MatSetBlockSizes(C,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));
2138:   MatGetType(A,&mtype);
2139:   MatSetType(C,mtype);
2140:   MatMPIAIJSetPreallocation(C,0,dnz,0,onz);
2141:   MatPreallocateFinalize(dnz,onz);
2142:   MatSetBlockSize(C,1);
2143:   for (i=0; i<pn; i++) {
2144:     row  = i + rstart;
2145:     nnz  = bi[i+1] - bi[i];
2146:     Jptr = bj + bi[i];
2147:     MatSetValues(C,1,&row,nnz,Jptr,NULL,INSERT_VALUES);
2148:   }
2149:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
2150:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
2151:   merge->bi        = bi;
2152:   merge->bj        = bj;
2153:   merge->coi       = coi;
2154:   merge->coj       = coj;
2155:   merge->buf_ri    = buf_ri;
2156:   merge->buf_rj    = buf_rj;
2157:   merge->owners_co = owners_co;

2159:   /* attach the supporting struct to C for reuse */
2160:   c = (Mat_MPIAIJ*)C->data;

2162:   c->ap       = ptap;
2163:   ptap->api   = NULL;
2164:   ptap->apj   = NULL;
2165:   ptap->merge = merge;
2166:   ptap->apa   = NULL;
2167:   ptap->destroy   = C->ops->destroy;
2168:   ptap->duplicate = C->ops->duplicate;

2170:   C->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
2171:   C->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
2172:   C->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;

2174: #if defined(PETSC_USE_INFO)
2175:   if (bi[pn] != 0) {
2176:     PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
2177:     PetscInfo1(C,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);
2178:   } else {
2179:     PetscInfo(C,"Empty matrix product\n");
2180:   }
2181: #endif
2182:   return(0);
2183: }

2185: /* ---------------------------------------------------------------- */
2186: static PetscErrorCode MatProductSymbolic_AtB_MPIAIJ_MPIAIJ(Mat C)
2187: {
2189:   Mat_Product    *product = C->product;
2190:   Mat            A=product->A,B=product->B;
2191:   PetscReal      fill=product->fill;
2192:   PetscBool      flg;

2195:   /* scalable */
2196:   PetscStrcmp(product->alg,"scalable",&flg);
2197:   if (flg) {
2198:     MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);
2199:     goto next;
2200:   }

2202:   /* nonscalable */
2203:   PetscStrcmp(product->alg,"nonscalable",&flg);
2204:   if (flg) {
2205:     MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);
2206:     goto next;
2207:   }

2209:   /* matmatmult */
2210:   PetscStrcmp(product->alg,"at*b",&flg);
2211:   if (flg) {
2212:     Mat         At;
2213:     Mat_APMPI   *ptap;
2214:     Mat_MPIAIJ  *c;
2215:     MatTranspose(A,MAT_INITIAL_MATRIX,&At);

2217:     MatMatMultSymbolic_MPIAIJ_MPIAIJ(At,B,fill,C);
2218:     c    = (Mat_MPIAIJ*)C->data;
2219:     ptap = c->ap;
2220:     if (ptap) {
2221:       ptap->Pt = At;
2222:       C->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;
2223:     }
2224:     C->ops->transposematmultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
2225:     goto next;
2226:   }

2228:   SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatProduct type is not supported");

2230: next:
2231:   C->ops->productnumeric = MatProductNumeric_AtB;

2233:   {
2234:     Mat_MPIAIJ *c  = (Mat_MPIAIJ*)C->data;
2235:     Mat_APMPI  *ap = c->ap;
2236:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatFreeIntermediateDataStructures","Mat");
2237:     ap->freestruct = PETSC_FALSE;
2238:     PetscOptionsBool("-mat_freeintermediatedatastructures","Free intermediate data structures", "MatFreeIntermediateDataStructures",ap->freestruct,&ap->freestruct, NULL);
2239:     PetscOptionsEnd();
2240:   }
2241:   return(0);
2242: }

2244: /* ---------------------------------------------------------------- */
2245: /* Set options for MatMatMultxxx_MPIAIJ_MPIAIJ */
2246: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_AB(Mat C)
2247: {
2249:   Mat_Product    *product = C->product;
2250:   Mat            A=product->A,B=product->B;
2251: #if defined(PETSC_HAVE_HYPRE)
2252:   const char     *algTypes[4] = {"scalable","nonscalable","seqmpi","hypre"};
2253:   PetscInt       nalg = 4;
2254: #else
2255:   const char     *algTypes[3] = {"scalable","nonscalable","seqmpi"};
2256:   PetscInt       nalg = 3;
2257: #endif
2258:   PetscInt       alg = 1; /* set nonscalable algorithm as default */
2259:   PetscBool      flg;
2260:   MPI_Comm       comm;

2263:   /* Check matrix local sizes */
2264:   PetscObjectGetComm((PetscObject)C,&comm);
2265:   if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,B->rmap->rstart,B->rmap->rend);

2267:   /* Set "nonscalable" as default algorithm */
2268:   PetscStrcmp(C->product->alg,"default",&flg);
2269:   if (flg) {
2270:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);

2272:     /* Set "scalable" as default if BN and local nonzeros of A and B are large */
2273:     if (B->cmap->N > 100000) { /* may switch to scalable algorithm as default */
2274:       MatInfo     Ainfo,Binfo;
2275:       PetscInt    nz_local;
2276:       PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;

2278:       MatGetInfo(A,MAT_LOCAL,&Ainfo);
2279:       MatGetInfo(B,MAT_LOCAL,&Binfo);
2280:       nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated);

2282:       if (B->cmap->N > product->fill*nz_local) alg_scalable_loc = PETSC_TRUE;
2283:       MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);

2285:       if (alg_scalable) {
2286:         alg  = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
2287:         MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2288:         PetscInfo2(B,"Use scalable algorithm, BN %D, fill*nz_allocated %g\n",B->cmap->N,product->fill*nz_local);
2289:       }
2290:     }
2291:   }

2293:   /* Get runtime option */
2294:   if (product->api_user) {
2295:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatMult","Mat");
2296:     PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2297:     PetscOptionsEnd();
2298:   } else {
2299:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_AB","Mat");
2300:     PetscOptionsEList("-matproduct_ab_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2301:     PetscOptionsEnd();
2302:   }
2303:   if (flg) {
2304:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2305:   }

2307:   C->ops->productsymbolic = MatProductSymbolic_AB_MPIAIJ_MPIAIJ;
2308:   return(0);
2309: }

2311: /* Set options for MatTransposeMatMultXXX_MPIAIJ_MPIAIJ */
2312: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_AtB(Mat C)
2313: {
2315:   Mat_Product    *product = C->product;
2316:   Mat            A=product->A,B=product->B;
2317:   const char     *algTypes[3] = {"scalable","nonscalable","at*b"};
2318:   PetscInt       nalg = 3;
2319:   PetscInt       alg = 1; /* set default algorithm  */
2320:   PetscBool      flg;
2321:   MPI_Comm       comm;

2324:   /* Check matrix local sizes */
2325:   PetscObjectGetComm((PetscObject)C,&comm);
2326:   if (A->rmap->rstart != B->rmap->rstart || A->rmap->rend != B->rmap->rend) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != B (%D,%D)",A->rmap->rstart,A->rmap->rend,B->rmap->rstart,B->rmap->rend);

2328:   /* Set default algorithm */
2329:   PetscStrcmp(C->product->alg,"default",&flg);
2330:   if (flg) {
2331:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2332:   }

2334:   /* Set "scalable" as default if BN and local nonzeros of A and B are large */
2335:   if (alg && B->cmap->N > 100000) { /* may switch to scalable algorithm as default */
2336:     MatInfo     Ainfo,Binfo;
2337:     PetscInt    nz_local;
2338:     PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;

2340:     MatGetInfo(A,MAT_LOCAL,&Ainfo);
2341:     MatGetInfo(B,MAT_LOCAL,&Binfo);
2342:     nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated);

2344:     if (B->cmap->N > product->fill*nz_local) alg_scalable_loc = PETSC_TRUE;
2345:     MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);

2347:     if (alg_scalable) {
2348:       alg  = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
2349:       MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2350:       PetscInfo2(B,"Use scalable algorithm, BN %D, fill*nz_allocated %g\n",B->cmap->N,product->fill*nz_local);
2351:     }
2352:   }

2354:   /* Get runtime option */
2355:   if (product->api_user) {
2356:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatTransposeMatMult","Mat");
2357:     PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2358:     PetscOptionsEnd();
2359:   } else {
2360:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_AtB","Mat");
2361:     PetscOptionsEList("-matproduct_atb_via","Algorithmic approach","MatTransposeMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2362:     PetscOptionsEnd();
2363:   }
2364:   if (flg) {
2365:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2366:   }

2368:   C->ops->productsymbolic = MatProductSymbolic_AtB_MPIAIJ_MPIAIJ;
2369:   return(0);
2370: }

2372: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_PtAP(Mat C)
2373: {
2375:   Mat_Product    *product = C->product;
2376:   Mat            A=product->A,P=product->B;
2377:   MPI_Comm       comm;
2378:   PetscBool      flg;
2379:   PetscInt       alg=1; /* set default algorithm */
2380: #if !defined(PETSC_HAVE_HYPRE)
2381:   const char     *algTypes[4] = {"scalable","nonscalable","allatonce","allatonce_merged"};
2382:   PetscInt       nalg=4;
2383: #else
2384:   const char     *algTypes[5] = {"scalable","nonscalable","allatonce","allatonce_merged","hypre"};
2385:   PetscInt       nalg=5;
2386: #endif
2387:   PetscInt       pN=P->cmap->N;

2390:   /* Check matrix local sizes */
2391:   PetscObjectGetComm((PetscObject)C,&comm);
2392:   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Arow (%D, %D) != Prow (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
2393:   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Acol (%D, %D) != Prow (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend);

2395:   /* Set "nonscalable" as default algorithm */
2396:   PetscStrcmp(C->product->alg,"default",&flg);
2397:   if (flg) {
2398:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);

2400:     /* Set "scalable" as default if BN and local nonzeros of A and B are large */
2401:     if (pN > 100000) {
2402:       MatInfo     Ainfo,Pinfo;
2403:       PetscInt    nz_local;
2404:       PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;

2406:       MatGetInfo(A,MAT_LOCAL,&Ainfo);
2407:       MatGetInfo(P,MAT_LOCAL,&Pinfo);
2408:       nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated);

2410:       if (pN > product->fill*nz_local) alg_scalable_loc = PETSC_TRUE;
2411:       MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);

2413:       if (alg_scalable) {
2414:         alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
2415:         MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2416:       }
2417:     }
2418:   }

2420:   /* Get runtime option */
2421:   if (product->api_user) {
2422:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatPtAP","Mat");
2423:     PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);
2424:     PetscOptionsEnd();
2425:   } else {
2426:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_PtAP","Mat");
2427:     PetscOptionsEList("-matproduct_ptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);
2428:     PetscOptionsEnd();
2429:   }
2430:   if (flg) {
2431:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2432:   }

2434:   C->ops->productsymbolic = MatProductSymbolic_PtAP_MPIAIJ_MPIAIJ;
2435:   return(0);
2436: }

2438: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_RARt(Mat C)
2439: {
2440:   Mat_Product *product = C->product;
2441:   Mat         A = product->A,R=product->B;

2444:   /* Check matrix local sizes */
2445:   if (A->cmap->n != R->cmap->n || A->rmap->n != R->cmap->n) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A local (%D, %D), R local (%D,%D)",A->rmap->n,A->rmap->n,R->rmap->n,R->cmap->n);

2447:   C->ops->productsymbolic = MatProductSymbolic_RARt_MPIAIJ_MPIAIJ;
2448:   return(0);
2449: }

2451: /*
2452:  Set options for ABC = A*B*C = A*(B*C); ABC's algorithm must be chosen from AB's algorithm
2453: */
2454: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_ABC(Mat C)
2455: {
2457:   Mat_Product    *product = C->product;
2458:   PetscBool      flg = PETSC_FALSE;
2459:   PetscInt       alg = 1; /* default algorithm */
2460:   const char     *algTypes[3] = {"scalable","nonscalable","seqmpi"};
2461:   PetscInt       nalg = 3;

2464:   /* Set default algorithm */
2465:   PetscStrcmp(C->product->alg,"default",&flg);
2466:   if (flg) {
2467:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2468:   }

2470:   /* Get runtime option */
2471:   if (product->api_user) {
2472:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatMatMult","Mat");
2473:     PetscOptionsEList("-matmatmatmult_via","Algorithmic approach","MatMatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2474:     PetscOptionsEnd();
2475:   } else {
2476:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_ABC","Mat");
2477:     PetscOptionsEList("-matproduct_abc_via","Algorithmic approach","MatProduct_ABC",algTypes,nalg,algTypes[alg],&alg,&flg);
2478:     PetscOptionsEnd();
2479:   }
2480:   if (flg) {
2481:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2482:   }

2484:   C->ops->matmatmultsymbolic = MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ;
2485:   C->ops->productsymbolic    = MatProductSymbolic_ABC;
2486:   return(0);
2487: }

2489: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ(Mat C)
2490: {
2492:   Mat_Product    *product = C->product;

2495:   switch (product->type) {
2496:   case MATPRODUCT_AB:
2497:     MatProductSetFromOptions_MPIAIJ_AB(C);
2498:     break;
2499:   case MATPRODUCT_AtB:
2500:     MatProductSetFromOptions_MPIAIJ_AtB(C);
2501:     break;
2502:   case MATPRODUCT_PtAP:
2503:     MatProductSetFromOptions_MPIAIJ_PtAP(C);
2504:     break;
2505:   case MATPRODUCT_RARt:
2506:     MatProductSetFromOptions_MPIAIJ_RARt(C);
2507:     break;
2508:   case MATPRODUCT_ABC:
2509:     MatProductSetFromOptions_MPIAIJ_ABC(C);
2510:     break;
2511:   default: SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_SUP,"MatProduct type is not supported");
2512:   }
2513:   return(0);
2514: }