Actual source code: mpimatmatmult.c

petsc-master 2019-12-07
Report Typos and Errors

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

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

 18: PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
 19: {
 21: #if defined(PETSC_HAVE_HYPRE)
 22:   const char     *algTypes[4] = {"scalable","nonscalable","seqmpi","hypre"};
 23:   PetscInt       nalg = 4;
 24: #else
 25:   const char     *algTypes[3] = {"scalable","nonscalable","seqmpi"};
 26:   PetscInt       nalg = 3;
 27: #endif
 28:   PetscInt       alg = 1; /* set nonscalable algorithm as default */
 29:   MPI_Comm       comm;
 30:   PetscBool      flg;

 33:   if (scall == MAT_INITIAL_MATRIX) {
 34:     PetscObjectGetComm((PetscObject)A,&comm);
 35:     if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend) SETERRQ4(comm,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);

 37:     PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatMatMult","Mat");
 38:     PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[1],&alg,&flg);
 39:     PetscOptionsEnd();

 41:     if (!flg && B->cmap->N > 100000) { /* may switch to scalable algorithm as default */
 42:       MatInfo     Ainfo,Binfo;
 43:       PetscInt    nz_local;
 44:       PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;

 46:       MatGetInfo(A,MAT_LOCAL,&Ainfo);
 47:       MatGetInfo(B,MAT_LOCAL,&Binfo);
 48:       nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated);

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

 53:       if (alg_scalable) {
 54:         alg  = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
 55:         PetscInfo2(B,"Use scalable algorithm, BN %D, fill*nz_allocated %g\n",B->cmap->N,fill*nz_local);
 56:       }
 57:     }

 59:     PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
 60:     switch (alg) {
 61:     case 1:
 62:       MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);
 63:       break;
 64:     case 2:
 65:       MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(A,B,fill,C);
 66:       break;
 67: #if defined(PETSC_HAVE_HYPRE)
 68:     case 3:
 69:       MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);
 70:       break;
 71: #endif
 72:     default:
 73:       MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);
 74:       break;
 75:     }
 76:     PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);

 78:     if (alg == 0 || alg == 1) {
 79:       Mat_MPIAIJ *c  = (Mat_MPIAIJ*)(*C)->data;
 80:       Mat_APMPI  *ap = c->ap;
 81:       PetscOptionsBegin(PetscObjectComm((PetscObject)(*C)),((PetscObject)(*C))->prefix,"MatFreeIntermediateDataStructures","Mat");
 82:       ap->freestruct = PETSC_FALSE;
 83:       PetscOptionsBool("-mat_freeintermediatedatastructures","Free intermediate data structures", "MatFreeIntermediateDataStructures",ap->freestruct,&ap->freestruct, NULL);
 84:       PetscOptionsEnd();
 85:     }
 86:   }

 88:   PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
 89:   (*(*C)->ops->matmultnumeric)(A,B,*C);
 90:   PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
 91:   return(0);
 92: }

 94: PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
 95: {
 97:   Mat_MPIAIJ     *a    = (Mat_MPIAIJ*)A->data;
 98:   Mat_APMPI      *ptap = a->ap;

101:   PetscFree2(ptap->startsj_s,ptap->startsj_r);
102:   PetscFree(ptap->bufa);
103:   MatDestroy(&ptap->P_loc);
104:   MatDestroy(&ptap->P_oth);
105:   MatDestroy(&ptap->Pt);
106:   PetscFree(ptap->api);
107:   PetscFree(ptap->apj);
108:   PetscFree(ptap->apa);
109:   ptap->destroy(A);
110:   PetscFree(ptap);
111:   return(0);
112: }

114: PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C)
115: {
117:   Mat_MPIAIJ     *a  =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
118:   Mat_SeqAIJ     *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
119:   Mat_SeqAIJ     *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
120:   PetscScalar    *cda=cd->a,*coa=co->a;
121:   Mat_SeqAIJ     *p_loc,*p_oth;
122:   PetscScalar    *apa,*ca;
123:   PetscInt       cm   =C->rmap->n;
124:   Mat_APMPI      *ptap=c->ap;
125:   PetscInt       *api,*apj,*apJ,i,k;
126:   PetscInt       cstart=C->cmap->rstart;
127:   PetscInt       cdnz,conz,k0,k1;
128:   MPI_Comm       comm;
129:   PetscMPIInt    size;

132:   PetscObjectGetComm((PetscObject)A,&comm);
133:   MPI_Comm_size(comm,&size);

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

137:   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
138:   /*-----------------------------------------------------*/
139:   /* update numerical values of P_oth and P_loc */
140:   MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);
141:   MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);

143:   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
144:   /*----------------------------------------------------------*/
145:   /* get data from symbolic products */
146:   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
147:   p_oth = NULL;
148:   if (size >1) {
149:     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
150:   }

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

155:   api = ptap->api;
156:   apj = ptap->apj;
157:   for (i=0; i<cm; i++) {
158:     /* compute apa = A[i,:]*P */
159:     AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);

161:     /* set values in C */
162:     apJ  = apj + api[i];
163:     cdnz = cd->i[i+1] - cd->i[i];
164:     conz = co->i[i+1] - co->i[i];

166:     /* 1st off-diagoanl part of C */
167:     ca = coa + co->i[i];
168:     k  = 0;
169:     for (k0=0; k0<conz; k0++) {
170:       if (apJ[k] >= cstart) break;
171:       ca[k0]      = apa[apJ[k]];
172:       apa[apJ[k++]] = 0.0;
173:     }

175:     /* diagonal part of C */
176:     ca = cda + cd->i[i];
177:     for (k1=0; k1<cdnz; k1++) {
178:       ca[k1]      = apa[apJ[k]];
179:       apa[apJ[k++]] = 0.0;
180:     }

182:     /* 2nd off-diagoanl part of C */
183:     ca = coa + co->i[i];
184:     for (; k0<conz; k0++) {
185:       ca[k0]      = apa[apJ[k]];
186:       apa[apJ[k++]] = 0.0;
187:     }
188:   }
189:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
190:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

192:   if (ptap->freestruct) {
193:     MatFreeIntermediateDataStructures(C);
194:   }
195:   return(0);
196: }

