Actual source code: mpimatmatmult.c


  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/sfimpl.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:     return(0);
 33:   }

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

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

 49:   /* backend general code */
 50:   PetscStrcmp(alg,"backend",&flg);
 51:   if (flg) {
 52:     MatProductSymbolic_MPIAIJBACKEND(C);
 53:     return(0);
 54:   }

 56: #if defined(PETSC_HAVE_HYPRE)
 57:   PetscStrcmp(alg,"hypre",&flg);
 58:   if (flg) {
 59:     MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);
 60:     return(0);
 61:   }
 62: #endif
 63:   SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_SUP,"Mat Product Algorithm is not supported");
 64: }

 66: PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(void *data)
 67: {
 69:   Mat_APMPI      *ptap = (Mat_APMPI*)data;

 72:   PetscFree2(ptap->startsj_s,ptap->startsj_r);
 73:   PetscFree(ptap->bufa);
 74:   MatDestroy(&ptap->P_loc);
 75:   MatDestroy(&ptap->P_oth);
 76:   MatDestroy(&ptap->Pt);
 77:   PetscFree(ptap->api);
 78:   PetscFree(ptap->apj);
 79:   PetscFree(ptap->apa);
 80:   PetscFree(ptap);
 81:   return(0);
 82: }

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

103:   MatCheckProduct(C,3);
104:   ptap = (Mat_APMPI*)C->product->data;
105:   if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data");
106:   PetscObjectGetComm((PetscObject)A,&comm);
107:   MPI_Comm_size(comm,&size);

109:   if (!ptap->P_oth && size>1) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"AP cannot be reused. Do not call MatProductClear()");

111:   /* flag CPU mask for C */
112: #if defined(PETSC_HAVE_DEVICE)
113:   if (C->offloadmask != PETSC_OFFLOAD_UNALLOCATED) C->offloadmask = PETSC_OFFLOAD_CPU;
114:   if (c->A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) c->A->offloadmask = PETSC_OFFLOAD_CPU;
115:   if (c->B->offloadmask != PETSC_OFFLOAD_UNALLOCATED) c->B->offloadmask = PETSC_OFFLOAD_CPU;
116: #endif

118:   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
119:   /*-----------------------------------------------------*/
120:   /* update numerical values of P_oth and P_loc */
121:   MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);
122:   MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);

124:   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
125:   /*----------------------------------------------------------*/
126:   /* get data from symbolic products */
127:   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
128:   p_oth = NULL;
129:   if (size >1) {
130:     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
131:   }

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

136:   api = ptap->api;
137:   apj = ptap->apj;
138:   /* trigger copy to CPU */
139:   MatSeqAIJGetArrayRead(a->A,&dummy);
140:   MatSeqAIJRestoreArrayRead(a->A,&dummy);
141:   MatSeqAIJGetArrayRead(a->B,&dummy);
142:   MatSeqAIJRestoreArrayRead(a->B,&dummy);
143:   for (i=0; i<cm; i++) {
144:     /* compute apa = A[i,:]*P */
145:     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);

147:     /* set values in C */
148:     apJ  = apj + api[i];
149:     cdnz = cd->i[i+1] - cd->i[i];
150:     conz = co->i[i+1] - co->i[i];

152:     /* 1st off-diagonal part of C */
153:     ca = coa + co->i[i];
154:     k  = 0;
155:     for (k0=0; k0<conz; k0++) {
156:       if (apJ[k] >= cstart) break;
157:       ca[k0]      = apa[apJ[k]];
158:       apa[apJ[k++]] = 0.0;
159:     }

161:     /* diagonal part of C */
162:     ca = cda + cd->i[i];
163:     for (k1=0; k1<cdnz; k1++) {
164:       ca[k1]      = apa[apJ[k]];
165:       apa[apJ[k++]] = 0.0;
166:     }

168:     /* 2nd off-diagonal part of C */
169:     ca = coa + co->i[i];
170:     for (; k0<conz; k0++) {
171:       ca[k0]      = apa[apJ[k]];
172:       apa[apJ[k++]] = 0.0;
173:     }
174:   }
175:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
176:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
177:   return(0);
178: }

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

198:   MatCheckProduct(C,4);
199:   if (C->product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty");
200:   PetscObjectGetComm((PetscObject)A,&comm);
201:   MPI_Comm_size(comm,&size);

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

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

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

212:   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
213:   pi_loc = p_loc->i; pj_loc = p_loc->j;
214:   if (size > 1) {
215:     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
216:     pi_oth = p_oth->i; pj_oth = p_oth->j;
217:   } else {
218:     p_oth = NULL;
219:     pi_oth = NULL; pj_oth = NULL;
220:   }

222:   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
223:   /*-------------------------------------------------------------------*/
224:   PetscMalloc1(am+2,&api);
225:   ptap->api = api;
226:   api[0]    = 0;

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

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

235:   MatPreallocateInitialize(comm,am,pn,dnz,onz);
236:   for (i=0; i<am; i++) {
237:     /* diagonal portion of A */
238:     nzi = adi[i+1] - adi[i];
239:     for (j=0; j<nzi; j++) {
240:       row  = *adj++;
241:       pnz  = pi_loc[row+1] - pi_loc[row];
242:       Jptr = pj_loc + pi_loc[row];
243:       /* add non-zero cols of P into the sorted linked list lnk */
244:       PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);
245:     }
246:     /* off-diagonal portion of A */
247:     nzi = aoi[i+1] - aoi[i];
248:     for (j=0; j<nzi; j++) {
249:       row  = *aoj++;
250:       pnz  = pi_oth[row+1] - pi_oth[row];
251:       Jptr = pj_oth + pi_oth[row];
252:       PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);
253:     }
254:     /* add possible missing diagonal entry */
255:     if (C->force_diagonals) {
256:       j = i + rstart; /* column index */
257:       PetscLLCondensedAddSorted(1,&j,lnk,lnkbt);
258:     }

260:     apnz     = lnk[0];
261:     api[i+1] = api[i] + apnz;

263:     /* if free space is not available, double the total space in the list */
264:     if (current_space->local_remaining<apnz) {
265:       PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);
266:       nspacedouble++;
267:     }

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

273:     current_space->array           += apnz;
274:     current_space->local_used      += apnz;
275:     current_space->local_remaining -= apnz;
276:   }

278:   /* Allocate space for apj, initialize apj, and */
279:   /* destroy list of free space and other temporary array(s) */
280:   PetscMalloc1(api[am]+1,&ptap->apj);
281:   apj  = ptap->apj;
282:   PetscFreeSpaceContiguous(&free_space,ptap->apj);
283:   PetscLLDestroy(lnk,lnkbt);

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

288:   /* set and assemble symbolic parallel matrix C */
289:   /*---------------------------------------------*/
290:   MatSetSizes(C,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);
291:   MatSetBlockSizesFromMats(C,A,P);

293:   MatGetType(A,&mtype);
294:   MatSetType(C,mtype);
295:   MatMPIAIJSetPreallocation(C,0,dnz,0,onz);
296:   MatPreallocateFinalize(dnz,onz);

298:   MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(C, apj, api);
299:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
300:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
301:   MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);

303:   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
304:   C->ops->productnumeric = MatProductNumeric_AB;

306:   /* attach the supporting struct to C for reuse */
307:   C->product->data    = ptap;
308:   C->product->destroy = MatDestroy_MPIAIJ_MatMatMult;

310:   /* set MatInfo */
311:   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
312:   if (afill < 1.0) afill = 1.0;
313:   C->info.mallocs           = nspacedouble;
314:   C->info.fill_ratio_given  = fill;
315:   C->info.fill_ratio_needed = afill;

317: #if defined(PETSC_USE_INFO)
318:   if (api[am]) {
319:     PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
320:     PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
321:   } else {
322:     PetscInfo(C,"Empty matrix product\n");
323:   }
324: #endif
325:   return(0);
326: }

328: /* ------------------------------------------------------- */
329: static PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat);
330: static PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat);

332: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense_AB(Mat C)
333: {
334:   Mat_Product *product = C->product;
335:   Mat         A = product->A,B=product->B;

338:   if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend)
339:     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);

341:   C->ops->matmultsymbolic = MatMatMultSymbolic_MPIAIJ_MPIDense;
342:   C->ops->productsymbolic = MatProductSymbolic_AB;
343:   return(0);
344: }
345: /* -------------------------------------------------------------------- */
346: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense_AtB(Mat C)
347: {
348:   Mat_Product *product = C->product;
349:   Mat         A = product->A,B=product->B;

352:   if (A->rmap->rstart != B->rmap->rstart || A->rmap->rend != B->rmap->rend)
353:     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);

355:   C->ops->transposematmultsymbolic = MatTransposeMatMultSymbolic_MPIAIJ_MPIDense;
356:   C->ops->productsymbolic          = MatProductSymbolic_AtB;
357:   return(0);
358: }

360: /* --------------------------------------------------------------------- */
361: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense(Mat C)
362: {
364:   Mat_Product    *product = C->product;

367:   switch (product->type) {
368:   case MATPRODUCT_AB:
369:     MatProductSetFromOptions_MPIAIJ_MPIDense_AB(C);
370:     break;
371:   case MATPRODUCT_AtB:
372:     MatProductSetFromOptions_MPIAIJ_MPIDense_AtB(C);
373:     break;
374:   default:
375:     break;
376:   }
377:   return(0);
378: }
379: /* ------------------------------------------------------- */

381: typedef struct {
382:   Mat          workB,workB1;
383:   MPI_Request  *rwaits,*swaits;
384:   PetscInt     nsends,nrecvs;
385:   MPI_Datatype *stype,*rtype;
386:   PetscInt     blda;
387: } MPIAIJ_MPIDense;

389: PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
390: {
391:   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*)ctx;
392:   PetscErrorCode  ierr;
393:   PetscInt        i;

396:   MatDestroy(&contents->workB);
397:   MatDestroy(&contents->workB1);
398:   for (i=0; i<contents->nsends; i++) {
399:     MPI_Type_free(&contents->stype[i]);
400:   }
401:   for (i=0; i<contents->nrecvs; i++) {
402:     MPI_Type_free(&contents->rtype[i]);
403:   }
404:   PetscFree4(contents->stype,contents->rtype,contents->rwaits,contents->swaits);
405:   PetscFree(contents);
406:   return(0);
407: }

