Actual source code: mpimatmatmult.c

petsc-3.11.4 2019-09-28
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) 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:   PetscScalar        *apa;
214:   PetscReal          afill;
215:   MatType            mtype;

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

301:   ptap->apa = apa;

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

309:   MatGetType(A,&mtype);
310:   MatSetType(Cmpi,mtype);
311:   MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);

313:   MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);
314:   MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);
315:   MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);
316:   MatPreallocateFinalize(dnz,onz);

318:   ptap->destroy        = Cmpi->ops->destroy;
319:   ptap->duplicate      = Cmpi->ops->duplicate;
320:   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
321:   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;
322:   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;

324:   /* attach the supporting struct to Cmpi for reuse */
325:   c       = (Mat_MPIAIJ*)Cmpi->data;
326:   c->ap = ptap;

328:   *C = Cmpi;

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

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

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

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

364: typedef struct {
365:   Mat         workB;
366:   PetscScalar *rvalues,*svalues;
367:   MPI_Request *rwaits,*swaits;
368: } MPIAIJ_MPIDense;

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

376:   MatDestroy(&contents->workB);
377:   PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);
378:   PetscFree(contents);
379:   return(0);
380: }

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

386:   It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C
387: */
388: PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C)
389: {
390:   PetscErrorCode         ierr;
391:   PetscBool              flg;
392:   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
393:   PetscInt               nz   = aij->B->cmap->n,to_n,to_entries,from_n,from_entries;
394:   PetscContainer         container;
395:   MPIAIJ_MPIDense        *contents;
396:   VecScatter             ctx   = aij->Mvctx;

399:   PetscObjectTypeCompare((PetscObject)B,MATMPIDENSE,&flg);
400:   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Second matrix must be mpidense");

402:   /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/
403:   PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);
404:   if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"First matrix must be MPIAIJ");

406:   C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;

408:   PetscNew(&contents);
409:   /* Create work matrix used to store off processor rows of B needed for local product */
410:   MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);
411:   /* Create work arrays needed */
412:   VecScatterGetRemoteCount_Private(ctx,PETSC_TRUE/*send*/,&to_n,&to_entries);
413:   VecScatterGetRemoteCount_Private(ctx,PETSC_FALSE/*recv*/,&from_n,&from_entries);
414:   PetscMalloc4(B->cmap->N*from_entries,&contents->rvalues,B->cmap->N*to_entries,&contents->svalues,from_n,&contents->rwaits,to_n,&contents->swaits);

416:   PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);
417:   PetscContainerSetPointer(container,contents);
418:   PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);
419:   PetscObjectCompose((PetscObject)C,"workB",(PetscObject)container);
420:   PetscContainerDestroy(&container);

422:   (*C->ops->matmultnumeric)(A,B,C);
423:   return(0);
424: }

426: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
427: {
428:   PetscErrorCode         ierr;
429:   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*) A->data;
430:   PetscInt               nz   = aij->B->cmap->n,to_n,to_entries,from_n,from_entries;
431:   PetscContainer         container;
432:   MPIAIJ_MPIDense        *contents;
433:   VecScatter             ctx   = aij->Mvctx;
434:   PetscInt               m     = A->rmap->n,n=B->cmap->n;

437:   MatCreate(PetscObjectComm((PetscObject)B),C);
438:   MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);
439:   MatSetBlockSizesFromMats(*C,A,B);
440:   MatSetType(*C,MATMPIDENSE);
441:   MatMPIDenseSetPreallocation(*C,NULL);
442:   MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);
443:   MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);

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

447:   PetscNew(&contents);
448:   /* Create work matrix used to store off processor rows of B needed for local product */
449:   MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);
450:   /* Create work arrays needed */
451:   VecScatterGetRemoteCount_Private(ctx,PETSC_TRUE/*send*/,&to_n,&to_entries);
452:   VecScatterGetRemoteCount_Private(ctx,PETSC_FALSE/*recv*/,&from_n,&from_entries);
453:   PetscMalloc4(B->cmap->N*from_entries,&contents->rvalues,B->cmap->N*to_entries,&contents->svalues,from_n,&contents->rwaits,to_n,&contents->swaits);

455:   PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);
456:   PetscContainerSetPointer(container,contents);
457:   PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);
458:   PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);
459:   PetscContainerDestroy(&container);
460:   return(0);
461: }

463: /*
464:     Performs an efficient scatter on the rows of B needed by this process; this is
465:     a modification of the VecScatterBegin_() routines.
466: */
467: PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
468: {
469:   Mat_MPIAIJ             *aij = (Mat_MPIAIJ*)A->data;
470:   PetscErrorCode         ierr;
471:   PetscScalar            *b,*w,*svalues,*rvalues;
472:   VecScatter             ctx   = aij->Mvctx;
473:   PetscInt               i,j,k;
474:   const PetscInt         *sindices,*sstarts,*rindices,*rstarts;
475:   const PetscMPIInt      *sprocs,*rprocs;
476:   PetscInt               nsends,nrecvs,nrecvs2;
477:   MPI_Request            *swaits,*rwaits;
478:   MPI_Comm               comm;
479:   PetscMPIInt            tag  = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n,nsends_mpi,nrecvs_mpi;
480:   MPI_Status             status;
481:   MPIAIJ_MPIDense        *contents;
482:   PetscContainer         container;
483:   Mat                    workB;

486:   PetscObjectGetComm((PetscObject)A,&comm);
487:   PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);
488:   if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
489:   PetscContainerGetPointer(container,(void**)&contents);

491:   workB = *outworkB = contents->workB;
492:   if (nrows != workB->rmap->n) SETERRQ2(comm,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",nrows,workB->cmap->n);
493:   VecScatterGetRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,&sprocs,NULL/*bs*/);
494:   VecScatterGetRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,&rindices,&rprocs,NULL/*bs*/);
495:   PetscMPIIntCast(nsends,&nsends_mpi);
496:   PetscMPIIntCast(nrecvs,&nrecvs_mpi);
497:   svalues = contents->svalues;
498:   rvalues = contents->rvalues;
499:   swaits  = contents->swaits;
500:   rwaits  = contents->rwaits;

502:   MatDenseGetArray(B,&b);
503:   MatDenseGetArray(workB,&w);

505:   for (i=0; i<nrecvs; i++) {
506:     MPI_Irecv(rvalues+ncols*(rstarts[i]-rstarts[0]),ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);
507:   }

509:   for (i=0; i<nsends; i++) {
510:     /* pack a message at a time */
511:     for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
512:       for (k=0; k<ncols; k++) {
513:         svalues[ncols*(sstarts[i]-sstarts[0]+j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
514:       }
515:     }
516:     MPI_Isend(svalues+ncols*(sstarts[i]-sstarts[0]),ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);
517:   }

519:   nrecvs2 = nrecvs;
520:   while (nrecvs2) {
521:     MPI_Waitany(nrecvs_mpi,rwaits,&imdex,&status);
522:     nrecvs2--;
523:     /* unpack a message at a time */
524:     for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
525:       for (k=0; k<ncols; k++) {
526:         w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex]-rstarts[0]+j) + k];
527:       }
528:     }
529:   }
530:   if (nsends) {MPI_Waitall(nsends_mpi,swaits,MPI_STATUSES_IGNORE);}

532:   VecScatterRestoreRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,&sprocs,NULL/*bs*/);
533:   VecScatterRestoreRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,&rindices,&rprocs,NULL/*bs*/);
534:   MatDenseRestoreArray(B,&b);
535:   MatDenseRestoreArray(workB,&w);
536:   MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);
537:   MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);
538:   return(0);