198: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C)
199: {
200:   PetscErrorCode     ierr;
201:   MPI_Comm           comm;
202:   PetscMPIInt        size;
203:   Mat                Cmpi;
204:   Mat_APMPI          *ptap;
205:   PetscFreeSpaceList free_space=NULL,current_space=NULL;
206:   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data,*c;
207:   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
208:   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
209:   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
210:   PetscInt           *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
211:   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
212:   PetscBT            lnkbt;
213:   PetscReal          afill;
214:   MatType            mtype;

217:   PetscObjectGetComm((PetscObject)A,&comm);
218:   MPI_Comm_size(comm,&size);

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

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

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

229:   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
230:   pi_loc = p_loc->i; pj_loc = p_loc->j;
231:   if (size > 1) {
232:     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
233:     pi_oth = p_oth->i; pj_oth = p_oth->j;
234:   } else {
235:     p_oth = NULL;
236:     pi_oth = NULL; pj_oth = NULL;
237:   }

239:   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
240:   /*-------------------------------------------------------------------*/
241:   PetscMalloc1(am+2,&api);
242:   ptap->api = api;
243:   api[0]    = 0;

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

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

252:   MatPreallocateInitialize(comm,am,pn,dnz,onz);
253:   for (i=0; i<am; i++) {
254:     /* diagonal portion of A */
255:     nzi = adi[i+1] - adi[i];
256:     for (j=0; j<nzi; j++) {
257:       row  = *adj++;
258:       pnz  = pi_loc[row+1] - pi_loc[row];
259:       Jptr = pj_loc + pi_loc[row];
260:       /* add non-zero cols of P into the sorted linked list lnk */
261:       PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);
262:     }
263:     /* off-diagonal portion of A */
264:     nzi = aoi[i+1] - aoi[i];
265:     for (j=0; j<nzi; j++) {
266:       row  = *aoj++;
267:       pnz  = pi_oth[row+1] - pi_oth[row];
268:       Jptr = pj_oth + pi_oth[row];
269:       PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);
270:     }

272:     apnz     = lnk[0];
273:     api[i+1] = api[i] + apnz;

275:     /* if free space is not available, double the total space in the list */
276:     if (current_space->local_remaining<apnz) {
277:       PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);
278:       nspacedouble++;
279:     }

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

285:     current_space->array           += apnz;
286:     current_space->local_used      += apnz;
287:     current_space->local_remaining -= apnz;
288:   }

290:   /* Allocate space for apj, initialize apj, and */
291:   /* destroy list of free space and other temporary array(s) */
292:   PetscMalloc1(api[am]+1,&ptap->apj);
293:   apj  = ptap->apj;
294:   PetscFreeSpaceContiguous(&free_space,ptap->apj);
295:   PetscLLDestroy(lnk,lnkbt);

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

300:   /* create and assemble symbolic parallel matrix Cmpi */
301:   /*----------------------------------------------------*/
302:   MatCreate(comm,&Cmpi);
303:   MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);
304:   MatSetBlockSizesFromMats(Cmpi,A,P);

306:   MatGetType(A,&mtype);
307:   MatSetType(Cmpi,mtype);
308:   MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);

310:   MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);
311:   MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);
312:   MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);
313:   MatPreallocateFinalize(dnz,onz);

315:   ptap->destroy        = Cmpi->ops->destroy;
316:   ptap->duplicate      = Cmpi->ops->duplicate;
317:   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
318:   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
319:   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;

321:   /* attach the supporting struct to Cmpi for reuse */
322:   c       = (Mat_MPIAIJ*)Cmpi->data;
323:   c->ap = ptap;

325:   *C = Cmpi;

327:   /* set MatInfo */
328:   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
329:   if (afill < 1.0) afill = 1.0;
330:   Cmpi->info.mallocs           = nspacedouble;
331:   Cmpi->info.fill_ratio_given  = fill;
332:   Cmpi->info.fill_ratio_needed = afill;

334: #if defined(PETSC_USE_INFO)
335:   if (api[am]) {
336:     PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
337:     PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
338:   } else {
339:     PetscInfo(Cmpi,"Empty matrix product\n");
340:   }
341: #endif
342:   return(0);
343: }

345: PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
346: {

350:   if (scall == MAT_INITIAL_MATRIX) {
351:     *C = NULL;
352:     PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
353:     MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);
354:     PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
355:   }
356:   PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
357:   MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);
358:   PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
359:   return(0);
360: }

362: typedef struct {
363:   Mat          workB,Bb,Cb,workB1,Bb1,Cb1;
364:   MPI_Request  *rwaits,*swaits;
365:   PetscInt     numBb;  /* num of Bb matrices */
366:   PetscInt     nsends,nrecvs;
367:   MPI_Datatype *stype,*rtype;
368: } MPIAIJ_MPIDense;

370: PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
371: {
372:   MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*)ctx;
373:   PetscErrorCode  ierr;
374:   PetscInt        i;

377:   MatDestroy(&contents->workB);

379:   if (contents->numBb) {
380:     MatDestroy(&contents->Bb);
381:     MatDestroy(&contents->Cb);

383:     MatDestroy(&contents->workB1);
384:     MatDestroy(&contents->Bb1);
385:     MatDestroy(&contents->Cb1);
386:   }
387:   for (i=0; i<contents->nsends; i++) {
388:     MPI_Type_free(&contents->stype[i]);
389:   }
390:   for (i=0; i<contents->nrecvs; i++) {
391:     MPI_Type_free(&contents->rtype[i]);
392:   }
393:   PetscFree4(contents->stype,contents->rtype,contents->rwaits,contents->swaits);
394:   PetscFree(contents);
395:   return(0);
396: }

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

402:   It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C
403: */
404: PETSC_INTERN PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C)
405: {
406:   PetscBool      flg;

410:   PetscObjectTypeCompare((PetscObject)A,MATNEST,&flg);
411:   if (flg) {
412:     PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
413:     MatMatMultSymbolic_Nest_Dense(A,B,PETSC_DEFAULT,&C);
414:     PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
415:     C->ops->matmultnumeric = MatMatMultNumeric_Nest_Dense;
416:   } else {
417:     MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,PETSC_DEFAULT,&C);
418:   }
419:   (*C->ops->matmultnumeric)(A,B,C);
420:   return(0);
421: }

423: /*
424:   Create Bb, Cb, Bb1 and Cb1 matrices to be used by MatMatMult_MPIAIJ_MPIDense().
425:   These matrices are used as wrappers for sub-columns of B and C, thus their own matrix operations are not used.
426:   Modified from MatCreateDense().
427: */
428: PETSC_STATIC_INLINE PetscErrorCode MatCreateSubMPIDense_private(MPI_Comm comm,PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt rbs,PetscInt cbs,PetscScalar *data,Mat *A)
429: {

433:   MatCreate(comm,A);
434:   MatSetSizes(*A,m,n,M,N);
435:   MatSetBlockSizes(*A,rbs,cbs);
436:   MatSetType(*A,MATMPIDENSE);
437:   MatMPIDenseSetPreallocation(*A,data);
438:   (*A)->assembled = PETSC_TRUE;
439:   return(0);
440: }

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