409: static PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat C)
410: {
411:   PetscErrorCode  ierr;
412:   Mat_MPIAIJ      *aij=(Mat_MPIAIJ*)A->data;
413:   PetscInt        nz=aij->B->cmap->n,nsends,nrecvs,i,nrows_to,j,blda,clda,m,M,n,N;
414:   MPIAIJ_MPIDense *contents;
415:   VecScatter      ctx=aij->Mvctx;
416:   PetscInt        Am=A->rmap->n,Bm=B->rmap->n,BN=B->cmap->N,Bbn,Bbn1,bs,nrows_from,numBb;
417:   MPI_Comm        comm;
418:   MPI_Datatype    type1,*stype,*rtype;
419:   const PetscInt  *sindices,*sstarts,*rstarts;
420:   PetscMPIInt     *disp;
421:   PetscBool       cisdense;

424:   MatCheckProduct(C,4);
425:   if (C->product->data) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty");
426:   PetscObjectGetComm((PetscObject)A,&comm);
427:   PetscObjectBaseTypeCompare((PetscObject)C,MATMPIDENSE,&cisdense);
428:   if (!cisdense) {
429:     MatSetType(C,((PetscObject)B)->type_name);
430:   }
431:   MatGetLocalSize(C,&m,&n);
432:   MatGetSize(C,&M,&N);
433:   if (m == PETSC_DECIDE || n == PETSC_DECIDE || M == PETSC_DECIDE || N == PETSC_DECIDE) {
434:     MatSetSizes(C,Am,B->cmap->n,A->rmap->N,BN);
435:   }
436:   MatSetBlockSizesFromMats(C,A,B);
437:   MatSetUp(C);
438:   MatDenseGetLDA(B,&blda);
439:   MatDenseGetLDA(C,&clda);
440:   PetscNew(&contents);

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

445:   /* Create column block of B and C for memory scalability when BN is too large */
446:   /* Estimate Bbn, column size of Bb */
447:   if (nz) {
448:     Bbn1 = 2*Am*BN/nz;
449:     if (!Bbn1) Bbn1 = 1;
450:   } else Bbn1 = BN;

452:   bs = PetscAbs(B->cmap->bs);
453:   Bbn1 = Bbn1/bs *bs; /* Bbn1 is a multiple of bs */
454:   if (Bbn1 > BN) Bbn1 = BN;
455:   MPI_Allreduce(&Bbn1,&Bbn,1,MPIU_INT,MPI_MAX,comm);

457:   /* Enable runtime option for Bbn */
458:   PetscOptionsBegin(comm,((PetscObject)C)->prefix,"MatMatMult","Mat");
459:   PetscOptionsInt("-matmatmult_Bbn","Number of columns in Bb","MatMatMult",Bbn,&Bbn,NULL);
460:   PetscOptionsEnd();
461:   Bbn  = PetscMin(Bbn,BN);

463:   if (Bbn > 0 && Bbn < BN) {
464:     numBb = BN/Bbn;
465:     Bbn1 = BN - numBb*Bbn;
466:   } else numBb = 0;

468:   if (numBb) {
469:     PetscInfo3(C,"use Bb, BN=%D, Bbn=%D; numBb=%D\n",BN,Bbn,numBb);
470:     if (Bbn1) { /* Create workB1 for the remaining columns */
471:       PetscInfo2(C,"use Bb1, BN=%D, Bbn1=%D\n",BN,Bbn1);
472:       /* Create work matrix used to store off processor rows of B needed for local product */
473:       MatCreateSeqDense(PETSC_COMM_SELF,nz,Bbn1,NULL,&contents->workB1);
474:     } else contents->workB1 = NULL;
475:   }

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

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

485:   contents->rtype  = rtype;
486:   contents->nrecvs = nrecvs;
487:   contents->blda   = blda;

489:   PetscMalloc1(Bm+1,&disp);
490:   for (i=0; i<nsends; i++) {
491:     nrows_to = sstarts[i+1]-sstarts[i];
492:     for (j=0; j<nrows_to; j++) {
493:       disp[j] = sindices[sstarts[i]+j]; /* rowB to be sent */
494:     }
495:     MPI_Type_create_indexed_block(nrows_to,1,(const PetscMPIInt *)disp,MPIU_SCALAR,&type1);

497:     MPI_Type_create_resized(type1,0,blda*sizeof(PetscScalar),&stype[i]);
498:     MPI_Type_commit(&stype[i]);
499:     MPI_Type_free(&type1);
500:   }

502:   for (i=0; i<nrecvs; i++) {
503:     /* received values from a process form a (nrows_from x Bbn) row block in workB (column-wise) */
504:     nrows_from = rstarts[i+1]-rstarts[i];
505:     disp[0] = 0;
506:     MPI_Type_create_indexed_block(1, nrows_from, (const PetscMPIInt *)disp, MPIU_SCALAR, &type1);
507:     MPI_Type_create_resized(type1, 0, nz*sizeof(PetscScalar), &rtype[i]);
508:     MPI_Type_commit(&rtype[i]);
509:     MPI_Type_free(&type1);
510:   }

512:   PetscFree(disp);
513:   VecScatterRestoreRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,NULL,NULL);
514:   VecScatterRestoreRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,NULL,NULL);
515:   MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);
516:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
517:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
518:   MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);

520:   C->product->data = contents;
521:   C->product->destroy = MatMPIAIJ_MPIDenseDestroy;
522:   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
523:   return(0);
524: }

526: PETSC_INTERN PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat,const PetscBool);
527: /*
528:     Performs an efficient scatter on the rows of B needed by this process; this is
529:     a modification of the VecScatterBegin_() routines.

531:     Input: Bbidx = 0: B = Bb
532:                  = 1: B = Bb1, see MatMatMultSymbolic_MPIAIJ_MPIDense()
533: */
534: PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,PetscInt Bbidx,Mat C,Mat *outworkB)
535: {
536:   Mat_MPIAIJ        *aij = (Mat_MPIAIJ*)A->data;
537:   PetscErrorCode    ierr;
538:   const PetscScalar *b;
539:   PetscScalar       *rvalues;
540:   VecScatter        ctx = aij->Mvctx;
541:   const PetscInt    *sindices,*sstarts,*rstarts;
542:   const PetscMPIInt *sprocs,*rprocs;
543:   PetscInt          i,nsends,nrecvs;
544:   MPI_Request       *swaits,*rwaits;
545:   MPI_Comm          comm;
546:   PetscMPIInt       tag=((PetscObject)ctx)->tag,ncols=B->cmap->N,nrows=aij->B->cmap->n,nsends_mpi,nrecvs_mpi;
547:   MPIAIJ_MPIDense   *contents;
548:   Mat               workB;
549:   MPI_Datatype      *stype,*rtype;
550:   PetscInt          blda;

553:   MatCheckProduct(C,4);
554:   if (!C->product->data) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty");
555:   contents = (MPIAIJ_MPIDense*)C->product->data;
556:   VecScatterGetRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,&sprocs,NULL/*bs*/);
557:   VecScatterGetRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,&rprocs,NULL/*bs*/);
558:   PetscMPIIntCast(nsends,&nsends_mpi);
559:   PetscMPIIntCast(nrecvs,&nrecvs_mpi);
560:   if (Bbidx == 0) {
561:     workB = *outworkB = contents->workB;
562:   } else {
563:     workB = *outworkB = contents->workB1;
564:   }
565:   if (nrows != workB->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",workB->cmap->n,nrows);
566:   swaits = contents->swaits;
567:   rwaits = contents->rwaits;

569:   MatDenseGetArrayRead(B,&b);
570:   MatDenseGetLDA(B,&blda);
571:   if (blda != contents->blda) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot reuse an input matrix with lda %D != %D",blda,contents->blda);
572:   MatDenseGetArray(workB,&rvalues);

574:   /* Post recv, use MPI derived data type to save memory */
575:   PetscObjectGetComm((PetscObject)C,&comm);
576:   rtype = contents->rtype;
577:   for (i=0; i<nrecvs; i++) {
578:     MPI_Irecv(rvalues+(rstarts[i]-rstarts[0]),ncols,rtype[i],rprocs[i],tag,comm,rwaits+i);
579:   }

581:   stype = contents->stype;
582:   for (i=0; i<nsends; i++) {
583:     MPI_Isend(b,ncols,stype[i],sprocs[i],tag,comm,swaits+i);
584:   }

586:   if (nrecvs) {MPI_Waitall(nrecvs_mpi,rwaits,MPI_STATUSES_IGNORE);}
587:   if (nsends) {MPI_Waitall(nsends_mpi,swaits,MPI_STATUSES_IGNORE);}

589:   VecScatterRestoreRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,&sprocs,NULL);
590:   VecScatterRestoreRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,&rprocs,NULL);
591:   MatDenseRestoreArrayRead(B,&b);
592:   MatDenseRestoreArray(workB,&rvalues);
593:   return(0);
594: }

596: static PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
597: {
598:   PetscErrorCode  ierr;
599:   Mat_MPIAIJ      *aij    = (Mat_MPIAIJ*)A->data;
600:   Mat_MPIDense    *bdense = (Mat_MPIDense*)B->data;
601:   Mat_MPIDense    *cdense = (Mat_MPIDense*)C->data;
602:   Mat             workB;
603:   MPIAIJ_MPIDense *contents;

606:   MatCheckProduct(C,3);
607:   if (!C->product->data) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty");
608:   contents = (MPIAIJ_MPIDense*)C->product->data;
609:   /* diagonal block of A times all local rows of B */
610:   /* TODO: this calls a symbolic multiplication every time, which could be avoided */
611:   MatMatMult(aij->A,bdense->A,MAT_REUSE_MATRIX,PETSC_DEFAULT,&cdense->A);
612:   if (contents->workB->cmap->n == B->cmap->N) {
613:     /* get off processor parts of B needed to complete C=A*B */
614:     MatMPIDenseScatter(A,B,0,C,&workB);

616:     /* off-diagonal block of A times nonlocal rows of B */
617:     MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A,PETSC_TRUE);
618:   } else {
619:     Mat      Bb,Cb;
620:     PetscInt BN=B->cmap->N,n=contents->workB->cmap->n,i;
621:     if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Column block size %D must be positive",n);

623:     for (i=0; i<BN; i+=n) {
624:       MatDenseGetSubMatrix(B,i,PetscMin(i+n,BN),&Bb);
625:       MatDenseGetSubMatrix(C,i,PetscMin(i+n,BN),&Cb);

627:       /* get off processor parts of B needed to complete C=A*B */
628:       MatMPIDenseScatter(A,Bb,i+n>BN,C,&workB);

630:       /* off-diagonal block of A times nonlocal rows of B */
631:       cdense = (Mat_MPIDense*)Cb->data;
632:       MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A,PETSC_TRUE);

634:       MatDenseRestoreSubMatrix(B,&Bb);
635:       MatDenseRestoreSubMatrix(C,&Cb);
636:     }
637:   }
638:   return(0);
639: }