540: }
541: extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);

543: PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
544: {
546:   Mat_MPIAIJ     *aij    = (Mat_MPIAIJ*)A->data;
547:   Mat_MPIDense   *bdense = (Mat_MPIDense*)B->data;
548:   Mat_MPIDense   *cdense = (Mat_MPIDense*)C->data;
549:   Mat            workB;

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

555:   /* get off processor parts of B needed to complete the product */
556:   MatMPIDenseScatter(A,B,C,&workB);

558:   /* off-diagonal block of A times nonlocal rows of B */
559:   MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);
560:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
561:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
562:   return(0);
563: }

565: PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
566: {
568:   Mat_MPIAIJ     *a   = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
569:   Mat_SeqAIJ     *ad  = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
570:   Mat_SeqAIJ     *cd  = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
571:   PetscInt       *adi = ad->i,*adj,*aoi=ao->i,*aoj;
572:   PetscScalar    *ada,*aoa,*cda=cd->a,*coa=co->a;
573:   Mat_SeqAIJ     *p_loc,*p_oth;
574:   PetscInt       *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
575:   PetscScalar    *pa_loc,*pa_oth,*pa,valtmp,*ca;
576:   PetscInt       cm    = C->rmap->n,anz,pnz;
577:   Mat_APMPI      *ptap = c->ap;
578:   PetscScalar    *apa_sparse;
579:   PetscInt       *api,*apj,*apJ,i,j,k,row;
580:   PetscInt       cstart = C->cmap->rstart;
581:   PetscInt       cdnz,conz,k0,k1,nextp;
582:   MPI_Comm       comm;
583:   PetscMPIInt    size;

586:   PetscObjectGetComm((PetscObject)C,&comm);
587:   MPI_Comm_size(comm,&size);

589:   if (!ptap) {
590:     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"AP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
591:   }
592:   apa_sparse = ptap->apa;

594:   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
595:   /*-----------------------------------------------------*/
596:   /* update numerical values of P_oth and P_loc */
597:   MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);
598:   MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);

600:   /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
601:   /*----------------------------------------------------------*/
602:   /* get data from symbolic products */
603:   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
604:   pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a;
605:   if (size >1) {
606:     p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
607:     pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a;
608:   } else {
609:     p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL;
610:   }

612:   api = ptap->api;
613:   apj = ptap->apj;
614:   for (i=0; i<cm; i++) {
615:     apJ = apj + api[i];

617:     /* diagonal portion of A */
618:     anz = adi[i+1] - adi[i];
619:     adj = ad->j + adi[i];
620:     ada = ad->a + adi[i];
621:     for (j=0; j<anz; j++) {
622:       row = adj[j];
623:       pnz = pi_loc[row+1] - pi_loc[row];
624:       pj  = pj_loc + pi_loc[row];
625:       pa  = pa_loc + pi_loc[row];
626:       /* perform sparse axpy */
627:       valtmp = ada[j];
628:       nextp  = 0;
629:       for (k=0; nextp<pnz; k++) {
630:         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
631:           apa_sparse[k] += valtmp*pa[nextp++];
632:         }
633:       }
634:       PetscLogFlops(2.0*pnz);
635:     }

637:     /* off-diagonal portion of A */
638:     anz = aoi[i+1] - aoi[i];
639:     aoj = ao->j + aoi[i];
640:     aoa = ao->a + aoi[i];
641:     for (j=0; j<anz; j++) {
642:       row = aoj[j];
643:       pnz = pi_oth[row+1] - pi_oth[row];
644:       pj  = pj_oth + pi_oth[row];
645:       pa  = pa_oth + pi_oth[row];
646:       /* perform sparse axpy */
647:       valtmp = aoa[j];
648:       nextp  = 0;
649:       for (k=0; nextp<pnz; k++) {
650:         if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
651:           apa_sparse[k] += valtmp*pa[nextp++];
652:         }
653:       }
654:       PetscLogFlops(2.0*pnz);
655:     }

657:     /* set values in C */
658:     cdnz = cd->i[i+1] - cd->i[i];
659:     conz = co->i[i+1] - co->i[i];

661:     /* 1st off-diagoanl part of C */
662:     ca = coa + co->i[i];
663:     k  = 0;
664:     for (k0=0; k0<conz; k0++) {
665:       if (apJ[k] >= cstart) break;
666:       ca[k0]        = apa_sparse[k];
667:       apa_sparse[k] = 0.0;
668:       k++;
669:     }

671:     /* diagonal part of C */
672:     ca = cda + cd->i[i];
673:     for (k1=0; k1<cdnz; k1++) {
674:       ca[k1]        = apa_sparse[k];
675:       apa_sparse[k] = 0.0;
676:       k++;
677:     }

679:     /* 2nd off-diagoanl part of C */
680:     ca = coa + co->i[i];
681:     for (; k0<conz; k0++) {
682:       ca[k0]        = apa_sparse[k];
683:       apa_sparse[k] = 0.0;
684:       k++;
685:     }
686:   }
687:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
688:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

690:   if (ptap->freestruct) {
691:     MatFreeIntermediateDataStructures(C);
692:   }
693:   return(0);
694: }

696: /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
697: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
698: {
699:   PetscErrorCode     ierr;
700:   MPI_Comm           comm;
701:   PetscMPIInt        size;
702:   Mat                Cmpi;
703:   Mat_APMPI          *ptap;
704:   PetscFreeSpaceList free_space = NULL,current_space=NULL;
705:   Mat_MPIAIJ         *a         = (Mat_MPIAIJ*)A->data,*c;
706:   Mat_SeqAIJ         *ad        = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
707:   PetscInt           *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
708:   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
709:   PetscInt           i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
710:   PetscInt           am=A->rmap->n,pn=P->cmap->n,pm=P->rmap->n,lsize=pn+20;
711:   PetscReal          afill;
712:   PetscScalar        *apa;
713:   MatType            mtype;

716:   PetscObjectGetComm((PetscObject)A,&comm);
717:   MPI_Comm_size(comm,&size);

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

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

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

728:   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
729:   pi_loc = p_loc->i; pj_loc = p_loc->j;
730:   if (size > 1) {
731:     p_oth  = (Mat_SeqAIJ*)(ptap->P_oth)->data;
732:     pi_oth = p_oth->i; pj_oth = p_oth->j;
733:   } else {
734:     p_oth  = NULL;
735:     pi_oth = NULL; pj_oth = NULL;
736:   }

738:   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
739:   /*-------------------------------------------------------------------*/
740:   PetscMalloc1(am+2,&api);
741:   ptap->api = api;
742:   api[0]    = 0;

744:   PetscLLCondensedCreate_Scalable(lsize,&lnk);

746:   /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
747:   PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);
748:   current_space = free_space;
749:   MatPreallocateInitialize(comm,am,pn,dnz,onz);
750:   for (i=0; i<am; i++) {
751:     /* diagonal portion of A */
752:     nzi = adi[i+1] - adi[i];
753:     for (j=0; j<nzi; j++) {
754:       row  = *adj++;
755:       pnz  = pi_loc[row+1] - pi_loc[row];
756:       Jptr = pj_loc + pi_loc[row];
757:       /* Expand list if it is not long enough */
758:       if (pnz+apnz_max > lsize) {
759:         lsize = pnz+apnz_max;
760:         PetscLLCondensedExpand_Scalable(lsize, &lnk);
761:       }
762:       /* add non-zero cols of P into the sorted linked list lnk */
763:       PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);
764:       apnz     = *lnk; /* The first element in the list is the number of items in the list */
765:       api[i+1] = api[i] + apnz;
766:       if (apnz > apnz_max) apnz_max = apnz;
767:     }
768:     /* off-diagonal portion of A */
769:     nzi = aoi[i+1] - aoi[i];
770:     for (j=0; j<nzi; j++) {
771:       row  = *aoj++;
772:       pnz  = pi_oth[row+1] - pi_oth[row];
773:       Jptr = pj_oth + pi_oth[row];
774:       /* Expand list if it is not long enough */
775:       if (pnz+apnz_max > lsize) {
776:         lsize = pnz + apnz_max;
777:         PetscLLCondensedExpand_Scalable(lsize, &lnk);
778:       }
779:       /* add non-zero cols of P into the sorted linked list lnk */
780:       PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);
781:       apnz     = *lnk;  /* The first element in the list is the number of items in the list */
782:       api[i+1] = api[i] + apnz;
783:       if (apnz > apnz_max) apnz_max = apnz;
784:     }
785:     apnz     = *lnk;
786:     api[i+1] = api[i] + apnz;
787:     if (apnz > apnz_max) apnz_max = apnz;