459:   PetscObjectGetComm((PetscObject)A,&comm);
460:   if (!(*C)) {
461:     MatCreate(comm,C);
462:     MatSetSizes(*C,Am,Bn,A->rmap->N,BN);
463:     MatSetBlockSizesFromMats(*C,A,B);
464:     MatSetType(*C,MATMPIDENSE);
465:     MatMPIDenseSetPreallocation(*C,NULL);
466:     MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);
467:     MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);
468:   } else {
469:     /* Check matrix size */
470:     if ((*C)->rmap->n != A->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local row dimensions are incompatible, %D != %D",(*C)->rmap->n,A->rmap->n);
471:     if ((*C)->cmap->n != B->cmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local column dimensions are incompatible, %D != %D",(*C)->cmap->n,B->cmap->n);
472:   }

474:   (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

558:   PetscContainerCreate(comm,&container);
559:   PetscContainerSetPointer(container,contents);
560:   PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);
561:   PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);
562:   PetscContainerDestroy(&container);
563:   return(0);
564: }

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

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

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

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

612:   MatDenseGetArrayRead(B,&b);
613:   MatDenseGetArray(workB,&rvalues);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

745:   PetscObjectGetComm((PetscObject)C,&comm);
746:   MPI_Comm_size(comm,&size);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

902:   PetscLLCondensedCreate_Scalable(lsize,&lnk);

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

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

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

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

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

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

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

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

986:   ptap->destroy             = Cmpi->ops->destroy;
987:   ptap->duplicate           = Cmpi->ops->duplicate;
988:   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
989:   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
990:   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;

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

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

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

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

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

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

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

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

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

1099:   *size4 = l;
1100: }

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

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

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

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

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


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

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

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

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

1162:   /* Symbolic calc of A_loc_diag * P_loc_diag */
1163:   MatGetOptionsPrefix(A,&prefix);
1164:   MatSetOptionsPrefix(a->A,prefix);
1165:   MatAppendOptionsPrefix(a->A,"inner_diag_");
1166:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->A, p->A, fill, &adpd);
1167:   adpd_seq = (Mat_SeqAIJ*)((adpd)->data);
1168:   adpdi = adpd_seq->i; adpdj = adpd_seq->j;
1169:   p_off = (Mat_SeqAIJ*)((p->B)->data);
1170:   poff_i = p_off->i; poff_j = p_off->j;

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


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

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

1195:     adponz     = lnk[0];
1196:     adpoi[i+1] = adpoi[i] + adponz;

1198:     /* if free space is not available, double the total space in the list */
1199:     if (current_space->local_remaining<adponz) {
1200:       PetscFreeSpaceGet(PetscIntSumTruncate(adponz,current_space->total_array_size),&current_space);
1201:       nspacedouble++;
1202:     }

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

1207:     current_space->array           += adponz;
1208:     current_space->local_used      += adponz;
1209:     current_space->local_remaining -= adponz;
1210:   }

1212:   /* Symbolic calc of A_off * P_oth */
1213:   MatSetOptionsPrefix(a->B,prefix);
1214:   MatAppendOptionsPrefix(a->B,"inner_offdiag_");
1215:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->B, ptap->P_oth, fill, &aopoth);
1216:   aopoth_seq = (Mat_SeqAIJ*)((aopoth)->data);
1217:   aopothi = aopoth_seq->i; aopothj = aopoth_seq->j;

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

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

1225:   /* Copy from linked list to j-array */
1226:   PetscFreeSpaceContiguous(&free_space_diag,adpoj);
1227:   PetscLLDestroy(lnk,lnkbt);

1229:   adpoJ = adpoj;
1230:   adpdJ = adpdj;
1231:   aopJ = aopothj;
1232:   apj  = ptap->apj;
1233:   apJ = apj; /* still empty */

1235:   /* Merge j-arrays of A_off * P, A_diag * P_loc_off, and */
1236:   /* A_diag * P_loc_diag to get A*P */
1237:   for (i = 0; i < am; i++) {
1238:     aopnz  =  aopothi[i+1] -  aopothi[i];
1239:     adponz = adpoi[i+1] - adpoi[i];
1240:     adpdnz = adpdi[i+1] - adpdi[i];

1242:     /* Correct indices from A_diag*P_diag */
1243:     for(i1 = 0; i1 < adpdnz; i1++) {
1244:       adpdJ[i1] += p_colstart;
1245:     }
1246:     /* Merge j-arrays of A_diag * P_loc_off and A_diag * P_loc_diag and A_off * P_oth */
1247:     Merge3SortedArrays(adponz, adpoJ, adpdnz, adpdJ, aopnz, aopJ, &apnz, apJ);
1248:     MatPreallocateSet(i+rstart, apnz, apJ, dnz, onz);

1250:     aopJ += aopnz;
1251:     adpoJ += adponz;
1252:     adpdJ += adpdnz;
1253:     apJ += apnz;
1254:     api[i+1] = api[i] + apnz;
1255:   }

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

1260:   /* create and assemble symbolic parallel matrix Cmpi */
1261:   MatCreate(comm,&Cmpi);
1262:   MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);
1263:   MatSetBlockSizesFromMats(Cmpi,A,P);
1264:   MatGetType(A,&mtype);
1265:   MatSetType(Cmpi,mtype);
1266:   MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);


1269:   MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);
1270:   MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);
1271:   MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);
1272:   MatPreallocateFinalize(dnz,onz);


1275:   ptap->destroy        = Cmpi->ops->destroy;
1276:   ptap->duplicate      = Cmpi->ops->duplicate;
1277:   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1278:   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;

1280:   /* attach the supporting struct to Cmpi for reuse */
1281:   c       = (Mat_MPIAIJ*)Cmpi->data;
1282:   c->ap = ptap;
1283:   *C = Cmpi;

1285:   /* set MatInfo */
1286:   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
1287:   if (afill < 1.0) afill = 1.0;
1288:   Cmpi->info.mallocs           = nspacedouble;
1289:   Cmpi->info.fill_ratio_given  = fill;
1290:   Cmpi->info.fill_ratio_needed = afill;

1292: #if defined(PETSC_USE_INFO)
1293:   if (api[am]) {
1294:     PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
1295:     PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
1296:   } else {
1297:     PetscInfo(Cmpi,"Empty matrix product\n");
1298:   }
1299: #endif

1301:   MatDestroy(&aopoth);
1302:   MatDestroy(&adpd);
1303:   PetscFree(j_temp);
1304:   PetscFree(adpoj);
1305:   PetscFree(adpoi);
1306:   return(0);
1307: }