641: PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
642: {
643:   PetscErrorCode    ierr;
644:   Mat_MPIAIJ        *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
645:   Mat_SeqAIJ        *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
646:   Mat_SeqAIJ        *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
647:   PetscInt          *adi = ad->i,*adj,*aoi=ao->i,*aoj;
648:   PetscScalar       *ada,*aoa,*cda=cd->a,*coa=co->a;
649:   Mat_SeqAIJ        *p_loc,*p_oth;
650:   PetscInt          *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
651:   PetscScalar       *pa_loc,*pa_oth,*pa,valtmp,*ca;
652:   PetscInt          cm    = C->rmap->n,anz,pnz;
653:   Mat_APMPI         *ptap;
654:   PetscScalar       *apa_sparse;
655:   const PetscScalar *dummy;
656:   PetscInt          *api,*apj,*apJ,i,j,k,row;
657:   PetscInt          cstart = C->cmap->rstart;
658:   PetscInt          cdnz,conz,k0,k1,nextp;
659:   MPI_Comm          comm;
660:   PetscMPIInt       size;

663:   MatCheckProduct(C,3);
664:   ptap = (Mat_APMPI*)C->product->data;
665:   if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data");
666:   PetscObjectGetComm((PetscObject)C,&comm);
667:   MPI_Comm_size(comm,&size);
668:   if (!ptap->P_oth && size>1) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"AP cannot be reused. Do not call MatProductClear()");

670:   /* flag CPU mask for C */
671: #if defined(PETSC_HAVE_DEVICE)
672:   if (C->offloadmask != PETSC_OFFLOAD_UNALLOCATED) C->offloadmask = PETSC_OFFLOAD_CPU;
673:   if (c->A->offloadmask != PETSC_OFFLOAD_UNALLOCATED) c->A->offloadmask = PETSC_OFFLOAD_CPU;
674:   if (c->B->offloadmask != PETSC_OFFLOAD_UNALLOCATED) c->B->offloadmask = PETSC_OFFLOAD_CPU;
675: #endif
676:   apa_sparse = ptap->apa;

678:   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
679:   /*-----------------------------------------------------*/
680:   /* update numerical values of P_oth and P_loc */
681:   MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);
682:   MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);

684:   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
685:   /*----------------------------------------------------------*/
686:   /* get data from symbolic products */
687:   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
688:   pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a;
689:   if (size >1) {
690:     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
691:     pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a;
692:   } else {
693:     p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL;
694:   }

696:   /* trigger copy to CPU */
697:   MatSeqAIJGetArrayRead(a->A,&dummy);
698:   MatSeqAIJRestoreArrayRead(a->A,&dummy);
699:   MatSeqAIJGetArrayRead(a->B,&dummy);
700:   MatSeqAIJRestoreArrayRead(a->B,&dummy);
701:   api = ptap->api;
702:   apj = ptap->apj;
703:   for (i=0; i<cm; i++) {
704:     apJ = apj + api[i];

706:     /* diagonal portion of A */
707:     anz = adi[i+1] - adi[i];
708:     adj = ad->j + adi[i];
709:     ada = ad->a + adi[i];
710:     for (j=0; j<anz; j++) {
711:       row = adj[j];
712:       pnz = pi_loc[row+1] - pi_loc[row];
713:       pj  = pj_loc + pi_loc[row];
714:       pa  = pa_loc + pi_loc[row];
715:       /* perform sparse axpy */
716:       valtmp = ada[j];
717:       nextp  = 0;
718:       for (k=0; nextp<pnz; k++) {
719:         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
720:           apa_sparse[k] += valtmp*pa[nextp++];
721:         }
722:       }
723:       PetscLogFlops(2.0*pnz);
724:     }

726:     /* off-diagonal portion of A */
727:     anz = aoi[i+1] - aoi[i];
728:     aoj = ao->j + aoi[i];
729:     aoa = ao->a + aoi[i];
730:     for (j=0; j<anz; j++) {
731:       row = aoj[j];
732:       pnz = pi_oth[row+1] - pi_oth[row];
733:       pj  = pj_oth + pi_oth[row];
734:       pa  = pa_oth + pi_oth[row];
735:       /* perform sparse axpy */
736:       valtmp = aoa[j];
737:       nextp  = 0;
738:       for (k=0; nextp<pnz; k++) {
739:         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
740:           apa_sparse[k] += valtmp*pa[nextp++];
741:         }
742:       }
743:       PetscLogFlops(2.0*pnz);
744:     }

746:     /* set values in C */
747:     cdnz = cd->i[i+1] - cd->i[i];
748:     conz = co->i[i+1] - co->i[i];

750:     /* 1st off-diagonal part of C */
751:     ca = coa + co->i[i];
752:     k  = 0;
753:     for (k0=0; k0<conz; k0++) {
754:       if (apJ[k] >= cstart) break;
755:       ca[k0]        = apa_sparse[k];
756:       apa_sparse[k] = 0.0;
757:       k++;
758:     }

760:     /* diagonal part of C */
761:     ca = cda + cd->i[i];
762:     for (k1=0; k1<cdnz; k1++) {
763:       ca[k1]        = apa_sparse[k];
764:       apa_sparse[k] = 0.0;
765:       k++;
766:     }

768:     /* 2nd off-diagonal part of C */
769:     ca = coa + co->i[i];
770:     for (; k0<conz; k0++) {
771:       ca[k0]        = apa_sparse[k];
772:       apa_sparse[k] = 0.0;
773:       k++;
774:     }
775:   }
776:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
777:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
778:   return(0);
779: }

781: /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
782: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat C)
783: {
784:   PetscErrorCode     ierr;
785:   MPI_Comm           comm;
786:   PetscMPIInt        size;
787:   Mat_APMPI          *ptap;
788:   PetscFreeSpaceList free_space = NULL,current_space=NULL;
789:   Mat_MPIAIJ         *a  = (Mat_MPIAIJ*)A->data;
790:   Mat_SeqAIJ         *ad = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
791:   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
792:   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
793:   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=1;
794:   PetscInt           am=A->rmap->n,pn=P->cmap->n,pm=P->rmap->n,lsize=pn+20;
795:   PetscReal          afill;
796:   MatType            mtype;

799:   MatCheckProduct(C,4);
800:   if (C->product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty");
801:   PetscObjectGetComm((PetscObject)A,&comm);
802:   MPI_Comm_size(comm,&size);

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

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

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

813:   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
814:   pi_loc = p_loc->i; pj_loc = p_loc->j;
815:   if (size > 1) {
816:     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
817:     pi_oth = p_oth->i; pj_oth = p_oth->j;
818:   } else {
819:     p_oth  = NULL;
820:     pi_oth = NULL; pj_oth = NULL;
821:   }

823:   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
824:   /*-------------------------------------------------------------------*/
825:   PetscMalloc1(am+2,&api);
826:   ptap->api = api;
827:   api[0]    = 0;

829:   PetscLLCondensedCreate_Scalable(lsize,&lnk);

831:   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
832:   PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);
833:   current_space = free_space;
834:   MatPreallocateInitialize(comm,am,pn,dnz,onz);
835:   for (i=0; i<am; i++) {
836:     /* diagonal portion of A */
837:     nzi = adi[i+1] - adi[i];
838:     for (j=0; j<nzi; j++) {
839:       row  = *adj++;
840:       pnz  = pi_loc[row+1] - pi_loc[row];
841:       Jptr = pj_loc + pi_loc[row];
842:       /* Expand list if it is not long enough */
843:       if (pnz+apnz_max > lsize) {
844:         lsize = pnz+apnz_max;
845:         PetscLLCondensedExpand_Scalable(lsize, &lnk);
846:       }
847:       /* add non-zero cols of P into the sorted linked list lnk */
848:       PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);
849:       apnz     = *lnk; /* The first element in the list is the number of items in the list */
850:       api[i+1] = api[i] + apnz;
851:       if (apnz > apnz_max) apnz_max = apnz + 1; /* '1' for diagonal entry */
852:     }
853:     /* off-diagonal portion of A */
854:     nzi = aoi[i+1] - aoi[i];
855:     for (j=0; j<nzi; j++) {
856:       row  = *aoj++;
857:       pnz  = pi_oth[row+1] - pi_oth[row];
858:       Jptr = pj_oth + pi_oth[row];
859:       /* Expand list if it is not long enough */
860:       if (pnz+apnz_max > lsize) {
861:         lsize = pnz + apnz_max;
862:         PetscLLCondensedExpand_Scalable(lsize, &lnk);
863:       }
864:       /* add non-zero cols of P into the sorted linked list lnk */
865:       PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);
866:       apnz     = *lnk;  /* The first element in the list is the number of items in the list */
867:       api[i+1] = api[i] + apnz;
868:       if (apnz > apnz_max) apnz_max = apnz + 1; /* '1' for diagonal entry */
869:     }

871:     /* add missing diagonal entry */
872:     if (C->force_diagonals) {
873:       j = i + rstart; /* column index */
874:       PetscLLCondensedAddSorted_Scalable(1,&j,lnk);
875:     }

877:     apnz     = *lnk;
878:     api[i+1] = api[i] + apnz;
879:     if (apnz > apnz_max) apnz_max = apnz;

881:     /* if free space is not available, double the total space in the list */
882:     if (current_space->local_remaining<apnz) {
883:       PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);
884:       nspacedouble++;
885:     }

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

891:     current_space->array           += apnz;
892:     current_space->local_used      += apnz;
893:     current_space->local_remaining -= apnz;
894:   }

896:   /* Allocate space for apj, initialize apj, and */
897:   /* destroy list of free space and other temporary array(s) */
898:   PetscMalloc1(api[am]+1,&ptap->apj);
899:   apj  = ptap->apj;
900:   PetscFreeSpaceContiguous(&free_space,ptap->apj);
901:   PetscLLCondensedDestroy_Scalable(lnk);

903:   /* create and assemble symbolic parallel matrix C */
904:   /*----------------------------------------------------*/
905:   MatSetSizes(C,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);
906:   MatSetBlockSizesFromMats(C,A,P);
907:   MatGetType(A,&mtype);
908:   MatSetType(C,mtype);
909:   MatMPIAIJSetPreallocation(C,0,dnz,0,onz);
910:   MatPreallocateFinalize(dnz,onz);

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

915:   MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(C, apj, api);
916:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
917:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
918:   MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);

920:   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
921:   C->ops->productnumeric = MatProductNumeric_AB;

923:   /* attach the supporting struct to C for reuse */
924:   C->product->data    = ptap;
925:   C->product->destroy = MatDestroy_MPIAIJ_MatMatMult;

927:   /* set MatInfo */
928:   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
929:   if (afill < 1.0) afill = 1.0;
930:   C->info.mallocs           = nspacedouble;
931:   C->info.fill_ratio_given  = fill;
932:   C->info.fill_ratio_needed = afill;