789:     /* if free space is not available, double the total space in the list */
790:     if (current_space->local_remaining<apnz) {
791:       PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),&current_space);
792:       nspacedouble++;
793:     }

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

799:     current_space->array           += apnz;
800:     current_space->local_used      += apnz;
801:     current_space->local_remaining -= apnz;
802:   }

804:   /* Allocate space for apj, initialize apj, and */
805:   /* destroy list of free space and other temporary array(s) */
806:   PetscMalloc1(api[am]+1,&ptap->apj);
807:   apj  = ptap->apj;
808:   PetscFreeSpaceContiguous(&free_space,ptap->apj);
809:   PetscLLCondensedDestroy_Scalable(lnk);

811:   /* create and assemble symbolic parallel matrix Cmpi */
812:   /*----------------------------------------------------*/
813:   MatCreate(comm,&Cmpi);
814:   MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);
815:   MatSetBlockSizesFromMats(Cmpi,A,P);
816:   MatGetType(A,&mtype);
817:   MatSetType(Cmpi,mtype);
818:   MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);

820:   /* malloc apa for assembly Cmpi */
821:   PetscCalloc1(apnz_max,&apa);
822:   ptap->apa = apa;

824:   MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);
825:   MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);
826:   MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);
827:   MatPreallocateFinalize(dnz,onz);

829:   ptap->destroy             = Cmpi->ops->destroy;
830:   ptap->duplicate           = Cmpi->ops->duplicate;
831:   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
832:   Cmpi->ops->destroy        = MatDestroy_MPIAIJ_MatMatMult;
833:   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;

835:   /* attach the supporting struct to Cmpi for reuse */
836:   c       = (Mat_MPIAIJ*)Cmpi->data;
837:   c->ap = ptap;
838:   *C = Cmpi;

840:   /* set MatInfo */
841:   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
842:   if (afill < 1.0) afill = 1.0;
843:   Cmpi->info.mallocs           = nspacedouble;
844:   Cmpi->info.fill_ratio_given  = fill;
845:   Cmpi->info.fill_ratio_needed = afill;

847: #if defined(PETSC_USE_INFO)
848:   if (api[am]) {
849:     PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
850:     PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
851:   } else {
852:     PetscInfo(Cmpi,"Empty matrix product\n");
853:   }
854: #endif
855:   return(0);
856: }

858: /* This function is needed for the seqMPI matrix-matrix multiplication.  */
859: /* Three input arrays are merged to one output array. The size of the    */
860: /* output array is also output. Duplicate entries only show up once.     */
861: static void Merge3SortedArrays(PetscInt  size1, PetscInt *in1,
862:                                PetscInt  size2, PetscInt *in2,
863:                                PetscInt  size3, PetscInt *in3,
864:                                PetscInt *size4, PetscInt *out)
865: {
866:   int i = 0, j = 0, k = 0, l = 0;

868:   /* Traverse all three arrays */
869:   while (i<size1 && j<size2 && k<size3) {
870:     if (in1[i] < in2[j] && in1[i] < in3[k]) {
871:       out[l++] = in1[i++];
872:     }
873:     else if(in2[j] < in1[i] && in2[j] < in3[k]) {
874:       out[l++] = in2[j++];
875:     }
876:     else if(in3[k] < in1[i] && in3[k] < in2[j]) {
877:       out[l++] = in3[k++];
878:     }
879:     else if(in1[i] == in2[j] && in1[i] < in3[k]) {
880:       out[l++] = in1[i];
881:       i++, j++;
882:     }
883:     else if(in1[i] == in3[k] && in1[i] < in2[j]) {
884:       out[l++] = in1[i];
885:       i++, k++;
886:     }
887:     else if(in3[k] == in2[j] && in2[j] < in1[i])  {
888:       out[l++] = in2[j];
889:       k++, j++;
890:     }
891:     else if(in1[i] == in2[j] && in1[i] == in3[k]) {
892:       out[l++] = in1[i];
893:       i++, j++, k++;
894:     }
895:   }

897:   /* Traverse two remaining arrays */
898:   while (i<size1 && j<size2) {
899:     if (in1[i] < in2[j]) {
900:       out[l++] = in1[i++];
901:     }
902:     else if(in1[i] > in2[j]) {
903:       out[l++] = in2[j++];
904:     }
905:     else {
906:       out[l++] = in1[i];
907:       i++, j++;
908:     }
909:   }

911:   while (i<size1 && k<size3) {
912:     if (in1[i] < in3[k]) {
913:       out[l++] = in1[i++];
914:     }
915:     else if(in1[i] > in3[k]) {
916:       out[l++] = in3[k++];
917:     }
918:     else {
919:       out[l++] = in1[i];
920:       i++, k++;
921:     }
922:   }

924:   while (k<size3 && j<size2)  {
925:     if (in3[k] < in2[j]) {
926:       out[l++] = in3[k++];
927:     }
928:     else if(in3[k] > in2[j]) {
929:       out[l++] = in2[j++];
930:     }
931:     else {
932:       out[l++] = in3[k];
933:       k++, j++;
934:     }
935:   }

937:   /* Traverse one remaining array */
938:   while (i<size1) out[l++] = in1[i++];
939:   while (j<size2) out[l++] = in2[j++];
940:   while (k<size3) out[l++] = in3[k++];

942:   *size4 = l;
943: }