1310: /*-------------------------------------------------------------------------*/
1311: PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
1312: {
1314:   const char     *algTypes[3] = {"scalable","nonscalable","matmatmult"};
1315:   PetscInt       aN=A->cmap->N,alg=1; /* set default algorithm */
1316:   PetscBool      flg;

1319:   if (scall == MAT_INITIAL_MATRIX) {
1320:     PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatTransposeMatMult","Mat");
1321:     PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[1],&alg,&flg);
1322:     PetscOptionsEnd();

1324:     PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);
1325:     switch (alg) {
1326:     case 1:
1327:       if (!flg && aN > 100000) { /* may switch to scalable algorithm as default */
1328:         MatInfo     Ainfo,Pinfo;
1329:         PetscInt    nz_local;
1330:         PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;
1331:         MPI_Comm    comm;

1333:         MatGetInfo(A,MAT_LOCAL,&Ainfo);
1334:         MatGetInfo(P,MAT_LOCAL,&Pinfo);
1335:         nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated); /* estimated local nonzero entries */

1337:         if (aN > fill*nz_local) alg_scalable_loc = PETSC_TRUE;
1338:         PetscObjectGetComm((PetscObject)A,&comm);
1339:         MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);

1341:         if (alg_scalable) {
1342:           alg  = 0; /* scalable algorithm would slower than nonscalable algorithm */
1343:           MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);
1344:           break;
1345:         }
1346:       }
1347:       MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);
1348:       break;
1349:     case 2:
1350:     {
1351:       Mat         Pt;
1352:       Mat_APMPI   *ptap;
1353:       Mat_MPIAIJ  *c;
1354:       MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);
1355:       MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);
1356:       c        = (Mat_MPIAIJ*)(*C)->data;
1357:       ptap     = c->ap;
1358:       if (ptap) {
1359:        ptap->Pt = Pt;
1360:        (*C)->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;
1361:       }
1362:       (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
1363:       return(0);
1364:     }
1365:       break;
1366:     default: /* scalable algorithm */
1367:       MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);
1368:       break;
1369:     }
1370:     PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);

1372:     {
1373:       Mat_MPIAIJ *c  = (Mat_MPIAIJ*)(*C)->data;
1374:       Mat_APMPI  *ap = c->ap;
1375:       PetscOptionsBegin(PetscObjectComm((PetscObject)(*C)),((PetscObject)(*C))->prefix,"MatFreeIntermediateDataStructures","Mat");
1376:       ap->freestruct = PETSC_FALSE;
1377:       PetscOptionsBool("-mat_freeintermediatedatastructures","Free intermediate data structures", "MatFreeIntermediateDataStructures",ap->freestruct,&ap->freestruct, NULL);
1378:       PetscOptionsEnd();
1379:     }
1380:   }

1382:   PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);
1383:   (*(*C)->ops->mattransposemultnumeric)(P,A,*C);
1384:   PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);
1385:   return(0);
1386: }

1388: /* This routine only works when scall=MAT_REUSE_MATRIX! */
1389: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
1390: {
1392:   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
1393:   Mat_APMPI      *ptap= c->ap;
1394:   Mat            Pt;

1397:   if (!ptap->Pt) {
1398:     MPI_Comm comm;
1399:     PetscObjectGetComm((PetscObject)C,&comm);
1400:     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
1401:   }

1403:   Pt=ptap->Pt;
1404:   MatTranspose(P,MAT_REUSE_MATRIX,&Pt);
1405:   MatMatMultNumeric(Pt,A,C);

1407:   /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */
1408:   if (ptap->freestruct) {
1409:     MatFreeIntermediateDataStructures(C);
1410:   }
1411:   return(0);
1412: }

1414: /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1415: PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C)
1416: {
1417:   PetscErrorCode      ierr;
1418:   Mat_APMPI           *ptap;
1419:   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*c;
1420:   MPI_Comm            comm;
1421:   PetscMPIInt         size,rank;
1422:   Mat                 Cmpi;
1423:   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1424:   PetscInt            pn=P->cmap->n,aN=A->cmap->N,an=A->cmap->n;
1425:   PetscInt            *lnk,i,k,nsend;
1426:   PetscBT             lnkbt;
1427:   PetscMPIInt         tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
1428:   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1429:   PetscInt            len,proc,*dnz,*onz,*owners,nzi;
1430:   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1431:   MPI_Request         *swaits,*rwaits;
1432:   MPI_Status          *sstatus,rstatus;
1433:   PetscLayout         rowmap;
1434:   PetscInt            *owners_co,*coi,*coj;    /* i and j array of (p->B)^T*A*P - used in the communication */
1435:   PetscMPIInt         *len_r,*id_r;    /* array of length of comm->size, store send/recv matrix values */
1436:   PetscInt            *Jptr,*prmap=p->garray,con,j,Crmax;
1437:   Mat_SeqAIJ          *a_loc,*c_loc,*c_oth;
1438:   PetscTable          ta;
1439:   MatType             mtype;
1440:   const char          *prefix;

1443:   PetscObjectGetComm((PetscObject)A,&comm);
1444:   MPI_Comm_size(comm,&size);
1445:   MPI_Comm_rank(comm,&rank);

1447:   /* create symbolic parallel matrix Cmpi */
1448:   MatCreate(comm,&Cmpi);
1449:   MatGetType(A,&mtype);
1450:   MatSetType(Cmpi,mtype);

1452:   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;

1454:   /* create struct Mat_APMPI and attached it to C later */
1455:   PetscNew(&ptap);
1456:   ptap->reuse = MAT_INITIAL_MATRIX;

1458:   /* (0) compute Rd = Pd^T, Ro = Po^T  */
1459:   /* --------------------------------- */
1460:   MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);
1461:   MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);

1463:   /* (1) compute symbolic A_loc */
1464:   /* ---------------------------*/
1465:   MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&ptap->A_loc);

1467:   /* (2-1) compute symbolic C_oth = Ro*A_loc  */
1468:   /* ------------------------------------ */
1469:   MatGetOptionsPrefix(A,&prefix);
1470:   MatSetOptionsPrefix(ptap->Ro,prefix);
1471:   MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");
1472:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->A_loc,fill,&ptap->C_oth);

1474:   /* (3) send coj of C_oth to other processors  */
1475:   /* ------------------------------------------ */
1476:   /* determine row ownership */
1477:   PetscLayoutCreate(comm,&rowmap);
1478:   rowmap->n  = pn;
1479:   rowmap->bs = 1;
1480:   PetscLayoutSetUp(rowmap);
1481:   owners = rowmap->range;

1483:   /* determine the number of messages to send, their lengths */
1484:   PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);
1485:   PetscArrayzero(len_s,size);
1486:   PetscArrayzero(len_si,size);