934: #if defined(PETSC_USE_INFO)
935:   if (api[am]) {
936:     PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
937:     PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
938:   } else {
939:     PetscInfo(C,"Empty matrix product\n");
940:   }
941: #endif
942:   return(0);
943: }

945: /* This function is needed for the seqMPI matrix-matrix multiplication.  */
946: /* Three input arrays are merged to one output array. The size of the    */
947: /* output array is also output. Duplicate entries only show up once.     */
948: static void Merge3SortedArrays(PetscInt  size1, PetscInt *in1,
949:                                PetscInt  size2, PetscInt *in2,
950:                                PetscInt  size3, PetscInt *in3,
951:                                PetscInt *size4, PetscInt *out)
952: {
953:   int i = 0, j = 0, k = 0, l = 0;

955:   /* Traverse all three arrays */
956:   while (i<size1 && j<size2 && k<size3) {
957:     if (in1[i] < in2[j] && in1[i] < in3[k]) {
958:       out[l++] = in1[i++];
959:     }
960:     else if (in2[j] < in1[i] && in2[j] < in3[k]) {
961:       out[l++] = in2[j++];
962:     }
963:     else if (in3[k] < in1[i] && in3[k] < in2[j]) {
964:       out[l++] = in3[k++];
965:     }
966:     else if (in1[i] == in2[j] && in1[i] < in3[k]) {
967:       out[l++] = in1[i];
968:       i++, j++;
969:     }
970:     else if (in1[i] == in3[k] && in1[i] < in2[j]) {
971:       out[l++] = in1[i];
972:       i++, k++;
973:     }
974:     else if (in3[k] == in2[j] && in2[j] < in1[i])  {
975:       out[l++] = in2[j];
976:       k++, j++;
977:     }
978:     else if (in1[i] == in2[j] && in1[i] == in3[k]) {
979:       out[l++] = in1[i];
980:       i++, j++, k++;
981:     }
982:   }

984:   /* Traverse two remaining arrays */
985:   while (i<size1 && j<size2) {
986:     if (in1[i] < in2[j]) {
987:       out[l++] = in1[i++];
988:     }
989:     else if (in1[i] > in2[j]) {
990:       out[l++] = in2[j++];
991:     }
992:     else {
993:       out[l++] = in1[i];
994:       i++, j++;
995:     }
996:   }

998:   while (i<size1 && k<size3) {
999:     if (in1[i] < in3[k]) {
1000:       out[l++] = in1[i++];
1001:     }
1002:     else if (in1[i] > in3[k]) {
1003:       out[l++] = in3[k++];
1004:     }
1005:     else {
1006:       out[l++] = in1[i];
1007:       i++, k++;
1008:     }
1009:   }

1011:   while (k<size3 && j<size2)  {
1012:     if (in3[k] < in2[j]) {
1013:       out[l++] = in3[k++];
1014:     }
1015:     else if (in3[k] > in2[j]) {
1016:       out[l++] = in2[j++];
1017:     }
1018:     else {
1019:       out[l++] = in3[k];
1020:       k++, j++;
1021:     }
1022:   }

1024:   /* Traverse one remaining array */
1025:   while (i<size1) out[l++] = in1[i++];
1026:   while (j<size2) out[l++] = in2[j++];
1027:   while (k<size3) out[l++] = in3[k++];

1029:   *size4 = l;
1030: }

1032: /* This matrix-matrix multiplication algorithm divides the multiplication into three multiplications and  */
1033: /* adds up the products. Two of these three multiplications are performed with existing (sequential)      */
1034: /* matrix-matrix multiplications.  */
1035: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat A, Mat P, PetscReal fill, Mat C)
1036: {
1037:   PetscErrorCode     ierr;
1038:   MPI_Comm           comm;
1039:   PetscMPIInt        size;
1040:   Mat_APMPI          *ptap;
1041:   PetscFreeSpaceList free_space_diag=NULL, current_space=NULL;
1042:   Mat_MPIAIJ         *a  =(Mat_MPIAIJ*)A->data;
1043:   Mat_SeqAIJ         *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc;
1044:   Mat_MPIAIJ         *p  =(Mat_MPIAIJ*)P->data;
1045:   Mat_SeqAIJ         *adpd_seq, *p_off, *aopoth_seq;
1046:   PetscInt           adponz, adpdnz;
1047:   PetscInt           *pi_loc,*dnz,*onz;
1048:   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,rstart=A->rmap->rstart;
1049:   PetscInt           *lnk,i, i1=0,pnz,row,*adpoi,*adpoj, *api, *adpoJ, *aopJ, *apJ,*Jptr, aopnz, nspacedouble=0,j,nzi,
1050:                      *apj,apnz, *adpdi, *adpdj, *adpdJ, *poff_i, *poff_j, *j_temp, *aopothi, *aopothj;
1051:   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n, p_colstart, p_colend;
1052:   PetscBT            lnkbt;
1053:   PetscReal          afill;
1054:   PetscMPIInt        rank;
1055:   Mat                adpd, aopoth;
1056:   MatType            mtype;
1057:   const char         *prefix;

1060:   MatCheckProduct(C,4);
1061:   if (C->product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty");
1062:   PetscObjectGetComm((PetscObject)A,&comm);
1063:   MPI_Comm_size(comm,&size);
1064:   MPI_Comm_rank(comm, &rank);
1065:   MatGetOwnershipRangeColumn(P, &p_colstart, &p_colend);

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

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

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

1076:   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1077:   pi_loc = p_loc->i;

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

1083:   adpoi[0]    = 0;
1084:   ptap->api = api;
1085:   api[0] = 0;

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

1091:   /* Symbolic calc of A_loc_diag * P_loc_diag */
1092:   MatGetOptionsPrefix(A,&prefix);
1093:   MatProductCreate(a->A,p->A,NULL,&adpd);
1094:   MatGetOptionsPrefix(A,&prefix);
1095:   MatSetOptionsPrefix(adpd,prefix);
1096:   MatAppendOptionsPrefix(adpd,"inner_diag_");

1098:   MatProductSetType(adpd,MATPRODUCT_AB);
1099:   MatProductSetAlgorithm(adpd,"sorted");
1100:   MatProductSetFill(adpd,fill);
1101:   MatProductSetFromOptions(adpd);

1103:   adpd->force_diagonals = C->force_diagonals;
1104:   MatProductSymbolic(adpd);

1106:   adpd_seq = (Mat_SeqAIJ*)((adpd)->data);
1107:   adpdi = adpd_seq->i; adpdj = adpd_seq->j;
1108:   p_off = (Mat_SeqAIJ*)((p->B)->data);
1109:   poff_i = p_off->i; poff_j = p_off->j;

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

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

1119:   for (i=0; i<am; i++) {
1120:     /* A_diag * P_loc_off */
1121:     nzi = adi[i+1] - adi[i];
1122:     for (j=0; j<nzi; j++) {
1123:       row  = *adj++;
1124:       pnz  = poff_i[row+1] - poff_i[row];
1125:       Jptr = poff_j + poff_i[row];
1126:       for (i1 = 0; i1 < pnz; i1++) {
1127:         j_temp[i1] = p->garray[Jptr[i1]];
1128:       }
1129:       /* add non-zero cols of P into the sorted linked list lnk */
1130:       PetscLLCondensedAddSorted(pnz,j_temp,lnk,lnkbt);
1131:     }

1133:     adponz     = lnk[0];
1134:     adpoi[i+1] = adpoi[i] + adponz;

1136:     /* if free space is not available, double the total space in the list */
1137:     if (current_space->local_remaining<adponz) {
1138:       PetscFreeSpaceGet(PetscIntSumTruncate(adponz,current_space->total_array_size),&current_space);
1139:       nspacedouble++;
1140:     }

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

1145:     current_space->array           += adponz;
1146:     current_space->local_used      += adponz;
1147:     current_space->local_remaining -= adponz;
1148:   }

1150:   /* Symbolic calc of A_off * P_oth */
1151:   MatSetOptionsPrefix(a->B,prefix);
1152:   MatAppendOptionsPrefix(a->B,"inner_offdiag_");
1153:   MatCreate(PETSC_COMM_SELF,&aopoth);
1154:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->B, ptap->P_oth, fill, aopoth);
1155:   aopoth_seq = (Mat_SeqAIJ*)((aopoth)->data);
1156:   aopothi = aopoth_seq->i; aopothj = aopoth_seq->j;

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

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

1164:   /* Copy from linked list to j-array */
1165:   PetscFreeSpaceContiguous(&free_space_diag,adpoj);
1166:   PetscLLDestroy(lnk,lnkbt);

1168:   adpoJ = adpoj;
1169:   adpdJ = adpdj;
1170:   aopJ = aopothj;
1171:   apj  = ptap->apj;
1172:   apJ = apj; /* still empty */

1174:   /* Merge j-arrays of A_off * P, A_diag * P_loc_off, and */
1175:   /* A_diag * P_loc_diag to get A*P */
1176:   for (i = 0; i < am; i++) {
1177:     aopnz  =  aopothi[i+1] -  aopothi[i];
1178:     adponz = adpoi[i+1] - adpoi[i];
1179:     adpdnz = adpdi[i+1] - adpdi[i];

1181:     /* Correct indices from A_diag*P_diag */
1182:     for (i1 = 0; i1 < adpdnz; i1++) {
1183:       adpdJ[i1] += p_colstart;
1184:     }
1185:     /* Merge j-arrays of A_diag * P_loc_off and A_diag * P_loc_diag and A_off * P_oth */
1186:     Merge3SortedArrays(adponz, adpoJ, adpdnz, adpdJ, aopnz, aopJ, &apnz, apJ);
1187:     MatPreallocateSet(i+rstart, apnz, apJ, dnz, onz);

1189:     aopJ += aopnz;
1190:     adpoJ += adponz;
1191:     adpdJ += adpdnz;
1192:     apJ += apnz;
1193:     api[i+1] = api[i] + apnz;
1194:   }

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

1199:   /* create and assemble symbolic parallel matrix C */
1200:   MatSetSizes(C,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);
1201:   MatSetBlockSizesFromMats(C,A,P);
1202:   MatGetType(A,&mtype);
1203:   MatSetType(C,mtype);
1204:   MatMPIAIJSetPreallocation(C,0,dnz,0,onz);
1205:   MatPreallocateFinalize(dnz,onz);

1207:   MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(C, apj, api);
1208:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1209:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1210:   MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);

1212:   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1213:   C->ops->productnumeric = MatProductNumeric_AB;

1215:   /* attach the supporting struct to C for reuse */
1216:   C->product->data    = ptap;
1217:   C->product->destroy = MatDestroy_MPIAIJ_MatMatMult;

