Actual source code: mpimatmatmult.c
petsc-3.10.5 2019-03-28
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>
13: #if defined(PETSC_HAVE_HYPRE)
14: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*);
15: #endif
17: PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIAIJ(Mat A,Mat B,MatReuse scall,PetscReal fill, Mat *C)
18: {
20: #if defined(PETSC_HAVE_HYPRE)
21: const char *algTypes[4] = {"scalable","nonscalable","seqmpi","hypre"};
22: PetscInt nalg = 4;
23: #else
24: const char *algTypes[3] = {"scalable","nonscalable","seqmpi"};
25: PetscInt nalg = 3;
26: #endif
27: PetscInt alg = 1; /* set nonscalable algorithm as default */
28: MPI_Comm comm;
29: PetscBool flg;
32: if (scall == MAT_INITIAL_MATRIX) {
33: PetscObjectGetComm((PetscObject)A,&comm);
34: 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);
36: PetscObjectOptionsBegin((PetscObject)A);
37: PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */
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);
77: }
78: PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
79: (*(*C)->ops->matmultnumeric)(A,B,*C);
80: PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
81: return(0);
82: }
84: PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(Mat A)
85: {
87: Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data;
88: Mat_PtAPMPI *ptap = a->ptap;
91: PetscFree2(ptap->startsj_s,ptap->startsj_r);
92: PetscFree(ptap->bufa);
93: MatDestroy(&ptap->P_loc);
94: MatDestroy(&ptap->P_oth);
95: MatDestroy(&ptap->Pt);
96: PetscFree(ptap->api);
97: PetscFree(ptap->apj);
98: PetscFree(ptap->apa);
99: ptap->destroy(A);
100: PetscFree(ptap);
101: return(0);
102: }
104: PetscErrorCode MatDuplicate_MPIAIJ_MatMatMult(Mat A, MatDuplicateOption op, Mat *M)
105: {
107: Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data;
108: Mat_PtAPMPI *ptap = a->ptap;
111: (*ptap->duplicate)(A,op,M);
113: (*M)->ops->destroy = ptap->destroy; /* = MatDestroy_MPIAIJ, *M doesn't duplicate A's special structure! */
114: (*M)->ops->duplicate = ptap->duplicate; /* = MatDuplicate_MPIAIJ */
115: return(0);
116: }
118: PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C)
119: {
121: Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
122: Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
123: Mat_SeqAIJ *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
124: PetscScalar *cda=cd->a,*coa=co->a;
125: Mat_SeqAIJ *p_loc,*p_oth;
126: PetscScalar *apa,*ca;
127: PetscInt cm =C->rmap->n;
128: Mat_PtAPMPI *ptap=c->ptap;
129: PetscInt *api,*apj,*apJ,i,k;
130: PetscInt cstart=C->cmap->rstart;
131: PetscInt cdnz,conz,k0,k1;
132: MPI_Comm comm;
133: PetscMPIInt size;
136: PetscObjectGetComm((PetscObject)A,&comm);
137: MPI_Comm_size(comm,&size);
139: /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */
140: /*-----------------------------------------------------*/
141: /* update numerical values of P_oth and P_loc */
142: MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);
143: MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);
145: /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
146: /*----------------------------------------------------------*/
147: /* get data from symbolic products */
148: p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
149: p_oth = NULL;
150: if (size >1) {
151: p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
152: }
154: /* get apa for storing dense row A[i,:]*P */
155: apa = ptap->apa;
157: api = ptap->api;
158: apj = ptap->apj;
159: for (i=0; i<cm; i++) {
160: /* compute apa = A[i,:]*P */
161: AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
163: /* set values in C */
164: apJ = apj + api[i];
165: cdnz = cd->i[i+1] - cd->i[i];
166: conz = co->i[i+1] - co->i[i];
168: /* 1st off-diagoanl part of C */
169: ca = coa + co->i[i];
170: k = 0;
171: for (k0=0; k0<conz; k0++) {
172: if (apJ[k] >= cstart) break;
173: ca[k0] = apa[apJ[k]];
174: apa[apJ[k++]] = 0.0;
175: }
177: /* diagonal part of C */
178: ca = cda + cd->i[i];
179: for (k1=0; k1<cdnz; k1++) {
180: ca[k1] = apa[apJ[k]];
181: apa[apJ[k++]] = 0.0;
182: }
184: /* 2nd off-diagoanl part of C */
185: ca = coa + co->i[i];
186: for (; k0<conz; k0++) {
187: ca[k0] = apa[apJ[k]];
188: apa[apJ[k++]] = 0.0;
189: }
190: }
191: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
192: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
193: return(0);
194: }
196: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat *C)
197: {
198: PetscErrorCode ierr;
199: MPI_Comm comm;
200: PetscMPIInt size;
201: Mat Cmpi;
202: Mat_PtAPMPI *ptap;
203: PetscFreeSpaceList free_space=NULL,current_space=NULL;
204: Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*c;
205: Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
206: PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
207: PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
208: PetscInt *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
209: PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
210: PetscBT lnkbt;
211: PetscScalar *apa;
212: PetscReal afill;
215: PetscObjectGetComm((PetscObject)A,&comm);
216: MPI_Comm_size(comm,&size);
218: /* create struct Mat_PtAPMPI and attached it to C later */
219: PetscNew(&ptap);
221: /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
222: MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);
224: /* get P_loc by taking all local rows of P */
225: MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);
227: p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
228: pi_loc = p_loc->i; pj_loc = p_loc->j;
229: if (size > 1) {
230: p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
231: pi_oth = p_oth->i; pj_oth = p_oth->j;
232: } else {
233: p_oth = NULL;
234: pi_oth = NULL; pj_oth = NULL;
235: }
237: /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
238: /*-------------------------------------------------------------------*/
239: PetscMalloc1(am+2,&api);
240: ptap->api = api;
241: api[0] = 0;
243: /* create and initialize a linked list */
244: PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);
246: /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
247: PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);
248: current_space = free_space;
250: MatPreallocateInitialize(comm,am,pn,dnz,onz);
251: for (i=0; i<am; i++) {
252: /* diagonal portion of A */
253: nzi = adi[i+1] - adi[i];
254: for (j=0; j<nzi; j++) {
255: row = *adj++;
256: pnz = pi_loc[row+1] - pi_loc[row];
257: Jptr = pj_loc + pi_loc[row];
258: /* add non-zero cols of P into the sorted linked list lnk */
259: PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);
260: }
261: /* off-diagonal portion of A */
262: nzi = aoi[i+1] - aoi[i];
263: for (j=0; j<nzi; j++) {
264: row = *aoj++;
265: pnz = pi_oth[row+1] - pi_oth[row];
266: Jptr = pj_oth + pi_oth[row];
267: PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);
268: }
270: apnz = lnk[0];
271: api[i+1] = api[i] + apnz;
273: /* if free space is not available, double the total space in the list */
274: if (current_space->local_remaining<apnz) {
275: PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);
276: nspacedouble++;
277: }
279: /* Copy data into free space, then initialize lnk */
280: PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);
281: MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);
283: current_space->array += apnz;
284: current_space->local_used += apnz;
285: current_space->local_remaining -= apnz;
286: }
288: /* Allocate space for apj, initialize apj, and */
289: /* destroy list of free space and other temporary array(s) */
290: PetscMalloc1(api[am]+1,&ptap->apj);
291: apj = ptap->apj;
292: PetscFreeSpaceContiguous(&free_space,ptap->apj);
293: PetscLLDestroy(lnk,lnkbt);
295: /* malloc apa to store dense row A[i,:]*P */
296: PetscCalloc1(pN,&apa);
298: ptap->apa = 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: MatSetType(Cmpi,MATMPIAIJ);
307: MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);
309: MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);
310: MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);
311: MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);
312: MatPreallocateFinalize(dnz,onz);
314: ptap->destroy = Cmpi->ops->destroy;
315: ptap->duplicate = Cmpi->ops->duplicate;
316: Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
317: Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult;
318: Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
320: /* attach the supporting struct to Cmpi for reuse */
321: c = (Mat_MPIAIJ*)Cmpi->data;
322: c->ptap = ptap;
324: *C = Cmpi;
326: /* set MatInfo */
327: afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
328: if (afill < 1.0) afill = 1.0;
329: Cmpi->info.mallocs = nspacedouble;
330: Cmpi->info.fill_ratio_given = fill;
331: Cmpi->info.fill_ratio_needed = afill;
333: #if defined(PETSC_USE_INFO)
334: if (api[am]) {
335: PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
336: PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
337: } else {
338: PetscInfo(Cmpi,"Empty matrix product\n");
339: }
340: #endif
341: return(0);
342: }
344: PETSC_INTERN PetscErrorCode MatMatMult_MPIAIJ_MPIDense(Mat A,Mat B,MatReuse scall,PetscReal fill,Mat *C)
345: {
349: if (scall == MAT_INITIAL_MATRIX) {
350: PetscLogEventBegin(MAT_MatMultSymbolic,A,B,0,0);
351: MatMatMultSymbolic_MPIAIJ_MPIDense(A,B,fill,C);
352: PetscLogEventEnd(MAT_MatMultSymbolic,A,B,0,0);
353: }
354: PetscLogEventBegin(MAT_MatMultNumeric,A,B,0,0);
355: MatMatMultNumeric_MPIAIJ_MPIDense(A,B,*C);
356: PetscLogEventEnd(MAT_MatMultNumeric,A,B,0,0);
357: return(0);
358: }
360: typedef struct {
361: Mat workB;
362: PetscScalar *rvalues,*svalues;
363: MPI_Request *rwaits,*swaits;
364: } MPIAIJ_MPIDense;
366: PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
367: {
368: MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*) ctx;
369: PetscErrorCode ierr;
372: MatDestroy(&contents->workB);
373: PetscFree4(contents->rvalues,contents->svalues,contents->rwaits,contents->swaits);
374: PetscFree(contents);
375: return(0);
376: }
378: #include <petsc/private/vecscatterimpl.h>
379: /*
380: This is a "dummy function" that handles the case where matrix C was created as a dense matrix
381: directly by the user and passed to MatMatMult() with the MAT_REUSE_MATRIX option
383: It is the same as MatMatMultSymbolic_MPIAIJ_MPIDense() except does not create C
385: Developer Notes: This directly accesses information inside the VecScatter associated with the matrix-vector product
386: for this matrix. This is not desirable..
388: */
389: PetscErrorCode MatMatMultNumeric_MPIDense(Mat A,Mat B,Mat C)
390: {
391: PetscErrorCode ierr;
392: PetscBool flg;
393: Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data;
394: PetscInt nz = aij->B->cmap->n;
395: PetscContainer container;
396: MPIAIJ_MPIDense *contents;
397: VecScatter ctx = aij->Mvctx;
398: VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
399: VecScatter_MPI_General *to = (VecScatter_MPI_General*) ctx->todata;
402: PetscObjectTypeCompare((PetscObject)B,MATMPIDENSE,&flg);
403: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Second matrix must be mpidense");
405: /* Handle case where where user provided the final C matrix rather than calling MatMatMult() with MAT_INITIAL_MATRIX*/
406: PetscObjectTypeCompare((PetscObject)A,MATMPIAIJ,&flg);
407: if (!flg) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"First matrix must be MPIAIJ");
409: C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
411: PetscNew(&contents);
412: /* Create work matrix used to store off processor rows of B needed for local product */
413: MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);
414: /* Create work arrays needed */
415: PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
416: B->cmap->N*to->starts[to->n],&contents->svalues,
417: from->n,&contents->rwaits,
418: to->n,&contents->swaits);
420: PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);
421: PetscContainerSetPointer(container,contents);
422: PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);
423: PetscObjectCompose((PetscObject)C,"workB",(PetscObject)container);
424: PetscContainerDestroy(&container);
426: (*C->ops->matmultnumeric)(A,B,C);
427: return(0);
428: }
430: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat *C)
431: {
432: PetscErrorCode ierr;
433: Mat_MPIAIJ *aij = (Mat_MPIAIJ*) A->data;
434: PetscInt nz = aij->B->cmap->n;
435: PetscContainer container;
436: MPIAIJ_MPIDense *contents;
437: VecScatter ctx = aij->Mvctx;
438: VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
439: VecScatter_MPI_General *to = (VecScatter_MPI_General*) ctx->todata;
440: PetscInt m = A->rmap->n,n=B->cmap->n;
443: MatCreate(PetscObjectComm((PetscObject)B),C);
444: MatSetSizes(*C,m,n,A->rmap->N,B->cmap->N);
445: MatSetBlockSizesFromMats(*C,A,B);
446: MatSetType(*C,MATMPIDENSE);
447: MatMPIDenseSetPreallocation(*C,NULL);
448: MatAssemblyBegin(*C,MAT_FINAL_ASSEMBLY);
449: MatAssemblyEnd(*C,MAT_FINAL_ASSEMBLY);
451: (*C)->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
453: PetscNew(&contents);
454: /* Create work matrix used to store off processor rows of B needed for local product */
455: MatCreateSeqDense(PETSC_COMM_SELF,nz,B->cmap->N,NULL,&contents->workB);
456: /* Create work arrays needed */
457: PetscMalloc4(B->cmap->N*from->starts[from->n],&contents->rvalues,
458: B->cmap->N*to->starts[to->n],&contents->svalues,
459: from->n,&contents->rwaits,
460: to->n,&contents->swaits);
462: PetscContainerCreate(PetscObjectComm((PetscObject)A),&container);
463: PetscContainerSetPointer(container,contents);
464: PetscContainerSetUserDestroy(container,MatMPIAIJ_MPIDenseDestroy);
465: PetscObjectCompose((PetscObject)(*C),"workB",(PetscObject)container);
466: PetscContainerDestroy(&container);
467: return(0);
468: }
470: /*
471: Performs an efficient scatter on the rows of B needed by this process; this is
472: a modification of the VecScatterBegin_() routines.
473: */
474: PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,Mat C,Mat *outworkB)
475: {
476: Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data;
477: PetscErrorCode ierr;
478: PetscScalar *b,*w,*svalues,*rvalues;
479: VecScatter ctx = aij->Mvctx;
480: VecScatter_MPI_General *from = (VecScatter_MPI_General*) ctx->fromdata;
481: VecScatter_MPI_General *to = (VecScatter_MPI_General*) ctx->todata;
482: PetscInt i,j,k;
483: PetscInt *sindices,*sstarts,*rindices,*rstarts;
484: PetscMPIInt *sprocs,*rprocs,nrecvs;
485: MPI_Request *swaits,*rwaits;
486: MPI_Comm comm;
487: PetscMPIInt tag = ((PetscObject)ctx)->tag,ncols = B->cmap->N, nrows = aij->B->cmap->n,imdex,nrowsB = B->rmap->n;
488: MPI_Status status;
489: MPIAIJ_MPIDense *contents;
490: PetscContainer container;
491: Mat workB;
494: PetscObjectGetComm((PetscObject)A,&comm);
495: PetscObjectQuery((PetscObject)C,"workB",(PetscObject*)&container);
496: if (!container) SETERRQ(comm,PETSC_ERR_PLIB,"Container does not exist");
497: PetscContainerGetPointer(container,(void**)&contents);
499: workB = *outworkB = contents->workB;
500: 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);
501: sindices = to->indices;
502: sstarts = to->starts;
503: sprocs = to->procs;
504: swaits = contents->swaits;
505: svalues = contents->svalues;
507: rindices = from->indices;
508: rstarts = from->starts;
509: rprocs = from->procs;
510: rwaits = contents->rwaits;
511: rvalues = contents->rvalues;
513: MatDenseGetArray(B,&b);
514: MatDenseGetArray(workB,&w);
516: for (i=0; i<from->n; i++) {
517: MPI_Irecv(rvalues+ncols*rstarts[i],ncols*(rstarts[i+1]-rstarts[i]),MPIU_SCALAR,rprocs[i],tag,comm,rwaits+i);
518: }
520: for (i=0; i<to->n; i++) {
521: /* pack a message at a time */
522: for (j=0; j<sstarts[i+1]-sstarts[i]; j++) {
523: for (k=0; k<ncols; k++) {
524: svalues[ncols*(sstarts[i] + j) + k] = b[sindices[sstarts[i]+j] + nrowsB*k];
525: }
526: }
527: MPI_Isend(svalues+ncols*sstarts[i],ncols*(sstarts[i+1]-sstarts[i]),MPIU_SCALAR,sprocs[i],tag,comm,swaits+i);
528: }
530: nrecvs = from->n;
531: while (nrecvs) {
532: MPI_Waitany(from->n,rwaits,&imdex,&status);
533: nrecvs--;
534: /* unpack a message at a time */
535: for (j=0; j<rstarts[imdex+1]-rstarts[imdex]; j++) {
536: for (k=0; k<ncols; k++) {
537: w[rindices[rstarts[imdex]+j] + nrows*k] = rvalues[ncols*(rstarts[imdex] + j) + k];
538: }
539: }
540: }
541: if (to->n) {MPI_Waitall(to->n,swaits,to->sstatus);}
543: MatDenseRestoreArray(B,&b);
544: MatDenseRestoreArray(workB,&w);
545: MatAssemblyBegin(workB,MAT_FINAL_ASSEMBLY);
546: MatAssemblyEnd(workB,MAT_FINAL_ASSEMBLY);
547: return(0);
548: }
549: extern PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat);
551: PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
552: {
554: Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data;
555: Mat_MPIDense *bdense = (Mat_MPIDense*)B->data;
556: Mat_MPIDense *cdense = (Mat_MPIDense*)C->data;
557: Mat workB;
560: /* diagonal block of A times all local rows of B*/
561: MatMatMultNumeric_SeqAIJ_SeqDense(aij->A,bdense->A,cdense->A);
563: /* get off processor parts of B needed to complete the product */
564: MatMPIDenseScatter(A,B,C,&workB);
566: /* off-diagonal block of A times nonlocal rows of B */
567: MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A);
568: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
569: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
570: return(0);
571: }
573: PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
574: {
576: Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
577: Mat_SeqAIJ *ad = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
578: Mat_SeqAIJ *cd = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
579: PetscInt *adi = ad->i,*adj,*aoi=ao->i,*aoj;
580: PetscScalar *ada,*aoa,*cda=cd->a,*coa=co->a;
581: Mat_SeqAIJ *p_loc,*p_oth;
582: PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
583: PetscScalar *pa_loc,*pa_oth,*pa,valtmp,*ca;
584: PetscInt cm = C->rmap->n,anz,pnz;
585: Mat_PtAPMPI *ptap = c->ptap;
586: PetscScalar *apa_sparse = ptap->apa;
587: PetscInt *api,*apj,*apJ,i,j,k,row;
588: PetscInt cstart = C->cmap->rstart;
589: PetscInt cdnz,conz,k0,k1,nextp;
590: MPI_Comm comm;
591: PetscMPIInt size;
594: PetscObjectGetComm((PetscObject)A,&comm);
595: MPI_Comm_size(comm,&size);
597: /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */
598: /*-----------------------------------------------------*/
599: /* update numerical values of P_oth and P_loc */
600: MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);
601: MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);
603: /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
604: /*----------------------------------------------------------*/
605: /* get data from symbolic products */
606: p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
607: pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a;
608: if (size >1) {
609: p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
610: pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a;
611: } else {
612: p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL;
613: }
615: api = ptap->api;
616: apj = ptap->apj;
617: for (i=0; i<cm; i++) {
618: apJ = apj + api[i];
620: /* diagonal portion of A */
621: anz = adi[i+1] - adi[i];
622: adj = ad->j + adi[i];
623: ada = ad->a + adi[i];
624: for (j=0; j<anz; j++) {
625: row = adj[j];
626: pnz = pi_loc[row+1] - pi_loc[row];
627: pj = pj_loc + pi_loc[row];
628: pa = pa_loc + pi_loc[row];
629: /* perform sparse axpy */
630: valtmp = ada[j];
631: nextp = 0;
632: for (k=0; nextp<pnz; k++) {
633: if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
634: apa_sparse[k] += valtmp*pa[nextp++];
635: }
636: }
637: PetscLogFlops(2.0*pnz);
638: }
640: /* off-diagonal portion of A */
641: anz = aoi[i+1] - aoi[i];
642: aoj = ao->j + aoi[i];
643: aoa = ao->a + aoi[i];
644: for (j=0; j<anz; j++) {
645: row = aoj[j];
646: pnz = pi_oth[row+1] - pi_oth[row];
647: pj = pj_oth + pi_oth[row];
648: pa = pa_oth + pi_oth[row];
649: /* perform sparse axpy */
650: valtmp = aoa[j];
651: nextp = 0;
652: for (k=0; nextp<pnz; k++) {
653: if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
654: apa_sparse[k] += valtmp*pa[nextp++];
655: }
656: }
657: PetscLogFlops(2.0*pnz);
658: }
660: /* set values in C */
661: cdnz = cd->i[i+1] - cd->i[i];
662: conz = co->i[i+1] - co->i[i];
664: /* 1st off-diagoanl part of C */
665: ca = coa + co->i[i];
666: k = 0;
667: for (k0=0; k0<conz; k0++) {
668: if (apJ[k] >= cstart) break;
669: ca[k0] = apa_sparse[k];
670: apa_sparse[k] = 0.0;
671: k++;
672: }
674: /* diagonal part of C */
675: ca = cda + cd->i[i];
676: for (k1=0; k1<cdnz; k1++) {
677: ca[k1] = apa_sparse[k];
678: apa_sparse[k] = 0.0;
679: k++;
680: }
682: /* 2nd off-diagoanl part of C */
683: ca = coa + co->i[i];
684: for (; k0<conz; k0++) {
685: ca[k0] = apa_sparse[k];
686: apa_sparse[k] = 0.0;
687: k++;
688: }
689: }
690: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
691: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
692: return(0);
693: }
695: /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
696: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
697: {
698: PetscErrorCode ierr;
699: MPI_Comm comm;
700: PetscMPIInt size;
701: Mat Cmpi;
702: Mat_PtAPMPI *ptap;
703: PetscFreeSpaceList free_space = NULL,current_space=NULL;
704: Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data,*c;
705: Mat_SeqAIJ *ad = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
706: PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
707: PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
708: PetscInt i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
709: PetscInt am=A->rmap->n,pn=P->cmap->n,pm=P->rmap->n,lsize=pn+20;
710: PetscReal afill;
711: PetscScalar *apa;
714: PetscObjectGetComm((PetscObject)A,&comm);
715: MPI_Comm_size(comm,&size);
717: /* create struct Mat_PtAPMPI and attached it to C later */
718: PetscNew(&ptap);
720: /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
721: MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);
723: /* get P_loc by taking all local rows of P */
724: MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);
726: p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
727: pi_loc = p_loc->i; pj_loc = p_loc->j;
728: if (size > 1) {
729: p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
730: pi_oth = p_oth->i; pj_oth = p_oth->j;
731: } else {
732: p_oth = NULL;
733: pi_oth = NULL; pj_oth = NULL;
734: }
736: /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
737: /*-------------------------------------------------------------------*/
738: PetscMalloc1(am+2,&api);
739: ptap->api = api;
740: api[0] = 0;
742: PetscLLCondensedCreate_Scalable(lsize,&lnk);
744: /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
745: PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);
746: current_space = free_space;
747: MatPreallocateInitialize(comm,am,pn,dnz,onz);
748: for (i=0; i<am; i++) {
749: /* diagonal portion of A */
750: nzi = adi[i+1] - adi[i];
751: for (j=0; j<nzi; j++) {
752: row = *adj++;
753: pnz = pi_loc[row+1] - pi_loc[row];
754: Jptr = pj_loc + pi_loc[row];
755: /* Expand list if it is not long enough */
756: if (pnz+apnz_max > lsize) {
757: lsize = pnz+apnz_max;
758: PetscLLCondensedExpand_Scalable(lsize, &lnk);
759: }
760: /* add non-zero cols of P into the sorted linked list lnk */
761: PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);
762: apnz = *lnk; /* The first element in the list is the number of items in the list */
763: api[i+1] = api[i] + apnz;
764: if (apnz > apnz_max) apnz_max = apnz;
765: }
766: /* off-diagonal portion of A */
767: nzi = aoi[i+1] - aoi[i];
768: for (j=0; j<nzi; j++) {
769: row = *aoj++;
770: pnz = pi_oth[row+1] - pi_oth[row];
771: Jptr = pj_oth + pi_oth[row];
772: /* Expand list if it is not long enough */
773: if (pnz+apnz_max > lsize) {
774: lsize = pnz + apnz_max;
775: PetscLLCondensedExpand_Scalable(lsize, &lnk);
776: }
777: /* add non-zero cols of P into the sorted linked list lnk */
778: PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);
779: apnz = *lnk; /* The first element in the list is the number of items in the list */
780: api[i+1] = api[i] + apnz;
781: if (apnz > apnz_max) apnz_max = apnz;
782: }
783: apnz = *lnk;
784: api[i+1] = api[i] + apnz;
785: if (apnz > apnz_max) apnz_max = apnz;
787: /* if free space is not available, double the total space in the list */
788: if (current_space->local_remaining<apnz) {
789: PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);
790: nspacedouble++;
791: }
793: /* Copy data into free space, then initialize lnk */
794: PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);
795: MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);
797: current_space->array += apnz;
798: current_space->local_used += apnz;
799: current_space->local_remaining -= apnz;
800: }
802: /* Allocate space for apj, initialize apj, and */
803: /* destroy list of free space and other temporary array(s) */
804: PetscMalloc1(api[am]+1,&ptap->apj);
805: apj = ptap->apj;
806: PetscFreeSpaceContiguous(&free_space,ptap->apj);
807: PetscLLCondensedDestroy_Scalable(lnk);
809: /* create and assemble symbolic parallel matrix Cmpi */
810: /*----------------------------------------------------*/
811: MatCreate(comm,&Cmpi);
812: MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);
813: MatSetBlockSizesFromMats(Cmpi,A,P);
814: MatSetType(Cmpi,MATMPIAIJ);
815: MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);
817: /* malloc apa for assembly Cmpi */
818: PetscCalloc1(apnz_max,&apa);
819: ptap->apa = apa;
821: MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);
822: MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);
823: MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);
824: MatPreallocateFinalize(dnz,onz);
826: ptap->destroy = Cmpi->ops->destroy;
827: ptap->duplicate = Cmpi->ops->duplicate;
828: Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
829: Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult;
830: Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
832: /* attach the supporting struct to Cmpi for reuse */
833: c = (Mat_MPIAIJ*)Cmpi->data;
834: c->ptap = ptap;
835: *C = Cmpi;
837: /* set MatInfo */
838: afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
839: if (afill < 1.0) afill = 1.0;
840: Cmpi->info.mallocs = nspacedouble;
841: Cmpi->info.fill_ratio_given = fill;
842: Cmpi->info.fill_ratio_needed = afill;
844: #if defined(PETSC_USE_INFO)
845: if (api[am]) {
846: PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
847: PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
848: } else {
849: PetscInfo(Cmpi,"Empty matrix product\n");
850: }
851: #endif
852: return(0);
853: }
855: /* This function is needed for the seqMPI matrix-matrix multiplication. */
856: /* Three input arrays are merged to one output array. The size of the */
857: /* output array is also output. Duplicate entries only show up once. */
858: static void Merge3SortedArrays(PetscInt size1, PetscInt *in1,
859: PetscInt size2, PetscInt *in2,
860: PetscInt size3, PetscInt *in3,
861: PetscInt *size4, PetscInt *out)
862: {
863: int i = 0, j = 0, k = 0, l = 0;
865: /* Traverse all three arrays */
866: while (i<size1 && j<size2 && k<size3) {
867: if (in1[i] < in2[j] && in1[i] < in3[k]) {
868: out[l++] = in1[i++];
869: }
870: else if(in2[j] < in1[i] && in2[j] < in3[k]) {
871: out[l++] = in2[j++];
872: }
873: else if(in3[k] < in1[i] && in3[k] < in2[j]) {
874: out[l++] = in3[k++];
875: }
876: else if(in1[i] == in2[j] && in1[i] < in3[k]) {
877: out[l++] = in1[i];
878: i++, j++;
879: }
880: else if(in1[i] == in3[k] && in1[i] < in2[j]) {
881: out[l++] = in1[i];
882: i++, k++;
883: }
884: else if(in3[k] == in2[j] && in2[j] < in1[i]) {
885: out[l++] = in2[j];
886: k++, j++;
887: }
888: else if(in1[i] == in2[j] && in1[i] == in3[k]) {
889: out[l++] = in1[i];
890: i++, j++, k++;
891: }
892: }
894: /* Traverse two remaining arrays */
895: while (i<size1 && j<size2) {
896: if (in1[i] < in2[j]) {
897: out[l++] = in1[i++];
898: }
899: else if(in1[i] > in2[j]) {
900: out[l++] = in2[j++];
901: }
902: else {
903: out[l++] = in1[i];
904: i++, j++;
905: }
906: }
908: while (i<size1 && k<size3) {
909: if (in1[i] < in3[k]) {
910: out[l++] = in1[i++];
911: }
912: else if(in1[i] > in3[k]) {
913: out[l++] = in3[k++];
914: }
915: else {
916: out[l++] = in1[i];
917: i++, k++;
918: }
919: }
921: while (k<size3 && j<size2) {
922: if (in3[k] < in2[j]) {
923: out[l++] = in3[k++];
924: }
925: else if(in3[k] > in2[j]) {
926: out[l++] = in2[j++];
927: }
928: else {
929: out[l++] = in3[k];
930: k++, j++;
931: }
932: }
934: /* Traverse one remaining array */
935: while (i<size1) out[l++] = in1[i++];
936: while (j<size2) out[l++] = in2[j++];
937: while (k<size3) out[l++] = in3[k++];
939: *size4 = l;
940: }
942: /* This matrix-matrix multiplication algorithm divides the multiplication into three multiplications and */
943: /* adds up the products. Two of these three multiplications are performed with existing (sequential) */
944: /* matrix-matrix multiplications. */
947: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat A, Mat P, PetscReal fill, Mat *C)
948: {
949: PetscErrorCode ierr;
950: MPI_Comm comm;
951: PetscMPIInt size;
952: Mat Cmpi;
953: Mat_PtAPMPI *ptap;
954: PetscFreeSpaceList free_space_diag=NULL, current_space=NULL;
955: Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data;
956: Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc;
957: Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data;
958: Mat_MPIAIJ *c;
959: Mat_SeqAIJ *adpd_seq, *p_off, *aopoth_seq;
960: PetscInt adponz, adpdnz;
961: PetscInt *pi_loc,*dnz,*onz;
962: PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,rstart=A->rmap->rstart;
963: PetscInt *lnk,i, i1=0,pnz,row,*adpoi,*adpoj, *api, *adpoJ, *aopJ, *apJ,*Jptr, aopnz, nspacedouble=0,j,nzi,
964: *apj,apnz, *adpdi, *adpdj, *adpdJ, *poff_i, *poff_j, *j_temp, *aopothi, *aopothj;
965: PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n, p_colstart, p_colend;
966: PetscBT lnkbt;
967: PetscScalar *apa;
968: PetscReal afill;
969: PetscMPIInt rank;
970: Mat adpd, aopoth;
973: PetscObjectGetComm((PetscObject)A,&comm);
974: MPI_Comm_size(comm,&size);
975: MPI_Comm_rank(comm, &rank);
976: MatGetOwnershipRangeColumn(P, &p_colstart, &p_colend);
978: /* create struct Mat_PtAPMPI and attached it to C later */
979: PetscNew(&ptap);
981: /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
982: MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);
984: /* get P_loc by taking all local rows of P */
985: MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);
988: p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
989: pi_loc = p_loc->i;
991: /* Allocate memory for the i arrays of the matrices A*P, A_diag*P_off and A_offd * P */
992: PetscMalloc1(am+2,&api);
993: PetscMalloc1(am+2,&adpoi);
995: adpoi[0] = 0;
996: ptap->api = api;
997: api[0] = 0;
999: /* create and initialize a linked list, will be used for both A_diag * P_loc_off and A_offd * P_oth */
1000: PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);
1001: MatPreallocateInitialize(comm,am,pn,dnz,onz);
1003: /* Symbolic calc of A_loc_diag * P_loc_diag */
1004: MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->A, p->A, fill, &adpd);
1005: adpd_seq = (Mat_SeqAIJ*)((adpd)->data);
1006: adpdi = adpd_seq->i; adpdj = adpd_seq->j;
1007: p_off = (Mat_SeqAIJ*)((p->B)->data);
1008: poff_i = p_off->i; poff_j = p_off->j;
1010: /* j_temp stores indices of a result row before they are added to the linked list */
1011: PetscMalloc1(pN+2,&j_temp);
1014: /* Symbolic calc of the A_diag * p_loc_off */
1015: /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
1016: PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space_diag);
1017: current_space = free_space_diag;
1019: for (i=0; i<am; i++) {
1020: /* A_diag * P_loc_off */
1021: nzi = adi[i+1] - adi[i];
1022: for (j=0; j<nzi; j++) {
1023: row = *adj++;
1024: pnz = poff_i[row+1] - poff_i[row];
1025: Jptr = poff_j + poff_i[row];
1026: for(i1 = 0; i1 < pnz; i1++) {
1027: j_temp[i1] = p->garray[Jptr[i1]];
1028: }
1029: /* add non-zero cols of P into the sorted linked list lnk */
1030: PetscLLCondensedAddSorted(pnz,j_temp,lnk,lnkbt);
1031: }
1033: adponz = lnk[0];
1034: adpoi[i+1] = adpoi[i] + adponz;
1036: /* if free space is not available, double the total space in the list */
1037: if (current_space->local_remaining<adponz) {
1038: PetscFreeSpaceGet(PetscIntSumTruncate(adponz,current_space->total_array_size),¤t_space);
1039: nspacedouble++;
1040: }
1042: /* Copy data into free space, then initialize lnk */
1043: PetscLLCondensedClean(pN,adponz,current_space->array,lnk,lnkbt);
1045: current_space->array += adponz;
1046: current_space->local_used += adponz;
1047: current_space->local_remaining -= adponz;
1048: }
1050: /* Symbolic calc of A_off * P_oth */
1051: MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->B, ptap->P_oth, fill, &aopoth);
1052: aopoth_seq = (Mat_SeqAIJ*)((aopoth)->data);
1053: aopothi = aopoth_seq->i; aopothj = aopoth_seq->j;
1055: /* Allocate space for apj, adpj, aopj, ... */
1056: /* destroy lists of free space and other temporary array(s) */
1058: PetscMalloc1(aopothi[am] + adpoi[am] + adpdi[am]+2, &ptap->apj);
1059: PetscMalloc1(adpoi[am]+2, &adpoj);
1061: /* Copy from linked list to j-array */
1062: PetscFreeSpaceContiguous(&free_space_diag,adpoj);
1063: PetscLLDestroy(lnk,lnkbt);
1065: adpoJ = adpoj;
1066: adpdJ = adpdj;
1067: aopJ = aopothj;
1068: apj = ptap->apj;
1069: apJ = apj; /* still empty */
1071: /* Merge j-arrays of A_off * P, A_diag * P_loc_off, and */
1072: /* A_diag * P_loc_diag to get A*P */
1073: for (i = 0; i < am; i++) {
1074: aopnz = aopothi[i+1] - aopothi[i];
1075: adponz = adpoi[i+1] - adpoi[i];
1076: adpdnz = adpdi[i+1] - adpdi[i];
1078: /* Correct indices from A_diag*P_diag */
1079: for(i1 = 0; i1 < adpdnz; i1++) {
1080: adpdJ[i1] += p_colstart;
1081: }
1082: /* Merge j-arrays of A_diag * P_loc_off and A_diag * P_loc_diag and A_off * P_oth */
1083: Merge3SortedArrays(adponz, adpoJ, adpdnz, adpdJ, aopnz, aopJ, &apnz, apJ);
1084: MatPreallocateSet(i+rstart, apnz, apJ, dnz, onz);
1086: aopJ += aopnz;
1087: adpoJ += adponz;
1088: adpdJ += adpdnz;
1089: apJ += apnz;
1090: api[i+1] = api[i] + apnz;
1091: }
1093: /* malloc apa to store dense row A[i,:]*P */
1094: PetscCalloc1(pN+2,&apa);
1096: ptap->apa = apa;
1097: /* create and assemble symbolic parallel matrix Cmpi */
1098: MatCreate(comm,&Cmpi);
1099: MatSetSizes(Cmpi,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);
1100: MatSetBlockSizesFromMats(Cmpi,A,P);
1102: MatSetType(Cmpi,MATMPIAIJ);
1103: MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);
1106: MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(Cmpi, apj, api);
1107: MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);
1108: MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);
1109: MatPreallocateFinalize(dnz,onz);
1112: ptap->destroy = Cmpi->ops->destroy;
1113: ptap->duplicate = Cmpi->ops->duplicate;
1114: Cmpi->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1115: Cmpi->ops->destroy = MatDestroy_MPIAIJ_MatMatMult;
1116: Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatMatMult;
1118: /* attach the supporting struct to Cmpi for reuse */
1119: c = (Mat_MPIAIJ*)Cmpi->data;
1120: c->ptap = ptap;
1121: *C = Cmpi;
1123: /* set MatInfo */
1124: afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
1125: if (afill < 1.0) afill = 1.0;
1126: Cmpi->info.mallocs = nspacedouble;
1127: Cmpi->info.fill_ratio_given = fill;
1128: Cmpi->info.fill_ratio_needed = afill;
1130: #if defined(PETSC_USE_INFO)
1131: if (api[am]) {
1132: PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
1133: PetscInfo1(Cmpi,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
1134: } else {
1135: PetscInfo(Cmpi,"Empty matrix product\n");
1136: }
1137: #endif
1139: MatDestroy(&aopoth);
1140: MatDestroy(&adpd);
1141: PetscFree(j_temp);
1142: PetscFree(adpoj);
1143: PetscFree(adpoi);
1144: return(0);
1145: }
1148: /*-------------------------------------------------------------------------*/
1149: PetscErrorCode MatTransposeMatMult_MPIAIJ_MPIAIJ(Mat P,Mat A,MatReuse scall,PetscReal fill,Mat *C)
1150: {
1152: const char *algTypes[3] = {"scalable","nonscalable","matmatmult"};
1153: PetscInt aN=A->cmap->N,alg=1; /* set default algorithm */
1154: PetscBool flg;
1157: if (scall == MAT_INITIAL_MATRIX) {
1158: PetscObjectOptionsBegin((PetscObject)A);
1159: PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */
1160: PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,3,algTypes[1],&alg,&flg);
1161: PetscOptionsEnd();
1163: PetscLogEventBegin(MAT_TransposeMatMultSymbolic,P,A,0,0);
1164: switch (alg) {
1165: case 1:
1166: if (!flg && aN > 100000) { /* may switch to scalable algorithm as default */
1167: MatInfo Ainfo,Pinfo;
1168: PetscInt nz_local;
1169: PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable;
1170: MPI_Comm comm;
1172: MatGetInfo(A,MAT_LOCAL,&Ainfo);
1173: MatGetInfo(P,MAT_LOCAL,&Pinfo);
1174: nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated); /* estimated local nonzero entries */
1176: if (aN > fill*nz_local) alg_scalable_loc = PETSC_TRUE;
1177: PetscObjectGetComm((PetscObject)A,&comm);
1178: MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);
1180: if (alg_scalable) {
1181: alg = 0; /* scalable algorithm would slower than nonscalable algorithm */
1182: MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);
1183: break;
1184: }
1185: }
1186: MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(P,A,fill,C);
1187: break;
1188: case 2:
1189: {
1190: Mat Pt;
1191: Mat_PtAPMPI *ptap;
1192: Mat_MPIAIJ *c;
1193: MatTranspose(P,MAT_INITIAL_MATRIX,&Pt);
1194: MatMatMult(Pt,A,MAT_INITIAL_MATRIX,fill,C);
1195: c = (Mat_MPIAIJ*)(*C)->data;
1196: ptap = c->ptap;
1197: ptap->Pt = Pt;
1198: (*C)->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
1199: return(0);
1200: }
1201: break;
1202: default: /* scalable algorithm */
1203: MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(P,A,fill,C);
1204: break;
1205: }
1206: PetscLogEventEnd(MAT_TransposeMatMultSymbolic,P,A,0,0);
1207: }
1208: PetscLogEventBegin(MAT_TransposeMatMultNumeric,P,A,0,0);
1209: (*(*C)->ops->mattransposemultnumeric)(P,A,*C);
1210: PetscLogEventEnd(MAT_TransposeMatMultNumeric,P,A,0,0);
1211: return(0);
1212: }
1214: /* This routine only works when scall=MAT_REUSE_MATRIX! */
1215: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
1216: {
1218: Mat_MPIAIJ *c=(Mat_MPIAIJ*)C->data;
1219: Mat_PtAPMPI *ptap= c->ptap;
1220: Mat Pt=ptap->Pt;
1223: MatTranspose(P,MAT_REUSE_MATRIX,&Pt);
1224: MatMatMultNumeric(Pt,A,C);
1225: return(0);
1226: }
1228: PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat,MatDuplicateOption,Mat*);
1230: /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1231: PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat *C)
1232: {
1233: PetscErrorCode ierr;
1234: Mat_PtAPMPI *ptap;
1235: Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c;
1236: MPI_Comm comm;
1237: PetscMPIInt size,rank;
1238: Mat Cmpi;
1239: PetscFreeSpaceList free_space=NULL,current_space=NULL;
1240: PetscInt pn=P->cmap->n,aN=A->cmap->N,an=A->cmap->n;
1241: PetscInt *lnk,i,k,nsend;
1242: PetscBT lnkbt;
1243: PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv;
1244: PetscInt **buf_rj,**buf_ri,**buf_ri_k;
1245: PetscInt len,proc,*dnz,*onz,*owners,nzi;
1246: PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1247: MPI_Request *swaits,*rwaits;
1248: MPI_Status *sstatus,rstatus;
1249: PetscLayout rowmap;
1250: PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */
1251: PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */
1252: PetscInt *Jptr,*prmap=p->garray,con,j,Crmax;
1253: Mat_SeqAIJ *a_loc,*c_loc,*c_oth;
1254: PetscTable ta;
1257: PetscObjectGetComm((PetscObject)A,&comm);
1258: MPI_Comm_size(comm,&size);
1259: MPI_Comm_rank(comm,&rank);
1261: /* create symbolic parallel matrix Cmpi */
1262: MatCreate(comm,&Cmpi);
1263: MatSetType(Cmpi,MATMPIAIJ);
1265: Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1267: /* create struct Mat_PtAPMPI and attached it to C later */
1268: PetscNew(&ptap);
1269: ptap->reuse = MAT_INITIAL_MATRIX;
1271: /* (0) compute Rd = Pd^T, Ro = Po^T */
1272: /* --------------------------------- */
1273: MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);
1274: MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);
1276: /* (1) compute symbolic A_loc */
1277: /* ---------------------------*/
1278: MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&ptap->A_loc);
1280: /* (2-1) compute symbolic C_oth = Ro*A_loc */
1281: /* ------------------------------------ */
1282: MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->A_loc,fill,&ptap->C_oth);
1284: /* (3) send coj of C_oth to other processors */
1285: /* ------------------------------------------ */
1286: /* determine row ownership */
1287: PetscLayoutCreate(comm,&rowmap);
1288: rowmap->n = pn;
1289: rowmap->bs = 1;
1290: PetscLayoutSetUp(rowmap);
1291: owners = rowmap->range;
1293: /* determine the number of messages to send, their lengths */
1294: PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);
1295: PetscMemzero(len_s,size*sizeof(PetscMPIInt));
1296: PetscMemzero(len_si,size*sizeof(PetscMPIInt));
1298: c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
1299: coi = c_oth->i; coj = c_oth->j;
1300: con = ptap->C_oth->rmap->n;
1301: proc = 0;
1302: for (i=0; i<con; i++) {
1303: while (prmap[i] >= owners[proc+1]) proc++;
1304: len_si[proc]++; /* num of rows in Co(=Pt*A) to be sent to [proc] */
1305: len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
1306: }
1308: len = 0; /* max length of buf_si[], see (4) */
1309: owners_co[0] = 0;
1310: nsend = 0;
1311: for (proc=0; proc<size; proc++) {
1312: owners_co[proc+1] = owners_co[proc] + len_si[proc];
1313: if (len_s[proc]) {
1314: nsend++;
1315: len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
1316: len += len_si[proc];
1317: }
1318: }
1320: /* determine the number and length of messages to receive for coi and coj */
1321: PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);
1322: PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);
1324: /* post the Irecv and Isend of coj */
1325: PetscCommGetNewTag(comm,&tagj);
1326: PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);
1327: PetscMalloc1(nsend+1,&swaits);
1328: for (proc=0, k=0; proc<size; proc++) {
1329: if (!len_s[proc]) continue;
1330: i = owners_co[proc];
1331: MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);
1332: k++;
1333: }
1335: /* (2-2) compute symbolic C_loc = Rd*A_loc */
1336: /* ---------------------------------------- */
1337: MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->A_loc,fill,&ptap->C_loc);
1338: c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
1340: /* receives coj are complete */
1341: for (i=0; i<nrecv; i++) {
1342: MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);
1343: }
1344: PetscFree(rwaits);
1345: if (nsend) {MPI_Waitall(nsend,swaits,sstatus);}
1347: /* add received column indices into ta to update Crmax */
1348: a_loc = (Mat_SeqAIJ*)(ptap->A_loc)->data;
1350: /* create and initialize a linked list */
1351: PetscTableCreate(an,aN,&ta); /* for compute Crmax */
1352: MatRowMergeMax_SeqAIJ(a_loc,ptap->A_loc->rmap->N,ta);
1354: for (k=0; k<nrecv; k++) {/* k-th received message */
1355: Jptr = buf_rj[k];
1356: for (j=0; j<len_r[k]; j++) {
1357: PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);
1358: }
1359: }
1360: PetscTableGetCount(ta,&Crmax);
1361: PetscTableDestroy(&ta);
1363: /* (4) send and recv coi */
1364: /*-----------------------*/
1365: PetscCommGetNewTag(comm,&tagi);
1366: PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);
1367: PetscMalloc1(len+1,&buf_s);
1368: buf_si = buf_s; /* points to the beginning of k-th msg to be sent */
1369: for (proc=0,k=0; proc<size; proc++) {
1370: if (!len_s[proc]) continue;
1371: /* form outgoing message for i-structure:
1372: buf_si[0]: nrows to be sent
1373: [1:nrows]: row index (global)
1374: [nrows+1:2*nrows+1]: i-structure index
1375: */
1376: /*-------------------------------------------*/
1377: nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
1378: buf_si_i = buf_si + nrows+1;
1379: buf_si[0] = nrows;
1380: buf_si_i[0] = 0;
1381: nrows = 0;
1382: for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1383: nzi = coi[i+1] - coi[i];
1384: buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1385: buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */
1386: nrows++;
1387: }
1388: MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);
1389: k++;
1390: buf_si += len_si[proc];
1391: }
1392: for (i=0; i<nrecv; i++) {
1393: MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);
1394: }
1395: PetscFree(rwaits);
1396: if (nsend) {MPI_Waitall(nsend,swaits,sstatus);}
1398: PetscFree4(len_s,len_si,sstatus,owners_co);
1399: PetscFree(len_ri);
1400: PetscFree(swaits);
1401: PetscFree(buf_s);
1403: /* (5) compute the local portion of Cmpi */
1404: /* ------------------------------------------ */
1405: /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */
1406: PetscFreeSpaceGet(Crmax,&free_space);
1407: current_space = free_space;
1409: PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);
1410: for (k=0; k<nrecv; k++) {
1411: buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1412: nrows = *buf_ri_k[k];
1413: nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */
1414: nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */
1415: }
1417: MatPreallocateInitialize(comm,pn,an,dnz,onz);
1418: PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);
1419: for (i=0; i<pn; i++) {
1420: /* add C_loc into Cmpi */
1421: nzi = c_loc->i[i+1] - c_loc->i[i];
1422: Jptr = c_loc->j + c_loc->i[i];
1423: PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);
1425: /* add received col data into lnk */
1426: for (k=0; k<nrecv; k++) { /* k-th received message */
1427: if (i == *nextrow[k]) { /* i-th row */
1428: nzi = *(nextci[k]+1) - *nextci[k];
1429: Jptr = buf_rj[k] + *nextci[k];
1430: PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);
1431: nextrow[k]++; nextci[k]++;
1432: }
1433: }
1434: nzi = lnk[0];
1436: /* copy data into free space, then initialize lnk */
1437: PetscLLCondensedClean(aN,nzi,current_space->array,lnk,lnkbt);
1438: MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);
1439: }
1440: PetscFree3(buf_ri_k,nextrow,nextci);
1441: PetscLLDestroy(lnk,lnkbt);
1442: PetscFreeSpaceDestroy(free_space);
1444: /* local sizes and preallocation */
1445: MatSetSizes(Cmpi,pn,an,PETSC_DETERMINE,PETSC_DETERMINE);
1446: MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));
1447: MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);
1448: MatPreallocateFinalize(dnz,onz);
1450: /* members in merge */
1451: PetscFree(id_r);
1452: PetscFree(len_r);
1453: PetscFree(buf_ri[0]);
1454: PetscFree(buf_ri);
1455: PetscFree(buf_rj[0]);
1456: PetscFree(buf_rj);
1457: PetscLayoutDestroy(&rowmap);
1459: /* attach the supporting struct to Cmpi for reuse */
1460: c = (Mat_MPIAIJ*)Cmpi->data;
1461: c->ptap = ptap;
1462: ptap->duplicate = Cmpi->ops->duplicate;
1463: ptap->destroy = Cmpi->ops->destroy;
1465: /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
1466: Cmpi->assembled = PETSC_FALSE;
1467: Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP;
1468: Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
1469: *C = Cmpi;
1470: return(0);
1471: }
1473: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
1474: {
1475: PetscErrorCode ierr;
1476: Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1477: Mat_SeqAIJ *c_seq;
1478: Mat_PtAPMPI *ptap = c->ptap;
1479: Mat A_loc,C_loc,C_oth;
1480: PetscInt i,rstart,rend,cm,ncols,row;
1481: const PetscInt *cols;
1482: const PetscScalar *vals;
1485: MatZeroEntries(C);
1487: if (ptap->reuse == MAT_REUSE_MATRIX) {
1488: /* These matrices are obtained in MatTransposeMatMultSymbolic() */
1489: /* 1) get R = Pd^T, Ro = Po^T */
1490: /*----------------------------*/
1491: MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);
1492: MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);
1494: /* 2) compute numeric A_loc */
1495: /*--------------------------*/
1496: MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&ptap->A_loc);
1497: }
1499: /* 3) C_loc = Rd*A_loc, C_oth = Ro*A_loc */
1500: A_loc = ptap->A_loc;
1501: ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,A_loc,ptap->C_loc);
1502: ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,A_loc,ptap->C_oth);
1503: C_loc = ptap->C_loc;
1504: C_oth = ptap->C_oth;
1506: /* add C_loc and Co to to C */
1507: MatGetOwnershipRange(C,&rstart,&rend);
1509: /* C_loc -> C */
1510: cm = C_loc->rmap->N;
1511: c_seq = (Mat_SeqAIJ*)C_loc->data;
1512: cols = c_seq->j;
1513: vals = c_seq->a;
1514: for (i=0; i<cm; i++) {
1515: ncols = c_seq->i[i+1] - c_seq->i[i];
1516: row = rstart + i;
1517: MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);
1518: cols += ncols; vals += ncols;
1519: }
1521: /* Co -> C, off-processor part */
1522: cm = C_oth->rmap->N;
1523: c_seq = (Mat_SeqAIJ*)C_oth->data;
1524: cols = c_seq->j;
1525: vals = c_seq->a;
1526: for (i=0; i<cm; i++) {
1527: ncols = c_seq->i[i+1] - c_seq->i[i];
1528: row = p->garray[i];
1529: MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);
1530: cols += ncols; vals += ncols;
1531: }
1532: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1533: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1535: ptap->reuse = MAT_REUSE_MATRIX;
1536: return(0);
1537: }
1539: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1540: {
1541: PetscErrorCode ierr;
1542: Mat_Merge_SeqsToMPI *merge;
1543: Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
1544: Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1545: Mat_PtAPMPI *ptap;
1546: PetscInt *adj;
1547: PetscInt i,j,k,anz,pnz,row,*cj,nexta;
1548: MatScalar *ada,*ca,valtmp;
1549: PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1550: MPI_Comm comm;
1551: PetscMPIInt size,rank,taga,*len_s;
1552: PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1553: PetscInt **buf_ri,**buf_rj;
1554: PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
1555: MPI_Request *s_waits,*r_waits;
1556: MPI_Status *status;
1557: MatScalar **abuf_r,*ba_i,*pA,*coa,*ba;
1558: PetscInt *ai,*aj,*coi,*coj,*poJ,*pdJ;
1559: Mat A_loc;
1560: Mat_SeqAIJ *a_loc;
1563: PetscObjectGetComm((PetscObject)C,&comm);
1564: MPI_Comm_size(comm,&size);
1565: MPI_Comm_rank(comm,&rank);
1567: ptap = c->ptap;
1568: merge = ptap->merge;
1570: /* 2) compute numeric C_seq = P_loc^T*A_loc */
1571: /*------------------------------------------*/
1572: /* get data from symbolic products */
1573: coi = merge->coi; coj = merge->coj;
1574: PetscCalloc1(coi[pon]+1,&coa);
1575: bi = merge->bi; bj = merge->bj;
1576: owners = merge->rowmap->range;
1577: PetscCalloc1(bi[cm]+1,&ba);
1579: /* get A_loc by taking all local rows of A */
1580: A_loc = ptap->A_loc;
1581: MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);
1582: a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1583: ai = a_loc->i;
1584: aj = a_loc->j;
1586: for (i=0; i<am; i++) {
1587: anz = ai[i+1] - ai[i];
1588: adj = aj + ai[i];
1589: ada = a_loc->a + ai[i];
1591: /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1592: /*-------------------------------------------------------------*/
1593: /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1594: pnz = po->i[i+1] - po->i[i];
1595: poJ = po->j + po->i[i];
1596: pA = po->a + po->i[i];
1597: for (j=0; j<pnz; j++) {
1598: row = poJ[j];
1599: cj = coj + coi[row];
1600: ca = coa + coi[row];
1601: /* perform sparse axpy */
1602: nexta = 0;
1603: valtmp = pA[j];
1604: for (k=0; nexta<anz; k++) {
1605: if (cj[k] == adj[nexta]) {
1606: ca[k] += valtmp*ada[nexta];
1607: nexta++;
1608: }
1609: }
1610: PetscLogFlops(2.0*anz);
1611: }
1613: /* put the value into Cd (diagonal part) */
1614: pnz = pd->i[i+1] - pd->i[i];
1615: pdJ = pd->j + pd->i[i];
1616: pA = pd->a + pd->i[i];
1617: for (j=0; j<pnz; j++) {
1618: row = pdJ[j];
1619: cj = bj + bi[row];
1620: ca = ba + bi[row];
1621: /* perform sparse axpy */
1622: nexta = 0;
1623: valtmp = pA[j];
1624: for (k=0; nexta<anz; k++) {
1625: if (cj[k] == adj[nexta]) {
1626: ca[k] += valtmp*ada[nexta];
1627: nexta++;
1628: }
1629: }
1630: PetscLogFlops(2.0*anz);
1631: }
1632: }
1634: /* 3) send and recv matrix values coa */
1635: /*------------------------------------*/
1636: buf_ri = merge->buf_ri;
1637: buf_rj = merge->buf_rj;
1638: len_s = merge->len_s;
1639: PetscCommGetNewTag(comm,&taga);
1640: PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);
1642: PetscMalloc2(merge->nsend+1,&s_waits,size,&status);
1643: for (proc=0,k=0; proc<size; proc++) {
1644: if (!len_s[proc]) continue;
1645: i = merge->owners_co[proc];
1646: MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);
1647: k++;
1648: }
1649: if (merge->nrecv) {MPI_Waitall(merge->nrecv,r_waits,status);}
1650: if (merge->nsend) {MPI_Waitall(merge->nsend,s_waits,status);}
1652: PetscFree2(s_waits,status);
1653: PetscFree(r_waits);
1654: PetscFree(coa);
1656: /* 4) insert local Cseq and received values into Cmpi */
1657: /*----------------------------------------------------*/
1658: PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);
1659: for (k=0; k<merge->nrecv; k++) {
1660: buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1661: nrows = *(buf_ri_k[k]);
1662: nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */
1663: nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */
1664: }
1666: for (i=0; i<cm; i++) {
1667: row = owners[rank] + i; /* global row index of C_seq */
1668: bj_i = bj + bi[i]; /* col indices of the i-th row of C */
1669: ba_i = ba + bi[i];
1670: bnz = bi[i+1] - bi[i];
1671: /* add received vals into ba */
1672: for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1673: /* i-th row */
1674: if (i == *nextrow[k]) {
1675: cnz = *(nextci[k]+1) - *nextci[k];
1676: cj = buf_rj[k] + *(nextci[k]);
1677: ca = abuf_r[k] + *(nextci[k]);
1678: nextcj = 0;
1679: for (j=0; nextcj<cnz; j++) {
1680: if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1681: ba_i[j] += ca[nextcj++];
1682: }
1683: }
1684: nextrow[k]++; nextci[k]++;
1685: PetscLogFlops(2.0*cnz);
1686: }
1687: }
1688: MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);
1689: }
1690: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1691: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1693: PetscFree(ba);
1694: PetscFree(abuf_r[0]);
1695: PetscFree(abuf_r);
1696: PetscFree3(buf_ri_k,nextrow,nextci);
1697: return(0);
1698: }
1700: PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat *C)
1701: {
1702: PetscErrorCode ierr;
1703: Mat Cmpi,A_loc,POt,PDt;
1704: Mat_PtAPMPI *ptap;
1705: PetscFreeSpaceList free_space=NULL,current_space=NULL;
1706: Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*a=(Mat_MPIAIJ*)A->data,*c;
1707: PetscInt *pdti,*pdtj,*poti,*potj,*ptJ;
1708: PetscInt nnz;
1709: PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1710: PetscInt am =A->rmap->n,pn=P->cmap->n;
1711: MPI_Comm comm;
1712: PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1713: PetscInt **buf_rj,**buf_ri,**buf_ri_k;
1714: PetscInt len,proc,*dnz,*onz,*owners;
1715: PetscInt nzi,*bi,*bj;
1716: PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1717: MPI_Request *swaits,*rwaits;
1718: MPI_Status *sstatus,rstatus;
1719: Mat_Merge_SeqsToMPI *merge;
1720: PetscInt *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1721: PetscReal afill =1.0,afill_tmp;
1722: PetscInt rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax;
1723: PetscScalar *vals;
1724: Mat_SeqAIJ *a_loc,*pdt,*pot;
1725: PetscTable ta;
1728: PetscObjectGetComm((PetscObject)A,&comm);
1729: /* check if matrix local sizes are compatible */
1730: 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);
1732: MPI_Comm_size(comm,&size);
1733: MPI_Comm_rank(comm,&rank);
1735: /* create struct Mat_PtAPMPI and attached it to C later */
1736: PetscNew(&ptap);
1738: /* get A_loc by taking all local rows of A */
1739: MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);
1741: ptap->A_loc = A_loc;
1742: a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1743: ai = a_loc->i;
1744: aj = a_loc->j;
1746: /* determine symbolic Co=(p->B)^T*A - send to others */
1747: /*----------------------------------------------------*/
1748: MatTransposeSymbolic_SeqAIJ(p->A,&PDt);
1749: pdt = (Mat_SeqAIJ*)PDt->data;
1750: pdti = pdt->i; pdtj = pdt->j;
1752: MatTransposeSymbolic_SeqAIJ(p->B,&POt);
1753: pot = (Mat_SeqAIJ*)POt->data;
1754: poti = pot->i; potj = pot->j;
1756: /* then, compute symbolic Co = (p->B)^T*A */
1757: pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1758: >= (num of nonzero rows of C_seq) - pn */
1759: PetscMalloc1(pon+1,&coi);
1760: coi[0] = 0;
1762: /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1763: nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am]));
1764: PetscFreeSpaceGet(nnz,&free_space);
1765: current_space = free_space;
1767: /* create and initialize a linked list */
1768: PetscTableCreate(A->cmap->n + a->B->cmap->N,aN,&ta);
1769: MatRowMergeMax_SeqAIJ(a_loc,am,ta);
1770: PetscTableGetCount(ta,&Armax);
1772: PetscLLCondensedCreate_Scalable(Armax,&lnk);
1774: for (i=0; i<pon; i++) {
1775: pnz = poti[i+1] - poti[i];
1776: ptJ = potj + poti[i];
1777: for (j=0; j<pnz; j++) {
1778: row = ptJ[j]; /* row of A_loc == col of Pot */
1779: anz = ai[row+1] - ai[row];
1780: Jptr = aj + ai[row];
1781: /* add non-zero cols of AP into the sorted linked list lnk */
1782: PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);
1783: }
1784: nnz = lnk[0];
1786: /* If free space is not available, double the total space in the list */
1787: if (current_space->local_remaining<nnz) {
1788: PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);
1789: nspacedouble++;
1790: }
1792: /* Copy data into free space, and zero out denserows */
1793: PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);
1795: current_space->array += nnz;
1796: current_space->local_used += nnz;
1797: current_space->local_remaining -= nnz;
1799: coi[i+1] = coi[i] + nnz;
1800: }
1802: PetscMalloc1(coi[pon]+1,&coj);
1803: PetscFreeSpaceContiguous(&free_space,coj);
1804: PetscLLCondensedDestroy_Scalable(lnk); /* must destroy to get a new one for C */
1806: afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1807: if (afill_tmp > afill) afill = afill_tmp;
1809: /* send j-array (coj) of Co to other processors */
1810: /*----------------------------------------------*/
1811: /* determine row ownership */
1812: PetscNew(&merge);
1813: PetscLayoutCreate(comm,&merge->rowmap);
1815: merge->rowmap->n = pn;
1816: merge->rowmap->bs = 1;
1818: PetscLayoutSetUp(merge->rowmap);
1819: owners = merge->rowmap->range;
1821: /* determine the number of messages to send, their lengths */
1822: PetscCalloc1(size,&len_si);
1823: PetscMalloc1(size,&merge->len_s);
1825: len_s = merge->len_s;
1826: merge->nsend = 0;
1828: PetscMalloc1(size+2,&owners_co);
1829: PetscMemzero(len_s,size*sizeof(PetscMPIInt));
1831: proc = 0;
1832: for (i=0; i<pon; i++) {
1833: while (prmap[i] >= owners[proc+1]) proc++;
1834: len_si[proc]++; /* num of rows in Co to be sent to [proc] */
1835: len_s[proc] += coi[i+1] - coi[i];
1836: }
1838: len = 0; /* max length of buf_si[] */
1839: owners_co[0] = 0;
1840: for (proc=0; proc<size; proc++) {
1841: owners_co[proc+1] = owners_co[proc] + len_si[proc];
1842: if (len_si[proc]) {
1843: merge->nsend++;
1844: len_si[proc] = 2*(len_si[proc] + 1);
1845: len += len_si[proc];
1846: }
1847: }
1849: /* determine the number and length of messages to receive for coi and coj */
1850: PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);
1851: PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);
1853: /* post the Irecv and Isend of coj */
1854: PetscCommGetNewTag(comm,&tagj);
1855: PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);
1856: PetscMalloc1(merge->nsend+1,&swaits);
1857: for (proc=0, k=0; proc<size; proc++) {
1858: if (!len_s[proc]) continue;
1859: i = owners_co[proc];
1860: MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);
1861: k++;
1862: }
1864: /* receives and sends of coj are complete */
1865: PetscMalloc1(size,&sstatus);
1866: for (i=0; i<merge->nrecv; i++) {
1867: PetscMPIInt icompleted;
1868: MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);
1869: }
1870: PetscFree(rwaits);
1871: if (merge->nsend) {MPI_Waitall(merge->nsend,swaits,sstatus);}
1873: /* add received column indices into table to update Armax */
1874: /* Armax can be as large as aN if a P[row,:] is dense, see src/ksp/ksp/examples/tutorials/ex56.c! */
1875: for (k=0; k<merge->nrecv; k++) {/* k-th received message */
1876: Jptr = buf_rj[k];
1877: for (j=0; j<merge->len_r[k]; j++) {
1878: PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);
1879: }
1880: }
1881: PetscTableGetCount(ta,&Armax);
1882: /* 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); */
1884: /* send and recv coi */
1885: /*-------------------*/
1886: PetscCommGetNewTag(comm,&tagi);
1887: PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);
1888: PetscMalloc1(len+1,&buf_s);
1889: buf_si = buf_s; /* points to the beginning of k-th msg to be sent */
1890: for (proc=0,k=0; proc<size; proc++) {
1891: if (!len_s[proc]) continue;
1892: /* form outgoing message for i-structure:
1893: buf_si[0]: nrows to be sent
1894: [1:nrows]: row index (global)
1895: [nrows+1:2*nrows+1]: i-structure index
1896: */
1897: /*-------------------------------------------*/
1898: nrows = len_si[proc]/2 - 1;
1899: buf_si_i = buf_si + nrows+1;
1900: buf_si[0] = nrows;
1901: buf_si_i[0] = 0;
1902: nrows = 0;
1903: for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1904: nzi = coi[i+1] - coi[i];
1905: buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1906: buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */
1907: nrows++;
1908: }
1909: MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);
1910: k++;
1911: buf_si += len_si[proc];
1912: }
1913: i = merge->nrecv;
1914: while (i--) {
1915: PetscMPIInt icompleted;
1916: MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);
1917: }
1918: PetscFree(rwaits);
1919: if (merge->nsend) {MPI_Waitall(merge->nsend,swaits,sstatus);}
1920: PetscFree(len_si);
1921: PetscFree(len_ri);
1922: PetscFree(swaits);
1923: PetscFree(sstatus);
1924: PetscFree(buf_s);
1926: /* compute the local portion of C (mpi mat) */
1927: /*------------------------------------------*/
1928: /* allocate bi array and free space for accumulating nonzero column info */
1929: PetscMalloc1(pn+1,&bi);
1930: bi[0] = 0;
1932: /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1933: nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am])));
1934: PetscFreeSpaceGet(nnz,&free_space);
1935: current_space = free_space;
1937: PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);
1938: for (k=0; k<merge->nrecv; k++) {
1939: buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1940: nrows = *buf_ri_k[k];
1941: nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */
1942: nextci[k] = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th recieved i-structure */
1943: }
1945: PetscLLCondensedCreate_Scalable(Armax,&lnk);
1946: MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);
1947: rmax = 0;
1948: for (i=0; i<pn; i++) {
1949: /* add pdt[i,:]*AP into lnk */
1950: pnz = pdti[i+1] - pdti[i];
1951: ptJ = pdtj + pdti[i];
1952: for (j=0; j<pnz; j++) {
1953: row = ptJ[j]; /* row of AP == col of Pt */
1954: anz = ai[row+1] - ai[row];
1955: Jptr = aj + ai[row];
1956: /* add non-zero cols of AP into the sorted linked list lnk */
1957: PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);
1958: }
1960: /* add received col data into lnk */
1961: for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1962: if (i == *nextrow[k]) { /* i-th row */
1963: nzi = *(nextci[k]+1) - *nextci[k];
1964: Jptr = buf_rj[k] + *nextci[k];
1965: PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);
1966: nextrow[k]++; nextci[k]++;
1967: }
1968: }
1969: nnz = lnk[0];
1971: /* if free space is not available, make more free space */
1972: if (current_space->local_remaining<nnz) {
1973: PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);
1974: nspacedouble++;
1975: }
1976: /* copy data into free space, then initialize lnk */
1977: PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);
1978: MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);
1980: current_space->array += nnz;
1981: current_space->local_used += nnz;
1982: current_space->local_remaining -= nnz;
1984: bi[i+1] = bi[i] + nnz;
1985: if (nnz > rmax) rmax = nnz;
1986: }
1987: PetscFree3(buf_ri_k,nextrow,nextci);
1989: PetscMalloc1(bi[pn]+1,&bj);
1990: PetscFreeSpaceContiguous(&free_space,bj);
1991: afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
1992: if (afill_tmp > afill) afill = afill_tmp;
1993: PetscLLCondensedDestroy_Scalable(lnk);
1994: PetscTableDestroy(&ta);
1996: MatDestroy(&POt);
1997: MatDestroy(&PDt);
1999: /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part */
2000: /*----------------------------------------------------------------------------------*/
2001: PetscCalloc1(rmax+1,&vals);
2003: MatCreate(comm,&Cmpi);
2004: MatSetSizes(Cmpi,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);
2005: MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));
2006: MatSetType(Cmpi,MATMPIAIJ);
2007: MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);
2008: MatPreallocateFinalize(dnz,onz);
2009: MatSetBlockSize(Cmpi,1);
2010: for (i=0; i<pn; i++) {
2011: row = i + rstart;
2012: nnz = bi[i+1] - bi[i];
2013: Jptr = bj + bi[i];
2014: MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);
2015: }
2016: MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);
2017: MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);
2018: PetscFree(vals);
2020: merge->bi = bi;
2021: merge->bj = bj;
2022: merge->coi = coi;
2023: merge->coj = coj;
2024: merge->buf_ri = buf_ri;
2025: merge->buf_rj = buf_rj;
2026: merge->owners_co = owners_co;
2028: /* attach the supporting struct to Cmpi for reuse */
2029: c = (Mat_MPIAIJ*)Cmpi->data;
2031: c->ptap = ptap;
2032: ptap->api = NULL;
2033: ptap->apj = NULL;
2034: ptap->merge = merge;
2035: ptap->apa = NULL;
2036: ptap->destroy = Cmpi->ops->destroy;
2037: ptap->duplicate = Cmpi->ops->duplicate;
2039: Cmpi->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
2040: Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP;
2041: Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
2043: *C = Cmpi;
2044: #if defined(PETSC_USE_INFO)
2045: if (bi[pn] != 0) {
2046: PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
2047: PetscInfo1(Cmpi,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);
2048: } else {
2049: PetscInfo(Cmpi,"Empty matrix product\n");
2050: }
2051: #endif
2052: return(0);
2053: }