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: PetscObjectContainer container;
23: PetscObjectQuery((PetscObject)A,"MatMergeSeqsToMPI",(PetscObject *)&container);
24: if (container) {
25: PetscObjectContainerGetPointer(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: PetscObjectContainerDestroy(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: PetscObjectContainer container;
54: PetscObjectQuery((PetscObject)A,"MatMergeSeqsToMPI",(PetscObject *)&container);
55: if (container) {
56: PetscObjectContainerGetPointer(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: PetscObjectContainer 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;
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: PetscObjectContainerDestroy(container);
131: PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",0);
132: }
134: /* create the container 'Mat_MatMatMultMPI' and attach it to P */
135: PetscNew(Mat_MatMatMultMPI,&ap);
136: ap->abi=PETSC_NULL; ap->abj=PETSC_NULL;
137: ap->abnz_max = 0;
139: PetscObjectContainerCreate(PETSC_COMM_SELF,&container);
140: PetscObjectContainerSetPointer(container,ap);
141: PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject)container);
142: ap->MatDestroy = P->ops->destroy;
143: P->ops->destroy = MatDestroy_MPIAIJ_MatMatMult;
144: ap->reuse = MAT_INITIAL_MATRIX;
145: PetscObjectContainerSetUserDestroy(container,PetscObjectContainerDestroy_Mat_MatMatMultMPI);
147: /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
148: MatGetBrowsOfAoCols(A,P,MAT_INITIAL_MATRIX,&ap->startsj,&ap->bufa,&ap->B_oth);
149: /* get P_loc by taking all local rows of P */
150: MatGetLocalMat(P,MAT_INITIAL_MATRIX,&ap->B_loc);
152: p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data;
153: p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data;
154: pi_loc = p_loc->i; pj_loc = p_loc->j;
155: pi_oth = p_oth->i; pj_oth = p_oth->j;
157: /* first, compute symbolic AP = A_loc*P */
158: /*--------------------------------------*/
159: PetscMalloc((am+2)*sizeof(PetscInt),&api);
160: ap->abi = api;
161: api[0] = 0;
163: /* create and initialize a linked list */
164: nlnk = pN+1;
165: PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);
167: /* Initial FreeSpace size is fill*nnz(A) */
168: PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am])),&free_space);
169: current_space = free_space;
171: for (i=0;i<am;i++) {
172: apnz = 0;
173: /* diagonal portion of A */
174: nzi = adi[i+1] - adi[i];
175: for (j=0; j<nzi; j++){
176: row = *adj++;
177: pnz = pi_loc[row+1] - pi_loc[row];
178: Jptr = pj_loc + pi_loc[row];
179: /* add non-zero cols of P into the sorted linked list lnk */
180: PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);
181: apnz += nlnk;
182: }
183: /* off-diagonal portion of A */
184: nzi = aoi[i+1] - aoi[i];
185: for (j=0; j<nzi; j++){
186: row = *aoj++;
187: pnz = pi_oth[row+1] - pi_oth[row];
188: Jptr = pj_oth + pi_oth[row];
189: PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);
190: apnz += nlnk;
191: }
193: api[i+1] = api[i] + apnz;
194: if (ap->abnz_max < apnz) ap->abnz_max = apnz;
196: /* if free space is not available, double the total space in the list */
197: if (current_space->local_remaining<apnz) {
198: PetscFreeSpaceGet(current_space->total_array_size,¤t_space);
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(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: merge->rowmap.n = pn;
264: merge->rowmap.N = PETSC_DECIDE;
265: merge->rowmap.bs = 1;
266: PetscMapInitialize(comm,&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);}
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: }
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: PetscMalloc((3*merge->nrecv+1)*sizeof(PetscInt**),&buf_ri_k);
384: nextrow = buf_ri_k + merge->nrecv;
385: nextci = nextrow + merge->nrecv;
386: for (k=0; k<merge->nrecv; k++){
387: buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
388: nrows = *buf_ri_k[k];
389: nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */
390: nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */
391: }
392: MatPreallocateInitialize(comm,pn,pn,dnz,onz);
393: for (i=0; i<pn; i++) {
394: /* add pdt[i,:]*AP into lnk */
395: nnz = 0;
396: pnz = pdti[i+1] - pdti[i];
397: j = pnz;
398: ptJ = pdtj + pdti[i+1];
399: while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
400: j--; ptJ--;
401: row = *ptJ; /* row of AP == col of Pt */
402: apnz = api[row+1] - api[row];
403: Jptr = apj + api[row];
404: /* add non-zero cols of AP into the sorted linked list lnk */
405: PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);
406: nnz += nlnk;
407: }
408: /* add received col data into lnk */
409: for (k=0; k<merge->nrecv; k++){ /* k-th received message */
410: if (i == *nextrow[k]) { /* i-th row */
411: nzi = *(nextci[k]+1) - *nextci[k];
412: Jptr = buf_rj[k] + *nextci[k];
413: PetscLLAdd(nzi,Jptr,pN,nlnk,lnk,lnkbt);
414: nnz += nlnk;
415: nextrow[k]++; nextci[k]++;
416: }
417: }
419: /* if free space is not available, make more free space */
420: if (current_space->local_remaining<nnz) {
421: PetscFreeSpaceGet(current_space->total_array_size,¤t_space);
422: }
423: /* copy data into free space, then initialize lnk */
424: PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);
425: MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);
426: current_space->array += nnz;
427: current_space->local_used += nnz;
428: current_space->local_remaining -= nnz;
429: bi[i+1] = bi[i] + nnz;
430: }
431: MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);
432: PetscFree(buf_ri_k);
434: PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);
435: PetscFreeSpaceContiguous(&free_space,bj);
436: PetscLLDestroy(lnk,lnkbt);
438: /* create symbolic parallel matrix B_mpi */
439: /*---------------------------------------*/
440: MatCreate(comm,&B_mpi);
441: MatSetSizes(B_mpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);
442: MatSetType(B_mpi,MATMPIAIJ);
443: MatMPIAIJSetPreallocation(B_mpi,0,dnz,0,onz);
444: MatPreallocateFinalize(dnz,onz);
446: merge->bi = bi;
447: merge->bj = bj;
448: merge->coi = coi;
449: merge->coj = coj;
450: merge->buf_ri = buf_ri;
451: merge->buf_rj = buf_rj;
452: merge->owners_co = owners_co;
453: merge->MatDestroy = B_mpi->ops->destroy;
454: merge->MatDuplicate = B_mpi->ops->duplicate;
456: /* B_mpi is not ready for use - assembly will be done by MatPtAPNumeric() */
457: B_mpi->assembled = PETSC_FALSE;
458: B_mpi->ops->destroy = MatDestroy_MPIAIJ_MatPtAP;
459: B_mpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
461: /* attach the supporting struct to B_mpi for reuse */
462: PetscObjectContainerCreate(PETSC_COMM_SELF,&container);
463: PetscObjectContainerSetPointer(container,merge);
464: PetscObjectCompose((PetscObject)B_mpi,"MatMergeSeqsToMPI",(PetscObject)container);
465: *C = B_mpi;
466: return(0);
467: }
471: PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
472: {
473: PetscErrorCode ierr;
474: Mat_Merge_SeqsToMPI *merge;
475: Mat_MatMatMultMPI *ap;
476: PetscObjectContainer cont_merge,cont_ptap;
477: Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data;
478: Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
479: Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
480: Mat_SeqAIJ *p_loc,*p_oth;
481: PetscInt *adi=ad->i,*aoi=ao->i,*adj=ad->j,*aoj=ao->j,*apJ,nextp,flops=0;
482: PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
483: PetscInt i,j,k,anz,pnz,apnz,nextap,row,*cj;
484: MatScalar *ada=ad->a,*aoa=ao->a,*apa,*pa,*ca,*pa_loc,*pa_oth;
485: PetscInt am=A->rmap.n,cm=C->rmap.n,pon=(p->B)->cmap.n;
486: MPI_Comm comm=C->comm;
487: PetscMPIInt size,rank,taga,*len_s;
488: PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
489: PetscInt **buf_ri,**buf_rj;
490: PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
491: MPI_Request *s_waits,*r_waits;
492: MPI_Status *status;
493: MatScalar **abuf_r,*ba_i,*pA,*coa,*ba;
494: PetscInt *api,*apj,*coi,*coj;
495: PetscInt *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap.rstart,pcend=P->cmap.rend;
498: MPI_Comm_size(comm,&size);
499: MPI_Comm_rank(comm,&rank);
501: PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&cont_merge);
502: if (cont_merge) {
503: PetscObjectContainerGetPointer(cont_merge,(void **)&merge);
504: } else {
505: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses an object container");
506: }
508: PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&cont_ptap);
509: if (cont_ptap) {
510: PetscObjectContainerGetPointer(cont_ptap,(void **)&ap);
511: if (ap->reuse == MAT_INITIAL_MATRIX){
512: ap->reuse = MAT_REUSE_MATRIX;
513: } else { /* update numerical values of P_oth and P_loc */
514: MatGetBrowsOfAoCols(A,P,MAT_REUSE_MATRIX,&ap->startsj,&ap->bufa,&ap->B_oth);
515: MatGetLocalMat(P,MAT_REUSE_MATRIX,&ap->B_loc);
516: }
517: } else {
518: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container");
519: }
521: /* get data from symbolic products */
522: p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data;
523: p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data;
524: pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a,pA=pa_loc;
525: pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
526:
527: coi = merge->coi; coj = merge->coj;
528: PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);
529: PetscMemzero(coa,coi[pon]*sizeof(MatScalar));
531: bi = merge->bi; bj = merge->bj;
532: owners = merge->rowmap.range;
533: PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);
534: PetscMemzero(ba,bi[cm]*sizeof(MatScalar));
536: /* get data from symbolic A*P */
537: PetscMalloc((ap->abnz_max+1)*sizeof(MatScalar),&apa);
538: PetscMemzero(apa,ap->abnz_max*sizeof(MatScalar));
540: /* compute numeric C_seq=P_loc^T*A_loc*P */
541: api = ap->abi; apj = ap->abj;
542: for (i=0;i<am;i++) {
543: /* form i-th sparse row of A*P */
544: apnz = api[i+1] - api[i];
545: apJ = apj + api[i];
546: /* diagonal portion of A */
547: anz = adi[i+1] - adi[i];
548: for (j=0;j<anz;j++) {
549: row = *adj++;
550: pnz = pi_loc[row+1] - pi_loc[row];
551: pj = pj_loc + pi_loc[row];
552: pa = pa_loc + pi_loc[row];
553: nextp = 0;
554: for (k=0; nextp<pnz; k++) {
555: if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
556: apa[k] += (*ada)*pa[nextp++];
557: }
558: }
559: flops += 2*pnz;
560: ada++;
561: }
562: /* off-diagonal portion of A */
563: anz = aoi[i+1] - aoi[i];
564: for (j=0; j<anz; j++) {
565: row = *aoj++;
566: pnz = pi_oth[row+1] - pi_oth[row];
567: pj = pj_oth + pi_oth[row];
568: pa = pa_oth + pi_oth[row];
569: nextp = 0;
570: for (k=0; nextp<pnz; k++) {
571: if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
572: apa[k] += (*aoa)*pa[nextp++];
573: }
574: }
575: flops += 2*pnz;
576: aoa++;
577: }
579: /* Compute P_loc[i,:]^T*AP[i,:] using outer product */
580: pnz = pi_loc[i+1] - pi_loc[i];
581: for (j=0; j<pnz; j++) {
582: nextap = 0;
583: row = *pJ++; /* global index */
584: if (row < pcstart || row >=pcend) { /* put the value into Co */
585: cj = coj + coi[*poJ];
586: ca = coa + coi[*poJ++];
587: } else { /* put the value into Cd */
588: cj = bj + bi[*pdJ];
589: ca = ba + bi[*pdJ++];
590: }
591: for (k=0; nextap<apnz; k++) {
592: if (cj[k]==apJ[nextap]) ca[k] += (*pA)*apa[nextap++];
593: }
594: flops += 2*apnz;
595: pA++;
596: }
598: /* zero the current row info for A*P */
599: PetscMemzero(apa,apnz*sizeof(MatScalar));
600: }
601: PetscFree(apa);
602:
603: /* send and recv matrix values */
604: /*-----------------------------*/
605: buf_ri = merge->buf_ri;
606: buf_rj = merge->buf_rj;
607: len_s = merge->len_s;
608: PetscCommGetNewTag(comm,&taga);
609: PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);
611: PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&s_waits);
612: for (proc=0,k=0; proc<size; proc++){
613: if (!len_s[proc]) continue;
614: i = merge->owners_co[proc];
615: MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);
616: k++;
617: }
618: PetscMalloc(size*sizeof(MPI_Status),&status);
619: if (merge->nrecv) {MPI_Waitall(merge->nrecv,r_waits,status);}
620: if (merge->nsend) {MPI_Waitall(merge->nsend,s_waits,status);}
621: PetscFree(status);
623: PetscFree(s_waits);
624: PetscFree(r_waits);
625: PetscFree(coa);
627: /* insert local and received values into C */
628: /*-----------------------------------------*/
629: PetscMalloc((3*merge->nrecv+1)*sizeof(PetscInt**),&buf_ri_k);
630: nextrow = buf_ri_k + merge->nrecv;
631: nextci = nextrow + merge->nrecv;
633: for (k=0; k<merge->nrecv; k++){
634: buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
635: nrows = *(buf_ri_k[k]);
636: nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */
637: nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */
638: }
640: for (i=0; i<cm; i++) {
641: row = owners[rank] + i; /* global row index of C_seq */
642: bj_i = bj + bi[i]; /* col indices of the i-th row of C */
643: ba_i = ba + bi[i];
644: bnz = bi[i+1] - bi[i];
645: /* add received vals into ba */
646: for (k=0; k<merge->nrecv; k++){ /* k-th received message */
647: /* i-th row */
648: if (i == *nextrow[k]) {
649: cnz = *(nextci[k]+1) - *nextci[k];
650: cj = buf_rj[k] + *(nextci[k]);
651: ca = abuf_r[k] + *(nextci[k]);
652: nextcj = 0;
653: for (j=0; nextcj<cnz; j++){
654: if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
655: ba_i[j] += ca[nextcj++];
656: }
657: }
658: nextrow[k]++; nextci[k]++;
659: }
660: }
661: MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);
662: flops += 2*cnz;
663: }
664: MatSetOption(C,MAT_COLUMNS_SORTED); /* -cause delay? */
665: MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
666: MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);
668: PetscFree(ba);
669: PetscFree(abuf_r);
670: PetscFree(buf_ri_k);
671: PetscLogFlops(flops);
672: return(0);
673: }