1219:   /* set MatInfo */
1220:   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
1221:   if (afill < 1.0) afill = 1.0;
1222:   C->info.mallocs           = nspacedouble;
1223:   C->info.fill_ratio_given  = fill;
1224:   C->info.fill_ratio_needed = afill;

1226: #if defined(PETSC_USE_INFO)
1227:   if (api[am]) {
1228:     PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
1229:     PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
1230:   } else {
1231:     PetscInfo(C,"Empty matrix product\n");
1232:   }
1233: #endif

1235:   MatDestroy(&aopoth);
1236:   MatDestroy(&adpd);
1237:   PetscFree(j_temp);
1238:   PetscFree(adpoj);
1239:   PetscFree(adpoi);
1240:   return(0);
1241: }

1243: /*-------------------------------------------------------------------------*/
1244: /* This routine only works when scall=MAT_REUSE_MATRIX! */
1245: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
1246: {
1248:   Mat_APMPI      *ptap;
1249:   Mat            Pt;

1252:   MatCheckProduct(C,3);
1253:   ptap = (Mat_APMPI*)C->product->data;
1254:   if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data");
1255:   if (!ptap->Pt) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatProductClear()");

1257:   Pt   = ptap->Pt;
1258:   MatTranspose(P,MAT_REUSE_MATRIX,&Pt);
1259:   MatMatMultNumeric_MPIAIJ_MPIAIJ(Pt,A,C);
1260:   return(0);
1261: }

1263: /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1264: PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat C)
1265: {
1266:   PetscErrorCode      ierr;
1267:   Mat_APMPI           *ptap;
1268:   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data;
1269:   MPI_Comm            comm;
1270:   PetscMPIInt         size,rank;
1271:   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1272:   PetscInt            pn=P->cmap->n,aN=A->cmap->N,an=A->cmap->n;
1273:   PetscInt            *lnk,i,k,nsend,rstart;
1274:   PetscBT             lnkbt;
1275:   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,nrecv;
1276:   PETSC_UNUSED PetscMPIInt icompleted=0;
1277:   PetscInt            **buf_rj,**buf_ri,**buf_ri_k,row,ncols,*cols;
1278:   PetscInt            len,proc,*dnz,*onz,*owners,nzi;
1279:   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1280:   MPI_Request         *swaits,*rwaits;
1281:   MPI_Status          *sstatus,rstatus;
1282:   PetscLayout         rowmap;
1283:   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
1284:   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
1285:   PetscInt            *Jptr,*prmap=p->garray,con,j,Crmax;
1286:   Mat_SeqAIJ          *a_loc,*c_loc,*c_oth;
1287:   PetscTable          ta;
1288:   MatType             mtype;
1289:   const char          *prefix;

1292:   PetscObjectGetComm((PetscObject)A,&comm);
1293:   MPI_Comm_size(comm,&size);
1294:   MPI_Comm_rank(comm,&rank);

1296:   /* create symbolic parallel matrix C */
1297:   MatGetType(A,&mtype);
1298:   MatSetType(C,mtype);

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

1302:   /* create struct Mat_APMPI and attached it to C later */
1303:   PetscNew(&ptap);
1304:   ptap->reuse = MAT_INITIAL_MATRIX;

1306:   /* (0) compute Rd = Pd^T, Ro = Po^T  */
1307:   /* --------------------------------- */
1308:   MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);
1309:   MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);

1311:   /* (1) compute symbolic A_loc */
1312:   /* ---------------------------*/
1313:   MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&ptap->A_loc);

1315:   /* (2-1) compute symbolic C_oth = Ro*A_loc  */
1316:   /* ------------------------------------ */
1317:   MatGetOptionsPrefix(A,&prefix);
1318:   MatSetOptionsPrefix(ptap->Ro,prefix);
1319:   MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");
1320:   MatCreate(PETSC_COMM_SELF,&ptap->C_oth);
1321:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->A_loc,fill,ptap->C_oth);

1323:   /* (3) send coj of C_oth to other processors  */
1324:   /* ------------------------------------------ */
1325:   /* determine row ownership */
1326:   PetscLayoutCreate(comm,&rowmap);
1327:   rowmap->n  = pn;
1328:   rowmap->bs = 1;
1329:   PetscLayoutSetUp(rowmap);
1330:   owners = rowmap->range;

1332:   /* determine the number of messages to send, their lengths */
1333:   PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);
1334:   PetscArrayzero(len_s,size);
1335:   PetscArrayzero(len_si,size);

1337:   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
1338:   coi   = c_oth->i; coj = c_oth->j;
1339:   con   = ptap->C_oth->rmap->n;
1340:   proc  = 0;
1341:   for (i=0; i<con; i++) {
1342:     while (prmap[i] >= owners[proc+1]) proc++;
1343:     len_si[proc]++;               /* num of rows in Co(=Pt*A) to be sent to [proc] */
1344:     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
1345:   }

1347:   len          = 0; /* max length of buf_si[], see (4) */
1348:   owners_co[0] = 0;
1349:   nsend        = 0;
1350:   for (proc=0; proc<size; proc++) {
1351:     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1352:     if (len_s[proc]) {
1353:       nsend++;
1354:       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
1355:       len         += len_si[proc];
1356:     }
1357:   }

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

1363:   /* post the Irecv and Isend of coj */
1364:   PetscCommGetNewTag(comm,&tagj);
1365:   PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);
1366:   PetscMalloc1(nsend+1,&swaits);
1367:   for (proc=0, k=0; proc<size; proc++) {
1368:     if (!len_s[proc]) continue;
1369:     i    = owners_co[proc];
1370:     MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);
1371:     k++;
1372:   }

1374:   /* (2-2) compute symbolic C_loc = Rd*A_loc */
1375:   /* ---------------------------------------- */
1376:   MatSetOptionsPrefix(ptap->Rd,prefix);
1377:   MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");
1378:   MatCreate(PETSC_COMM_SELF,&ptap->C_loc);
1379:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->A_loc,fill,ptap->C_loc);
1380:   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;

1382:   /* receives coj are complete */
1383:   for (i=0; i<nrecv; i++) {
1384:     MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);
1385:   }
1386:   PetscFree(rwaits);
1387:   if (nsend) {MPI_Waitall(nsend,swaits,sstatus);}

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

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

1396:   for (k=0; k<nrecv; k++) {/* k-th received message */
1397:     Jptr = buf_rj[k];
1398:     for (j=0; j<len_r[k]; j++) {
1399:       PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);
1400:     }
1401:   }
1402:   PetscTableGetCount(ta,&Crmax);
1403:   PetscTableDestroy(&ta);

1405:   /* (4) send and recv coi */
1406:   /*-----------------------*/
1407:   PetscCommGetNewTag(comm,&tagi);
1408:   PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);
1409:   PetscMalloc1(len+1,&buf_s);
1410:   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1411:   for (proc=0,k=0; proc<size; proc++) {
1412:     if (!len_s[proc]) continue;
1413:     /* form outgoing message for i-structure:
1414:          buf_si[0]:                 nrows to be sent
1415:                [1:nrows]:           row index (global)
1416:                [nrows+1:2*nrows+1]: i-structure index
1417:     */
1418:     /*-------------------------------------------*/
1419:     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
1420:     buf_si_i    = buf_si + nrows+1;
1421:     buf_si[0]   = nrows;
1422:     buf_si_i[0] = 0;
1423:     nrows       = 0;
1424:     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1425:       nzi = coi[i+1] - coi[i];
1426:       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1427:       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1428:       nrows++;
1429:     }
1430:     MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);
1431:     k++;
1432:     buf_si += len_si[proc];
1433:   }
1434:   for (i=0; i<nrecv; i++) {
1435:     MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);
1436:   }
1437:   PetscFree(rwaits);
1438:   if (nsend) {MPI_Waitall(nsend,swaits,sstatus);}

1440:   PetscFree4(len_s,len_si,sstatus,owners_co);
1441:   PetscFree(len_ri);
1442:   PetscFree(swaits);
1443:   PetscFree(buf_s);

1445:   /* (5) compute the local portion of C      */
1446:   /* ------------------------------------------ */
1447:   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of C */
1448:   PetscFreeSpaceGet(Crmax,&free_space);
1449:   current_space = free_space;

1451:   PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);
1452:   for (k=0; k<nrecv; k++) {
1453:     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1454:     nrows       = *buf_ri_k[k];
1455:     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1456:     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recved i-structure  */
1457:   }

1459:   MatPreallocateInitialize(comm,pn,an,dnz,onz);
1460:   PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);
1461:   for (i=0; i<pn; i++) { /* for each local row of C */
1462:     /* add C_loc into C */
1463:     nzi  = c_loc->i[i+1] - c_loc->i[i];
1464:     Jptr = c_loc->j + c_loc->i[i];
1465:     PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);

1467:     /* add received col data into lnk */
1468:     for (k=0; k<nrecv; k++) { /* k-th received message */
1469:       if (i == *nextrow[k]) { /* i-th row */
1470:         nzi  = *(nextci[k]+1) - *nextci[k];
1471:         Jptr = buf_rj[k] + *nextci[k];
1472:         PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);
1473:         nextrow[k]++; nextci[k]++;
1474:       }
1475:     }

1477:     /* add missing diagonal entry */
1478:     if (C->force_diagonals) {
1479:       k = i + owners[rank]; /* column index */
1480:       PetscLLCondensedAddSorted(1,&k,lnk,lnkbt);
1481:     }

1483:     nzi = lnk[0];

1485:     /* copy data into free space, then initialize lnk */
1486:     PetscLLCondensedClean(aN,nzi,current_space->array,lnk,lnkbt);
1487:     MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);
1488:   }
1489:   PetscFree3(buf_ri_k,nextrow,nextci);
1490:   PetscLLDestroy(lnk,lnkbt);
1491:   PetscFreeSpaceDestroy(free_space);

1493:   /* local sizes and preallocation */
1494:   MatSetSizes(C,pn,an,PETSC_DETERMINE,PETSC_DETERMINE);
1495:   if (P->cmap->bs > 0) {PetscLayoutSetBlockSize(C->rmap,P->cmap->bs);}
1496:   if (A->cmap->bs > 0) {PetscLayoutSetBlockSize(C->cmap,A->cmap->bs);}
1497:   MatMPIAIJSetPreallocation(C,0,dnz,0,onz);
1498:   MatPreallocateFinalize(dnz,onz);

