Actual source code: mpiptap.c
2: /*
3: Defines projective product routines where A is a MPIAIJ matrix
4: C = P^T * A * P
5: */
7: #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/
8: #include <../src/mat/utils/freespace.h>
9: #include <../src/mat/impls/aij/mpi/mpiaij.h>
10: #include <petscbt.h>
15: PetscErrorCode MatDestroy_MPIAIJ_MatPtAP(Mat A)
16: {
17: PetscErrorCode ierr;
18: PetscContainer container;
21: PetscObjectQuery((PetscObject)A,"MatMergeSeqsToMPI",(PetscObject *)&container);
22: if (container) {
23: Mat_Merge_SeqsToMPI *merge;
24: PetscContainerGetPointer(container,(void **)&merge);
25: PetscFree(merge->id_r);
26: PetscFree(merge->len_s);
27: PetscFree(merge->len_r);
28: PetscFree(merge->bi);
29: PetscFree(merge->bj);
30: PetscFree(merge->buf_ri[0]);
31: PetscFree(merge->buf_ri);
32: PetscFree(merge->buf_rj[0]);
33: PetscFree(merge->buf_rj);
34: PetscFree(merge->coi);
35: PetscFree(merge->coj);
36: PetscFree(merge->owners_co);
37: PetscLayoutDestroy(&merge->rowmap);
38: merge->destroy(A);
39: PetscFree(merge);
40: PetscObjectCompose((PetscObject)A,"MatMergeSeqsToMPI",0);
41: }
42: return(0);
43: }
47: PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
48: {
49: PetscErrorCode ierr;
50: Mat_Merge_SeqsToMPI *merge;
51: PetscContainer container;
54: PetscObjectQuery((PetscObject)A,"MatMergeSeqsToMPI",(PetscObject *)&container);
55: if (container) {
56: PetscContainerGetPointer(container,(void **)&merge);
57: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
58: (*merge->duplicate)(A,op,M);
59: (*M)->ops->destroy = merge->destroy; /* =MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */
60: (*M)->ops->duplicate = merge->duplicate; /* =MatDuplicate_ MPIAIJ */
61: return(0);
62: }
66: PetscErrorCode MatPtAPSymbolic_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
67: {
71: if (!P->ops->ptapsymbolic_mpiaij) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",((PetscObject)A)->type_name,((PetscObject)P)->type_name);
72: (*P->ops->ptapsymbolic_mpiaij)(A,P,fill,C);
73: return(0);
74: }
78: PetscErrorCode MatPtAPNumeric_MPIAIJ(Mat A,Mat P,Mat C)
79: {
83: if (!P->ops->ptapnumeric_mpiaij) SETERRQ2(((PetscObject)A)->comm,PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",((PetscObject)A)->type_name,((PetscObject)P)->type_name);
84: (*P->ops->ptapnumeric_mpiaij)(A,P,C);
85: return(0);
86: }
90: PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
91: {
92: PetscErrorCode ierr;
93: Mat B_mpi;
94: Mat_MatMatMultMPI *ap;
95: PetscContainer container;
96: PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL;
97: Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data;
98: Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
99: Mat_SeqAIJ *p_loc,*p_oth;
100: PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
101: PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz;
102: PetscInt nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row;
103: PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
104: PetscBT lnkbt;
105: MPI_Comm comm=((PetscObject)A)->comm;
106: PetscMPIInt size,rank,tag,*len_si,*len_s,*len_ri;
107: PetscInt **buf_rj,**buf_ri,**buf_ri_k;
108: PetscInt len,proc,*dnz,*onz,*owners;
109: PetscInt nzi,*bi,*bj;
110: PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
111: MPI_Request *swaits,*rwaits;
112: MPI_Status *sstatus,rstatus;
113: Mat_Merge_SeqsToMPI *merge;
114: PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0;
115: PetscMPIInt j;
118: MPI_Comm_size(comm,&size);
119: MPI_Comm_rank(comm,&rank);
121: /* destroy the container 'Mat_MatMatMultMPI' in case that P is attached to it */
122: PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&container);
123: if (container) {
124: /* reset functions */
125: PetscContainerGetPointer(container,(void **)&ap);
126: P->ops->destroy = ap->destroy;
127: P->ops->duplicate = ap->duplicate;
128: /* destroy container and contents */
129: PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",0);
130: }
132: /* create the container 'Mat_MatMatMultMPI' and attach it to P */
133: PetscNew(Mat_MatMatMultMPI,&ap);
134: ap->abi=PETSC_NULL; ap->abj=PETSC_NULL;
135: ap->abnz_max = 0;
137: PetscContainerCreate(PETSC_COMM_SELF,&container);
138: PetscContainerSetPointer(container,ap);
139: PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);
140: PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject)container);
141: PetscContainerDestroy(&container);
142: ap->destroy = P->ops->destroy;
143: P->ops->destroy = MatDestroy_MPIAIJ_MatMatMult;
144: ap->reuse = MAT_INITIAL_MATRIX;
146: /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
147: MatGetBrowsOfAoCols(A,P,MAT_INITIAL_MATRIX,&ap->startsj,&ap->startsj_r,&ap->bufa,&ap->B_oth);
148: /* get P_loc by taking all local rows of P */
149: MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ap->B_loc);
151: p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data;
152: p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data;
153: pi_loc = p_loc->i; pj_loc = p_loc->j;
154: pi_oth = p_oth->i; pj_oth = p_oth->j;
156: /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
157: /*-------------------------------------------------------------------*/
158: PetscMalloc((am+2)*sizeof(PetscInt),&api);
159: ap->abi = api;
160: api[0] = 0;
162: /* create and initialize a linked list */
163: nlnk = pN+1;
164: PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);
166: /* Initial FreeSpace size is fill*nnz(A) */
167: PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am])),&free_space);
168: current_space = free_space;
170: for (i=0;i<am;i++) {
171: apnz = 0;
172: /* diagonal portion of A */
173: nzi = adi[i+1] - adi[i];
174: for (j=0; j<nzi; j++){
175: row = *adj++;
176: pnz = pi_loc[row+1] - pi_loc[row];
177: Jptr = pj_loc + pi_loc[row];
178: /* add non-zero cols of P into the sorted linked list lnk */
179: PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);
180: apnz += nlnk;
181: }
182: /* off-diagonal portion of A */
183: nzi = aoi[i+1] - aoi[i];
184: for (j=0; j<nzi; j++){
185: row = *aoj++;
186: pnz = pi_oth[row+1] - pi_oth[row];
187: Jptr = pj_oth + pi_oth[row];
188: PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);
189: apnz += nlnk;
190: }
192: api[i+1] = api[i] + apnz;
193: if (ap->abnz_max < apnz) ap->abnz_max = apnz;
195: /* if free space is not available, double the total space in the list */
196: if (current_space->local_remaining<apnz) {
197: PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);
198: nspacedouble++;
199: }
201: /* Copy data into free space, then initialize lnk */
202: PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);
203: current_space->array += apnz;
204: current_space->local_used += apnz;
205: current_space->local_remaining -= apnz;
206: }
207: /* Allocate space for apj, initialize apj, and */
208: /* destroy list of free space and other temporary array(s) */
209: PetscMalloc((api[am]+1)*sizeof(PetscInt),&ap->abj);
210: apj = ap->abj;
211: PetscFreeSpaceContiguous(&free_space,ap->abj);
213: /* determine symbolic Co=(p->B)^T*AP - send to others */
214: /*----------------------------------------------------*/
215: MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);
217: /* then, compute symbolic Co = (p->B)^T*AP */
218: pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
219: >= (num of nonzero rows of C_seq) - pn */
220: PetscMalloc((pon+1)*sizeof(PetscInt),&coi);
221: coi[0] = 0;
223: /* set initial free space to be 3*pon*max( nnz(AP) per row) */
224: nnz = 3*pon*ap->abnz_max + 1;
225: PetscFreeSpaceGet(nnz,&free_space);
226: current_space = free_space;
228: for (i=0; i<pon; i++) {
229: nnz = 0;
230: pnz = poti[i+1] - poti[i];
231: j = pnz;
232: ptJ = potj + poti[i+1];
233: while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
234: j--; ptJ--;
235: row = *ptJ; /* row of AP == col of Pot */
236: apnz = api[row+1] - api[row];
237: Jptr = apj + api[row];
238: /* add non-zero cols of AP into the sorted linked list lnk */
239: PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);
240: nnz += nlnk;
241: }
243: /* If free space is not available, double the total space in the list */
244: if (current_space->local_remaining<nnz) {
245: PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);
246: }
248: /* Copy data into free space, and zero out denserows */
249: PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);
250: current_space->array += nnz;
251: current_space->local_used += nnz;
252: current_space->local_remaining -= nnz;
253: coi[i+1] = coi[i] + nnz;
254: }
255: PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);
256: PetscFreeSpaceContiguous(&free_space,coj);
257: MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);
259: /* send j-array (coj) of Co to other processors */
260: /*----------------------------------------------*/
261: /* determine row ownership */
262: PetscNew(Mat_Merge_SeqsToMPI,&merge);
263: PetscLayoutCreate(comm,&merge->rowmap);
264: merge->rowmap->n = pn;
265: merge->rowmap->bs = 1;
266: PetscLayoutSetUp(merge->rowmap);
267: owners = merge->rowmap->range;
269: /* determine the number of messages to send, their lengths */
270: PetscMalloc(size*sizeof(PetscMPIInt),&len_si);
271: PetscMemzero(len_si,size*sizeof(PetscMPIInt));
272: PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);
273: len_s = merge->len_s;
274: merge->nsend = 0;
275:
276: PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);
277: PetscMemzero(len_s,size*sizeof(PetscMPIInt));
279: proc = 0;
280: for (i=0; i<pon; i++){
281: while (prmap[i] >= owners[proc+1]) proc++;
282: len_si[proc]++; /* num of rows in Co to be sent to [proc] */
283: len_s[proc] += coi[i+1] - coi[i];
284: }
286: len = 0; /* max length of buf_si[] */
287: owners_co[0] = 0;
288: for (proc=0; proc<size; proc++){
289: owners_co[proc+1] = owners_co[proc] + len_si[proc];
290: if (len_si[proc]){
291: merge->nsend++;
292: len_si[proc] = 2*(len_si[proc] + 1);
293: len += len_si[proc];
294: }
295: }
297: /* determine the number and length of messages to receive for coi and coj */
298: PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);
299: PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);
301: /* post the Irecv and Isend of coj */
302: PetscCommGetNewTag(comm,&tag);
303: PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);
305: PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);
307: for (proc=0, k=0; proc<size; proc++){
308: if (!len_s[proc]) continue;
309: i = owners_co[proc];
310: MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tag,comm,swaits+k);
311: k++;
312: }
314: /* receives and sends of coj are complete */
315: PetscMalloc(size*sizeof(MPI_Status),&sstatus);
316: i = merge->nrecv;
317: while (i--) {
318: MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);
319: }
320: PetscFree(rwaits);
321: if (merge->nsend) {MPI_Waitall(merge->nsend,swaits,sstatus);}
322:
323: /* send and recv coi */
324: /*-------------------*/
325: PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);
326:
327: PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);
328: buf_si = buf_s; /* points to the beginning of k-th msg to be sent */
329: for (proc=0,k=0; proc<size; proc++){
330: if (!len_s[proc]) continue;
331: /* form outgoing message for i-structure:
332: buf_si[0]: nrows to be sent
333: [1:nrows]: row index (global)
334: [nrows+1:2*nrows+1]: i-structure index
335: */
336: /*-------------------------------------------*/
337: nrows = len_si[proc]/2 - 1;
338: buf_si_i = buf_si + nrows+1;
339: buf_si[0] = nrows;
340: buf_si_i[0] = 0;
341: nrows = 0;
342: for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
343: nzi = coi[i+1] - coi[i];
344: buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
345: buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
346: nrows++;
347: }
348: MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tag,comm,swaits+k);
349: k++;
350: buf_si += len_si[proc];
351: }
352: i = merge->nrecv;
353: while (i--) {
354: MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);
355: }
356: PetscFree(rwaits);
357: if (merge->nsend) {MPI_Waitall(merge->nsend,swaits,sstatus);}
358: /*
359: PetscInfo2(A,"nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);
360: for (i=0; i<merge->nrecv; i++){
361: PetscInfo3(A,"recv len_ri=%d, len_rj=%d from [%d]\n",len_ri[i],merge->len_r[i],merge->id_r[i]);
362: }
363: */
364: PetscFree(len_si);
365: PetscFree(len_ri);
366: PetscFree(swaits);
367: PetscFree(sstatus);
368: PetscFree(buf_s);
370: /* compute the local portion of C (mpi mat) */
371: /*------------------------------------------*/
372: MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);
374: /* allocate bi array and free space for accumulating nonzero column info */
375: PetscMalloc((pn+1)*sizeof(PetscInt),&bi);
376: bi[0] = 0;
377:
378: /* set initial free space to be 3*pn*max( nnz(AP) per row) */
379: nnz = 3*pn*ap->abnz_max + 1;
380: PetscFreeSpaceGet(nnz,&free_space);
381: current_space = free_space;
383: PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);
384: for (k=0; k<merge->nrecv; k++){
385: buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
386: nrows = *buf_ri_k[k];
387: nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */
388: nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */
389: }
390: MatPreallocateInitialize(comm,pn,pn,dnz,onz);
391: for (i=0; i<pn; i++) {
392: /* add pdt[i,:]*AP into lnk */
393: nnz = 0;
394: pnz = pdti[i+1] - pdti[i];
395: j = pnz;
396: ptJ = pdtj + pdti[i+1];
397: while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
398: j--; ptJ--;
399: row = *ptJ; /* row of AP == col of Pt */
400: apnz = api[row+1] - api[row];
401: Jptr = apj + api[row];
402: /* add non-zero cols of AP into the sorted linked list lnk */
403: PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);
404: nnz += nlnk;
405: }
406: /* add received col data into lnk */
407: for (k=0; k<merge->nrecv; k++){ /* k-th received message */
408: if (i == *nextrow[k]) { /* i-th row */
409: nzi = *(nextci[k]+1) - *nextci[k];
410: Jptr = buf_rj[k] + *nextci[k];
411: PetscLLAdd(nzi,Jptr,pN,nlnk,lnk,lnkbt);
412: nnz += nlnk;
413: nextrow[k]++; nextci[k]++;
414: }
415: }
417: /* if free space is not available, make more free space */
418: if (current_space->local_remaining<nnz) {
419: PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);
420: }
421: /* copy data into free space, then initialize lnk */
422: PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);
423: MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);
424: current_space->array += nnz;
425: current_space->local_used += nnz;
426: current_space->local_remaining -= nnz;
427: bi[i+1] = bi[i] + nnz;
428: }
429: MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);
430: PetscFree3(buf_ri_k,nextrow,nextci);
432: PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);
433: PetscFreeSpaceContiguous(&free_space,bj);
434: PetscLLDestroy(lnk,lnkbt);
436: /* create symbolic parallel matrix B_mpi */
437: /*---------------------------------------*/
438: MatCreate(comm,&B_mpi);
439: MatSetSizes(B_mpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);
440: MatSetType(B_mpi,MATMPIAIJ);
441: MatMPIAIJSetPreallocation(B_mpi,0,dnz,0,onz);
442: MatPreallocateFinalize(dnz,onz);
444: merge->bi = bi;
445: merge->bj = bj;
446: merge->coi = coi;
447: merge->coj = coj;
448: merge->buf_ri = buf_ri;
449: merge->buf_rj = buf_rj;
450: merge->owners_co = owners_co;
451: merge->destroy = B_mpi->ops->destroy;
452: merge->duplicate = B_mpi->ops->duplicate;
454: /* B_mpi is not ready for use - assembly will be done by MatPtAPNumeric() */
455: B_mpi->assembled = PETSC_FALSE;
456: B_mpi->ops->destroy = MatDestroy_MPIAIJ_MatPtAP;
457: B_mpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
459: /* attach the supporting struct to B_mpi for reuse */
460: PetscContainerCreate(PETSC_COMM_SELF,&container);
461: PetscContainerSetPointer(container,merge);
462: PetscObjectCompose((PetscObject)B_mpi,"MatMergeSeqsToMPI",(PetscObject)container);
463: PetscContainerDestroy(&container);
464: *C = B_mpi;
465: #if defined(PETSC_USE_INFO)
466: if (bi[pn] != 0) {
467: PetscReal afill = ((PetscReal)bi[pn])/(adi[am]+aoi[am]);
468: if (afill < 1.0) afill = 1.0;
469: PetscInfo3(B_mpi,"Reallocs %D; Fill ratio: given %G needed %G when computing A*P.\n",nspacedouble,fill,afill);
470: PetscInfo1(B_mpi,"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);
471: } else {
472: PetscInfo(B_mpi,"Empty matrix product\n");
473: }
474: #endif
475: return(0);
476: }
480: PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
481: {
482: PetscErrorCode ierr;
483: Mat_Merge_SeqsToMPI *merge;
484: Mat_MatMatMultMPI *ap;
485: PetscContainer cont_merge,cont_ptap;
486: Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data;
487: Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
488: Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
489: Mat_SeqAIJ *p_loc,*p_oth;
490: PetscInt *adi=ad->i,*aoi=ao->i,*adj=ad->j,*aoj=ao->j,*apJ,nextp;
491: PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
492: PetscInt i,j,k,anz,pnz,apnz,nextap,row,*cj;
493: MatScalar *ada=ad->a,*aoa=ao->a,*apa,*pa,*ca,*pa_loc,*pa_oth;
494: PetscInt am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
495: MPI_Comm comm=((PetscObject)C)->comm;
496: PetscMPIInt size,rank,taga,*len_s;
497: PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
498: PetscInt **buf_ri,**buf_rj;
499: PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
500: MPI_Request *s_waits,*r_waits;
501: MPI_Status *status;
502: MatScalar **abuf_r,*ba_i,*pA,*coa,*ba;
503: PetscInt *api,*apj,*coi,*coj;
504: PetscInt *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
507: MPI_Comm_size(comm,&size);
508: MPI_Comm_rank(comm,&rank);
510: PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&cont_merge);
511: if (cont_merge) {
512: PetscContainerGetPointer(cont_merge,(void **)&merge);
513: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses an object container");
515: PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&cont_ptap);
516: if (cont_ptap) {
517: PetscContainerGetPointer(cont_ptap,(void **)&ap);
518: if (ap->reuse == MAT_INITIAL_MATRIX){
519: ap->reuse = MAT_REUSE_MATRIX;
520: } else { /* update numerical values of P_oth and P_loc */
521: MatGetBrowsOfAoCols(A,P,MAT_REUSE_MATRIX,&ap->startsj,&ap->startsj_r,&ap->bufa,&ap->B_oth);
522: MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ap->B_loc);
523: }
524: } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container");
526: /* get data from symbolic products */
527: p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data;
528: p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data;
529: pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a,pA=pa_loc;
530: pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
531:
532: coi = merge->coi; coj = merge->coj;
533: PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);
534: PetscMemzero(coa,coi[pon]*sizeof(MatScalar));
536: bi = merge->bi; bj = merge->bj;
537: owners = merge->rowmap->range;
538: PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);
539: PetscMemzero(ba,bi[cm]*sizeof(MatScalar));
541: /* get data from symbolic A*P */
542: PetscMalloc((ap->abnz_max+1)*sizeof(MatScalar),&apa);
543: PetscMemzero(apa,ap->abnz_max*sizeof(MatScalar));
545: /* compute numeric C_seq=P_loc^T*A_loc*P */
546: api = ap->abi; apj = ap->abj;
547: for (i=0;i<am;i++) {
548: /* form i-th sparse row of A*P */
549: apnz = api[i+1] - api[i];
550: apJ = apj + api[i];
551: /* diagonal portion of A */
552: anz = adi[i+1] - adi[i];
553: for (j=0;j<anz;j++) {
554: row = *adj++;
555: pnz = pi_loc[row+1] - pi_loc[row];
556: pj = pj_loc + pi_loc[row];
557: pa = pa_loc + pi_loc[row];
558: nextp = 0;
559: for (k=0; nextp<pnz; k++) {
560: if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
561: apa[k] += (*ada)*pa[nextp++];
562: }
563: }
564: PetscLogFlops(2.0*pnz);
565: ada++;
566: }
567: /* off-diagonal portion of A */
568: anz = aoi[i+1] - aoi[i];
569: for (j=0; j<anz; j++) {
570: row = *aoj++;
571: pnz = pi_oth[row+1] - pi_oth[row];
572: pj = pj_oth + pi_oth[row];
573: pa = pa_oth + pi_oth[row];
574: nextp = 0;
575: for (k=0; nextp<pnz; k++) {
576: if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
577: apa[k] += (*aoa)*pa[nextp++];
578: }
579: }
580: PetscLogFlops(2.0*pnz);
581: aoa++;
582: }
584: /* Compute P_loc[i,:]^T*AP[i,:] using outer product */
585: pnz = pi_loc[i+1] - pi_loc[i];
586: for (j=0; j<pnz; j++) {
587: nextap = 0;
588: row = *pJ++; /* global index */
589: if (row < pcstart || row >=pcend) { /* put the value into Co */
590: cj = coj + coi[*poJ];
591: ca = coa + coi[*poJ++];
592: } else { /* put the value into Cd */
593: cj = bj + bi[*pdJ];
594: ca = ba + bi[*pdJ++];
595: }
596: for (k=0; nextap<apnz; k++) {
597: if (cj[k]==apJ[nextap]) ca[k] += (*pA)*apa[nextap++];
598: }
599: PetscLogFlops(2.0*apnz);
600: pA++;
601: }
603: /* zero the current row info for A*P */
604: PetscMemzero(apa,apnz*sizeof(MatScalar));
605: }
606: PetscFree(apa);
607:
608: /* send and recv matrix values */
609: /*-----------------------------*/
610: buf_ri = merge->buf_ri;
611: buf_rj = merge->buf_rj;
612: len_s = merge->len_s;
613: PetscCommGetNewTag(comm,&taga);
614: PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);
616: PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&s_waits);
617: for (proc=0,k=0; proc<size; proc++){
618: if (!len_s[proc]) continue;
619: i = merge->owners_co[proc];
620: MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);
621: k++;
622: }
623: PetscMalloc(size*sizeof(MPI_Status),&status);
624: if (merge->nrecv) {MPI_Waitall(merge->nrecv,r_waits,status);}
625: if (merge->nsend) {MPI_Waitall(merge->nsend,s_waits,status);}
626: PetscFree(status);
628: PetscFree(s_waits);
629: PetscFree(r_waits);
630: PetscFree(coa);
632: /* insert local and received values into C */
633: /*-----------------------------------------*/
634: PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);
636: for (k=0; k<merge->nrecv; k++){
637: buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
638: nrows = *(buf_ri_k[k]);
639: nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */
640: nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */
641: }
643: for (i=0; i<cm; i++) {
644: row = owners[rank] + i; /* global row index of C_seq */
645: bj_i = bj + bi[i]; /* col indices of the i-th row of C */
646: ba_i = ba + bi[i];
647: bnz = bi[i+1] - bi[i];
648: /* add received vals into ba */
649: for (k=0; k<merge->nrecv; k++){ /* k-th received message */
650: /* i-th row */
651: if (i == *nextrow[k]) {
652: cnz = *(nextci[k]+1) - *nextci[k];
653: cj = buf_rj[k] + *(nextci[k]);
654: ca = abuf_r[k] + *(nextci[k]);
655: nextcj = 0;
656: for (j=0; nextcj<cnz; j++){
657: if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
658: ba_i[j] += ca[nextcj++];
659: }
660: }
661: nextrow[k]++; nextci[k]++;
662: }
663: }
664: MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);
665: PetscLogFlops(2.0*cnz);
666: }
667: MatSetBlockSize(C,1);
668: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
669: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
671: PetscFree(ba);
672: PetscFree(abuf_r[0]);
673: PetscFree(abuf_r);
674: PetscFree3(buf_ri_k,nextrow,nextci);
675: return(0);
676: }