1488:   c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
1489:   coi   = c_oth->i; coj = c_oth->j;
1490:   con   = ptap->C_oth->rmap->n;
1491:   proc  = 0;
1492:   for (i=0; i<con; i++) {
1493:     while (prmap[i] >= owners[proc+1]) proc++;
1494:     len_si[proc]++;               /* num of rows in Co(=Pt*A) to be sent to [proc] */
1495:     len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
1496:   }

1498:   len          = 0; /* max length of buf_si[], see (4) */
1499:   owners_co[0] = 0;
1500:   nsend        = 0;
1501:   for (proc=0; proc<size; proc++) {
1502:     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1503:     if (len_s[proc]) {
1504:       nsend++;
1505:       len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
1506:       len         += len_si[proc];
1507:     }
1508:   }

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

1514:   /* post the Irecv and Isend of coj */
1515:   PetscCommGetNewTag(comm,&tagj);
1516:   PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);
1517:   PetscMalloc1(nsend+1,&swaits);
1518:   for (proc=0, k=0; proc<size; proc++) {
1519:     if (!len_s[proc]) continue;
1520:     i    = owners_co[proc];
1521:     MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);
1522:     k++;
1523:   }

1525:   /* (2-2) compute symbolic C_loc = Rd*A_loc */
1526:   /* ---------------------------------------- */
1527:   MatSetOptionsPrefix(ptap->Rd,prefix);
1528:   MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");
1529:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->A_loc,fill,&ptap->C_loc);
1530:   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;

1532:   /* receives coj are complete */
1533:   for (i=0; i<nrecv; i++) {
1534:     MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);
1535:   }
1536:   PetscFree(rwaits);
1537:   if (nsend) {MPI_Waitall(nsend,swaits,sstatus);}

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

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

1546:   for (k=0; k<nrecv; k++) {/* k-th received message */
1547:     Jptr = buf_rj[k];
1548:     for (j=0; j<len_r[k]; j++) {
1549:       PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);
1550:     }
1551:   }
1552:   PetscTableGetCount(ta,&Crmax);
1553:   PetscTableDestroy(&ta);

1555:   /* (4) send and recv coi */
1556:   /*-----------------------*/
1557:   PetscCommGetNewTag(comm,&tagi);
1558:   PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);
1559:   PetscMalloc1(len+1,&buf_s);
1560:   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1561:   for (proc=0,k=0; proc<size; proc++) {
1562:     if (!len_s[proc]) continue;
1563:     /* form outgoing message for i-structure:
1564:          buf_si[0]:                 nrows to be sent
1565:                [1:nrows]:           row index (global)
1566:                [nrows+1:2*nrows+1]: i-structure index
1567:     */
1568:     /*-------------------------------------------*/
1569:     nrows       = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
1570:     buf_si_i    = buf_si + nrows+1;
1571:     buf_si[0]   = nrows;
1572:     buf_si_i[0] = 0;
1573:     nrows       = 0;
1574:     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1575:       nzi = coi[i+1] - coi[i];
1576:       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1577:       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1578:       nrows++;
1579:     }
1580:     MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);
1581:     k++;
1582:     buf_si += len_si[proc];
1583:   }
1584:   for (i=0; i<nrecv; i++) {
1585:     MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);
1586:   }
1587:   PetscFree(rwaits);
1588:   if (nsend) {MPI_Waitall(nsend,swaits,sstatus);}

1590:   PetscFree4(len_s,len_si,sstatus,owners_co);
1591:   PetscFree(len_ri);
1592:   PetscFree(swaits);
1593:   PetscFree(buf_s);

1595:   /* (5) compute the local portion of Cmpi      */
1596:   /* ------------------------------------------ */
1597:   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
1598:   PetscFreeSpaceGet(Crmax,&free_space);
1599:   current_space = free_space;

1601:   PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);
1602:   for (k=0; k<nrecv; k++) {
1603:     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1604:     nrows       = *buf_ri_k[k];
1605:     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1606:     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1607:   }

1609:   MatPreallocateInitialize(comm,pn,an,dnz,onz);
1610:   PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);
1611:   for (i=0; i<pn; i++) {
1612:     /* add C_loc into Cmpi */
1613:     nzi  = c_loc->i[i+1] - c_loc->i[i];
1614:     Jptr = c_loc->j + c_loc->i[i];
1615:     PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);

1617:     /* add received col data into lnk */
1618:     for (k=0; k<nrecv; k++) { /* k-th received message */
1619:       if (i == *nextrow[k]) { /* i-th row */
1620:         nzi  = *(nextci[k]+1) - *nextci[k];
1621:         Jptr = buf_rj[k] + *nextci[k];
1622:         PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);
1623:         nextrow[k]++; nextci[k]++;
1624:       }
1625:     }
1626:     nzi = lnk[0];

1628:     /* copy data into free space, then initialize lnk */
1629:     PetscLLCondensedClean(aN,nzi,current_space->array,lnk,lnkbt);
1630:     MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);
1631:   }
1632:   PetscFree3(buf_ri_k,nextrow,nextci);
1633:   PetscLLDestroy(lnk,lnkbt);
1634:   PetscFreeSpaceDestroy(free_space);

1636:   /* local sizes and preallocation */
1637:   MatSetSizes(Cmpi,pn,an,PETSC_DETERMINE,PETSC_DETERMINE);
1638:   if (P->cmap->bs > 0) {PetscLayoutSetBlockSize(Cmpi->rmap,P->cmap->bs);}
1639:   if (A->cmap->bs > 0) {PetscLayoutSetBlockSize(Cmpi->cmap,A->cmap->bs);}
1640:   MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);
1641:   MatPreallocateFinalize(dnz,onz);

1643:   /* members in merge */
1644:   PetscFree(id_r);
1645:   PetscFree(len_r);
1646:   PetscFree(buf_ri[0]);
1647:   PetscFree(buf_ri);
1648:   PetscFree(buf_rj[0]);
1649:   PetscFree(buf_rj);
1650:   PetscLayoutDestroy(&rowmap);

1652:   /* attach the supporting struct to Cmpi for reuse */
1653:   c = (Mat_MPIAIJ*)Cmpi->data;
1654:   c->ap         = ptap;
1655:   ptap->destroy = Cmpi->ops->destroy;

1657:   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
1658:   Cmpi->assembled        = PETSC_FALSE;
1659:   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
1660:   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;

1662:   *C                     = Cmpi;
1663:   return(0);
1664: }