1500:   /* add C_loc and C_oth to C */
1501:   MatGetOwnershipRange(C,&rstart,NULL);
1502:   for (i=0; i<pn; i++) {
1503:     ncols = c_loc->i[i+1] - c_loc->i[i];
1504:     cols  = c_loc->j + c_loc->i[i];
1505:     row   = rstart + i;
1506:     MatSetValues(C,1,(const PetscInt*)&row,ncols,(const PetscInt*)cols,NULL,INSERT_VALUES);

1508:     if (C->force_diagonals) {
1509:       MatSetValues(C,1,(const PetscInt*)&row,1,(const PetscInt*)&row,NULL,INSERT_VALUES);
1510:     }
1511:   }
1512:   for (i=0; i<con; i++) {
1513:     ncols = c_oth->i[i+1] - c_oth->i[i];
1514:     cols  = c_oth->j + c_oth->i[i];
1515:     row   = prmap[i];
1516:     MatSetValues(C,1,(const PetscInt*)&row,ncols,(const PetscInt*)cols,NULL,INSERT_VALUES);
1517:   }
1518:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1519:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1520:   MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);

1522:   /* members in merge */
1523:   PetscFree(id_r);
1524:   PetscFree(len_r);
1525:   PetscFree(buf_ri[0]);
1526:   PetscFree(buf_ri);
1527:   PetscFree(buf_rj[0]);
1528:   PetscFree(buf_rj);
1529:   PetscLayoutDestroy(&rowmap);

1531:   /* attach the supporting struct to C for reuse */
1532:   C->product->data    = ptap;
1533:   C->product->destroy = MatDestroy_MPIAIJ_PtAP;
1534:   return(0);
1535: }

1537: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
1538: {
1539:   PetscErrorCode    ierr;
1540:   Mat_MPIAIJ        *p=(Mat_MPIAIJ*)P->data;
1541:   Mat_SeqAIJ        *c_seq;
1542:   Mat_APMPI         *ptap;
1543:   Mat               A_loc,C_loc,C_oth;
1544:   PetscInt          i,rstart,rend,cm,ncols,row;
1545:   const PetscInt    *cols;
1546:   const PetscScalar *vals;

1549:   MatCheckProduct(C,3);
1550:   ptap = (Mat_APMPI*)C->product->data;
1551:   if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data");
1552:   if (!ptap->A_loc) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatProductClear()");
1553:   MatZeroEntries(C);

1555:   if (ptap->reuse == MAT_REUSE_MATRIX) {
1556:     /* These matrices are obtained in MatTransposeMatMultSymbolic() */
1557:     /* 1) get R = Pd^T, Ro = Po^T */
1558:     /*----------------------------*/
1559:     MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);
1560:     MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);

1562:     /* 2) compute numeric A_loc */
1563:     /*--------------------------*/
1564:     MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&ptap->A_loc);
1565:   }

1567:   /* 3) C_loc = Rd*A_loc, C_oth = Ro*A_loc */
1568:   A_loc = ptap->A_loc;
1569:   ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,A_loc,ptap->C_loc);
1570:   ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,A_loc,ptap->C_oth);
1571:   C_loc = ptap->C_loc;
1572:   C_oth = ptap->C_oth;

1574:   /* add C_loc and C_oth to C */
1575:   MatGetOwnershipRange(C,&rstart,&rend);

1577:   /* C_loc -> C */
1578:   cm    = C_loc->rmap->N;
1579:   c_seq = (Mat_SeqAIJ*)C_loc->data;
1580:   cols = c_seq->j;
1581:   vals = c_seq->a;
1582:   for (i=0; i<cm; i++) {
1583:     ncols = c_seq->i[i+1] - c_seq->i[i];
1584:     row = rstart + i;
1585:     MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);
1586:     cols += ncols; vals += ncols;
1587:   }

1589:   /* Co -> C, off-processor part */
1590:   cm    = C_oth->rmap->N;
1591:   c_seq = (Mat_SeqAIJ*)C_oth->data;
1592:   cols  = c_seq->j;
1593:   vals  = c_seq->a;
1594:   for (i=0; i<cm; i++) {
1595:     ncols = c_seq->i[i+1] - c_seq->i[i];
1596:     row = p->garray[i];
1597:     MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);
1598:     cols += ncols; vals += ncols;
1599:   }
1600:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1601:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1602:   MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);

1604:   ptap->reuse = MAT_REUSE_MATRIX;
1605:   return(0);
1606: }

1608: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1609: {
1610:   PetscErrorCode      ierr;
1611:   Mat_Merge_SeqsToMPI *merge;
1612:   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data;
1613:   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1614:   Mat_APMPI           *ptap;
1615:   PetscInt            *adj;
1616:   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1617:   MatScalar           *ada,*ca,valtmp;
1618:   PetscInt            am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1619:   MPI_Comm            comm;
1620:   PetscMPIInt         size,rank,taga,*len_s;
1621:   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1622:   PetscInt            **buf_ri,**buf_rj;
1623:   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1624:   MPI_Request         *s_waits,*r_waits;
1625:   MPI_Status          *status;
1626:   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1627:   const PetscScalar   *dummy;
1628:   PetscInt            *ai,*aj,*coi,*coj,*poJ,*pdJ;
1629:   Mat                 A_loc;
1630:   Mat_SeqAIJ          *a_loc;

1633:   MatCheckProduct(C,3);
1634:   ptap = (Mat_APMPI*)C->product->data;
1635:   if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data");
1636:   if (!ptap->A_loc) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatProductClear()");
1637:   PetscObjectGetComm((PetscObject)C,&comm);
1638:   MPI_Comm_size(comm,&size);
1639:   MPI_Comm_rank(comm,&rank);

1641:   merge = ptap->merge;

1643:   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1644:   /*------------------------------------------*/
1645:   /* get data from symbolic products */
1646:   coi    = merge->coi; coj = merge->coj;
1647:   PetscCalloc1(coi[pon]+1,&coa);
1648:   bi     = merge->bi; bj = merge->bj;
1649:   owners = merge->rowmap->range;
1650:   PetscCalloc1(bi[cm]+1,&ba);

1652:   /* get A_loc by taking all local rows of A */
1653:   A_loc = ptap->A_loc;
1654:   MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);
1655:   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1656:   ai    = a_loc->i;
1657:   aj    = a_loc->j;

1659:   /* trigger copy to CPU */
1660:   MatSeqAIJGetArrayRead(p->A,&dummy);
1661:   MatSeqAIJRestoreArrayRead(p->A,&dummy);
1662:   MatSeqAIJGetArrayRead(p->B,&dummy);
1663:   MatSeqAIJRestoreArrayRead(p->B,&dummy);
1664:   for (i=0; i<am; i++) {
1665:     anz = ai[i+1] - ai[i];
1666:     adj = aj + ai[i];
1667:     ada = a_loc->a + ai[i];

1669:     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1670:     /*-------------------------------------------------------------*/
1671:     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1672:     pnz = po->i[i+1] - po->i[i];
1673:     poJ = po->j + po->i[i];
1674:     pA  = po->a + po->i[i];
1675:     for (j=0; j<pnz; j++) {
1676:       row = poJ[j];
1677:       cj  = coj + coi[row];
1678:       ca  = coa + coi[row];
1679:       /* perform sparse axpy */
1680:       nexta  = 0;
1681:       valtmp = pA[j];
1682:       for (k=0; nexta<anz; k++) {
1683:         if (cj[k] == adj[nexta]) {
1684:           ca[k] += valtmp*ada[nexta];
1685:           nexta++;
1686:         }
1687:       }
1688:       PetscLogFlops(2.0*anz);
1689:     }

1691:     /* put the value into Cd (diagonal part) */
1692:     pnz = pd->i[i+1] - pd->i[i];
1693:     pdJ = pd->j + pd->i[i];
1694:     pA  = pd->a + pd->i[i];
1695:     for (j=0; j<pnz; j++) {
1696:       row = pdJ[j];
1697:       cj  = bj + bi[row];
1698:       ca  = ba + bi[row];
1699:       /* perform sparse axpy */
1700:       nexta  = 0;
1701:       valtmp = pA[j];
1702:       for (k=0; nexta<anz; k++) {
1703:         if (cj[k] == adj[nexta]) {
1704:           ca[k] += valtmp*ada[nexta];
1705:           nexta++;
1706:         }
1707:       }
1708:       PetscLogFlops(2.0*anz);
1709:     }
1710:   }

1712:   /* 3) send and recv matrix values coa */
1713:   /*------------------------------------*/
1714:   buf_ri = merge->buf_ri;
1715:   buf_rj = merge->buf_rj;
1716:   len_s  = merge->len_s;
1717:   PetscCommGetNewTag(comm,&taga);
1718:   PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);

1720:   PetscMalloc2(merge->nsend+1,&s_waits,size,&status);
1721:   for (proc=0,k=0; proc<size; proc++) {
1722:     if (!len_s[proc]) continue;
1723:     i    = merge->owners_co[proc];
1724:     MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);
1725:     k++;
1726:   }
1727:   if (merge->nrecv) {MPI_Waitall(merge->nrecv,r_waits,status);}
1728:   if (merge->nsend) {MPI_Waitall(merge->nsend,s_waits,status);}

1730:   PetscFree2(s_waits,status);
1731:   PetscFree(r_waits);
1732:   PetscFree(coa);

1734:   /* 4) insert local Cseq and received values into Cmpi */
1735:   /*----------------------------------------------------*/
1736:   PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);
1737:   for (k=0; k<merge->nrecv; k++) {
1738:     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1739:     nrows       = *(buf_ri_k[k]);
1740:     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1741:     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recved i-structure  */
1742:   }

1744:   for (i=0; i<cm; i++) {
1745:     row  = owners[rank] + i; /* global row index of C_seq */
1746:     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1747:     ba_i = ba + bi[i];
1748:     bnz  = bi[i+1] - bi[i];
1749:     /* add received vals into ba */
1750:     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1751:       /* i-th row */
1752:       if (i == *nextrow[k]) {
1753:         cnz    = *(nextci[k]+1) - *nextci[k];
1754:         cj     = buf_rj[k] + *(nextci[k]);
1755:         ca     = abuf_r[k] + *(nextci[k]);
1756:         nextcj = 0;
1757:         for (j=0; nextcj<cnz; j++) {
1758:           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1759:             ba_i[j] += ca[nextcj++];
1760:           }
1761:         }
1762:         nextrow[k]++; nextci[k]++;
1763:         PetscLogFlops(2.0*cnz);
1764:       }
1765:     }
1766:     MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);
1767:   }
1768:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1769:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

1771:   PetscFree(ba);
1772:   PetscFree(abuf_r[0]);
1773:   PetscFree(abuf_r);
1774:   PetscFree3(buf_ri_k,nextrow,nextci);
1775:   return(0);
1776: }