945: /* This matrix-matrix multiplication algorithm divides the multiplication into three multiplications and  */
946: /* adds up the products. Two of these three multiplications are performed with existing (sequential)      */
947: /* matrix-matrix multiplications.  */
948: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat A, Mat P, PetscReal fill, Mat *C)
949: {
950:   PetscErrorCode     ierr;
951:   MPI_Comm           comm;
952:   PetscMPIInt        size;
953:   Mat                Cmpi;
954:   Mat_APMPI          *ptap;
955:   PetscFreeSpaceList free_space_diag=NULL, current_space=NULL;
956:   Mat_MPIAIJ         *a        =(Mat_MPIAIJ*)A->data;
957:   Mat_SeqAIJ         *ad       =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc;
958:   Mat_MPIAIJ         *p        =(Mat_MPIAIJ*)P->data;
959:   Mat_MPIAIJ         *c;
960:   Mat_SeqAIJ         *adpd_seq, *p_off, *aopoth_seq;
961:   PetscInt           adponz, adpdnz;
962:   PetscInt           *pi_loc,*dnz,*onz;
963:   PetscInt           *adi=ad->i,*adj=ad->j,*aoi=ao->i,rstart=A->rmap->rstart;
964:   PetscInt           *lnk,i, i1=0,pnz,row,*adpoi,*adpoj, *api, *adpoJ, *aopJ, *apJ,*Jptr, aopnz, nspacedouble=0,j,nzi,
965:                      *apj,apnz, *adpdi, *adpdj, *adpdJ, *poff_i, *poff_j, *j_temp, *aopothi, *aopothj;
966:   PetscInt           am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n, p_colstart, p_colend;
967:   PetscBT            lnkbt;
968:   PetscScalar        *apa;
969:   PetscReal          afill;
970:   PetscMPIInt        rank;
971:   Mat                adpd, aopoth;
972:   MatType            mtype;
973:   const char         *prefix;

976:   PetscObjectGetComm((PetscObject)A,&comm);
977:   MPI_Comm_size(comm,&size);
978:   MPI_Comm_rank(comm, &rank);
979:   MatGetOwnershipRangeColumn(P, &p_colstart, &p_colend);

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

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

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


991:   p_loc  = (Mat_SeqAIJ*)(ptap->P_loc)->data;
992:   pi_loc = p_loc->i;

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

998:   adpoi[0]    = 0;
999:   ptap->api = api;
1000:   api[0] = 0;

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

1006:   /* Symbolic calc of A_loc_diag * P_loc_diag */
1007:   MatGetOptionsPrefix(A,&prefix);
1008:   MatSetOptionsPrefix(a->A,prefix);
1009:   MatAppendOptionsPrefix(a->A,"inner_diag_");
1010:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->A, p->A, fill, &adpd);
1011:   adpd_seq = (Mat_SeqAIJ*)((adpd)->data);
1012:   adpdi = adpd_seq->i; adpdj = adpd_seq->j;
1013:   p_off = (Mat_SeqAIJ*)((p->B)->data);
1014:   poff_i = p_off->i; poff_j = p_off->j;

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


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

1025:   for (i=0; i<am; i++) {
1026:     /* A_diag * P_loc_off */
1027:     nzi = adi[i+1] - adi[i];
1028:     for (j=0; j<nzi; j++) {
1029:       row  = *adj++;
1030:       pnz  = poff_i[row+1] - poff_i[row];
1031:       Jptr = poff_j + poff_i[row];
1032:       for(i1 = 0; i1 < pnz; i1++) {
1033:         j_temp[i1] = p->garray[Jptr[i1]];
1034:       }
1035:       /* add non-zero cols of P into the sorted linked list lnk */
1036:       PetscLLCondensedAddSorted(pnz,j_temp,lnk,lnkbt);
1037:     }

1039:     adponz     = lnk[0];
1040:     adpoi[i+1] = adpoi[i] + adponz;

1042:     /* if free space is not available, double the total space in the list */
1043:     if (current_space->local_remaining<adponz) {
1044:       PetscFreeSpaceGet(PetscIntSumTruncate(adponz,current_space->total_array_size),&current_space);
1045:       nspacedouble++;
1046:     }

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

1051:     current_space->array           += adponz;
1052:     current_space->local_used      += adponz;
1053:     current_space->local_remaining -= adponz;
1054:   }

1056:   /* Symbolic calc of A_off * P_oth */
1057:   MatSetOptionsPrefix(a->B,prefix);
1058:   MatAppendOptionsPrefix(a->B,"inner_offdiag_");
1059:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->B, ptap->P_oth, fill, &aopoth);
1060:   aopoth_seq = (Mat_SeqAIJ*)((aopoth)->data);
1061:   aopothi = aopoth_seq->i; aopothj = aopoth_seq->j;

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

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

1069:   /* Copy from linked list to j-array */
1070:   PetscFreeSpaceContiguous(&free_space_diag,adpoj);
1071:   PetscLLDestroy(lnk,lnkbt);

1073:   adpoJ = adpoj;
1074:   adpdJ = adpdj;
1075:   aopJ = aopothj;
1076:   apj  = ptap->apj;
1077:   apJ = apj; /* still empty */

1079:   /* Merge j-arrays of A_off * P, A_diag * P_loc_off, and */
1080:   /* A_diag * P_loc_diag to get A*P */
1081:   for (i = 0; i < am; i++) {
1082:     aopnz  =  aopothi[i+1] -  aopothi[i];
1083:     adponz = adpoi[i+1] - adpoi[i];
1084:     adpdnz = adpdi[i+1] - adpdi[i];

1086:     /* Correct indices from A_diag*P_diag */
1087:     for(i1 = 0; i1 < adpdnz; i1++) {
1088:       adpdJ[i1] += p_colstart;
1089:     }
1090:     /* Merge j-arrays of A_diag * P_loc_off and A_diag * P_loc_diag and A_off * P_oth */
1091:     Merge3SortedArrays(adponz, adpoJ, adpdnz, adpdJ, aopnz, aopJ, &apnz, apJ);
1092:     MatPreallocateSet(i+rstart, apnz, apJ, dnz, onz);

1094:     aopJ += aopnz;
1095:     adpoJ += adponz;
1096:     adpdJ += adpdnz;
1097:     apJ += apnz;
1098:     api[i+1] = api[i] + apnz;
1099:   }

1101:   /* malloc apa to store dense row A[i,:]*P */
1102:   PetscCalloc1(pN+2,&apa);

1104:   ptap->apa = apa;
1105:   /* create and assemble symbolic parallel matrix Cmpi */
1106:   MatCreate(comm,&Cmpi);
1107:   MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);
1108:   MatSetBlockSizesFromMats(Cmpi,A,P);
1109:   MatGetType(A,&mtype);
1110:   MatSetType(Cmpi,mtype);
1111:   MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);


1114:   MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);
1115:   MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);
1116:   MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);
1117:   MatPreallocateFinalize(dnz,onz);


1120:   ptap->destroy        = Cmpi->ops->destroy;
1121:   ptap->duplicate      = Cmpi->ops->duplicate;
1122:   Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1123:   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_MatMatMult;

1125:   /* attach the supporting struct to Cmpi for reuse */
1126:   c       = (Mat_MPIAIJ*)Cmpi->data;
1127:   c->ap = ptap;
1128:   *C = Cmpi;

1130:   /* set MatInfo */
1131:   afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
1132:   if (afill < 1.0) afill = 1.0;
1133:   Cmpi->info.mallocs           = nspacedouble;
1134:   Cmpi->info.fill_ratio_given  = fill;
1135:   Cmpi->info.fill_ratio_needed = afill;

1137: #if defined(PETSC_USE_INFO)
1138:   if (api[am]) {
1139:     PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
1140:     PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
1141:   } else {
1142:     PetscInfo(Cmpi,"Empty matrix product\n");
1143:   }
1144: #endif

1146:   MatDestroy(&aopoth);
1147:   MatDestroy(&adpd);
1148:   PetscFree(j_temp);
1149:   PetscFree(adpoj);
1150:   PetscFree(adpoi);
1151:   return(0);
1152: }