1666: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
1667: {
1668:   PetscErrorCode    ierr;
1669:   Mat_MPIAIJ        *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1670:   Mat_SeqAIJ        *c_seq;
1671:   Mat_APMPI         *ptap = c->ap;
1672:   Mat               A_loc,C_loc,C_oth;
1673:   PetscInt          i,rstart,rend,cm,ncols,row;
1674:   const PetscInt    *cols;
1675:   const PetscScalar *vals;

1678:   if (!ptap->A_loc) {
1679:     MPI_Comm comm;
1680:     PetscObjectGetComm((PetscObject)C,&comm);
1681:     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
1682:   }

1684:   MatZeroEntries(C);

1686:   if (ptap->reuse == MAT_REUSE_MATRIX) {
1687:     /* These matrices are obtained in MatTransposeMatMultSymbolic() */
1688:     /* 1) get R = Pd^T, Ro = Po^T */
1689:     /*----------------------------*/
1690:     MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);
1691:     MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);

1693:     /* 2) compute numeric A_loc */
1694:     /*--------------------------*/
1695:     MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&ptap->A_loc);
1696:   }

1698:   /* 3) C_loc = Rd*A_loc, C_oth = Ro*A_loc */
1699:   A_loc = ptap->A_loc;
1700:   ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,A_loc,ptap->C_loc);
1701:   ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,A_loc,ptap->C_oth);
1702:   C_loc = ptap->C_loc;
1703:   C_oth = ptap->C_oth;

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

1708:   /* C_loc -> C */
1709:   cm    = C_loc->rmap->N;
1710:   c_seq = (Mat_SeqAIJ*)C_loc->data;
1711:   cols = c_seq->j;
1712:   vals = c_seq->a;
1713:   for (i=0; i<cm; i++) {
1714:     ncols = c_seq->i[i+1] - c_seq->i[i];
1715:     row = rstart + i;
1716:     MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);
1717:     cols += ncols; vals += ncols;
1718:   }

1720:   /* Co -> C, off-processor part */
1721:   cm    = C_oth->rmap->N;
1722:   c_seq = (Mat_SeqAIJ*)C_oth->data;
1723:   cols  = c_seq->j;
1724:   vals  = c_seq->a;
1725:   for (i=0; i<cm; i++) {
1726:     ncols = c_seq->i[i+1] - c_seq->i[i];
1727:     row = p->garray[i];
1728:     MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);
1729:     cols += ncols; vals += ncols;
1730:   }
1731:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1732:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

1734:   ptap->reuse = MAT_REUSE_MATRIX;

1736:   /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */
1737:   if (ptap->freestruct) {
1738:     MatFreeIntermediateDataStructures(C);
1739:   }
1740:   return(0);
1741: }

1743: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1744: {
1745:   PetscErrorCode      ierr;
1746:   Mat_Merge_SeqsToMPI *merge;
1747:   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1748:   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1749:   Mat_APMPI           *ptap;
1750:   PetscInt            *adj;
1751:   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1752:   MatScalar           *ada,*ca,valtmp;
1753:   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1754:   MPI_Comm            comm;
1755:   PetscMPIInt         size,rank,taga,*len_s;
1756:   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1757:   PetscInt            **buf_ri,**buf_rj;
1758:   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1759:   MPI_Request         *s_waits,*r_waits;
1760:   MPI_Status          *status;
1761:   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1762:   PetscInt            *ai,*aj,*coi,*coj,*poJ,*pdJ;
1763:   Mat                 A_loc;
1764:   Mat_SeqAIJ          *a_loc;

1767:   PetscObjectGetComm((PetscObject)C,&comm);
1768:   MPI_Comm_size(comm,&size);
1769:   MPI_Comm_rank(comm,&rank);

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

1775:   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1776:   /*------------------------------------------*/
1777:   /* get data from symbolic products */
1778:   coi    = merge->coi; coj = merge->coj;
1779:   PetscCalloc1(coi[pon]+1,&coa);
1780:   bi     = merge->bi; bj = merge->bj;
1781:   owners = merge->rowmap->range;
1782:   PetscCalloc1(bi[cm]+1,&ba);

1784:   /* get A_loc by taking all local rows of A */
1785:   A_loc = ptap->A_loc;
1786:   MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);
1787:   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1788:   ai    = a_loc->i;
1789:   aj    = a_loc->j;

1791:   for (i=0; i<am; i++) {
1792:     anz = ai[i+1] - ai[i];
1793:     adj = aj + ai[i];
1794:     ada = a_loc->a + ai[i];

1796:     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1797:     /*-------------------------------------------------------------*/
1798:     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1799:     pnz = po->i[i+1] - po->i[i];
1800:     poJ = po->j + po->i[i];
1801:     pA  = po->a + po->i[i];
1802:     for (j=0; j<pnz; j++) {
1803:       row = poJ[j];
1804:       cj  = coj + coi[row];
1805:       ca  = coa + coi[row];
1806:       /* perform sparse axpy */
1807:       nexta  = 0;
1808:       valtmp = pA[j];
1809:       for (k=0; nexta<anz; k++) {
1810:         if (cj[k] == adj[nexta]) {
1811:           ca[k] += valtmp*ada[nexta];
1812:           nexta++;
1813:         }
1814:       }
1815:       PetscLogFlops(2.0*anz);
1816:     }

1818:     /* put the value into Cd (diagonal part) */
1819:     pnz = pd->i[i+1] - pd->i[i];
1820:     pdJ = pd->j + pd->i[i];
1821:     pA  = pd->a + pd->i[i];
1822:     for (j=0; j<pnz; j++) {
1823:       row = pdJ[j];
1824:       cj  = bj + bi[row];
1825:       ca  = ba + bi[row];
1826:       /* perform sparse axpy */
1827:       nexta  = 0;
1828:       valtmp = pA[j];
1829:       for (k=0; nexta<anz; k++) {
1830:         if (cj[k] == adj[nexta]) {
1831:           ca[k] += valtmp*ada[nexta];
1832:           nexta++;
1833:         }
1834:       }
1835:       PetscLogFlops(2.0*anz);
1836:     }
1837:   }

1839:   /* 3) send and recv matrix values coa */
1840:   /*------------------------------------*/
1841:   buf_ri = merge->buf_ri;
1842:   buf_rj = merge->buf_rj;
1843:   len_s  = merge->len_s;
1844:   PetscCommGetNewTag(comm,&taga);
1845:   PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);

1847:   PetscMalloc2(merge->nsend+1,&s_waits,size,&status);
1848:   for (proc=0,k=0; proc<size; proc++) {
1849:     if (!len_s[proc]) continue;
1850:     i    = merge->owners_co[proc];
1851:     MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);
1852:     k++;
1853:   }
1854:   if (merge->nrecv) {MPI_Waitall(merge->nrecv,r_waits,status);}
1855:   if (merge->nsend) {MPI_Waitall(merge->nsend,s_waits,status);}

