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