1155: /*-------------------------------------------------------------------------*/
1156: PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
1157: {
1159:   const char     *algTypes[3] = {"scalable","nonscalable","matmatmult"};
1160:   PetscInt       aN=A->cmap->N,alg=1; /* set default algorithm */
1161:   PetscBool      flg;

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

1169:     PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);
1170:     switch (alg) {
1171:     case 1:
1172:       if (!flg && aN > 100000) { /* may switch to scalable algorithm as default */
1173:         MatInfo     Ainfo,Pinfo;
1174:         PetscInt    nz_local;
1175:         PetscBool   alg_scalable_loc=PETSC_FALSE,alg_scalable;
1176:         MPI_Comm    comm;

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

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

1186:         if (alg_scalable) {
1187:           alg  = 0; /* scalable algorithm would slower than nonscalable algorithm */
1188:           MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);
1189:           break;
1190:         }
1191:       }
1192:       MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);
1193:       break;
1194:     case 2:
1195:     {
1196:       Mat         Pt;
1197:       Mat_APMPI   *ptap;
1198:       Mat_MPIAIJ  *c;
1199:       MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);
1200:       MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);
1201:       c        = (Mat_MPIAIJ*)(*C)->data;
1202:       ptap     = c->ap;
1203:       if (ptap) {
1204:        ptap->Pt = Pt;
1205:        (*C)->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;
1206:       }
1207:       (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
1208:       return(0);
1209:     }
1210:       break;
1211:     default: /* scalable algorithm */
1212:       MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);
1213:       break;
1214:     }
1215:     PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);

1217:     {
1218:       Mat_MPIAIJ *c  = (Mat_MPIAIJ*)(*C)->data;
1219:       Mat_APMPI  *ap = c->ap;
1220:       PetscOptionsBegin(PetscObjectComm((PetscObject)(*C)),((PetscObject)(*C))->prefix,"MatFreeIntermediateDataStructures","Mat");
1221:       ap->freestruct = PETSC_FALSE;
1222:       PetscOptionsBool("-mat_freeintermediatedatastructures","Free intermediate data structures", "MatFreeIntermediateDataStructures",ap->freestruct,&ap->freestruct, NULL);
1223:       PetscOptionsEnd();
1224:     }
1225:   }

1227:   PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);
1228:   (*(*C)->ops->mattransposemultnumeric)(P,A,*C);
1229:   PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);
1230:   return(0);
1231: }

1233: /* This routine only works when scall=MAT_REUSE_MATRIX! */
1234: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
1235: {
1237:   Mat_MPIAIJ     *c=(Mat_MPIAIJ*)C->data;
1238:   Mat_APMPI      *ptap= c->ap;
1239:   Mat            Pt;

1242:   if (!ptap) {
1243:     MPI_Comm comm;
1244:     PetscObjectGetComm((PetscObject)C,&comm);
1245:     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
1246:   }

1248:   Pt=ptap->Pt;
1249:   MatTranspose(P,MAT_REUSE_MATRIX,&Pt);
1250:   MatMatMultNumeric(Pt,A,C);

1252:   /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */
1253:   if (ptap->freestruct) {
1254:     MatFreeIntermediateDataStructures(C);
1255:   }
1256:   return(0);
1257: }

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

1288:   PetscObjectGetComm((PetscObject)A,&comm);
1289:   MPI_Comm_size(comm,&size);
1290:   MPI_Comm_rank(comm,&rank);

1292:   /* create symbolic parallel matrix Cmpi */
1293:   MatCreate(comm,&Cmpi);
1294:   MatGetType(A,&mtype);
1295:   MatSetType(Cmpi,mtype);

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

1299:   /* create struct Mat_APMPI and attached it to C later */
1300:   PetscNew(&ptap);
1301:   ptap->reuse = MAT_INITIAL_MATRIX;

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

1308:   /* (1) compute symbolic A_loc */
1309:   /* ---------------------------*/
1310:   MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&ptap->A_loc);

1312:   /* (2-1) compute symbolic C_oth = Ro*A_loc  */
1313:   /* ------------------------------------ */
1314:   MatGetOptionsPrefix(A,&prefix);
1315:   MatSetOptionsPrefix(ptap->Ro,prefix);
1316:   MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");
1317:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->A_loc,fill,&ptap->C_oth);

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

1328:   /* determine the number of messages to send, their lengths */
1329:   PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);
1330:   PetscMemzero(len_s,size*sizeof(PetscMPIInt));
1331:   PetscMemzero(len_si,size*sizeof(PetscMPIInt));

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

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

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

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

1370:   /* (2-2) compute symbolic C_loc = Rd*A_loc */
1371:   /* ---------------------------------------- */
1372:   MatSetOptionsPrefix(ptap->Rd,prefix);
1373:   MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");
1374:   MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->A_loc,fill,&ptap->C_loc);
1375:   c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;

1377:   /* receives coj are complete */
1378:   for (i=0; i<nrecv; i++) {
1379:     MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);
1380:   }
1381:   PetscFree(rwaits);
1382:   if (nsend) {MPI_Waitall(nsend,swaits,sstatus);}

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

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

1391:   for (k=0; k<nrecv; k++) {/* k-th received message */
1392:     Jptr = buf_rj[k];
1393:     for (j=0; j<len_r[k]; j++) {
1394:       PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);
1395:     }
1396:   }
1397:   PetscTableGetCount(ta,&Crmax);
1398:   PetscTableDestroy(&ta);

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

1435:   PetscFree4(len_s,len_si,sstatus,owners_co);
1436:   PetscFree(len_ri);
1437:   PetscFree(swaits);
1438:   PetscFree(buf_s);

1440:   /* (5) compute the local portion of Cmpi      */
1441:   /* ------------------------------------------ */
1442:   /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
1443:   PetscFreeSpaceGet(Crmax,&free_space);
1444:   current_space = free_space;

1446:   PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);
1447:   for (k=0; k<nrecv; k++) {
1448:     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1449:     nrows       = *buf_ri_k[k];
1450:     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1451:     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1452:   }

1454:   MatPreallocateInitialize(comm,pn,an,dnz,onz);
1455:   PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);
1456:   for (i=0; i<pn; i++) {
1457:     /* add C_loc into Cmpi */
1458:     nzi  = c_loc->i[i+1] - c_loc->i[i];
1459:     Jptr = c_loc->j + c_loc->i[i];
1460:     PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);

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

1473:     /* copy data into free space, then initialize lnk */
1474:     PetscLLCondensedClean(aN,nzi,current_space->array,lnk,lnkbt);
1475:     MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);
1476:   }
1477:   PetscFree3(buf_ri_k,nextrow,nextci);
1478:   PetscLLDestroy(lnk,lnkbt);
1479:   PetscFreeSpaceDestroy(free_space);

1481:   /* local sizes and preallocation */
1482:   MatSetSizes(Cmpi,pn,an,PETSC_DETERMINE,PETSC_DETERMINE);
1483:   if (P->cmap->bs > 0) {PetscLayoutSetBlockSize(Cmpi->rmap,P->cmap->bs);}
1484:   if (A->cmap->bs > 0) {PetscLayoutSetBlockSize(Cmpi->cmap,A->cmap->bs);}
1485:   MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);
1486:   MatPreallocateFinalize(dnz,onz);

1488:   /* members in merge */
1489:   PetscFree(id_r);
1490:   PetscFree(len_r);
1491:   PetscFree(buf_ri[0]);
1492:   PetscFree(buf_ri);
1493:   PetscFree(buf_rj[0]);
1494:   PetscFree(buf_rj);
1495:   PetscLayoutDestroy(&rowmap);