1857:   PetscFree2(s_waits,status);
1858:   PetscFree(r_waits);
1859:   PetscFree(coa);

1861:   /* 4) insert local Cseq and received values into Cmpi */
1862:   /*----------------------------------------------------*/
1863:   PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);
1864:   for (k=0; k<merge->nrecv; k++) {
1865:     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1866:     nrows       = *(buf_ri_k[k]);
1867:     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1868:     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1869:   }

1871:   for (i=0; i<cm; i++) {
1872:     row  = owners[rank] + i; /* global row index of C_seq */
1873:     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1874:     ba_i = ba + bi[i];
1875:     bnz  = bi[i+1] - bi[i];
1876:     /* add received vals into ba */
1877:     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1878:       /* i-th row */
1879:       if (i == *nextrow[k]) {
1880:         cnz    = *(nextci[k]+1) - *nextci[k];
1881:         cj     = buf_rj[k] + *(nextci[k]);
1882:         ca     = abuf_r[k] + *(nextci[k]);
1883:         nextcj = 0;
1884:         for (j=0; nextcj<cnz; j++) {
1885:           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1886:             ba_i[j] += ca[nextcj++];
1887:           }
1888:         }
1889:         nextrow[k]++; nextci[k]++;
1890:         PetscLogFlops(2.0*cnz);
1891:       }
1892:     }
1893:     MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);
1894:   }
1895:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1896:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

1898:   PetscFree(ba);
1899:   PetscFree(abuf_r[0]);
1900:   PetscFree(abuf_r);
1901:   PetscFree3(buf_ri_k,nextrow,nextci);

1903:   if (ptap->freestruct) {
1904:     MatFreeIntermediateDataStructures(C);
1905:   }
1906:   return(0);
1907: }

1909: PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1910: {
1911:   PetscErrorCode      ierr;
1912:   Mat                 Cmpi,A_loc,POt,PDt;
1913:   Mat_APMPI           *ptap;
1914:   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1915:   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*a=(Mat_MPIAIJ*)A->data,*c;
1916:   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1917:   PetscInt            nnz;
1918:   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1919:   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1920:   MPI_Comm            comm;
1921:   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1922:   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1923:   PetscInt            len,proc,*dnz,*onz,*owners;
1924:   PetscInt            nzi,*bi,*bj;
1925:   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1926:   MPI_Request         *swaits,*rwaits;
1927:   MPI_Status          *sstatus,rstatus;
1928:   Mat_Merge_SeqsToMPI *merge;
1929:   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1930:   PetscReal           afill  =1.0,afill_tmp;
1931:   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax;
1932:   Mat_SeqAIJ          *a_loc,*pdt,*pot;
1933:   PetscTable          ta;
1934:   MatType             mtype;

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

1941:   MPI_Comm_size(comm,&size);
1942:   MPI_Comm_rank(comm,&rank);

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

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

1950:   ptap->A_loc = A_loc;
1951:   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1952:   ai          = a_loc->i;
1953:   aj          = a_loc->j;

1955:   /* determine symbolic Co=(p->B)^T*A - send to others */
1956:   /*----------------------------------------------------*/
1957:   MatTransposeSymbolic_SeqAIJ(p->A,&PDt);
1958:   pdt  = (Mat_SeqAIJ*)PDt->data;
1959:   pdti = pdt->i; pdtj = pdt->j;

1961:   MatTransposeSymbolic_SeqAIJ(p->B,&POt);
1962:   pot  = (Mat_SeqAIJ*)POt->data;
1963:   poti = pot->i; potj = pot->j;

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

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

1976:   /* create and initialize a linked list */
1977:   PetscTableCreate(A->cmap->n + a->B->cmap->N,aN,&ta);
1978:   MatRowMergeMax_SeqAIJ(a_loc,am,ta);
1979:   PetscTableGetCount(ta,&Armax);

1981:   PetscLLCondensedCreate_Scalable(Armax,&lnk);

1983:   for (i=0; i<pon; i++) {
1984:     pnz = poti[i+1] - poti[i];
1985:     ptJ = potj + poti[i];
1986:     for (j=0; j<pnz; j++) {
1987:       row  = ptJ[j]; /* row of A_loc == col of Pot */
1988:       anz  = ai[row+1] - ai[row];
1989:       Jptr = aj + ai[row];
1990:       /* add non-zero cols of AP into the sorted linked list lnk */
1991:       PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);
1992:     }
1993:     nnz = lnk[0];

1995:     /* If free space is not available, double the total space in the list */
1996:     if (current_space->local_remaining<nnz) {
1997:       PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);
1998:       nspacedouble++;
1999:     }

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

2004:     current_space->array           += nnz;
2005:     current_space->local_used      += nnz;
2006:     current_space->local_remaining -= nnz;

2008:     coi[i+1] = coi[i] + nnz;
2009:   }

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

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

2018:   /* send j-array (coj) of Co to other processors */
2019:   /*----------------------------------------------*/
2020:   /* determine row ownership */
2021:   PetscNew(&merge);
2022:   PetscLayoutCreate(comm,&merge->rowmap);

2024:   merge->rowmap->n  = pn;
2025:   merge->rowmap->bs = 1;

2027:   PetscLayoutSetUp(merge->rowmap);
2028:   owners = merge->rowmap->range;

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

2034:   len_s        = merge->len_s;
2035:   merge->nsend = 0;

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

2039:   proc = 0;
2040:   for (i=0; i<pon; i++) {
2041:     while (prmap[i] >= owners[proc+1]) proc++;
2042:     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
2043:     len_s[proc] += coi[i+1] - coi[i];
2044:   }

2046:   len          = 0; /* max length of buf_si[] */
2047:   owners_co[0] = 0;
2048:   for (proc=0; proc<size; proc++) {
2049:     owners_co[proc+1] = owners_co[proc] + len_si[proc];
2050:     if (len_si[proc]) {
2051:       merge->nsend++;
2052:       len_si[proc] = 2*(len_si[proc] + 1);
2053:       len         += len_si[proc];
2054:     }
2055:   }

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

2061:   /* post the Irecv and Isend of coj */
2062:   PetscCommGetNewTag(comm,&tagj);
2063:   PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);
2064:   PetscMalloc1(merge->nsend+1,&swaits);
2065:   for (proc=0, k=0; proc<size; proc++) {
2066:     if (!len_s[proc]) continue;
2067:     i    = owners_co[proc];
2068:     MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);
2069:     k++;
2070:   }

2072:   /* receives and sends of coj are complete */
2073:   PetscMalloc1(size,&sstatus);
2074:   for (i=0; i<merge->nrecv; i++) {
2075:     PetscMPIInt icompleted;
2076:     MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);
2077:   }
2078:   PetscFree(rwaits);
2079:   if (merge->nsend) {MPI_Waitall(merge->nsend,swaits,sstatus);}