1778: PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat C)
1779: {
1780:   PetscErrorCode      ierr;
1781:   Mat                 A_loc;
1782:   Mat_APMPI           *ptap;
1783:   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1784:   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*a=(Mat_MPIAIJ*)A->data;
1785:   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1786:   PetscInt            nnz;
1787:   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1788:   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1789:   MPI_Comm            comm;
1790:   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1791:   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1792:   PetscInt            len,proc,*dnz,*onz,*owners;
1793:   PetscInt            nzi,*bi,*bj;
1794:   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1795:   MPI_Request         *swaits,*rwaits;
1796:   MPI_Status          *sstatus,rstatus;
1797:   Mat_Merge_SeqsToMPI *merge;
1798:   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1799:   PetscReal           afill  =1.0,afill_tmp;
1800:   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax;
1801:   Mat_SeqAIJ          *a_loc;
1802:   PetscTable          ta;
1803:   MatType             mtype;

1806:   PetscObjectGetComm((PetscObject)A,&comm);
1807:   /* check if matrix local sizes are compatible */
1808:   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);

1810:   MPI_Comm_size(comm,&size);
1811:   MPI_Comm_rank(comm,&rank);

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

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

1819:   ptap->A_loc = A_loc;
1820:   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1821:   ai          = a_loc->i;
1822:   aj          = a_loc->j;

1824:   /* determine symbolic Co=(p->B)^T*A - send to others */
1825:   /*----------------------------------------------------*/
1826:   MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);
1827:   MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);
1828:   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1829:                          >= (num of nonzero rows of C_seq) - pn */
1830:   PetscMalloc1(pon+1,&coi);
1831:   coi[0] = 0;

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

1838:   /* create and initialize a linked list */
1839:   PetscTableCreate(A->cmap->n + a->B->cmap->N,aN,&ta);
1840:   MatRowMergeMax_SeqAIJ(a_loc,am,ta);
1841:   PetscTableGetCount(ta,&Armax);

1843:   PetscLLCondensedCreate_Scalable(Armax,&lnk);

1845:   for (i=0; i<pon; i++) {
1846:     pnz = poti[i+1] - poti[i];
1847:     ptJ = potj + poti[i];
1848:     for (j=0; j<pnz; j++) {
1849:       row  = ptJ[j]; /* row of A_loc == col of Pot */
1850:       anz  = ai[row+1] - ai[row];
1851:       Jptr = aj + ai[row];
1852:       /* add non-zero cols of AP into the sorted linked list lnk */
1853:       PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);
1854:     }
1855:     nnz = lnk[0];

1857:     /* If free space is not available, double the total space in the list */
1858:     if (current_space->local_remaining<nnz) {
1859:       PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);
1860:       nspacedouble++;
1861:     }

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

1866:     current_space->array           += nnz;
1867:     current_space->local_used      += nnz;
1868:     current_space->local_remaining -= nnz;

1870:     coi[i+1] = coi[i] + nnz;
1871:   }

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

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

1880:   /* send j-array (coj) of Co to other processors */
1881:   /*----------------------------------------------*/
1882:   /* determine row ownership */
1883:   PetscNew(&merge);
1884:   PetscLayoutCreate(comm,&merge->rowmap);

1886:   merge->rowmap->n  = pn;
1887:   merge->rowmap->bs = 1;

1889:   PetscLayoutSetUp(merge->rowmap);
1890:   owners = merge->rowmap->range;

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

1896:   len_s        = merge->len_s;
1897:   merge->nsend = 0;

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

1901:   proc = 0;
1902:   for (i=0; i<pon; i++) {
1903:     while (prmap[i] >= owners[proc+1]) proc++;
1904:     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1905:     len_s[proc] += coi[i+1] - coi[i];
1906:   }

1908:   len          = 0; /* max length of buf_si[] */
1909:   owners_co[0] = 0;
1910:   for (proc=0; proc<size; proc++) {
1911:     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1912:     if (len_si[proc]) {
1913:       merge->nsend++;
1914:       len_si[proc] = 2*(len_si[proc] + 1);
1915:       len         += len_si[proc];
1916:     }
1917:   }

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

1923:   /* post the Irecv and Isend of coj */
1924:   PetscCommGetNewTag(comm,&tagj);
1925:   PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);
1926:   PetscMalloc1(merge->nsend+1,&swaits);
1927:   for (proc=0, k=0; proc<size; proc++) {
1928:     if (!len_s[proc]) continue;
1929:     i    = owners_co[proc];
1930:     MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);
1931:     k++;
1932:   }

1934:   /* receives and sends of coj are complete */
1935:   PetscMalloc1(size,&sstatus);
1936:   for (i=0; i<merge->nrecv; i++) {
1937:     PETSC_UNUSED PetscMPIInt icompleted;
1938:     MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);
1939:   }
1940:   PetscFree(rwaits);
1941:   if (merge->nsend) {MPI_Waitall(merge->nsend,swaits,sstatus);}

1943:   /* add received column indices into table to update Armax */
1944:   /* Armax can be as large as aN if a P[row,:] is dense, see src/ksp/ksp/tutorials/ex56.c! */
1945:   for (k=0; k<merge->nrecv; k++) {/* k-th received message */
1946:     Jptr = buf_rj[k];
1947:     for (j=0; j<merge->len_r[k]; j++) {
1948:       PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);
1949:     }
1950:   }
1951:   PetscTableGetCount(ta,&Armax);
1952:   /* 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); */

1954:   /* send and recv coi */
1955:   /*-------------------*/
1956:   PetscCommGetNewTag(comm,&tagi);
1957:   PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);
1958:   PetscMalloc1(len+1,&buf_s);
1959:   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1960:   for (proc=0,k=0; proc<size; proc++) {
1961:     if (!len_s[proc]) continue;
1962:     /* form outgoing message for i-structure:
1963:          buf_si[0]:                 nrows to be sent
1964:                [1:nrows]:           row index (global)
1965:                [nrows+1:2*nrows+1]: i-structure index
1966:     */
1967:     /*-------------------------------------------*/
1968:     nrows       = len_si[proc]/2 - 1;
1969:     buf_si_i    = buf_si + nrows+1;
1970:     buf_si[0]   = nrows;
1971:     buf_si_i[0] = 0;
1972:     nrows       = 0;
1973:     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1974:       nzi               = coi[i+1] - coi[i];
1975:       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1976:       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1977:       nrows++;
1978:     }
1979:     MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);
1980:     k++;
1981:     buf_si += len_si[proc];
1982:   }
1983:   i = merge->nrecv;
1984:   while (i--) {
1985:     PETSC_UNUSED PetscMPIInt icompleted;
1986:     MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);
1987:   }
1988:   PetscFree(rwaits);
1989:   if (merge->nsend) {MPI_Waitall(merge->nsend,swaits,sstatus);}
1990:   PetscFree(len_si);
1991:   PetscFree(len_ri);
1992:   PetscFree(swaits);
1993:   PetscFree(sstatus);
1994:   PetscFree(buf_s);

1996:   /* compute the local portion of C (mpi mat) */
1997:   /*------------------------------------------*/
1998:   /* allocate bi array and free space for accumulating nonzero column info */
1999:   PetscMalloc1(pn+1,&bi);
2000:   bi[0] = 0;

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

2007:   PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);
2008:   for (k=0; k<merge->nrecv; k++) {
2009:     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
2010:     nrows       = *buf_ri_k[k];
2011:     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
2012:     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th received i-structure  */
2013:   }

2015:   PetscLLCondensedCreate_Scalable(Armax,&lnk);
2016:   MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);
2017:   rmax = 0;
2018:   for (i=0; i<pn; i++) {
2019:     /* add pdt[i,:]*AP into lnk */
2020:     pnz = pdti[i+1] - pdti[i];
2021:     ptJ = pdtj + pdti[i];
2022:     for (j=0; j<pnz; j++) {
2023:       row  = ptJ[j];  /* row of AP == col of Pt */
2024:       anz  = ai[row+1] - ai[row];
2025:       Jptr = aj + ai[row];
2026:       /* add non-zero cols of AP into the sorted linked list lnk */
2027:       PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);
2028:     }

2030:     /* add received col data into lnk */
2031:     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
2032:       if (i == *nextrow[k]) { /* i-th row */
2033:         nzi  = *(nextci[k]+1) - *nextci[k];
2034:         Jptr = buf_rj[k] + *nextci[k];
2035:         PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);
2036:         nextrow[k]++; nextci[k]++;
2037:       }
2038:     }

2040:     /* add missing diagonal entry */
2041:     if (C->force_diagonals) {
2042:       k = i + owners[rank]; /* column index */
2043:       PetscLLCondensedAddSorted_Scalable(1,&k,lnk);
2044:     }

2046:     nnz = lnk[0];

2048:     /* if free space is not available, make more free space */
2049:     if (current_space->local_remaining<nnz) {
2050:       PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);
2051:       nspacedouble++;
2052:     }
2053:     /* copy data into free space, then initialize lnk */
2054:     PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);
2055:     MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);

2057:     current_space->array           += nnz;
2058:     current_space->local_used      += nnz;
2059:     current_space->local_remaining -= nnz;

2061:     bi[i+1] = bi[i] + nnz;
2062:     if (nnz > rmax) rmax = nnz;
2063:   }
2064:   PetscFree3(buf_ri_k,nextrow,nextci);

2066:   PetscMalloc1(bi[pn]+1,&bj);
2067:   PetscFreeSpaceContiguous(&free_space,bj);
2068:   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
2069:   if (afill_tmp > afill) afill = afill_tmp;
2070:   PetscLLCondensedDestroy_Scalable(lnk);
2071:   PetscTableDestroy(&ta);
2072:   MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);
2073:   MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);

2075:   /* create symbolic parallel matrix C - why cannot be assembled in Numeric part   */
2076:   /*-------------------------------------------------------------------------------*/
2077:   MatSetSizes(C,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);
2078:   MatSetBlockSizes(C,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));
2079:   MatGetType(A,&mtype);
2080:   MatSetType(C,mtype);
2081:   MatMPIAIJSetPreallocation(C,0,dnz,0,onz);
2082:   MatPreallocateFinalize(dnz,onz);
2083:   MatSetBlockSize(C,1);
2084:   MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);
2085:   for (i=0; i<pn; i++) {
2086:     row  = i + rstart;
2087:     nnz  = bi[i+1] - bi[i];
2088:     Jptr = bj + bi[i];
2089:     MatSetValues(C,1,&row,nnz,Jptr,NULL,INSERT_VALUES);
2090:   }
2091:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
2092:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
2093:   MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2094:   merge->bi        = bi;
2095:   merge->bj        = bj;
2096:   merge->coi       = coi;
2097:   merge->coj       = coj;
2098:   merge->buf_ri    = buf_ri;
2099:   merge->buf_rj    = buf_rj;
2100:   merge->owners_co = owners_co;

2102:   /* attach the supporting struct to C for reuse */
2103:   C->product->data    = ptap;
2104:   C->product->destroy = MatDestroy_MPIAIJ_PtAP;
2105:   ptap->merge         = merge;

