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