2081:   /* add received column indices into table to update Armax */
2082:   /* Armax can be as large as aN if a P[row,:] is dense, see src/ksp/ksp/examples/tutorials/ex56.c! */
2083:   for (k=0; k<merge->nrecv; k++) {/* k-th received message */
2084:     Jptr = buf_rj[k];
2085:     for (j=0; j<merge->len_r[k]; j++) {
2086:       PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);
2087:     }
2088:   }
2089:   PetscTableGetCount(ta,&Armax);
2090:   /* 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); */

2092:   /* send and recv coi */
2093:   /*-------------------*/
2094:   PetscCommGetNewTag(comm,&tagi);
2095:   PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);
2096:   PetscMalloc1(len+1,&buf_s);
2097:   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
2098:   for (proc=0,k=0; proc<size; proc++) {
2099:     if (!len_s[proc]) continue;
2100:     /* form outgoing message for i-structure:
2101:          buf_si[0]:                 nrows to be sent
2102:                [1:nrows]:           row index (global)
2103:                [nrows+1:2*nrows+1]: i-structure index
2104:     */
2105:     /*-------------------------------------------*/
2106:     nrows       = len_si[proc]/2 - 1;
2107:     buf_si_i    = buf_si + nrows+1;
2108:     buf_si[0]   = nrows;
2109:     buf_si_i[0] = 0;
2110:     nrows       = 0;
2111:     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
2112:       nzi               = coi[i+1] - coi[i];
2113:       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
2114:       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
2115:       nrows++;
2116:     }
2117:     MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);
2118:     k++;
2119:     buf_si += len_si[proc];
2120:   }
2121:   i = merge->nrecv;
2122:   while (i--) {
2123:     PetscMPIInt icompleted;
2124:     MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);
2125:   }
2126:   PetscFree(rwaits);
2127:   if (merge->nsend) {MPI_Waitall(merge->nsend,swaits,sstatus);}
2128:   PetscFree(len_si);
2129:   PetscFree(len_ri);
2130:   PetscFree(swaits);
2131:   PetscFree(sstatus);
2132:   PetscFree(buf_s);

2134:   /* compute the local portion of C (mpi mat) */
2135:   /*------------------------------------------*/
2136:   /* allocate bi array and free space for accumulating nonzero column info */
2137:   PetscMalloc1(pn+1,&bi);
2138:   bi[0] = 0;

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

2145:   PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);
2146:   for (k=0; k<merge->nrecv; k++) {
2147:     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
2148:     nrows       = *buf_ri_k[k];
2149:     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
2150:     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
2151:   }

2153:   PetscLLCondensedCreate_Scalable(Armax,&lnk);
2154:   MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);
2155:   rmax = 0;
2156:   for (i=0; i<pn; i++) {
2157:     /* add pdt[i,:]*AP into lnk */
2158:     pnz = pdti[i+1] - pdti[i];
2159:     ptJ = pdtj + pdti[i];
2160:     for (j=0; j<pnz; j++) {
2161:       row  = ptJ[j];  /* row of AP == col of Pt */
2162:       anz  = ai[row+1] - ai[row];
2163:       Jptr = aj + ai[row];
2164:       /* add non-zero cols of AP into the sorted linked list lnk */
2165:       PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);
2166:     }

2168:     /* add received col data into lnk */
2169:     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
2170:       if (i == *nextrow[k]) { /* i-th row */
2171:         nzi  = *(nextci[k]+1) - *nextci[k];
2172:         Jptr = buf_rj[k] + *nextci[k];
2173:         PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);
2174:         nextrow[k]++; nextci[k]++;
2175:       }
2176:     }
2177:     nnz = lnk[0];

2179:     /* if free space is not available, make more free space */
2180:     if (current_space->local_remaining<nnz) {
2181:       PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);
2182:       nspacedouble++;
2183:     }
2184:     /* copy data into free space, then initialize lnk */
2185:     PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);
2186:     MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);

2188:     current_space->array           += nnz;
2189:     current_space->local_used      += nnz;
2190:     current_space->local_remaining -= nnz;

2192:     bi[i+1] = bi[i] + nnz;
2193:     if (nnz > rmax) rmax = nnz;
2194:   }
2195:   PetscFree3(buf_ri_k,nextrow,nextci);

2197:   PetscMalloc1(bi[pn]+1,&bj);
2198:   PetscFreeSpaceContiguous(&free_space,bj);
2199:   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
2200:   if (afill_tmp > afill) afill = afill_tmp;
2201:   PetscLLCondensedDestroy_Scalable(lnk);
2202:   PetscTableDestroy(&ta);

2204:   MatDestroy(&POt);
2205:   MatDestroy(&PDt);

2207:   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
2208:   /*----------------------------------------------------------------------------------*/
2209:   MatCreate(comm,&Cmpi);
2210:   MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);
2211:   MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));
2212:   MatGetType(A,&mtype);
2213:   MatSetType(Cmpi,mtype);
2214:   MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);
2215:   MatPreallocateFinalize(dnz,onz);
2216:   MatSetBlockSize(Cmpi,1);
2217:   for (i=0; i<pn; i++) {
2218:     row  = i + rstart;
2219:     nnz  = bi[i+1] - bi[i];
2220:     Jptr = bj + bi[i];
2221:     MatSetValues(Cmpi,1,&row,nnz,Jptr,NULL,INSERT_VALUES);
2222:   }
2223:   MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);
2224:   MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);
2225:   merge->bi        = bi;
2226:   merge->bj        = bj;
2227:   merge->coi       = coi;
2228:   merge->coj       = coj;
2229:   merge->buf_ri    = buf_ri;
2230:   merge->buf_rj    = buf_rj;
2231:   merge->owners_co = owners_co;

2233:   /* attach the supporting struct to Cmpi for reuse */
2234:   c = (Mat_MPIAIJ*)Cmpi->data;

2236:   c->ap       = ptap;
2237:   ptap->api   = NULL;
2238:   ptap->apj   = NULL;
2239:   ptap->merge = merge;
2240:   ptap->apa   = NULL;
2241:   ptap->destroy   = Cmpi->ops->destroy;
2242:   ptap->duplicate = Cmpi->ops->duplicate;

2244:   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
2245:   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
2246:   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;

2248:   *C = Cmpi;
2249: #if defined(PETSC_USE_INFO)
2250:   if (bi[pn] != 0) {
2251:     PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
2252:     PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);
2253:   } else {
2254:     PetscInfo(Cmpi,"Empty matrix product\n");
2255:   }
2256: #endif
2257:   return(0);
2258: }