1497:   /* attach the supporting struct to Cmpi for reuse */
1498:   c = (Mat_MPIAIJ*)Cmpi->data;
1499:   c->ap         = ptap;
1500:   ptap->destroy = Cmpi->ops->destroy;

1502:   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
1503:   Cmpi->assembled        = PETSC_FALSE;
1504:   Cmpi->ops->destroy     = MatDestroy_MPIAIJ_PtAP;
1505:   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;

1507:   *C                     = Cmpi;
1508:   return(0);
1509: }

1511: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
1512: {
1513:   PetscErrorCode    ierr;
1514:   Mat_MPIAIJ        *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1515:   Mat_SeqAIJ        *c_seq;
1516:   Mat_APMPI         *ptap = c->ap;
1517:   Mat               A_loc,C_loc,C_oth;
1518:   PetscInt          i,rstart,rend,cm,ncols,row;
1519:   const PetscInt    *cols;
1520:   const PetscScalar *vals;

1523:   if (!ptap) {
1524:     MPI_Comm comm;
1525:     PetscObjectGetComm((PetscObject)C,&comm);
1526:     SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'");
1527:   }

1529:   MatZeroEntries(C);

1531:   if (ptap->reuse == MAT_REUSE_MATRIX) {
1532:     /* These matrices are obtained in MatTransposeMatMultSymbolic() */
1533:     /* 1) get R = Pd^T, Ro = Po^T */
1534:     /*----------------------------*/
1535:     MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);
1536:     MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);

1538:     /* 2) compute numeric A_loc */
1539:     /*--------------------------*/
1540:     MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&ptap->A_loc);
1541:   }

1543:   /* 3) C_loc = Rd*A_loc, C_oth = Ro*A_loc */
1544:   A_loc = ptap->A_loc;
1545:   ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,A_loc,ptap->C_loc);
1546:   ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,A_loc,ptap->C_oth);
1547:   C_loc = ptap->C_loc;
1548:   C_oth = ptap->C_oth;

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

1553:   /* C_loc -> C */
1554:   cm    = C_loc->rmap->N;
1555:   c_seq = (Mat_SeqAIJ*)C_loc->data;
1556:   cols = c_seq->j;
1557:   vals = c_seq->a;
1558:   for (i=0; i<cm; i++) {
1559:     ncols = c_seq->i[i+1] - c_seq->i[i];
1560:     row = rstart + i;
1561:     MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);
1562:     cols += ncols; vals += ncols;
1563:   }

1565:   /* Co -> C, off-processor part */
1566:   cm    = C_oth->rmap->N;
1567:   c_seq = (Mat_SeqAIJ*)C_oth->data;
1568:   cols  = c_seq->j;
1569:   vals  = c_seq->a;
1570:   for (i=0; i<cm; i++) {
1571:     ncols = c_seq->i[i+1] - c_seq->i[i];
1572:     row = p->garray[i];
1573:     MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);
1574:     cols += ncols; vals += ncols;
1575:   }
1576:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1577:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

1579:   ptap->reuse = MAT_REUSE_MATRIX;

1581:   /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */
1582:   if (ptap->freestruct) {
1583:     MatFreeIntermediateDataStructures(C);
1584:   }
1585:   return(0);
1586: }

1588: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1589: {
1590:   PetscErrorCode      ierr;
1591:   Mat_Merge_SeqsToMPI *merge;
1592:   Mat_MPIAIJ          *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1593:   Mat_SeqAIJ          *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1594:   Mat_APMPI           *ptap;
1595:   PetscInt            *adj;
1596:   PetscInt            i,j,k,anz,pnz,row,*cj,nexta;
1597:   MatScalar           *ada,*ca,valtmp;
1598:   PetscInt            am  =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1599:   MPI_Comm            comm;
1600:   PetscMPIInt         size,rank,taga,*len_s;
1601:   PetscInt            *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1602:   PetscInt            **buf_ri,**buf_rj;
1603:   PetscInt            cnz=0,*bj_i,*bi,*bj,bnz,nextcj;  /* bi,bj,ba: local array of C(mpi mat) */
1604:   MPI_Request         *s_waits,*r_waits;
1605:   MPI_Status          *status;
1606:   MatScalar           **abuf_r,*ba_i,*pA,*coa,*ba;
1607:   PetscInt            *ai,*aj,*coi,*coj,*poJ,*pdJ;
1608:   Mat                 A_loc;
1609:   Mat_SeqAIJ          *a_loc;

1612:   PetscObjectGetComm((PetscObject)C,&comm);
1613:   MPI_Comm_size(comm,&size);
1614:   MPI_Comm_rank(comm,&rank);

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

1620:   /* 2) compute numeric C_seq = P_loc^T*A_loc */
1621:   /*------------------------------------------*/
1622:   /* get data from symbolic products */
1623:   coi    = merge->coi; coj = merge->coj;
1624:   PetscCalloc1(coi[pon]+1,&coa);
1625:   bi     = merge->bi; bj = merge->bj;
1626:   owners = merge->rowmap->range;
1627:   PetscCalloc1(bi[cm]+1,&ba);

1629:   /* get A_loc by taking all local rows of A */
1630:   A_loc = ptap->A_loc;
1631:   MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);
1632:   a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1633:   ai    = a_loc->i;
1634:   aj    = a_loc->j;

1636:   for (i=0; i<am; i++) {
1637:     anz = ai[i+1] - ai[i];
1638:     adj = aj + ai[i];
1639:     ada = a_loc->a + ai[i];

1641:     /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1642:     /*-------------------------------------------------------------*/
1643:     /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1644:     pnz = po->i[i+1] - po->i[i];
1645:     poJ = po->j + po->i[i];
1646:     pA  = po->a + po->i[i];
1647:     for (j=0; j<pnz; j++) {
1648:       row = poJ[j];
1649:       cj  = coj + coi[row];
1650:       ca  = coa + coi[row];
1651:       /* perform sparse axpy */
1652:       nexta  = 0;
1653:       valtmp = pA[j];
1654:       for (k=0; nexta<anz; k++) {
1655:         if (cj[k] == adj[nexta]) {
1656:           ca[k] += valtmp*ada[nexta];
1657:           nexta++;
1658:         }
1659:       }
1660:       PetscLogFlops(2.0*anz);
1661:     }

1663:     /* put the value into Cd (diagonal part) */
1664:     pnz = pd->i[i+1] - pd->i[i];
1665:     pdJ = pd->j + pd->i[i];
1666:     pA  = pd->a + pd->i[i];
1667:     for (j=0; j<pnz; j++) {
1668:       row = pdJ[j];
1669:       cj  = bj + bi[row];
1670:       ca  = ba + bi[row];
1671:       /* perform sparse axpy */
1672:       nexta  = 0;
1673:       valtmp = pA[j];
1674:       for (k=0; nexta<anz; k++) {
1675:         if (cj[k] == adj[nexta]) {
1676:           ca[k] += valtmp*ada[nexta];
1677:           nexta++;
1678:         }
1679:       }
1680:       PetscLogFlops(2.0*anz);
1681:     }
1682:   }

1684:   /* 3) send and recv matrix values coa */
1685:   /*------------------------------------*/
1686:   buf_ri = merge->buf_ri;
1687:   buf_rj = merge->buf_rj;
1688:   len_s  = merge->len_s;
1689:   PetscCommGetNewTag(comm,&taga);
1690:   PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);

