Actual source code: mpimatmatmult.c
petsc-3.14.3 2021-01-09
2: /*
3: Defines matrix-matrix product routines for pairs of MPIAIJ matrices
4: C = A * B
5: */
6: #include <../src/mat/impls/aij/seq/aij.h>
7: #include <../src/mat/utils/freespace.h>
8: #include <../src/mat/impls/aij/mpi/mpiaij.h>
9: #include <petscbt.h>
10: #include <../src/mat/impls/dense/mpi/mpidense.h>
11: #include <petsc/private/vecimpl.h>
12: #include <petsc/private/vecscatterimpl.h>
14: #if defined(PETSC_HAVE_HYPRE)
15: PETSC_INTERN PetscErrorCode MatMatMultSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat);
16: #endif
18: PETSC_INTERN PetscErrorCode MatProductSymbolic_AB_MPIAIJ_MPIAIJ(Mat C)
19: {
20: PetscErrorCode ierr;
21: Mat_Product *product = C->product;
22: Mat A=product->A,B=product->B;
23: MatProductAlgorithm alg=product->alg;
24: PetscReal fill=product->fill;
25: PetscBool flg;
28: /* scalable */
29: PetscStrcmp(alg,"scalable",&flg);
30: if (flg) {
31: MatMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);
32: return(0);
33: }
35: /* nonscalable */
36: PetscStrcmp(alg,"nonscalable",&flg);
37: if (flg) {
38: MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);
39: return(0);
40: }
42: /* seqmpi */
43: PetscStrcmp(alg,"seqmpi",&flg);
44: if (flg) {
45: MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(A,B,fill,C);
46: return(0);
47: }
49: #if defined(PETSC_HAVE_HYPRE)
50: PetscStrcmp(alg,"hypre",&flg);
51: if (flg) {
52: MatMatMultSymbolic_AIJ_AIJ_wHYPRE(A,B,fill,C);
53: return(0);
54: }
55: #endif
56: SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_SUP,"Mat Product Algorithm is not supported");
57: }
59: PetscErrorCode MatDestroy_MPIAIJ_MatMatMult(void *data)
60: {
62: Mat_APMPI *ptap = (Mat_APMPI*)data;
65: PetscFree2(ptap->startsj_s,ptap->startsj_r);
66: PetscFree(ptap->bufa);
67: MatDestroy(&ptap->P_loc);
68: MatDestroy(&ptap->P_oth);
69: MatDestroy(&ptap->Pt);
70: PetscFree(ptap->api);
71: PetscFree(ptap->apj);
72: PetscFree(ptap->apa);
73: PetscFree(ptap);
74: return(0);
75: }
77: PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,Mat C)
78: {
80: Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
81: Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
82: Mat_SeqAIJ *cd =(Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
83: PetscScalar *cda=cd->a,*coa=co->a;
84: Mat_SeqAIJ *p_loc,*p_oth;
85: PetscScalar *apa,*ca;
86: PetscInt cm =C->rmap->n;
87: Mat_APMPI *ptap;
88: PetscInt *api,*apj,*apJ,i,k;
89: PetscInt cstart=C->cmap->rstart;
90: PetscInt cdnz,conz,k0,k1;
91: MPI_Comm comm;
92: PetscMPIInt size;
95: MatCheckProduct(C,3);
96: ptap = (Mat_APMPI*)C->product->data;
97: if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data");
98: PetscObjectGetComm((PetscObject)A,&comm);
99: MPI_Comm_size(comm,&size);
101: if (!ptap->P_oth && size>1) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"AP cannot be reused. Do not call MatProductClear()");
103: /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */
104: /*-----------------------------------------------------*/
105: /* update numerical values of P_oth and P_loc */
106: MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);
107: MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);
109: /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
110: /*----------------------------------------------------------*/
111: /* get data from symbolic products */
112: p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
113: p_oth = NULL;
114: if (size >1) {
115: p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
116: }
118: /* get apa for storing dense row A[i,:]*P */
119: apa = ptap->apa;
121: api = ptap->api;
122: apj = ptap->apj;
123: for (i=0; i<cm; i++) {
124: /* compute apa = A[i,:]*P */
125: AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa);
127: /* set values in C */
128: apJ = apj + api[i];
129: cdnz = cd->i[i+1] - cd->i[i];
130: conz = co->i[i+1] - co->i[i];
132: /* 1st off-diagonal part of C */
133: ca = coa + co->i[i];
134: k = 0;
135: for (k0=0; k0<conz; k0++) {
136: if (apJ[k] >= cstart) break;
137: ca[k0] = apa[apJ[k]];
138: apa[apJ[k++]] = 0.0;
139: }
141: /* diagonal part of C */
142: ca = cda + cd->i[i];
143: for (k1=0; k1<cdnz; k1++) {
144: ca[k1] = apa[apJ[k]];
145: apa[apJ[k++]] = 0.0;
146: }
148: /* 2nd off-diagonal part of C */
149: ca = coa + co->i[i];
150: for (; k0<conz; k0++) {
151: ca[k0] = apa[apJ[k]];
152: apa[apJ[k++]] = 0.0;
153: }
154: }
155: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
156: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
157: return(0);
158: }
160: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat A,Mat P,PetscReal fill,Mat C)
161: {
162: PetscErrorCode ierr;
163: MPI_Comm comm;
164: PetscMPIInt size;
165: Mat_APMPI *ptap;
166: PetscFreeSpaceList free_space=NULL,current_space=NULL;
167: Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data;
168: Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
169: PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
170: PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
171: PetscInt *lnk,i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi;
172: PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n;
173: PetscBT lnkbt;
174: PetscReal afill;
175: MatType mtype;
178: MatCheckProduct(C,4);
179: if (C->product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty");
180: PetscObjectGetComm((PetscObject)A,&comm);
181: MPI_Comm_size(comm,&size);
183: /* create struct Mat_APMPI and attached it to C later */
184: PetscNew(&ptap);
186: /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
187: MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);
189: /* get P_loc by taking all local rows of P */
190: MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);
192: p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
193: pi_loc = p_loc->i; pj_loc = p_loc->j;
194: if (size > 1) {
195: p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
196: pi_oth = p_oth->i; pj_oth = p_oth->j;
197: } else {
198: p_oth = NULL;
199: pi_oth = NULL; pj_oth = NULL;
200: }
202: /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
203: /*-------------------------------------------------------------------*/
204: PetscMalloc1(am+2,&api);
205: ptap->api = api;
206: api[0] = 0;
208: /* create and initialize a linked list */
209: PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);
211: /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
212: PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);
213: current_space = free_space;
215: MatPreallocateInitialize(comm,am,pn,dnz,onz);
216: for (i=0; i<am; i++) {
217: /* diagonal portion of A */
218: nzi = adi[i+1] - adi[i];
219: for (j=0; j<nzi; j++) {
220: row = *adj++;
221: pnz = pi_loc[row+1] - pi_loc[row];
222: Jptr = pj_loc + pi_loc[row];
223: /* add non-zero cols of P into the sorted linked list lnk */
224: PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);
225: }
226: /* off-diagonal portion of A */
227: nzi = aoi[i+1] - aoi[i];
228: for (j=0; j<nzi; j++) {
229: row = *aoj++;
230: pnz = pi_oth[row+1] - pi_oth[row];
231: Jptr = pj_oth + pi_oth[row];
232: PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);
233: }
235: apnz = lnk[0];
236: api[i+1] = api[i] + apnz;
238: /* if free space is not available, double the total space in the list */
239: if (current_space->local_remaining<apnz) {
240: PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);
241: nspacedouble++;
242: }
244: /* Copy data into free space, then initialize lnk */
245: PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);
246: MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);
248: current_space->array += apnz;
249: current_space->local_used += apnz;
250: current_space->local_remaining -= apnz;
251: }
253: /* Allocate space for apj, initialize apj, and */
254: /* destroy list of free space and other temporary array(s) */
255: PetscMalloc1(api[am]+1,&ptap->apj);
256: apj = ptap->apj;
257: PetscFreeSpaceContiguous(&free_space,ptap->apj);
258: PetscLLDestroy(lnk,lnkbt);
260: /* malloc apa to store dense row A[i,:]*P */
261: PetscCalloc1(pN,&ptap->apa);
263: /* set and assemble symbolic parallel matrix C */
264: /*---------------------------------------------*/
265: MatSetSizes(C,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);
266: MatSetBlockSizesFromMats(C,A,P);
268: MatGetType(A,&mtype);
269: MatSetType(C,mtype);
270: MatMPIAIJSetPreallocation(C,0,dnz,0,onz);
271: MatPreallocateFinalize(dnz,onz);
273: MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(C, apj, api);
274: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
275: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
276: MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
278: C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
279: C->ops->productnumeric = MatProductNumeric_AB;
281: /* attach the supporting struct to C for reuse */
282: C->product->data = ptap;
283: C->product->destroy = MatDestroy_MPIAIJ_MatMatMult;
285: /* set MatInfo */
286: afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
287: if (afill < 1.0) afill = 1.0;
288: C->info.mallocs = nspacedouble;
289: C->info.fill_ratio_given = fill;
290: C->info.fill_ratio_needed = afill;
292: #if defined(PETSC_USE_INFO)
293: if (api[am]) {
294: PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
295: PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
296: } else {
297: PetscInfo(C,"Empty matrix product\n");
298: }
299: #endif
300: return(0);
301: }
303: /* ------------------------------------------------------- */
304: static PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat,Mat,PetscReal,Mat);
305: static PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat,Mat,Mat);
307: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense_AB(Mat C)
308: {
309: Mat_Product *product = C->product;
310: Mat A = product->A,B=product->B;
313: if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend)
314: SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,B->rmap->rstart,B->rmap->rend);
316: C->ops->matmultsymbolic = MatMatMultSymbolic_MPIAIJ_MPIDense;
317: C->ops->productsymbolic = MatProductSymbolic_AB;
318: return(0);
319: }
320: /* -------------------------------------------------------------------- */
321: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense_AtB(Mat C)
322: {
323: Mat_Product *product = C->product;
324: Mat A = product->A,B=product->B;
327: if (A->rmap->rstart != B->rmap->rstart || A->rmap->rend != B->rmap->rend)
328: SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->rmap->rstart,A->rmap->rend,B->rmap->rstart,B->rmap->rend);
330: C->ops->transposematmultsymbolic = MatTransposeMatMultSymbolic_MPIAIJ_MPIDense;
331: C->ops->productsymbolic = MatProductSymbolic_AtB;
332: return(0);
333: }
335: /* --------------------------------------------------------------------- */
336: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ_MPIDense(Mat C)
337: {
339: Mat_Product *product = C->product;
342: switch (product->type) {
343: case MATPRODUCT_AB:
344: MatProductSetFromOptions_MPIAIJ_MPIDense_AB(C);
345: break;
346: case MATPRODUCT_AtB:
347: MatProductSetFromOptions_MPIAIJ_MPIDense_AtB(C);
348: break;
349: default:
350: break;
351: }
352: return(0);
353: }
354: /* ------------------------------------------------------- */
356: typedef struct {
357: Mat workB,workB1;
358: MPI_Request *rwaits,*swaits;
359: PetscInt nsends,nrecvs;
360: MPI_Datatype *stype,*rtype;
361: PetscInt blda;
362: } MPIAIJ_MPIDense;
364: PetscErrorCode MatMPIAIJ_MPIDenseDestroy(void *ctx)
365: {
366: MPIAIJ_MPIDense *contents = (MPIAIJ_MPIDense*)ctx;
367: PetscErrorCode ierr;
368: PetscInt i;
371: MatDestroy(&contents->workB);
372: MatDestroy(&contents->workB1);
373: for (i=0; i<contents->nsends; i++) {
374: MPI_Type_free(&contents->stype[i]);
375: }
376: for (i=0; i<contents->nrecvs; i++) {
377: MPI_Type_free(&contents->rtype[i]);
378: }
379: PetscFree4(contents->stype,contents->rtype,contents->rwaits,contents->swaits);
380: PetscFree(contents);
381: return(0);
382: }
384: static PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIDense(Mat A,Mat B,PetscReal fill,Mat C)
385: {
386: PetscErrorCode ierr;
387: Mat_MPIAIJ *aij=(Mat_MPIAIJ*)A->data;
388: PetscInt nz=aij->B->cmap->n,nsends,nrecvs,i,nrows_to,j,blda,clda;
389: MPIAIJ_MPIDense *contents;
390: VecScatter ctx=aij->Mvctx;
391: PetscInt Am=A->rmap->n,Bm=B->rmap->n,BN=B->cmap->N,Bbn,Bbn1,bs,nrows_from,numBb;
392: MPI_Comm comm;
393: MPI_Datatype type1,*stype,*rtype;
394: const PetscInt *sindices,*sstarts,*rstarts;
395: PetscMPIInt *disp;
396: PetscBool cisdense;
399: MatCheckProduct(C,4);
400: if (C->product->data) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data not empty");
401: PetscObjectGetComm((PetscObject)A,&comm);
402: PetscObjectBaseTypeCompare((PetscObject)C,MATMPIDENSE,&cisdense);
403: if (!cisdense) {
404: MatSetType(C,((PetscObject)B)->type_name);
405: }
406: MatSetSizes(C,Am,B->cmap->n,A->rmap->N,BN);
407: MatSetBlockSizesFromMats(C,A,B);
408: MatSetUp(C);
409: MatDenseGetLDA(B,&blda);
410: MatDenseGetLDA(C,&clda);
411: PetscNew(&contents);
413: VecScatterGetRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,NULL,NULL);
414: VecScatterGetRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,NULL,NULL);
416: /* Create column block of B and C for memory scalability when BN is too large */
417: /* Estimate Bbn, column size of Bb */
418: if (nz) {
419: Bbn1 = 2*Am*BN/nz;
420: } else Bbn1 = BN;
422: bs = PetscAbs(B->cmap->bs);
423: Bbn1 = Bbn1/bs *bs; /* Bbn1 is a multiple of bs */
424: if (Bbn1 > BN) Bbn1 = BN;
425: MPI_Allreduce(&Bbn1,&Bbn,1,MPIU_INT,MPI_MAX,comm);
427: /* Enable runtime option for Bbn */
428: PetscOptionsBegin(comm,((PetscObject)C)->prefix,"MatMatMult","Mat");
429: PetscOptionsInt("-matmatmult_Bbn","Number of columns in Bb","MatMatMult",Bbn,&Bbn,NULL);
430: PetscOptionsEnd();
431: Bbn = PetscMin(Bbn,BN);
433: if (Bbn > 0 && Bbn < BN) {
434: numBb = BN/Bbn;
435: Bbn1 = BN - numBb*Bbn;
436: } else numBb = 0;
438: if (numBb) {
439: PetscInfo3(C,"use Bb, BN=%D, Bbn=%D; numBb=%D\n",BN,Bbn,numBb);
440: if (Bbn1) { /* Create workB1 for the remaining columns */
441: PetscInfo2(C,"use Bb1, BN=%D, Bbn1=%D\n",BN,Bbn1);
442: /* Create work matrix used to store off processor rows of B needed for local product */
443: MatCreateSeqDense(PETSC_COMM_SELF,nz,Bbn1,NULL,&contents->workB1);
444: } else contents->workB1 = NULL;
445: }
447: /* Create work matrix used to store off processor rows of B needed for local product */
448: MatCreateSeqDense(PETSC_COMM_SELF,nz,Bbn,NULL,&contents->workB);
450: /* Use MPI derived data type to reduce memory required by the send/recv buffers */
451: PetscMalloc4(nsends,&stype,nrecvs,&rtype,nrecvs,&contents->rwaits,nsends,&contents->swaits);
452: contents->stype = stype;
453: contents->nsends = nsends;
455: contents->rtype = rtype;
456: contents->nrecvs = nrecvs;
457: contents->blda = blda;
459: PetscMalloc1(Bm+1,&disp);
460: for (i=0; i<nsends; i++) {
461: nrows_to = sstarts[i+1]-sstarts[i];
462: for (j=0; j<nrows_to; j++){
463: disp[j] = sindices[sstarts[i]+j]; /* rowB to be sent */
464: }
465: MPI_Type_create_indexed_block(nrows_to,1,(const PetscMPIInt *)disp,MPIU_SCALAR,&type1);
467: MPI_Type_create_resized(type1,0,blda*sizeof(PetscScalar),&stype[i]);
468: MPI_Type_commit(&stype[i]);
469: MPI_Type_free(&type1);
470: }
472: for (i=0; i<nrecvs; i++) {
473: /* received values from a process form a (nrows_from x Bbn) row block in workB (column-wise) */
474: nrows_from = rstarts[i+1]-rstarts[i];
475: disp[0] = 0;
476: MPI_Type_create_indexed_block(1, nrows_from, (const PetscMPIInt *)disp, MPIU_SCALAR, &type1);
477: MPI_Type_create_resized(type1, 0, nz*sizeof(PetscScalar), &rtype[i]);
478: MPI_Type_commit(&rtype[i]);
479: MPI_Type_free(&type1);
480: }
482: PetscFree(disp);
483: VecScatterRestoreRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,NULL,NULL);
484: VecScatterRestoreRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,NULL,NULL);
485: MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);
486: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
487: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
488: MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
490: C->product->data = contents;
491: C->product->destroy = MatMPIAIJ_MPIDenseDestroy;
492: C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIDense;
493: return(0);
494: }
496: PETSC_INTERN PetscErrorCode MatMatMultNumericAdd_SeqAIJ_SeqDense(Mat,Mat,Mat,const PetscBool);
497: /*
498: Performs an efficient scatter on the rows of B needed by this process; this is
499: a modification of the VecScatterBegin_() routines.
501: Input: Bbidx = 0: B = Bb
502: = 1: B = Bb1, see MatMatMultSymbolic_MPIAIJ_MPIDense()
503: */
504: PetscErrorCode MatMPIDenseScatter(Mat A,Mat B,PetscInt Bbidx,Mat C,Mat *outworkB)
505: {
506: Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data;
507: PetscErrorCode ierr;
508: const PetscScalar *b;
509: PetscScalar *rvalues;
510: VecScatter ctx = aij->Mvctx;
511: const PetscInt *sindices,*sstarts,*rstarts;
512: const PetscMPIInt *sprocs,*rprocs;
513: PetscInt i,nsends,nrecvs;
514: MPI_Request *swaits,*rwaits;
515: MPI_Comm comm;
516: PetscMPIInt tag=((PetscObject)ctx)->tag,ncols=B->cmap->N,nrows=aij->B->cmap->n,nsends_mpi,nrecvs_mpi;
517: MPIAIJ_MPIDense *contents;
518: Mat workB;
519: MPI_Datatype *stype,*rtype;
520: PetscInt blda;
523: MatCheckProduct(C,4);
524: if (!C->product->data) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty");
525: contents = (MPIAIJ_MPIDense*)C->product->data;
526: VecScatterGetRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,&sprocs,NULL/*bs*/);
527: VecScatterGetRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,&rprocs,NULL/*bs*/);
528: PetscMPIIntCast(nsends,&nsends_mpi);
529: PetscMPIIntCast(nrecvs,&nrecvs_mpi);
530: if (Bbidx == 0) {
531: workB = *outworkB = contents->workB;
532: } else {
533: workB = *outworkB = contents->workB1;
534: }
535: if (nrows != workB->rmap->n) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Number of rows of workB %D not equal to columns of aij->B %D",workB->cmap->n,nrows);
536: swaits = contents->swaits;
537: rwaits = contents->rwaits;
539: MatDenseGetArrayRead(B,&b);
540: MatDenseGetLDA(B,&blda);
541: if (blda != contents->blda) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Cannot reuse an input matrix with lda %D != %D",blda,contents->blda);
542: MatDenseGetArray(workB,&rvalues);
544: /* Post recv, use MPI derived data type to save memory */
545: PetscObjectGetComm((PetscObject)C,&comm);
546: rtype = contents->rtype;
547: for (i=0; i<nrecvs; i++) {
548: MPI_Irecv(rvalues+(rstarts[i]-rstarts[0]),ncols,rtype[i],rprocs[i],tag,comm,rwaits+i);
549: }
551: stype = contents->stype;
552: for (i=0; i<nsends; i++) {
553: MPI_Isend(b,ncols,stype[i],sprocs[i],tag,comm,swaits+i);
554: }
556: if (nrecvs) {MPI_Waitall(nrecvs_mpi,rwaits,MPI_STATUSES_IGNORE);}
557: if (nsends) {MPI_Waitall(nsends_mpi,swaits,MPI_STATUSES_IGNORE);}
559: VecScatterRestoreRemote_Private(ctx,PETSC_TRUE/*send*/,&nsends,&sstarts,&sindices,&sprocs,NULL);
560: VecScatterRestoreRemoteOrdered_Private(ctx,PETSC_FALSE/*recv*/,&nrecvs,&rstarts,NULL,&rprocs,NULL);
561: MatDenseRestoreArrayRead(B,&b);
562: MatDenseRestoreArray(workB,&rvalues);
563: return(0);
564: }
566: static PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIDense(Mat A,Mat B,Mat C)
567: {
568: PetscErrorCode ierr;
569: Mat_MPIAIJ *aij = (Mat_MPIAIJ*)A->data;
570: Mat_MPIDense *bdense = (Mat_MPIDense*)B->data;
571: Mat_MPIDense *cdense = (Mat_MPIDense*)C->data;
572: Mat workB;
573: MPIAIJ_MPIDense *contents;
576: MatCheckProduct(C,3);
577: if (!C->product->data) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_PLIB,"Product data empty");
578: contents = (MPIAIJ_MPIDense*)C->product->data;
579: /* diagonal block of A times all local rows of B */
580: /* TODO: this calls a symbolic multiplication every time, which could be avoided */
581: MatMatMult(aij->A,bdense->A,MAT_REUSE_MATRIX,PETSC_DEFAULT,&cdense->A);
582: if (contents->workB->cmap->n == B->cmap->N) {
583: /* get off processor parts of B needed to complete C=A*B */
584: MatMPIDenseScatter(A,B,0,C,&workB);
586: /* off-diagonal block of A times nonlocal rows of B */
587: MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A,PETSC_TRUE);
588: } else {
589: Mat Bb,Cb;
590: PetscInt BN=B->cmap->N,n=contents->workB->cmap->n,i;
592: for (i=0; i<BN; i+=n) {
593: MatDenseGetSubMatrix(B,i,PetscMin(i+n,BN),&Bb);
594: MatDenseGetSubMatrix(C,i,PetscMin(i+n,BN),&Cb);
596: /* get off processor parts of B needed to complete C=A*B */
597: MatMPIDenseScatter(A,Bb,i+n>BN,C,&workB);
599: /* off-diagonal block of A times nonlocal rows of B */
600: cdense = (Mat_MPIDense*)Cb->data;
601: MatMatMultNumericAdd_SeqAIJ_SeqDense(aij->B,workB,cdense->A,PETSC_TRUE);
603: MatDenseRestoreSubMatrix(B,&Bb);
604: MatDenseRestoreSubMatrix(C,&Cb);
605: }
606: }
607: return(0);
608: }
610: PetscErrorCode MatMatMultNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
611: {
613: Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data,*c=(Mat_MPIAIJ*)C->data;
614: Mat_SeqAIJ *ad = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
615: Mat_SeqAIJ *cd = (Mat_SeqAIJ*)(c->A)->data,*co=(Mat_SeqAIJ*)(c->B)->data;
616: PetscInt *adi = ad->i,*adj,*aoi=ao->i,*aoj;
617: PetscScalar *ada,*aoa,*cda=cd->a,*coa=co->a;
618: Mat_SeqAIJ *p_loc,*p_oth;
619: PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pj;
620: PetscScalar *pa_loc,*pa_oth,*pa,valtmp,*ca;
621: PetscInt cm = C->rmap->n,anz,pnz;
622: Mat_APMPI *ptap;
623: PetscScalar *apa_sparse;
624: PetscInt *api,*apj,*apJ,i,j,k,row;
625: PetscInt cstart = C->cmap->rstart;
626: PetscInt cdnz,conz,k0,k1,nextp;
627: MPI_Comm comm;
628: PetscMPIInt size;
631: MatCheckProduct(C,3);
632: ptap = (Mat_APMPI*)C->product->data;
633: if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data");
634: PetscObjectGetComm((PetscObject)C,&comm);
635: MPI_Comm_size(comm,&size);
636: if (!ptap->P_oth && size>1) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"AP cannot be reused. Do not call MatProductClear()");
638: apa_sparse = ptap->apa;
640: /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */
641: /*-----------------------------------------------------*/
642: /* update numerical values of P_oth and P_loc */
643: MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);
644: MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);
646: /* 2) compute numeric C_loc = A_loc*P = Ad*P_loc + Ao*P_oth */
647: /*----------------------------------------------------------*/
648: /* get data from symbolic products */
649: p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
650: pi_loc = p_loc->i; pj_loc = p_loc->j; pa_loc = p_loc->a;
651: if (size >1) {
652: p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
653: pi_oth = p_oth->i; pj_oth = p_oth->j; pa_oth = p_oth->a;
654: } else {
655: p_oth = NULL; pi_oth = NULL; pj_oth = NULL; pa_oth = NULL;
656: }
658: api = ptap->api;
659: apj = ptap->apj;
660: for (i=0; i<cm; i++) {
661: apJ = apj + api[i];
663: /* diagonal portion of A */
664: anz = adi[i+1] - adi[i];
665: adj = ad->j + adi[i];
666: ada = ad->a + adi[i];
667: for (j=0; j<anz; j++) {
668: row = adj[j];
669: pnz = pi_loc[row+1] - pi_loc[row];
670: pj = pj_loc + pi_loc[row];
671: pa = pa_loc + pi_loc[row];
672: /* perform sparse axpy */
673: valtmp = ada[j];
674: nextp = 0;
675: for (k=0; nextp<pnz; k++) {
676: if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
677: apa_sparse[k] += valtmp*pa[nextp++];
678: }
679: }
680: PetscLogFlops(2.0*pnz);
681: }
683: /* off-diagonal portion of A */
684: anz = aoi[i+1] - aoi[i];
685: aoj = ao->j + aoi[i];
686: aoa = ao->a + aoi[i];
687: for (j=0; j<anz; j++) {
688: row = aoj[j];
689: pnz = pi_oth[row+1] - pi_oth[row];
690: pj = pj_oth + pi_oth[row];
691: pa = pa_oth + pi_oth[row];
692: /* perform sparse axpy */
693: valtmp = aoa[j];
694: nextp = 0;
695: for (k=0; nextp<pnz; k++) {
696: if (apJ[k] == pj[nextp]) { /* column of AP == column of P */
697: apa_sparse[k] += valtmp*pa[nextp++];
698: }
699: }
700: PetscLogFlops(2.0*pnz);
701: }
703: /* set values in C */
704: cdnz = cd->i[i+1] - cd->i[i];
705: conz = co->i[i+1] - co->i[i];
707: /* 1st off-diagonal part of C */
708: ca = coa + co->i[i];
709: k = 0;
710: for (k0=0; k0<conz; k0++) {
711: if (apJ[k] >= cstart) break;
712: ca[k0] = apa_sparse[k];
713: apa_sparse[k] = 0.0;
714: k++;
715: }
717: /* diagonal part of C */
718: ca = cda + cd->i[i];
719: for (k1=0; k1<cdnz; k1++) {
720: ca[k1] = apa_sparse[k];
721: apa_sparse[k] = 0.0;
722: k++;
723: }
725: /* 2nd off-diagonal part of C */
726: ca = coa + co->i[i];
727: for (; k0<conz; k0++) {
728: ca[k0] = apa_sparse[k];
729: apa_sparse[k] = 0.0;
730: k++;
731: }
732: }
733: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
734: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
735: return(0);
736: }
738: /* same as MatMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(), except using LLCondensed to avoid O(BN) memory requirement */
739: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat C)
740: {
741: PetscErrorCode ierr;
742: MPI_Comm comm;
743: PetscMPIInt size;
744: Mat_APMPI *ptap;
745: PetscFreeSpaceList free_space = NULL,current_space=NULL;
746: Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data;
747: Mat_SeqAIJ *ad = (Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc,*p_oth;
748: PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*dnz,*onz;
749: PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,rstart=A->rmap->rstart;
750: PetscInt i,pnz,row,*api,*apj,*Jptr,apnz,nspacedouble=0,j,nzi,*lnk,apnz_max=0;
751: PetscInt am=A->rmap->n,pn=P->cmap->n,pm=P->rmap->n,lsize=pn+20;
752: PetscReal afill;
753: MatType mtype;
756: MatCheckProduct(C,4);
757: if (C->product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty");
758: PetscObjectGetComm((PetscObject)A,&comm);
759: MPI_Comm_size(comm,&size);
761: /* create struct Mat_APMPI and attached it to C later */
762: PetscNew(&ptap);
764: /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
765: MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);
767: /* get P_loc by taking all local rows of P */
768: MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);
770: p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
771: pi_loc = p_loc->i; pj_loc = p_loc->j;
772: if (size > 1) {
773: p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
774: pi_oth = p_oth->i; pj_oth = p_oth->j;
775: } else {
776: p_oth = NULL;
777: pi_oth = NULL; pj_oth = NULL;
778: }
780: /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
781: /*-------------------------------------------------------------------*/
782: PetscMalloc1(am+2,&api);
783: ptap->api = api;
784: api[0] = 0;
786: PetscLLCondensedCreate_Scalable(lsize,&lnk);
788: /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
789: PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);
790: current_space = free_space;
791: MatPreallocateInitialize(comm,am,pn,dnz,onz);
792: for (i=0; i<am; i++) {
793: /* diagonal portion of A */
794: nzi = adi[i+1] - adi[i];
795: for (j=0; j<nzi; j++) {
796: row = *adj++;
797: pnz = pi_loc[row+1] - pi_loc[row];
798: Jptr = pj_loc + pi_loc[row];
799: /* Expand list if it is not long enough */
800: if (pnz+apnz_max > lsize) {
801: lsize = pnz+apnz_max;
802: PetscLLCondensedExpand_Scalable(lsize, &lnk);
803: }
804: /* add non-zero cols of P into the sorted linked list lnk */
805: PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);
806: apnz = *lnk; /* The first element in the list is the number of items in the list */
807: api[i+1] = api[i] + apnz;
808: if (apnz > apnz_max) apnz_max = apnz;
809: }
810: /* off-diagonal portion of A */
811: nzi = aoi[i+1] - aoi[i];
812: for (j=0; j<nzi; j++) {
813: row = *aoj++;
814: pnz = pi_oth[row+1] - pi_oth[row];
815: Jptr = pj_oth + pi_oth[row];
816: /* Expand list if it is not long enough */
817: if (pnz+apnz_max > lsize) {
818: lsize = pnz + apnz_max;
819: PetscLLCondensedExpand_Scalable(lsize, &lnk);
820: }
821: /* add non-zero cols of P into the sorted linked list lnk */
822: PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);
823: apnz = *lnk; /* The first element in the list is the number of items in the list */
824: api[i+1] = api[i] + apnz;
825: if (apnz > apnz_max) apnz_max = apnz;
826: }
827: apnz = *lnk;
828: api[i+1] = api[i] + apnz;
829: if (apnz > apnz_max) apnz_max = apnz;
831: /* if free space is not available, double the total space in the list */
832: if (current_space->local_remaining<apnz) {
833: PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);
834: nspacedouble++;
835: }
837: /* Copy data into free space, then initialize lnk */
838: PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);
839: MatPreallocateSet(i+rstart,apnz,current_space->array,dnz,onz);
841: current_space->array += apnz;
842: current_space->local_used += apnz;
843: current_space->local_remaining -= apnz;
844: }
846: /* Allocate space for apj, initialize apj, and */
847: /* destroy list of free space and other temporary array(s) */
848: PetscMalloc1(api[am]+1,&ptap->apj);
849: apj = ptap->apj;
850: PetscFreeSpaceContiguous(&free_space,ptap->apj);
851: PetscLLCondensedDestroy_Scalable(lnk);
853: /* create and assemble symbolic parallel matrix C */
854: /*----------------------------------------------------*/
855: MatSetSizes(C,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);
856: MatSetBlockSizesFromMats(C,A,P);
857: MatGetType(A,&mtype);
858: MatSetType(C,mtype);
859: MatMPIAIJSetPreallocation(C,0,dnz,0,onz);
860: MatPreallocateFinalize(dnz,onz);
862: /* malloc apa for assembly C */
863: PetscCalloc1(apnz_max,&ptap->apa);
865: MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(C, apj, api);
866: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
867: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
868: MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
870: C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ;
871: C->ops->productnumeric = MatProductNumeric_AB;
873: /* attach the supporting struct to C for reuse */
874: C->product->data = ptap;
875: C->product->destroy = MatDestroy_MPIAIJ_MatMatMult;
877: /* set MatInfo */
878: afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
879: if (afill < 1.0) afill = 1.0;
880: C->info.mallocs = nspacedouble;
881: C->info.fill_ratio_given = fill;
882: C->info.fill_ratio_needed = afill;
884: #if defined(PETSC_USE_INFO)
885: if (api[am]) {
886: PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
887: PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
888: } else {
889: PetscInfo(C,"Empty matrix product\n");
890: }
891: #endif
892: return(0);
893: }
895: /* This function is needed for the seqMPI matrix-matrix multiplication. */
896: /* Three input arrays are merged to one output array. The size of the */
897: /* output array is also output. Duplicate entries only show up once. */
898: static void Merge3SortedArrays(PetscInt size1, PetscInt *in1,
899: PetscInt size2, PetscInt *in2,
900: PetscInt size3, PetscInt *in3,
901: PetscInt *size4, PetscInt *out)
902: {
903: int i = 0, j = 0, k = 0, l = 0;
905: /* Traverse all three arrays */
906: while (i<size1 && j<size2 && k<size3) {
907: if (in1[i] < in2[j] && in1[i] < in3[k]) {
908: out[l++] = in1[i++];
909: }
910: else if (in2[j] < in1[i] && in2[j] < in3[k]) {
911: out[l++] = in2[j++];
912: }
913: else if (in3[k] < in1[i] && in3[k] < in2[j]) {
914: out[l++] = in3[k++];
915: }
916: else if (in1[i] == in2[j] && in1[i] < in3[k]) {
917: out[l++] = in1[i];
918: i++, j++;
919: }
920: else if (in1[i] == in3[k] && in1[i] < in2[j]) {
921: out[l++] = in1[i];
922: i++, k++;
923: }
924: else if (in3[k] == in2[j] && in2[j] < in1[i]) {
925: out[l++] = in2[j];
926: k++, j++;
927: }
928: else if (in1[i] == in2[j] && in1[i] == in3[k]) {
929: out[l++] = in1[i];
930: i++, j++, k++;
931: }
932: }
934: /* Traverse two remaining arrays */
935: while (i<size1 && j<size2) {
936: if (in1[i] < in2[j]) {
937: out[l++] = in1[i++];
938: }
939: else if (in1[i] > in2[j]) {
940: out[l++] = in2[j++];
941: }
942: else {
943: out[l++] = in1[i];
944: i++, j++;
945: }
946: }
948: while (i<size1 && k<size3) {
949: if (in1[i] < in3[k]) {
950: out[l++] = in1[i++];
951: }
952: else if (in1[i] > in3[k]) {
953: out[l++] = in3[k++];
954: }
955: else {
956: out[l++] = in1[i];
957: i++, k++;
958: }
959: }
961: while (k<size3 && j<size2) {
962: if (in3[k] < in2[j]) {
963: out[l++] = in3[k++];
964: }
965: else if (in3[k] > in2[j]) {
966: out[l++] = in2[j++];
967: }
968: else {
969: out[l++] = in3[k];
970: k++, j++;
971: }
972: }
974: /* Traverse one remaining array */
975: while (i<size1) out[l++] = in1[i++];
976: while (j<size2) out[l++] = in2[j++];
977: while (k<size3) out[l++] = in3[k++];
979: *size4 = l;
980: }
982: /* This matrix-matrix multiplication algorithm divides the multiplication into three multiplications and */
983: /* adds up the products. Two of these three multiplications are performed with existing (sequential) */
984: /* matrix-matrix multiplications. */
985: PetscErrorCode MatMatMultSymbolic_MPIAIJ_MPIAIJ_seqMPI(Mat A, Mat P, PetscReal fill, Mat C)
986: {
987: PetscErrorCode ierr;
988: MPI_Comm comm;
989: PetscMPIInt size;
990: Mat_APMPI *ptap;
991: PetscFreeSpaceList free_space_diag=NULL, current_space=NULL;
992: Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data;
993: Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_loc;
994: Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data;
995: Mat_SeqAIJ *adpd_seq, *p_off, *aopoth_seq;
996: PetscInt adponz, adpdnz;
997: PetscInt *pi_loc,*dnz,*onz;
998: PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,rstart=A->rmap->rstart;
999: PetscInt *lnk,i, i1=0,pnz,row,*adpoi,*adpoj, *api, *adpoJ, *aopJ, *apJ,*Jptr, aopnz, nspacedouble=0,j,nzi,
1000: *apj,apnz, *adpdi, *adpdj, *adpdJ, *poff_i, *poff_j, *j_temp, *aopothi, *aopothj;
1001: PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n,pm=P->rmap->n, p_colstart, p_colend;
1002: PetscBT lnkbt;
1003: PetscReal afill;
1004: PetscMPIInt rank;
1005: Mat adpd, aopoth;
1006: MatType mtype;
1007: const char *prefix;
1010: MatCheckProduct(C,4);
1011: if (C->product->data) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Extra product struct not empty");
1012: PetscObjectGetComm((PetscObject)A,&comm);
1013: MPI_Comm_size(comm,&size);
1014: MPI_Comm_rank(comm, &rank);
1015: MatGetOwnershipRangeColumn(P, &p_colstart, &p_colend);
1017: /* create struct Mat_APMPI and attached it to C later */
1018: PetscNew(&ptap);
1020: /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
1021: MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);
1023: /* get P_loc by taking all local rows of P */
1024: MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);
1027: p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
1028: pi_loc = p_loc->i;
1030: /* Allocate memory for the i arrays of the matrices A*P, A_diag*P_off and A_offd * P */
1031: PetscMalloc1(am+2,&api);
1032: PetscMalloc1(am+2,&adpoi);
1034: adpoi[0] = 0;
1035: ptap->api = api;
1036: api[0] = 0;
1038: /* create and initialize a linked list, will be used for both A_diag * P_loc_off and A_offd * P_oth */
1039: PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);
1040: MatPreallocateInitialize(comm,am,pn,dnz,onz);
1042: /* Symbolic calc of A_loc_diag * P_loc_diag */
1043: MatGetOptionsPrefix(A,&prefix);
1044: MatProductCreate(a->A,p->A,NULL,&adpd);
1045: MatGetOptionsPrefix(A,&prefix);
1046: MatSetOptionsPrefix(adpd,prefix);
1047: MatAppendOptionsPrefix(adpd,"inner_diag_");
1049: MatProductSetType(adpd,MATPRODUCT_AB);
1050: MatProductSetAlgorithm(adpd,"sorted");
1051: MatProductSetFill(adpd,fill);
1052: MatProductSetFromOptions(adpd);
1053: MatProductSymbolic(adpd);
1055: adpd_seq = (Mat_SeqAIJ*)((adpd)->data);
1056: adpdi = adpd_seq->i; adpdj = adpd_seq->j;
1057: p_off = (Mat_SeqAIJ*)((p->B)->data);
1058: poff_i = p_off->i; poff_j = p_off->j;
1060: /* j_temp stores indices of a result row before they are added to the linked list */
1061: PetscMalloc1(pN+2,&j_temp);
1064: /* Symbolic calc of the A_diag * p_loc_off */
1065: /* Initial FreeSpace size is fill*(nnz(A)+nnz(P)) */
1066: PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space_diag);
1067: current_space = free_space_diag;
1069: for (i=0; i<am; i++) {
1070: /* A_diag * P_loc_off */
1071: nzi = adi[i+1] - adi[i];
1072: for (j=0; j<nzi; j++) {
1073: row = *adj++;
1074: pnz = poff_i[row+1] - poff_i[row];
1075: Jptr = poff_j + poff_i[row];
1076: for (i1 = 0; i1 < pnz; i1++) {
1077: j_temp[i1] = p->garray[Jptr[i1]];
1078: }
1079: /* add non-zero cols of P into the sorted linked list lnk */
1080: PetscLLCondensedAddSorted(pnz,j_temp,lnk,lnkbt);
1081: }
1083: adponz = lnk[0];
1084: adpoi[i+1] = adpoi[i] + adponz;
1086: /* if free space is not available, double the total space in the list */
1087: if (current_space->local_remaining<adponz) {
1088: PetscFreeSpaceGet(PetscIntSumTruncate(adponz,current_space->total_array_size),¤t_space);
1089: nspacedouble++;
1090: }
1092: /* Copy data into free space, then initialize lnk */
1093: PetscLLCondensedClean(pN,adponz,current_space->array,lnk,lnkbt);
1095: current_space->array += adponz;
1096: current_space->local_used += adponz;
1097: current_space->local_remaining -= adponz;
1098: }
1100: /* Symbolic calc of A_off * P_oth */
1101: MatSetOptionsPrefix(a->B,prefix);
1102: MatAppendOptionsPrefix(a->B,"inner_offdiag_");
1103: MatCreate(PETSC_COMM_SELF,&aopoth);
1104: MatMatMultSymbolic_SeqAIJ_SeqAIJ(a->B, ptap->P_oth, fill, aopoth);
1105: aopoth_seq = (Mat_SeqAIJ*)((aopoth)->data);
1106: aopothi = aopoth_seq->i; aopothj = aopoth_seq->j;
1108: /* Allocate space for apj, adpj, aopj, ... */
1109: /* destroy lists of free space and other temporary array(s) */
1111: PetscMalloc1(aopothi[am] + adpoi[am] + adpdi[am]+2, &ptap->apj);
1112: PetscMalloc1(adpoi[am]+2, &adpoj);
1114: /* Copy from linked list to j-array */
1115: PetscFreeSpaceContiguous(&free_space_diag,adpoj);
1116: PetscLLDestroy(lnk,lnkbt);
1118: adpoJ = adpoj;
1119: adpdJ = adpdj;
1120: aopJ = aopothj;
1121: apj = ptap->apj;
1122: apJ = apj; /* still empty */
1124: /* Merge j-arrays of A_off * P, A_diag * P_loc_off, and */
1125: /* A_diag * P_loc_diag to get A*P */
1126: for (i = 0; i < am; i++) {
1127: aopnz = aopothi[i+1] - aopothi[i];
1128: adponz = adpoi[i+1] - adpoi[i];
1129: adpdnz = adpdi[i+1] - adpdi[i];
1131: /* Correct indices from A_diag*P_diag */
1132: for (i1 = 0; i1 < adpdnz; i1++) {
1133: adpdJ[i1] += p_colstart;
1134: }
1135: /* Merge j-arrays of A_diag * P_loc_off and A_diag * P_loc_diag and A_off * P_oth */
1136: Merge3SortedArrays(adponz, adpoJ, adpdnz, adpdJ, aopnz, aopJ, &apnz, apJ);
1137: MatPreallocateSet(i+rstart, apnz, apJ, dnz, onz);
1139: aopJ += aopnz;
1140: adpoJ += adponz;
1141: adpdJ += adpdnz;
1142: apJ += apnz;
1143: api[i+1] = api[i] + apnz;
1144: }
1146: /* malloc apa to store dense row A[i,:]*P */
1147: PetscCalloc1(pN+2,&ptap->apa);
1149: /* create and assemble symbolic parallel matrix C */
1150: MatSetSizes(C,am,pn,PETSC_DETERMINE,PETSC_DETERMINE);
1151: MatSetBlockSizesFromMats(C,A,P);
1152: MatGetType(A,&mtype);
1153: MatSetType(C,mtype);
1154: MatMPIAIJSetPreallocation(C,0,dnz,0,onz);
1155: MatPreallocateFinalize(dnz,onz);
1157: MatSetValues_MPIAIJ_CopyFromCSRFormat_Symbolic(C, apj, api);
1158: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1159: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1160: MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
1162: C->ops->matmultnumeric = MatMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1163: C->ops->productnumeric = MatProductNumeric_AB;
1165: /* attach the supporting struct to C for reuse */
1166: C->product->data = ptap;
1167: C->product->destroy = MatDestroy_MPIAIJ_MatMatMult;
1169: /* set MatInfo */
1170: afill = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1) + 1.e-5;
1171: if (afill < 1.0) afill = 1.0;
1172: C->info.mallocs = nspacedouble;
1173: C->info.fill_ratio_given = fill;
1174: C->info.fill_ratio_needed = afill;
1176: #if defined(PETSC_USE_INFO)
1177: if (api[am]) {
1178: PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
1179: PetscInfo1(C,"Use MatMatMult(A,B,MatReuse,%g,&C) for best performance.;\n",(double)afill);
1180: } else {
1181: PetscInfo(C,"Empty matrix product\n");
1182: }
1183: #endif
1185: MatDestroy(&aopoth);
1186: MatDestroy(&adpd);
1187: PetscFree(j_temp);
1188: PetscFree(adpoj);
1189: PetscFree(adpoi);
1190: return(0);
1191: }
1193: /*-------------------------------------------------------------------------*/
1194: /* This routine only works when scall=MAT_REUSE_MATRIX! */
1195: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult(Mat P,Mat A,Mat C)
1196: {
1198: Mat_APMPI *ptap;
1199: Mat Pt;
1202: MatCheckProduct(C,3);
1203: ptap = (Mat_APMPI*)C->product->data;
1204: if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data");
1205: if (!ptap->Pt) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatProductClear()");
1207: Pt = ptap->Pt;
1208: MatTranspose(P,MAT_REUSE_MATRIX,&Pt);
1209: MatMatMultNumeric_MPIAIJ_MPIAIJ(Pt,A,C);
1210: return(0);
1211: }
1213: /* This routine is modified from MatPtAPSymbolic_MPIAIJ_MPIAIJ() */
1214: PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,PetscReal fill,Mat C)
1215: {
1216: PetscErrorCode ierr;
1217: Mat_APMPI *ptap;
1218: Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data;
1219: MPI_Comm comm;
1220: PetscMPIInt size,rank;
1221: PetscFreeSpaceList free_space=NULL,current_space=NULL;
1222: PetscInt pn=P->cmap->n,aN=A->cmap->N,an=A->cmap->n;
1223: PetscInt *lnk,i,k,nsend,rstart;
1224: PetscBT lnkbt;
1225: PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,nrecv;
1226: PETSC_UNUSED PetscMPIInt icompleted=0;
1227: PetscInt **buf_rj,**buf_ri,**buf_ri_k;
1228: PetscInt len,proc,*dnz,*onz,*owners,nzi;
1229: PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1230: MPI_Request *swaits,*rwaits;
1231: MPI_Status *sstatus,rstatus;
1232: PetscLayout rowmap;
1233: PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */
1234: PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */
1235: PetscInt *Jptr,*prmap=p->garray,con,j,Crmax;
1236: Mat_SeqAIJ *a_loc,*c_loc,*c_oth;
1237: PetscTable ta;
1238: MatType mtype;
1239: const char *prefix;
1242: PetscObjectGetComm((PetscObject)A,&comm);
1243: MPI_Comm_size(comm,&size);
1244: MPI_Comm_rank(comm,&rank);
1246: /* create symbolic parallel matrix C */
1247: MatGetType(A,&mtype);
1248: MatSetType(C,mtype);
1250: C->ops->transposematmultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable;
1252: /* create struct Mat_APMPI and attached it to C later */
1253: PetscNew(&ptap);
1254: ptap->reuse = MAT_INITIAL_MATRIX;
1256: /* (0) compute Rd = Pd^T, Ro = Po^T */
1257: /* --------------------------------- */
1258: MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);
1259: MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);
1261: /* (1) compute symbolic A_loc */
1262: /* ---------------------------*/
1263: MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&ptap->A_loc);
1265: /* (2-1) compute symbolic C_oth = Ro*A_loc */
1266: /* ------------------------------------ */
1267: MatGetOptionsPrefix(A,&prefix);
1268: MatSetOptionsPrefix(ptap->Ro,prefix);
1269: MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");
1270: MatCreate(PETSC_COMM_SELF,&ptap->C_oth);
1271: MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->A_loc,fill,ptap->C_oth);
1273: /* (3) send coj of C_oth to other processors */
1274: /* ------------------------------------------ */
1275: /* determine row ownership */
1276: PetscLayoutCreate(comm,&rowmap);
1277: rowmap->n = pn;
1278: rowmap->bs = 1;
1279: PetscLayoutSetUp(rowmap);
1280: owners = rowmap->range;
1282: /* determine the number of messages to send, their lengths */
1283: PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);
1284: PetscArrayzero(len_s,size);
1285: PetscArrayzero(len_si,size);
1287: c_oth = (Mat_SeqAIJ*)ptap->C_oth->data;
1288: coi = c_oth->i; coj = c_oth->j;
1289: con = ptap->C_oth->rmap->n;
1290: proc = 0;
1291: for (i=0; i<con; i++) {
1292: while (prmap[i] >= owners[proc+1]) proc++;
1293: len_si[proc]++; /* num of rows in Co(=Pt*A) to be sent to [proc] */
1294: len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */
1295: }
1297: len = 0; /* max length of buf_si[], see (4) */
1298: owners_co[0] = 0;
1299: nsend = 0;
1300: for (proc=0; proc<size; proc++) {
1301: owners_co[proc+1] = owners_co[proc] + len_si[proc];
1302: if (len_s[proc]) {
1303: nsend++;
1304: len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */
1305: len += len_si[proc];
1306: }
1307: }
1309: /* determine the number and length of messages to receive for coi and coj */
1310: PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);
1311: PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);
1313: /* post the Irecv and Isend of coj */
1314: PetscCommGetNewTag(comm,&tagj);
1315: PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);
1316: PetscMalloc1(nsend+1,&swaits);
1317: for (proc=0, k=0; proc<size; proc++) {
1318: if (!len_s[proc]) continue;
1319: i = owners_co[proc];
1320: MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);
1321: k++;
1322: }
1324: /* (2-2) compute symbolic C_loc = Rd*A_loc */
1325: /* ---------------------------------------- */
1326: MatSetOptionsPrefix(ptap->Rd,prefix);
1327: MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");
1328: MatCreate(PETSC_COMM_SELF,&ptap->C_loc);
1329: MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->A_loc,fill,ptap->C_loc);
1330: c_loc = (Mat_SeqAIJ*)ptap->C_loc->data;
1332: /* receives coj are complete */
1333: for (i=0; i<nrecv; i++) {
1334: MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);
1335: }
1336: PetscFree(rwaits);
1337: if (nsend) {MPI_Waitall(nsend,swaits,sstatus);}
1339: /* add received column indices into ta to update Crmax */
1340: a_loc = (Mat_SeqAIJ*)(ptap->A_loc)->data;
1342: /* create and initialize a linked list */
1343: PetscTableCreate(an,aN,&ta); /* for compute Crmax */
1344: MatRowMergeMax_SeqAIJ(a_loc,ptap->A_loc->rmap->N,ta);
1346: for (k=0; k<nrecv; k++) {/* k-th received message */
1347: Jptr = buf_rj[k];
1348: for (j=0; j<len_r[k]; j++) {
1349: PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);
1350: }
1351: }
1352: PetscTableGetCount(ta,&Crmax);
1353: PetscTableDestroy(&ta);
1355: /* (4) send and recv coi */
1356: /*-----------------------*/
1357: PetscCommGetNewTag(comm,&tagi);
1358: PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);
1359: PetscMalloc1(len+1,&buf_s);
1360: buf_si = buf_s; /* points to the beginning of k-th msg to be sent */
1361: for (proc=0,k=0; proc<size; proc++) {
1362: if (!len_s[proc]) continue;
1363: /* form outgoing message for i-structure:
1364: buf_si[0]: nrows to be sent
1365: [1:nrows]: row index (global)
1366: [nrows+1:2*nrows+1]: i-structure index
1367: */
1368: /*-------------------------------------------*/
1369: nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */
1370: buf_si_i = buf_si + nrows+1;
1371: buf_si[0] = nrows;
1372: buf_si_i[0] = 0;
1373: nrows = 0;
1374: for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1375: nzi = coi[i+1] - coi[i];
1376: buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1377: buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */
1378: nrows++;
1379: }
1380: MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);
1381: k++;
1382: buf_si += len_si[proc];
1383: }
1384: for (i=0; i<nrecv; i++) {
1385: MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);
1386: }
1387: PetscFree(rwaits);
1388: if (nsend) {MPI_Waitall(nsend,swaits,sstatus);}
1390: PetscFree4(len_s,len_si,sstatus,owners_co);
1391: PetscFree(len_ri);
1392: PetscFree(swaits);
1393: PetscFree(buf_s);
1395: /* (5) compute the local portion of C */
1396: /* ------------------------------------------ */
1397: /* set initial free space to be Crmax, sufficient for holding nozeros in each row of C */
1398: PetscFreeSpaceGet(Crmax,&free_space);
1399: current_space = free_space;
1401: PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);
1402: for (k=0; k<nrecv; k++) {
1403: buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1404: nrows = *buf_ri_k[k];
1405: nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */
1406: nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */
1407: }
1409: MatPreallocateInitialize(comm,pn,an,dnz,onz);
1410: PetscLLCondensedCreate(Crmax,aN,&lnk,&lnkbt);
1411: for (i=0; i<pn; i++) {
1412: /* add C_loc into C */
1413: nzi = c_loc->i[i+1] - c_loc->i[i];
1414: Jptr = c_loc->j + c_loc->i[i];
1415: PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);
1417: /* add received col data into lnk */
1418: for (k=0; k<nrecv; k++) { /* k-th received message */
1419: if (i == *nextrow[k]) { /* i-th row */
1420: nzi = *(nextci[k]+1) - *nextci[k];
1421: Jptr = buf_rj[k] + *nextci[k];
1422: PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);
1423: nextrow[k]++; nextci[k]++;
1424: }
1425: }
1426: nzi = lnk[0];
1428: /* copy data into free space, then initialize lnk */
1429: PetscLLCondensedClean(aN,nzi,current_space->array,lnk,lnkbt);
1430: MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);
1431: }
1432: PetscFree3(buf_ri_k,nextrow,nextci);
1433: PetscLLDestroy(lnk,lnkbt);
1434: PetscFreeSpaceDestroy(free_space);
1436: /* local sizes and preallocation */
1437: MatSetSizes(C,pn,an,PETSC_DETERMINE,PETSC_DETERMINE);
1438: if (P->cmap->bs > 0) {PetscLayoutSetBlockSize(C->rmap,P->cmap->bs);}
1439: if (A->cmap->bs > 0) {PetscLayoutSetBlockSize(C->cmap,A->cmap->bs);}
1440: MatMPIAIJSetPreallocation(C,0,dnz,0,onz);
1441: MatPreallocateFinalize(dnz,onz);
1443: /* add C_loc and C_oth to C */
1444: MatGetOwnershipRange(C,&rstart,NULL);
1445: for (i=0; i<pn; i++) {
1446: const PetscInt ncols = c_loc->i[i+1] - c_loc->i[i];
1447: const PetscInt *cols = c_loc->j + c_loc->i[i];
1448: const PetscInt row = rstart + i;
1449: MatSetValues(C,1,&row,ncols,cols,NULL,INSERT_VALUES);
1450: }
1451: for (i=0; i<con; i++) {
1452: const PetscInt ncols = c_oth->i[i+1] - c_oth->i[i];
1453: const PetscInt *cols = c_oth->j + c_oth->i[i];
1454: const PetscInt row = prmap[i];
1455: MatSetValues(C,1,&row,ncols,cols,NULL,INSERT_VALUES);
1456: }
1457: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1458: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1459: MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
1461: /* members in merge */
1462: PetscFree(id_r);
1463: PetscFree(len_r);
1464: PetscFree(buf_ri[0]);
1465: PetscFree(buf_ri);
1466: PetscFree(buf_rj[0]);
1467: PetscFree(buf_rj);
1468: PetscLayoutDestroy(&rowmap);
1470: /* attach the supporting struct to C for reuse */
1471: C->product->data = ptap;
1472: C->product->destroy = MatDestroy_MPIAIJ_PtAP;
1473: return(0);
1474: }
1476: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_nonscalable(Mat P,Mat A,Mat C)
1477: {
1478: PetscErrorCode ierr;
1479: Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data;
1480: Mat_SeqAIJ *c_seq;
1481: Mat_APMPI *ptap;
1482: Mat A_loc,C_loc,C_oth;
1483: PetscInt i,rstart,rend,cm,ncols,row;
1484: const PetscInt *cols;
1485: const PetscScalar *vals;
1488: MatCheckProduct(C,3);
1489: ptap = (Mat_APMPI*)C->product->data;
1490: if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data");
1491: if (!ptap->A_loc) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatProductClear()");
1492: MatZeroEntries(C);
1494: if (ptap->reuse == MAT_REUSE_MATRIX) {
1495: /* These matrices are obtained in MatTransposeMatMultSymbolic() */
1496: /* 1) get R = Pd^T, Ro = Po^T */
1497: /*----------------------------*/
1498: MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);
1499: MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);
1501: /* 2) compute numeric A_loc */
1502: /*--------------------------*/
1503: MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&ptap->A_loc);
1504: }
1506: /* 3) C_loc = Rd*A_loc, C_oth = Ro*A_loc */
1507: A_loc = ptap->A_loc;
1508: ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,A_loc,ptap->C_loc);
1509: ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,A_loc,ptap->C_oth);
1510: C_loc = ptap->C_loc;
1511: C_oth = ptap->C_oth;
1513: /* add C_loc and C_oth to C */
1514: MatGetOwnershipRange(C,&rstart,&rend);
1516: /* C_loc -> C */
1517: cm = C_loc->rmap->N;
1518: c_seq = (Mat_SeqAIJ*)C_loc->data;
1519: cols = c_seq->j;
1520: vals = c_seq->a;
1521: for (i=0; i<cm; i++) {
1522: ncols = c_seq->i[i+1] - c_seq->i[i];
1523: row = rstart + i;
1524: MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);
1525: cols += ncols; vals += ncols;
1526: }
1528: /* Co -> C, off-processor part */
1529: cm = C_oth->rmap->N;
1530: c_seq = (Mat_SeqAIJ*)C_oth->data;
1531: cols = c_seq->j;
1532: vals = c_seq->a;
1533: for (i=0; i<cm; i++) {
1534: ncols = c_seq->i[i+1] - c_seq->i[i];
1535: row = p->garray[i];
1536: MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);
1537: cols += ncols; vals += ncols;
1538: }
1539: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1540: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1541: MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
1543: ptap->reuse = MAT_REUSE_MATRIX;
1544: return(0);
1545: }
1547: PetscErrorCode MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ(Mat P,Mat A,Mat C)
1548: {
1549: PetscErrorCode ierr;
1550: Mat_Merge_SeqsToMPI *merge;
1551: Mat_MPIAIJ *p =(Mat_MPIAIJ*)P->data;
1552: Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
1553: Mat_APMPI *ptap;
1554: PetscInt *adj;
1555: PetscInt i,j,k,anz,pnz,row,*cj,nexta;
1556: MatScalar *ada,*ca,valtmp;
1557: PetscInt am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
1558: MPI_Comm comm;
1559: PetscMPIInt size,rank,taga,*len_s;
1560: PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
1561: PetscInt **buf_ri,**buf_rj;
1562: PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
1563: MPI_Request *s_waits,*r_waits;
1564: MPI_Status *status;
1565: MatScalar **abuf_r,*ba_i,*pA,*coa,*ba;
1566: PetscInt *ai,*aj,*coi,*coj,*poJ,*pdJ;
1567: Mat A_loc;
1568: Mat_SeqAIJ *a_loc;
1571: MatCheckProduct(C,3);
1572: ptap = (Mat_APMPI*)C->product->data;
1573: if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data");
1574: if (!ptap->A_loc) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtA cannot be reused. Do not call MatProductClear()");
1575: PetscObjectGetComm((PetscObject)C,&comm);
1576: MPI_Comm_size(comm,&size);
1577: MPI_Comm_rank(comm,&rank);
1579: merge = ptap->merge;
1581: /* 2) compute numeric C_seq = P_loc^T*A_loc */
1582: /*------------------------------------------*/
1583: /* get data from symbolic products */
1584: coi = merge->coi; coj = merge->coj;
1585: PetscCalloc1(coi[pon]+1,&coa);
1586: bi = merge->bi; bj = merge->bj;
1587: owners = merge->rowmap->range;
1588: PetscCalloc1(bi[cm]+1,&ba);
1590: /* get A_loc by taking all local rows of A */
1591: A_loc = ptap->A_loc;
1592: MatMPIAIJGetLocalMat(A,MAT_REUSE_MATRIX,&A_loc);
1593: a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1594: ai = a_loc->i;
1595: aj = a_loc->j;
1597: for (i=0; i<am; i++) {
1598: anz = ai[i+1] - ai[i];
1599: adj = aj + ai[i];
1600: ada = a_loc->a + ai[i];
1602: /* 2-b) Compute Cseq = P_loc[i,:]^T*A[i,:] using outer product */
1603: /*-------------------------------------------------------------*/
1604: /* put the value into Co=(p->B)^T*A (off-diagonal part, send to others) */
1605: pnz = po->i[i+1] - po->i[i];
1606: poJ = po->j + po->i[i];
1607: pA = po->a + po->i[i];
1608: for (j=0; j<pnz; j++) {
1609: row = poJ[j];
1610: cj = coj + coi[row];
1611: ca = coa + coi[row];
1612: /* perform sparse axpy */
1613: nexta = 0;
1614: valtmp = pA[j];
1615: for (k=0; nexta<anz; k++) {
1616: if (cj[k] == adj[nexta]) {
1617: ca[k] += valtmp*ada[nexta];
1618: nexta++;
1619: }
1620: }
1621: PetscLogFlops(2.0*anz);
1622: }
1624: /* put the value into Cd (diagonal part) */
1625: pnz = pd->i[i+1] - pd->i[i];
1626: pdJ = pd->j + pd->i[i];
1627: pA = pd->a + pd->i[i];
1628: for (j=0; j<pnz; j++) {
1629: row = pdJ[j];
1630: cj = bj + bi[row];
1631: ca = ba + bi[row];
1632: /* perform sparse axpy */
1633: nexta = 0;
1634: valtmp = pA[j];
1635: for (k=0; nexta<anz; k++) {
1636: if (cj[k] == adj[nexta]) {
1637: ca[k] += valtmp*ada[nexta];
1638: nexta++;
1639: }
1640: }
1641: PetscLogFlops(2.0*anz);
1642: }
1643: }
1645: /* 3) send and recv matrix values coa */
1646: /*------------------------------------*/
1647: buf_ri = merge->buf_ri;
1648: buf_rj = merge->buf_rj;
1649: len_s = merge->len_s;
1650: PetscCommGetNewTag(comm,&taga);
1651: PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);
1653: PetscMalloc2(merge->nsend+1,&s_waits,size,&status);
1654: for (proc=0,k=0; proc<size; proc++) {
1655: if (!len_s[proc]) continue;
1656: i = merge->owners_co[proc];
1657: MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);
1658: k++;
1659: }
1660: if (merge->nrecv) {MPI_Waitall(merge->nrecv,r_waits,status);}
1661: if (merge->nsend) {MPI_Waitall(merge->nsend,s_waits,status);}
1663: PetscFree2(s_waits,status);
1664: PetscFree(r_waits);
1665: PetscFree(coa);
1667: /* 4) insert local Cseq and received values into Cmpi */
1668: /*----------------------------------------------------*/
1669: PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);
1670: for (k=0; k<merge->nrecv; k++) {
1671: buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1672: nrows = *(buf_ri_k[k]);
1673: nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */
1674: nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */
1675: }
1677: for (i=0; i<cm; i++) {
1678: row = owners[rank] + i; /* global row index of C_seq */
1679: bj_i = bj + bi[i]; /* col indices of the i-th row of C */
1680: ba_i = ba + bi[i];
1681: bnz = bi[i+1] - bi[i];
1682: /* add received vals into ba */
1683: for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1684: /* i-th row */
1685: if (i == *nextrow[k]) {
1686: cnz = *(nextci[k]+1) - *nextci[k];
1687: cj = buf_rj[k] + *(nextci[k]);
1688: ca = abuf_r[k] + *(nextci[k]);
1689: nextcj = 0;
1690: for (j=0; nextcj<cnz; j++) {
1691: if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */
1692: ba_i[j] += ca[nextcj++];
1693: }
1694: }
1695: nextrow[k]++; nextci[k]++;
1696: PetscLogFlops(2.0*cnz);
1697: }
1698: }
1699: MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);
1700: }
1701: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
1702: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
1704: PetscFree(ba);
1705: PetscFree(abuf_r[0]);
1706: PetscFree(abuf_r);
1707: PetscFree3(buf_ri_k,nextrow,nextci);
1708: return(0);
1709: }
1711: PetscErrorCode MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(Mat P,Mat A,PetscReal fill,Mat C)
1712: {
1713: PetscErrorCode ierr;
1714: Mat A_loc,POt,PDt;
1715: Mat_APMPI *ptap;
1716: PetscFreeSpaceList free_space=NULL,current_space=NULL;
1717: Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*a=(Mat_MPIAIJ*)A->data;
1718: PetscInt *pdti,*pdtj,*poti,*potj,*ptJ;
1719: PetscInt nnz;
1720: PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row;
1721: PetscInt am =A->rmap->n,pn=P->cmap->n;
1722: MPI_Comm comm;
1723: PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
1724: PetscInt **buf_rj,**buf_ri,**buf_ri_k;
1725: PetscInt len,proc,*dnz,*onz,*owners;
1726: PetscInt nzi,*bi,*bj;
1727: PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
1728: MPI_Request *swaits,*rwaits;
1729: MPI_Status *sstatus,rstatus;
1730: Mat_Merge_SeqsToMPI *merge;
1731: PetscInt *ai,*aj,*Jptr,anz,*prmap=p->garray,pon,nspacedouble=0,j;
1732: PetscReal afill =1.0,afill_tmp;
1733: PetscInt rstart = P->cmap->rstart,rmax,aN=A->cmap->N,Armax;
1734: Mat_SeqAIJ *a_loc,*pdt,*pot;
1735: PetscTable ta;
1736: MatType mtype;
1739: PetscObjectGetComm((PetscObject)A,&comm);
1740: /* check if matrix local sizes are compatible */
1741: 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);
1743: MPI_Comm_size(comm,&size);
1744: MPI_Comm_rank(comm,&rank);
1746: /* create struct Mat_APMPI and attached it to C later */
1747: PetscNew(&ptap);
1749: /* get A_loc by taking all local rows of A */
1750: MatMPIAIJGetLocalMat(A,MAT_INITIAL_MATRIX,&A_loc);
1752: ptap->A_loc = A_loc;
1753: a_loc = (Mat_SeqAIJ*)(A_loc)->data;
1754: ai = a_loc->i;
1755: aj = a_loc->j;
1757: /* determine symbolic Co=(p->B)^T*A - send to others */
1758: /*----------------------------------------------------*/
1759: MatTransposeSymbolic_SeqAIJ(p->A,&PDt);
1760: pdt = (Mat_SeqAIJ*)PDt->data;
1761: pdti = pdt->i; pdtj = pdt->j;
1763: MatTransposeSymbolic_SeqAIJ(p->B,&POt);
1764: pot = (Mat_SeqAIJ*)POt->data;
1765: poti = pot->i; potj = pot->j;
1767: /* then, compute symbolic Co = (p->B)^T*A */
1768: pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
1769: >= (num of nonzero rows of C_seq) - pn */
1770: PetscMalloc1(pon+1,&coi);
1771: coi[0] = 0;
1773: /* set initial free space to be fill*(nnz(p->B) + nnz(A)) */
1774: nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],ai[am]));
1775: PetscFreeSpaceGet(nnz,&free_space);
1776: current_space = free_space;
1778: /* create and initialize a linked list */
1779: PetscTableCreate(A->cmap->n + a->B->cmap->N,aN,&ta);
1780: MatRowMergeMax_SeqAIJ(a_loc,am,ta);
1781: PetscTableGetCount(ta,&Armax);
1783: PetscLLCondensedCreate_Scalable(Armax,&lnk);
1785: for (i=0; i<pon; i++) {
1786: pnz = poti[i+1] - poti[i];
1787: ptJ = potj + poti[i];
1788: for (j=0; j<pnz; j++) {
1789: row = ptJ[j]; /* row of A_loc == col of Pot */
1790: anz = ai[row+1] - ai[row];
1791: Jptr = aj + ai[row];
1792: /* add non-zero cols of AP into the sorted linked list lnk */
1793: PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);
1794: }
1795: nnz = lnk[0];
1797: /* If free space is not available, double the total space in the list */
1798: if (current_space->local_remaining<nnz) {
1799: PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);
1800: nspacedouble++;
1801: }
1803: /* Copy data into free space, and zero out denserows */
1804: PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);
1806: current_space->array += nnz;
1807: current_space->local_used += nnz;
1808: current_space->local_remaining -= nnz;
1810: coi[i+1] = coi[i] + nnz;
1811: }
1813: PetscMalloc1(coi[pon]+1,&coj);
1814: PetscFreeSpaceContiguous(&free_space,coj);
1815: PetscLLCondensedDestroy_Scalable(lnk); /* must destroy to get a new one for C */
1817: afill_tmp = (PetscReal)coi[pon]/(poti[pon] + ai[am]+1);
1818: if (afill_tmp > afill) afill = afill_tmp;
1820: /* send j-array (coj) of Co to other processors */
1821: /*----------------------------------------------*/
1822: /* determine row ownership */
1823: PetscNew(&merge);
1824: PetscLayoutCreate(comm,&merge->rowmap);
1826: merge->rowmap->n = pn;
1827: merge->rowmap->bs = 1;
1829: PetscLayoutSetUp(merge->rowmap);
1830: owners = merge->rowmap->range;
1832: /* determine the number of messages to send, their lengths */
1833: PetscCalloc1(size,&len_si);
1834: PetscCalloc1(size,&merge->len_s);
1836: len_s = merge->len_s;
1837: merge->nsend = 0;
1839: PetscMalloc1(size+2,&owners_co);
1841: proc = 0;
1842: for (i=0; i<pon; i++) {
1843: while (prmap[i] >= owners[proc+1]) proc++;
1844: len_si[proc]++; /* num of rows in Co to be sent to [proc] */
1845: len_s[proc] += coi[i+1] - coi[i];
1846: }
1848: len = 0; /* max length of buf_si[] */
1849: owners_co[0] = 0;
1850: for (proc=0; proc<size; proc++) {
1851: owners_co[proc+1] = owners_co[proc] + len_si[proc];
1852: if (len_si[proc]) {
1853: merge->nsend++;
1854: len_si[proc] = 2*(len_si[proc] + 1);
1855: len += len_si[proc];
1856: }
1857: }
1859: /* determine the number and length of messages to receive for coi and coj */
1860: PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);
1861: PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);
1863: /* post the Irecv and Isend of coj */
1864: PetscCommGetNewTag(comm,&tagj);
1865: PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);
1866: PetscMalloc1(merge->nsend+1,&swaits);
1867: for (proc=0, k=0; proc<size; proc++) {
1868: if (!len_s[proc]) continue;
1869: i = owners_co[proc];
1870: MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);
1871: k++;
1872: }
1874: /* receives and sends of coj are complete */
1875: PetscMalloc1(size,&sstatus);
1876: for (i=0; i<merge->nrecv; i++) {
1877: PETSC_UNUSED PetscMPIInt icompleted;
1878: MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);
1879: }
1880: PetscFree(rwaits);
1881: if (merge->nsend) {MPI_Waitall(merge->nsend,swaits,sstatus);}
1883: /* add received column indices into table to update Armax */
1884: /* Armax can be as large as aN if a P[row,:] is dense, see src/ksp/ksp/tutorials/ex56.c! */
1885: for (k=0; k<merge->nrecv; k++) {/* k-th received message */
1886: Jptr = buf_rj[k];
1887: for (j=0; j<merge->len_r[k]; j++) {
1888: PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);
1889: }
1890: }
1891: PetscTableGetCount(ta,&Armax);
1892: /* 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); */
1894: /* send and recv coi */
1895: /*-------------------*/
1896: PetscCommGetNewTag(comm,&tagi);
1897: PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);
1898: PetscMalloc1(len+1,&buf_s);
1899: buf_si = buf_s; /* points to the beginning of k-th msg to be sent */
1900: for (proc=0,k=0; proc<size; proc++) {
1901: if (!len_s[proc]) continue;
1902: /* form outgoing message for i-structure:
1903: buf_si[0]: nrows to be sent
1904: [1:nrows]: row index (global)
1905: [nrows+1:2*nrows+1]: i-structure index
1906: */
1907: /*-------------------------------------------*/
1908: nrows = len_si[proc]/2 - 1;
1909: buf_si_i = buf_si + nrows+1;
1910: buf_si[0] = nrows;
1911: buf_si_i[0] = 0;
1912: nrows = 0;
1913: for (i=owners_co[proc]; i<owners_co[proc+1]; i++) {
1914: nzi = coi[i+1] - coi[i];
1915: buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
1916: buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */
1917: nrows++;
1918: }
1919: MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);
1920: k++;
1921: buf_si += len_si[proc];
1922: }
1923: i = merge->nrecv;
1924: while (i--) {
1925: PETSC_UNUSED PetscMPIInt icompleted;
1926: MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);
1927: }
1928: PetscFree(rwaits);
1929: if (merge->nsend) {MPI_Waitall(merge->nsend,swaits,sstatus);}
1930: PetscFree(len_si);
1931: PetscFree(len_ri);
1932: PetscFree(swaits);
1933: PetscFree(sstatus);
1934: PetscFree(buf_s);
1936: /* compute the local portion of C (mpi mat) */
1937: /*------------------------------------------*/
1938: /* allocate bi array and free space for accumulating nonzero column info */
1939: PetscMalloc1(pn+1,&bi);
1940: bi[0] = 0;
1942: /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
1943: nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pdti[pn],PetscIntSumTruncate(poti[pon],ai[am])));
1944: PetscFreeSpaceGet(nnz,&free_space);
1945: current_space = free_space;
1947: PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);
1948: for (k=0; k<merge->nrecv; k++) {
1949: buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
1950: nrows = *buf_ri_k[k];
1951: nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */
1952: nextci[k] = buf_ri_k[k] + (nrows + 1); /* points to the next i-structure of k-th received i-structure */
1953: }
1955: PetscLLCondensedCreate_Scalable(Armax,&lnk);
1956: MatPreallocateInitialize(comm,pn,A->cmap->n,dnz,onz);
1957: rmax = 0;
1958: for (i=0; i<pn; i++) {
1959: /* add pdt[i,:]*AP into lnk */
1960: pnz = pdti[i+1] - pdti[i];
1961: ptJ = pdtj + pdti[i];
1962: for (j=0; j<pnz; j++) {
1963: row = ptJ[j]; /* row of AP == col of Pt */
1964: anz = ai[row+1] - ai[row];
1965: Jptr = aj + ai[row];
1966: /* add non-zero cols of AP into the sorted linked list lnk */
1967: PetscLLCondensedAddSorted_Scalable(anz,Jptr,lnk);
1968: }
1970: /* add received col data into lnk */
1971: for (k=0; k<merge->nrecv; k++) { /* k-th received message */
1972: if (i == *nextrow[k]) { /* i-th row */
1973: nzi = *(nextci[k]+1) - *nextci[k];
1974: Jptr = buf_rj[k] + *nextci[k];
1975: PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);
1976: nextrow[k]++; nextci[k]++;
1977: }
1978: }
1979: nnz = lnk[0];
1981: /* if free space is not available, make more free space */
1982: if (current_space->local_remaining<nnz) {
1983: PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);
1984: nspacedouble++;
1985: }
1986: /* copy data into free space, then initialize lnk */
1987: PetscLLCondensedClean_Scalable(nnz,current_space->array,lnk);
1988: MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);
1990: current_space->array += nnz;
1991: current_space->local_used += nnz;
1992: current_space->local_remaining -= nnz;
1994: bi[i+1] = bi[i] + nnz;
1995: if (nnz > rmax) rmax = nnz;
1996: }
1997: PetscFree3(buf_ri_k,nextrow,nextci);
1999: PetscMalloc1(bi[pn]+1,&bj);
2000: PetscFreeSpaceContiguous(&free_space,bj);
2001: afill_tmp = (PetscReal)bi[pn]/(pdti[pn] + poti[pon] + ai[am]+1);
2002: if (afill_tmp > afill) afill = afill_tmp;
2003: PetscLLCondensedDestroy_Scalable(lnk);
2004: PetscTableDestroy(&ta);
2006: MatDestroy(&POt);
2007: MatDestroy(&PDt);
2009: /* create symbolic parallel matrix C - why cannot be assembled in Numeric part */
2010: /*-------------------------------------------------------------------------------*/
2011: MatSetSizes(C,pn,A->cmap->n,PETSC_DETERMINE,PETSC_DETERMINE);
2012: MatSetBlockSizes(C,PetscAbs(P->cmap->bs),PetscAbs(A->cmap->bs));
2013: MatGetType(A,&mtype);
2014: MatSetType(C,mtype);
2015: MatMPIAIJSetPreallocation(C,0,dnz,0,onz);
2016: MatPreallocateFinalize(dnz,onz);
2017: MatSetBlockSize(C,1);
2018: MatSetOption(C,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE);
2019: for (i=0; i<pn; i++) {
2020: row = i + rstart;
2021: nnz = bi[i+1] - bi[i];
2022: Jptr = bj + bi[i];
2023: MatSetValues(C,1,&row,nnz,Jptr,NULL,INSERT_VALUES);
2024: }
2025: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
2026: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
2027: MatSetOption(C,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
2028: merge->bi = bi;
2029: merge->bj = bj;
2030: merge->coi = coi;
2031: merge->coj = coj;
2032: merge->buf_ri = buf_ri;
2033: merge->buf_rj = buf_rj;
2034: merge->owners_co = owners_co;
2036: /* attach the supporting struct to C for reuse */
2037: C->product->data = ptap;
2038: C->product->destroy = MatDestroy_MPIAIJ_PtAP;
2039: ptap->merge = merge;
2041: C->ops->mattransposemultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ;
2043: #if defined(PETSC_USE_INFO)
2044: if (bi[pn] != 0) {
2045: PetscInfo3(C,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);
2046: PetscInfo1(C,"Use MatTransposeMatMult(A,B,MatReuse,%g,&C) for best performance.\n",(double)afill);
2047: } else {
2048: PetscInfo(C,"Empty matrix product\n");
2049: }
2050: #endif
2051: return(0);
2052: }
2054: /* ---------------------------------------------------------------- */
2055: static PetscErrorCode MatProductSymbolic_AtB_MPIAIJ_MPIAIJ(Mat C)
2056: {
2058: Mat_Product *product = C->product;
2059: Mat A=product->A,B=product->B;
2060: PetscReal fill=product->fill;
2061: PetscBool flg;
2064: /* scalable */
2065: PetscStrcmp(product->alg,"scalable",&flg);
2066: if (flg) {
2067: MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ(A,B,fill,C);
2068: goto next;
2069: }
2071: /* nonscalable */
2072: PetscStrcmp(product->alg,"nonscalable",&flg);
2073: if (flg) {
2074: MatTransposeMatMultSymbolic_MPIAIJ_MPIAIJ_nonscalable(A,B,fill,C);
2075: goto next;
2076: }
2078: /* matmatmult */
2079: PetscStrcmp(product->alg,"at*b",&flg);
2080: if (flg) {
2081: Mat At;
2082: Mat_APMPI *ptap;
2084: MatTranspose(A,MAT_INITIAL_MATRIX,&At);
2085: MatMatMultSymbolic_MPIAIJ_MPIAIJ(At,B,fill,C);
2086: ptap = (Mat_APMPI*)C->product->data;
2087: if (ptap) {
2088: ptap->Pt = At;
2089: C->product->destroy = MatDestroy_MPIAIJ_PtAP;
2090: }
2091: C->ops->transposematmultnumeric = MatTransposeMatMultNumeric_MPIAIJ_MPIAIJ_matmatmult;
2092: goto next;
2093: }
2095: SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"MatProduct type is not supported");
2097: next:
2098: C->ops->productnumeric = MatProductNumeric_AtB;
2099: return(0);
2100: }
2102: /* ---------------------------------------------------------------- */
2103: /* Set options for MatMatMultxxx_MPIAIJ_MPIAIJ */
2104: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_AB(Mat C)
2105: {
2107: Mat_Product *product = C->product;
2108: Mat A=product->A,B=product->B;
2109: #if defined(PETSC_HAVE_HYPRE)
2110: const char *algTypes[4] = {"scalable","nonscalable","seqmpi","hypre"};
2111: PetscInt nalg = 4;
2112: #else
2113: const char *algTypes[3] = {"scalable","nonscalable","seqmpi"};
2114: PetscInt nalg = 3;
2115: #endif
2116: PetscInt alg = 1; /* set nonscalable algorithm as default */
2117: PetscBool flg;
2118: MPI_Comm comm;
2121: /* Check matrix local sizes */
2122: PetscObjectGetComm((PetscObject)C,&comm);
2123: if (A->cmap->rstart != B->rmap->rstart || A->cmap->rend != B->rmap->rend) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, (%D, %D) != (%D,%D)",A->cmap->rstart,A->cmap->rend,B->rmap->rstart,B->rmap->rend);
2125: /* Set "nonscalable" as default algorithm */
2126: PetscStrcmp(C->product->alg,"default",&flg);
2127: if (flg) {
2128: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2130: /* Set "scalable" as default if BN and local nonzeros of A and B are large */
2131: if (B->cmap->N > 100000) { /* may switch to scalable algorithm as default */
2132: MatInfo Ainfo,Binfo;
2133: PetscInt nz_local;
2134: PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable;
2136: MatGetInfo(A,MAT_LOCAL,&Ainfo);
2137: MatGetInfo(B,MAT_LOCAL,&Binfo);
2138: nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated);
2140: if (B->cmap->N > product->fill*nz_local) alg_scalable_loc = PETSC_TRUE;
2141: MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);
2143: if (alg_scalable) {
2144: alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
2145: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2146: PetscInfo2(B,"Use scalable algorithm, BN %D, fill*nz_allocated %g\n",B->cmap->N,product->fill*nz_local);
2147: }
2148: }
2149: }
2151: /* Get runtime option */
2152: if (product->api_user) {
2153: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatMult","Mat");
2154: PetscOptionsEList("-matmatmult_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2155: PetscOptionsEnd();
2156: } else {
2157: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_AB","Mat");
2158: PetscOptionsEList("-matproduct_ab_via","Algorithmic approach","MatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2159: PetscOptionsEnd();
2160: }
2161: if (flg) {
2162: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2163: }
2165: C->ops->productsymbolic = MatProductSymbolic_AB_MPIAIJ_MPIAIJ;
2166: return(0);
2167: }
2169: /* Set options for MatTransposeMatMultXXX_MPIAIJ_MPIAIJ */
2170: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_AtB(Mat C)
2171: {
2173: Mat_Product *product = C->product;
2174: Mat A=product->A,B=product->B;
2175: const char *algTypes[3] = {"scalable","nonscalable","at*b"};
2176: PetscInt nalg = 3;
2177: PetscInt alg = 1; /* set default algorithm */
2178: PetscBool flg;
2179: MPI_Comm comm;
2182: /* Check matrix local sizes */
2183: PetscObjectGetComm((PetscObject)C,&comm);
2184: if (A->rmap->rstart != B->rmap->rstart || A->rmap->rend != B->rmap->rend) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A (%D, %D) != B (%D,%D)",A->rmap->rstart,A->rmap->rend,B->rmap->rstart,B->rmap->rend);
2186: /* Set default algorithm */
2187: PetscStrcmp(C->product->alg,"default",&flg);
2188: if (flg) {
2189: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2190: }
2192: /* Set "scalable" as default if BN and local nonzeros of A and B are large */
2193: if (alg && B->cmap->N > 100000) { /* may switch to scalable algorithm as default */
2194: MatInfo Ainfo,Binfo;
2195: PetscInt nz_local;
2196: PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable;
2198: MatGetInfo(A,MAT_LOCAL,&Ainfo);
2199: MatGetInfo(B,MAT_LOCAL,&Binfo);
2200: nz_local = (PetscInt)(Ainfo.nz_allocated + Binfo.nz_allocated);
2202: if (B->cmap->N > product->fill*nz_local) alg_scalable_loc = PETSC_TRUE;
2203: MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);
2205: if (alg_scalable) {
2206: alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
2207: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2208: PetscInfo2(B,"Use scalable algorithm, BN %D, fill*nz_allocated %g\n",B->cmap->N,product->fill*nz_local);
2209: }
2210: }
2212: /* Get runtime option */
2213: if (product->api_user) {
2214: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatTransposeMatMult","Mat");
2215: PetscOptionsEList("-mattransposematmult_via","Algorithmic approach","MatTransposeMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2216: PetscOptionsEnd();
2217: } else {
2218: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_AtB","Mat");
2219: PetscOptionsEList("-matproduct_atb_via","Algorithmic approach","MatTransposeMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2220: PetscOptionsEnd();
2221: }
2222: if (flg) {
2223: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2224: }
2226: C->ops->productsymbolic = MatProductSymbolic_AtB_MPIAIJ_MPIAIJ;
2227: return(0);
2228: }
2230: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_PtAP(Mat C)
2231: {
2233: Mat_Product *product = C->product;
2234: Mat A=product->A,P=product->B;
2235: MPI_Comm comm;
2236: PetscBool flg;
2237: PetscInt alg=1; /* set default algorithm */
2238: #if !defined(PETSC_HAVE_HYPRE)
2239: const char *algTypes[4] = {"scalable","nonscalable","allatonce","allatonce_merged"};
2240: PetscInt nalg=4;
2241: #else
2242: const char *algTypes[5] = {"scalable","nonscalable","allatonce","allatonce_merged","hypre"};
2243: PetscInt nalg=5;
2244: #endif
2245: PetscInt pN=P->cmap->N;
2248: /* Check matrix local sizes */
2249: PetscObjectGetComm((PetscObject)C,&comm);
2250: if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Arow (%D, %D) != Prow (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
2251: if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Acol (%D, %D) != Prow (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend);
2253: /* Set "nonscalable" as default algorithm */
2254: PetscStrcmp(C->product->alg,"default",&flg);
2255: if (flg) {
2256: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2258: /* Set "scalable" as default if BN and local nonzeros of A and B are large */
2259: if (pN > 100000) {
2260: MatInfo Ainfo,Pinfo;
2261: PetscInt nz_local;
2262: PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable;
2264: MatGetInfo(A,MAT_LOCAL,&Ainfo);
2265: MatGetInfo(P,MAT_LOCAL,&Pinfo);
2266: nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated);
2268: if (pN > product->fill*nz_local) alg_scalable_loc = PETSC_TRUE;
2269: MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);
2271: if (alg_scalable) {
2272: alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */
2273: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2274: }
2275: }
2276: }
2278: /* Get runtime option */
2279: if (product->api_user) {
2280: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatPtAP","Mat");
2281: PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);
2282: PetscOptionsEnd();
2283: } else {
2284: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_PtAP","Mat");
2285: PetscOptionsEList("-matproduct_ptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);
2286: PetscOptionsEnd();
2287: }
2288: if (flg) {
2289: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2290: }
2292: C->ops->productsymbolic = MatProductSymbolic_PtAP_MPIAIJ_MPIAIJ;
2293: return(0);
2294: }
2296: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_RARt(Mat C)
2297: {
2298: Mat_Product *product = C->product;
2299: Mat A = product->A,R=product->B;
2302: /* Check matrix local sizes */
2303: if (A->cmap->n != R->cmap->n || A->rmap->n != R->cmap->n) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, A local (%D, %D), R local (%D,%D)",A->rmap->n,A->rmap->n,R->rmap->n,R->cmap->n);
2305: C->ops->productsymbolic = MatProductSymbolic_RARt_MPIAIJ_MPIAIJ;
2306: return(0);
2307: }
2309: /*
2310: Set options for ABC = A*B*C = A*(B*C); ABC's algorithm must be chosen from AB's algorithm
2311: */
2312: static PetscErrorCode MatProductSetFromOptions_MPIAIJ_ABC(Mat C)
2313: {
2315: Mat_Product *product = C->product;
2316: PetscBool flg = PETSC_FALSE;
2317: PetscInt alg = 1; /* default algorithm */
2318: const char *algTypes[3] = {"scalable","nonscalable","seqmpi"};
2319: PetscInt nalg = 3;
2322: /* Set default algorithm */
2323: PetscStrcmp(C->product->alg,"default",&flg);
2324: if (flg) {
2325: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2326: }
2328: /* Get runtime option */
2329: if (product->api_user) {
2330: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatMatMatMult","Mat");
2331: PetscOptionsEList("-matmatmatmult_via","Algorithmic approach","MatMatMatMult",algTypes,nalg,algTypes[alg],&alg,&flg);
2332: PetscOptionsEnd();
2333: } else {
2334: PetscOptionsBegin(PetscObjectComm((PetscObject)C),((PetscObject)C)->prefix,"MatProduct_ABC","Mat");
2335: PetscOptionsEList("-matproduct_abc_via","Algorithmic approach","MatProduct_ABC",algTypes,nalg,algTypes[alg],&alg,&flg);
2336: PetscOptionsEnd();
2337: }
2338: if (flg) {
2339: MatProductSetAlgorithm(C,(MatProductAlgorithm)algTypes[alg]);
2340: }
2342: C->ops->matmatmultsymbolic = MatMatMatMultSymbolic_MPIAIJ_MPIAIJ_MPIAIJ;
2343: C->ops->productsymbolic = MatProductSymbolic_ABC;
2344: return(0);
2345: }
2347: PETSC_INTERN PetscErrorCode MatProductSetFromOptions_MPIAIJ(Mat C)
2348: {
2350: Mat_Product *product = C->product;
2353: switch (product->type) {
2354: case MATPRODUCT_AB:
2355: MatProductSetFromOptions_MPIAIJ_AB(C);
2356: break;
2357: case MATPRODUCT_AtB:
2358: MatProductSetFromOptions_MPIAIJ_AtB(C);
2359: break;
2360: case MATPRODUCT_PtAP:
2361: MatProductSetFromOptions_MPIAIJ_PtAP(C);
2362: break;
2363: case MATPRODUCT_RARt:
2364: MatProductSetFromOptions_MPIAIJ_RARt(C);
2365: break;
2366: case MATPRODUCT_ABC:
2367: MatProductSetFromOptions_MPIAIJ_ABC(C);
2368: break;
2369: default:
2370: break;
2371: }
2372: return(0);
2373: }