2107:   C->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;

2109: #if defined(PETSC_USE_INFO)
2110:   if (bi[pn] != 0) {
2111:     PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
2112:     PetscInfo1(C,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);
2113:   } else {
2114:     PetscInfo(C,"Empty matrix product\n");
2115:   }
2116: #endif
2117:   return(0);
2118: }

2120: /* ---------------------------------------------------------------- */
2121: static PetscErrorCode MatProductSymbolic_AtB_MPIAIJ_MPIAIJ(Mat C)
2122: {
2124:   Mat_Product    *product = C->product;
2125:   Mat            A=product->A,B=product->B;
2126:   PetscReal      fill=product->fill;
2127:   PetscBool      flg;

2130:   /* scalable */
2131:   PetscStrcmp(product->alg,"scalable",&flg);
2132:   if (flg) {
2133:     MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);
2134:     goto next;
2135:   }

2137:   /* nonscalable */
2138:   PetscStrcmp(product->alg,"nonscalable",&flg);
2139:   if (flg) {
2140:     MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);
2141:     goto next;
2142:   }

2144:   /* matmatmult */
2145:   PetscStrcmp(product->alg,"at*b",&flg);
2146:   if (flg) {
2147:     Mat       At;
2148:     Mat_APMPI *ptap;

2150:     MatTranspose(A,MAT_INITIAL_MATRIX,&At);
2151:     MatMatMultSymbolic_MPIAIJ_MPIAIJ(At,B,fill,C);
2152:     ptap = (Mat_APMPI*)C->product->data;
2153:     if (ptap) {
2154:       ptap->Pt = At;
2155:       C->product->destroy = MatDestroy_MPIAIJ_PtAP;
2156:     }
2157:     C->ops->transposematmultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
2158:     goto next;
2159:   }

2161:   /* backend general code */
2162:   PetscStrcmp(product->alg,"backend",&flg);
2163:   if (flg) {
2164:     MatProductSymbolic_MPIAIJBACKEND(C);
2165:     return(0);
2166:   }

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

2170: next:
2171:   C->ops->productnumeric = MatProductNumeric_AtB;
2172:   return(0);
2173: }

2175: /* ---------------------------------------------------------------- */
2176: /* Set options for MatMatMultxxx_MPIAIJ_MPIAIJ */
2177: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_AB(Mat C)
2178: {
2180:   Mat_Product    *product = C->product;
2181:   Mat            A=product->A,B=product->B;
2182: #if defined(PETSC_HAVE_HYPRE)
2183:   const char     *algTypes[5] = {"scalable","nonscalable","seqmpi","backend","hypre"};
2184:   PetscInt       nalg = 5;
2185: #else
2186:   const char     *algTypes[4] = {"scalable","nonscalable","seqmpi","backend",};
2187:   PetscInt       nalg = 4;
2188: #endif
2189:   PetscInt       alg = 1; /* set nonscalable algorithm as default */
2190:   PetscBool      flg;
2191:   MPI_Comm       comm;

2194:   /* Check matrix local sizes */
2195:   PetscObjectGetComm((PetscObject)C,&comm);
2196:   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);

2198:   /* Set "nonscalable" as default algorithm */
2199:   PetscStrcmp(C->product->alg,"default",&flg);
2200:   if (flg) {
2201:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);

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

2209:       MatGetInfo(A,MAT_LOCAL,&Ainfo);
2210:       MatGetInfo(B,MAT_LOCAL,&Binfo);
2211:       nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated);

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

2216:       if (alg_scalable) {
2217:         alg  = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
2218:         MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2219:         PetscInfo2(B,"Use scalable algorithm, BN %D, fill*nz_allocated %g\n",B->cmap->N,product->fill*nz_local);
2220:       }
2221:     }
2222:   }

2224:   /* Get runtime option */
2225:   if (product->api_user) {
2226:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatMult","Mat");
2227:     PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2228:     PetscOptionsEnd();
2229:   } else {
2230:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_AB","Mat");
2231:     PetscOptionsEList("-matproduct_ab_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2232:     PetscOptionsEnd();
2233:   }
2234:   if (flg) {
2235:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2236:   }

2238:   C->ops->productsymbolic = MatProductSymbolic_AB_MPIAIJ_MPIAIJ;
2239:   return(0);
2240: }

2242: /* Set options for MatTransposeMatMultXXX_MPIAIJ_MPIAIJ */
2243: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_AtB(Mat C)
2244: {
2246:   Mat_Product    *product = C->product;
2247:   Mat            A=product->A,B=product->B;
2248:   const char     *algTypes[4] = {"scalable","nonscalable","at*b","backend"};
2249:   PetscInt       nalg = 4;
2250:   PetscInt       alg = 1; /* set default algorithm  */
2251:   PetscBool      flg;
2252:   MPI_Comm       comm;

2255:   /* Check matrix local sizes */
2256:   PetscObjectGetComm((PetscObject)C,&comm);
2257:   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);

2259:   /* Set default algorithm */
2260:   PetscStrcmp(C->product->alg,"default",&flg);
2261:   if (flg) {
2262:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2263:   }

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

2271:     MatGetInfo(A,MAT_LOCAL,&Ainfo);
2272:     MatGetInfo(B,MAT_LOCAL,&Binfo);
2273:     nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated);

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

2278:     if (alg_scalable) {
2279:       alg  = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
2280:       MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2281:       PetscInfo2(B,"Use scalable algorithm, BN %D, fill*nz_allocated %g\n",B->cmap->N,product->fill*nz_local);
2282:     }
2283:   }

2285:   /* Get runtime option */
2286:   if (product->api_user) {
2287:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatTransposeMatMult","Mat");
2288:     PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2289:     PetscOptionsEnd();
2290:   } else {
2291:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_AtB","Mat");
2292:     PetscOptionsEList("-matproduct_atb_via","Algorithmic approach","MatTransposeMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2293:     PetscOptionsEnd();
2294:   }
2295:   if (flg) {
2296:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2297:   }

2299:   C->ops->productsymbolic = MatProductSymbolic_AtB_MPIAIJ_MPIAIJ;
2300:   return(0);
2301: }

2303: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_PtAP(Mat C)
2304: {
2306:   Mat_Product    *product = C->product;
2307:   Mat            A=product->A,P=product->B;
2308:   MPI_Comm       comm;
2309:   PetscBool      flg;
2310:   PetscInt       alg=1; /* set default algorithm */
2311: #if !defined(PETSC_HAVE_HYPRE)
2312:   const char     *algTypes[5] = {"scalable","nonscalable","allatonce","allatonce_merged","backend"};
2313:   PetscInt       nalg=5;
2314: #else
2315:   const char     *algTypes[6] = {"scalable","nonscalable","allatonce","allatonce_merged","backend","hypre"};
2316:   PetscInt       nalg=6;
2317: #endif
2318:   PetscInt       pN=P->cmap->N;

2321:   /* Check matrix local sizes */
2322:   PetscObjectGetComm((PetscObject)C,&comm);
2323:   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);
2324:   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);

2326:   /* Set "nonscalable" as default algorithm */
2327:   PetscStrcmp(C->product->alg,"default",&flg);
2328:   if (flg) {
2329:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);

2331:     /* Set "scalable" as default if BN and local nonzeros of A and B are large */
2332:     if (pN > 100000) {
2333:       MatInfo     Ainfo,Pinfo;
2334:       PetscInt    nz_local;
2335:       PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;

2337:       MatGetInfo(A,MAT_LOCAL,&Ainfo);
2338:       MatGetInfo(P,MAT_LOCAL,&Pinfo);
2339:       nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated);

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

2344:       if (alg_scalable) {
2345:         alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
2346:         MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2347:       }
2348:     }
2349:   }

2351:   /* Get runtime option */
2352:   if (product->api_user) {
2353:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatPtAP","Mat");
2354:     PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);
2355:     PetscOptionsEnd();
2356:   } else {
2357:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_PtAP","Mat");
2358:     PetscOptionsEList("-matproduct_ptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);
2359:     PetscOptionsEnd();
2360:   }
2361:   if (flg) {
2362:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2363:   }

2365:   C->ops->productsymbolic = MatProductSymbolic_PtAP_MPIAIJ_MPIAIJ;
2366:   return(0);
2367: }

2369: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_RARt(Mat C)
2370: {
2371:   Mat_Product *product = C->product;
2372:   Mat         A = product->A,R=product->B;

2375:   /* Check matrix local sizes */
2376:   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);

2378:   C->ops->productsymbolic = MatProductSymbolic_RARt_MPIAIJ_MPIAIJ;
2379:   return(0);
2380: }

2382: /*
2383:  Set options for ABC = A*B*C = A*(B*C); ABC's algorithm must be chosen from AB's algorithm
2384: */
2385: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_ABC(Mat C)
2386: {
2388:   Mat_Product    *product = C->product;
2389:   PetscBool      flg = PETSC_FALSE;
2390:   PetscInt       alg = 1; /* default algorithm */
2391:   const char     *algTypes[3] = {"scalable","nonscalable","seqmpi"};
2392:   PetscInt       nalg = 3;

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

2401:   /* Get runtime option */
2402:   if (product->api_user) {
2403:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatMatMult","Mat");
2404:     PetscOptionsEList("-matmatmatmult_via","Algorithmic approach","MatMatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2405:     PetscOptionsEnd();
2406:   } else {
2407:     PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_ABC","Mat");
2408:     PetscOptionsEList("-matproduct_abc_via","Algorithmic approach","MatProduct_ABC",algTypes,nalg,algTypes[alg],&alg,&flg);
2409:     PetscOptionsEnd();
2410:   }
2411:   if (flg) {
2412:     MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2413:   }

2415:   C->ops->matmatmultsymbolic = MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ;
2416:   C->ops->productsymbolic    = MatProductSymbolic_ABC;
2417:   return(0);
2418: }

2420: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ(Mat C)
2421: {
2423:   Mat_Product    *product = C->product;

2426:   switch (product->type) {
2427:   case MATPRODUCT_AB:
2428:     MatProductSetFromOptions_MPIAIJ_AB(C);
2429:     break;
2430:   case MATPRODUCT_AtB:
2431:     MatProductSetFromOptions_MPIAIJ_AtB(C);
2432:     break;
2433:   case MATPRODUCT_PtAP:
2434:     MatProductSetFromOptions_MPIAIJ_PtAP(C);
2435:     break;
2436:   case MATPRODUCT_RARt:
2437:     MatProductSetFromOptions_MPIAIJ_RARt(C);
2438:     break;
2439:   case MATPRODUCT_ABC:
2440:     MatProductSetFromOptions_MPIAIJ_ABC(C);
2441:     break;
2442:   default:
2443:     break;
2444:   }
2445:   return(0);
2446: }