1692:   PetscMalloc2(merge->nsend+1,&s_waits,size,&status);
1693:   for (proc=0,k=0; proc<size; proc++) {
1694:     if (!len_s[proc]) continue;
1695:     i    = merge->owners_co[proc];
1696:     MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);
1697:     k++;
1698:   }
1699:   if (merge->nrecv) {MPI_Waitall(merge->nrecv,r_waits,status);}
1700:   if (merge->nsend) {MPI_Waitall(merge->nsend,s_waits,status);}

1702:   PetscFree2(s_waits,status);
1703:   PetscFree(r_waits);
1704:   PetscFree(coa);

1706:   /* 4) insert local Cseq and received values into Cmpi */
1707:   /*----------------------------------------------------*/
1708:   PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);
1709:   for (k=0; k<merge->nrecv; k++) {
1710:     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1711:     nrows       = *(buf_ri_k[k]);
1712:     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
1713:     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure  */
1714:   }

1716:   for (i=0; i<cm; i++) {
1717:     row  = owners[rank] + i; /* global row index of C_seq */
1718:     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
1719:     ba_i = ba + bi[i];
1720:     bnz  = bi[i+1] - bi[i];
1721:     /* add received vals into ba */
1722:     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1723:       /* i-th row */
1724:       if (i == *nextrow[k]) {
1725:         cnz    = *(nextci[k]+1) - *nextci[k];
1726:         cj     = buf_rj[k] + *(nextci[k]);
1727:         ca     = abuf_r[k] + *(nextci[k]);
1728:         nextcj = 0;
1729:         for (j=0; nextcj<cnz; j++) {
1730:           if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1731:             ba_i[j] += ca[nextcj++];
1732:           }
1733:         }
1734:         nextrow[k]++; nextci[k]++;
1735:         PetscLogFlops(2.0*cnz);
1736:       }
1737:     }
1738:     MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);
1739:   }
1740:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1741:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

1743:   PetscFree(ba);
1744:   PetscFree(abuf_r[0]);
1745:   PetscFree(abuf_r);
1746:   PetscFree3(buf_ri_k,nextrow,nextci);

1748:   if (ptap->freestruct) {
1749:     MatFreeIntermediateDataStructures(C);
1750:   }
1751:   return(0);
1752: }

1754: PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1755: {
1756:   PetscErrorCode      ierr;
1757:   Mat                 Cmpi,A_loc,POt,PDt;
1758:   Mat_APMPI           *ptap;
1759:   PetscFreeSpaceList  free_space=NULL,current_space=NULL;
1760:   Mat_MPIAIJ          *p=(Mat_MPIAIJ*)P->data,*a=(Mat_MPIAIJ*)A->data,*c;
1761:   PetscInt            *pdti,*pdtj,*poti,*potj,*ptJ;
1762:   PetscInt            nnz;
1763:   PetscInt            *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1764:   PetscInt            am  =A->rmap->n,pn=P->cmap->n;
1765:   MPI_Comm            comm;
1766:   PetscMPIInt         size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1767:   PetscInt            **buf_rj,**buf_ri,**buf_ri_k;
1768:   PetscInt            len,proc,*dnz,*onz,*owners;
1769:   PetscInt            nzi,*bi,*bj;
1770:   PetscInt            nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1771:   MPI_Request         *swaits,*rwaits;
1772:   MPI_Status          *sstatus,rstatus;
1773:   Mat_Merge_SeqsToMPI *merge;
1774:   PetscInt            *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1775:   PetscReal           afill  =1.0,afill_tmp;
1776:   PetscInt            rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax;
1777:   PetscScalar         *vals;
1778:   Mat_SeqAIJ          *a_loc,*pdt,*pot;
1779:   PetscTable          ta;
1780:   MatType             mtype;

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

1787:   MPI_Comm_size(comm,&size);
1788:   MPI_Comm_rank(comm,&rank);

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

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

1796:   ptap->A_loc = A_loc;
1797:   a_loc       = (Mat_SeqAIJ*)(A_loc)->data;
1798:   ai          = a_loc->i;
1799:   aj          = a_loc->j;

1801:   /* determine symbolic Co=(p->B)^T*A - send to others */
1802:   /*----------------------------------------------------*/
1803:   MatTransposeSymbolic_SeqAIJ(p->A,&PDt);
1804:   pdt  = (Mat_SeqAIJ*)PDt->data;
1805:   pdti = pdt->i; pdtj = pdt->j;

1807:   MatTransposeSymbolic_SeqAIJ(p->B,&POt);
1808:   pot  = (Mat_SeqAIJ*)POt->data;
1809:   poti = pot->i; potj = pot->j;

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

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

1822:   /* create and initialize a linked list */
1823:   PetscTableCreate(A->cmap->n + a->B->cmap->N,aN,&ta);
1824:   MatRowMergeMax_SeqAIJ(a_loc,am,ta);
1825:   PetscTableGetCount(ta,&Armax);

1827:   PetscLLCondensedCreate_Scalable(Armax,&lnk);

1829:   for (i=0; i<pon; i++) {
1830:     pnz = poti[i+1] - poti[i];
1831:     ptJ = potj + poti[i];
1832:     for (j=0; j<pnz; j++) {
1833:       row  = ptJ[j]; /* row of A_loc == col of Pot */
1834:       anz  = ai[row+1] - ai[row];
1835:       Jptr = aj + ai[row];
1836:       /* add non-zero cols of AP into the sorted linked list lnk */
1837:       PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);
1838:     }
1839:     nnz = lnk[0];

1841:     /* If free space is not available, double the total space in the list */
1842:     if (current_space->local_remaining<nnz) {
1843:       PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);
1844:       nspacedouble++;
1845:     }

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

1850:     current_space->array           += nnz;
1851:     current_space->local_used      += nnz;
1852:     current_space->local_remaining -= nnz;

1854:     coi[i+1] = coi[i] + nnz;
1855:   }

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

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

1864:   /* send j-array (coj) of Co to other processors */
1865:   /*----------------------------------------------*/
1866:   /* determine row ownership */
1867:   PetscNew(&merge);
1868:   PetscLayoutCreate(comm,&merge->rowmap);

1870:   merge->rowmap->n  = pn;
1871:   merge->rowmap->bs = 1;

1873:   PetscLayoutSetUp(merge->rowmap);
1874:   owners = merge->rowmap->range;

1876:   /* determine the number of messages to send, their lengths */
1877:   PetscCalloc1(size,&len_si);
1878:   PetscMalloc1(size,&merge->len_s);

1880:   len_s        = merge->len_s;
1881:   merge->nsend = 0;

1883:   PetscMalloc1(size+2,&owners_co);
1884:   PetscMemzero(len_s,size*sizeof(PetscMPIInt));

1886:   proc = 0;
1887:   for (i=0; i<pon; i++) {
1888:     while (prmap[i] >= owners[proc+1]) proc++;
1889:     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
1890:     len_s[proc] += coi[i+1] - coi[i];
1891:   }

1893:   len          = 0; /* max length of buf_si[] */
1894:   owners_co[0] = 0;
1895:   for (proc=0; proc<size; proc++) {
1896:     owners_co[proc+1] = owners_co[proc] + len_si[proc];
1897:     if (len_si[proc]) {
1898:       merge->nsend++;
1899:       len_si[proc] = 2*(len_si[proc] + 1);
1900:       len         += len_si[proc];
1901:     }
1902:   }

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

1908:   /* post the Irecv and Isend of coj */
1909:   PetscCommGetNewTag(comm,&tagj);
1910:   PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);
1911:   PetscMalloc1(merge->nsend+1,&swaits);
1912:   for (proc=0, k=0; proc<size; proc++) {
1913:     if (!len_s[proc]) continue;
1914:     i    = owners_co[proc];
1915:     MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);
1916:     k++;
1917:   }

1919:   /* receives and sends of coj are complete */
1920:   PetscMalloc1(size,&sstatus);
1921:   for (i=0; i<merge->nrecv; i++) {
1922:     PetscMPIInt icompleted;
1923:     MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);
1924:   }
1925:   PetscFree(rwaits);
1926:   if (merge->nsend) {MPI_Waitall(merge->nsend,swaits,sstatus);}

1928:   /* add received column indices into table to update Armax */
1929:   /* Armax can be as large as aN if a P[row,:] is dense, see src/ksp/ksp/examples/tutorials/ex56.c! */
1930:   for (k=0; k<merge->nrecv; k++) {/* k-th received message */
1931:     Jptr = buf_rj[k];
1932:     for (j=0; j<merge->len_r[k]; j++) {
1933:       PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);
1934:     }
1935:   }
1936:   PetscTableGetCount(ta,&Armax);
1937:   /* 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); */

1939:   /* send and recv coi */
1940:   /*-------------------*/
1941:   PetscCommGetNewTag(comm,&tagi);
1942:   PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);
1943:   PetscMalloc1(len+1,&buf_s);
1944:   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
1945:   for (proc=0,k=0; proc<size; proc++) {
1946:     if (!len_s[proc]) continue;
1947:     /* form outgoing message for i-structure:
1948:          buf_si[0]:                 nrows to be sent
1949:                [1:nrows]:           row index (global)
1950:                [nrows+1:2*nrows+1]: i-structure index
1951:     */
1952:     /*-------------------------------------------*/
1953:     nrows       = len_si[proc]/2 - 1;
1954:     buf_si_i    = buf_si + nrows+1;
1955:     buf_si[0]   = nrows;
1956:     buf_si_i[0] = 0;
1957:     nrows       = 0;
1958:     for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1959:       nzi               = coi[i+1] - coi[i];
1960:       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi;  /* i-structure */
1961:       buf_si[nrows+1]   = prmap[i] -owners[proc]; /* local row index */
1962:       nrows++;
1963:     }
1964:     MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);
1965:     k++;
1966:     buf_si += len_si[proc];
1967:   }
1968:   i = merge->nrecv;
1969:   while (i--) {
1970:     PetscMPIInt icompleted;
1971:     MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);
1972:   }
1973:   PetscFree(rwaits);
1974:   if (merge->nsend) {MPI_Waitall(merge->nsend,swaits,sstatus);}
1975:   PetscFree(len_si);
1976:   PetscFree(len_ri);
1977:   PetscFree(swaits);
1978:   PetscFree(sstatus);
1979:   PetscFree(buf_s);

1981:   /* compute the local portion of C (mpi mat) */
1982:   /*------------------------------------------*/
1983:   /* allocate bi array and free space for accumulating nonzero column info */
1984:   PetscMalloc1(pn+1,&bi);
1985:   bi[0] = 0;

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

1992:   PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);
1993:   for (k=0; k<merge->nrecv; k++) {
1994:     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1995:     nrows       = *buf_ri_k[k];
1996:     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
1997:     nextci[k]   = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure  */
1998:   }

2000:   PetscLLCondensedCreate_Scalable(Armax,&lnk);
2001:   MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);
2002:   rmax = 0;
2003:   for (i=0; i<pn; i++) {
2004:     /* add pdt[i,:]*AP into lnk */
2005:     pnz = pdti[i+1] - pdti[i];
2006:     ptJ = pdtj + pdti[i];
2007:     for (j=0; j<pnz; j++) {
2008:       row  = ptJ[j];  /* row of AP == col of Pt */
2009:       anz  = ai[row+1] - ai[row];
2010:       Jptr = aj + ai[row];
2011:       /* add non-zero cols of AP into the sorted linked list lnk */
2012:       PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);
2013:     }

2015:     /* add received col data into lnk */
2016:     for (k=0; k<merge->nrecv; k++) { /* k-th received message */
2017:       if (i == *nextrow[k]) { /* i-th row */
2018:         nzi  = *(nextci[k]+1) - *nextci[k];
2019:         Jptr = buf_rj[k] + *nextci[k];
2020:         PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);
2021:         nextrow[k]++; nextci[k]++;
2022:       }
2023:     }
2024:     nnz = lnk[0];

2026:     /* if free space is not available, make more free space */
2027:     if (current_space->local_remaining<nnz) {
2028:       PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),&current_space);
2029:       nspacedouble++;
2030:     }
2031:     /* copy data into free space, then initialize lnk */
2032:     PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);
2033:     MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);

2035:     current_space->array           += nnz;
2036:     current_space->local_used      += nnz;
2037:     current_space->local_remaining -= nnz;

2039:     bi[i+1] = bi[i] + nnz;
2040:     if (nnz > rmax) rmax = nnz;
2041:   }
2042:   PetscFree3(buf_ri_k,nextrow,nextci);

2044:   PetscMalloc1(bi[pn]+1,&bj);
2045:   PetscFreeSpaceContiguous(&free_space,bj);
2046:   afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
2047:   if (afill_tmp > afill) afill = afill_tmp;
2048:   PetscLLCondensedDestroy_Scalable(lnk);
2049:   PetscTableDestroy(&ta);

2051:   MatDestroy(&POt);
2052:   MatDestroy(&PDt);

2054:   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part   */
2055:   /*----------------------------------------------------------------------------------*/
2056:   PetscCalloc1(rmax+1,&vals);

2058:   MatCreate(comm,&Cmpi);
2059:   MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);
2060:   MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));
2061:   MatGetType(A,&mtype);
2062:   MatSetType(Cmpi,mtype);
2063:   MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);
2064:   MatPreallocateFinalize(dnz,onz);
2065:   MatSetBlockSize(Cmpi,1);
2066:   for (i=0; i<pn; i++) {
2067:     row  = i + rstart;
2068:     nnz  = bi[i+1] - bi[i];
2069:     Jptr = bj + bi[i];
2070:     MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);
2071:   }
2072:   MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);
2073:   MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);
2074:   PetscFree(vals);

2076:   merge->bi        = bi;
2077:   merge->bj        = bj;
2078:   merge->coi       = coi;
2079:   merge->coj       = coj;
2080:   merge->buf_ri    = buf_ri;
2081:   merge->buf_rj    = buf_rj;
2082:   merge->owners_co = owners_co;

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

2087:   c->ap       = ptap;
2088:   ptap->api   = NULL;
2089:   ptap->apj   = NULL;
2090:   ptap->merge = merge;
2091:   ptap->apa   = NULL;
2092:   ptap->destroy   = Cmpi->ops->destroy;
2093:   ptap->duplicate = Cmpi->ops->duplicate;

2095:   Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
2096:   Cmpi->ops->destroy                 = MatDestroy_MPIAIJ_PtAP;
2097:   Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP;

2099:   *C = Cmpi;
2100: #if defined(PETSC_USE_INFO)
2101:   if (bi[pn] != 0) {
2102:     PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
2103:     PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);
2104:   } else {
2105:     PetscInfo(Cmpi,"Empty matrix product\n");
2106:   }
2107: #endif
2108:   return(0);
2109: }