xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision e32f2f54e699d0aa6e733466c00da7e34666fe5e)
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 {
61     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
62   }
63   ierr = (*merge->MatDuplicate)(A,op,M);CHKERRQ(ierr);
64   (*M)->ops->destroy   = merge->MatDestroy;   /* =MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */
65   (*M)->ops->duplicate = merge->MatDuplicate; /* =MatDuplicate_ MPIAIJ */
66   PetscFunctionReturn(0);
67 }
68 
69 #undef __FUNCT__
70 #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ"
71 PetscErrorCode MatPtAPSymbolic_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
72 {
73   PetscErrorCode ierr;
74 
75   PetscFunctionBegin;
76   if (!P->ops->ptapsymbolic_mpiaij) {
77     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",((PetscObject)A)->type_name,((PetscObject)P)->type_name);
78   }
79   ierr = (*P->ops->ptapsymbolic_mpiaij)(A,P,fill,C);CHKERRQ(ierr);
80   PetscFunctionReturn(0);
81 }
82 
83 #undef __FUNCT__
84 #define __FUNCT__ "MatPtAPNumeric_MPIAIJ"
85 PetscErrorCode MatPtAPNumeric_MPIAIJ(Mat A,Mat P,Mat C)
86 {
87   PetscErrorCode ierr;
88 
89   PetscFunctionBegin;
90   if (!P->ops->ptapnumeric_mpiaij) {
91     SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",((PetscObject)A)->type_name,((PetscObject)P)->type_name);
92   }
93   ierr = (*P->ops->ptapnumeric_mpiaij)(A,P,C);CHKERRQ(ierr);
94   PetscFunctionReturn(0);
95 }
96 
97 #undef __FUNCT__
98 #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
99 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
100 {
101   PetscErrorCode       ierr;
102   Mat                  B_mpi;
103   Mat_MatMatMultMPI    *ap;
104   PetscContainer       container;
105   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
106   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data;
107   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
108   Mat_SeqAIJ           *p_loc,*p_oth;
109   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
110   PetscInt             *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz;
111   PetscInt             nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row;
112   PetscInt             am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
113   PetscBT              lnkbt;
114   MPI_Comm             comm=((PetscObject)A)->comm;
115   PetscMPIInt          size,rank,tag,*len_si,*len_s,*len_ri;
116   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
117   PetscInt             len,proc,*dnz,*onz,*owners;
118   PetscInt             nzi,*bi,*bj;
119   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
120   MPI_Request          *swaits,*rwaits;
121   MPI_Status           *sstatus,rstatus;
122   Mat_Merge_SeqsToMPI  *merge;
123   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0;
124   PetscMPIInt          j;
125 
126   PetscFunctionBegin;
127   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
128   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
129 
130   /* destroy the container 'Mat_MatMatMultMPI' in case that P is attached to it */
131   ierr = PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr);
132   if (container) {
133     /* reset functions */
134     ierr = PetscContainerGetPointer(container,(void **)&ap);CHKERRQ(ierr);
135     P->ops->destroy = ap->MatDestroy;
136     P->ops->duplicate = ap->MatDuplicate;
137     /* destroy container and contents */
138     ierr = PetscContainerDestroy(container);CHKERRQ(ierr);
139     ierr = PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",0);CHKERRQ(ierr);
140   }
141 
142   /* create the container 'Mat_MatMatMultMPI' and attach it to P */
143   ierr = PetscNew(Mat_MatMatMultMPI,&ap);CHKERRQ(ierr);
144   ap->abi=PETSC_NULL; ap->abj=PETSC_NULL;
145   ap->abnz_max = 0;
146 
147   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
148   ierr = PetscContainerSetPointer(container,ap);CHKERRQ(ierr);
149   ierr = PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr);
150   ap->MatDestroy  = P->ops->destroy;
151   P->ops->destroy = MatDestroy_MPIAIJ_MatMatMult;
152   ap->reuse       = MAT_INITIAL_MATRIX;
153   ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr);
154 
155   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
156   ierr = MatGetBrowsOfAoCols(A,P,MAT_INITIAL_MATRIX,&ap->startsj,&ap->startsj_r,&ap->bufa,&ap->B_oth);CHKERRQ(ierr);
157   /* get P_loc by taking all local rows of P */
158   ierr = MatGetLocalMat(P,MAT_INITIAL_MATRIX,&ap->B_loc);CHKERRQ(ierr);
159 
160   p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data;
161   p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data;
162   pi_loc = p_loc->i; pj_loc = p_loc->j;
163   pi_oth = p_oth->i; pj_oth = p_oth->j;
164 
165   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
166   /*-------------------------------------------------------------------*/
167   ierr  = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
168   ap->abi = api;
169   api[0] = 0;
170 
171   /* create and initialize a linked list */
172   nlnk = pN+1;
173   ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
174 
175   /* Initial FreeSpace size is fill*nnz(A) */
176   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am])),&free_space);CHKERRQ(ierr);
177   current_space = free_space;
178 
179   for (i=0;i<am;i++) {
180     apnz = 0;
181     /* diagonal portion of A */
182     nzi = adi[i+1] - adi[i];
183     for (j=0; j<nzi; j++){
184       row = *adj++;
185       pnz = pi_loc[row+1] - pi_loc[row];
186       Jptr  = pj_loc + pi_loc[row];
187       /* add non-zero cols of P into the sorted linked list lnk */
188       ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
189       apnz += nlnk;
190     }
191     /* off-diagonal portion of A */
192     nzi = aoi[i+1] - aoi[i];
193     for (j=0; j<nzi; j++){
194       row = *aoj++;
195       pnz = pi_oth[row+1] - pi_oth[row];
196       Jptr  = pj_oth + pi_oth[row];
197       ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
198       apnz += nlnk;
199     }
200 
201     api[i+1] = api[i] + apnz;
202     if (ap->abnz_max < apnz) ap->abnz_max = apnz;
203 
204     /* if free space is not available, double the total space in the list */
205     if (current_space->local_remaining<apnz) {
206       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
207       nspacedouble++;
208     }
209 
210     /* Copy data into free space, then initialize lnk */
211     ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
212     current_space->array           += apnz;
213     current_space->local_used      += apnz;
214     current_space->local_remaining -= apnz;
215   }
216   /* Allocate space for apj, initialize apj, and */
217   /* destroy list of free space and other temporary array(s) */
218   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ap->abj);CHKERRQ(ierr);
219   apj = ap->abj;
220   ierr = PetscFreeSpaceContiguous(&free_space,ap->abj);CHKERRQ(ierr);
221 
222   /* determine symbolic Co=(p->B)^T*AP - send to others */
223   /*----------------------------------------------------*/
224   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
225 
226   /* then, compute symbolic Co = (p->B)^T*AP */
227   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
228                          >= (num of nonzero rows of C_seq) - pn */
229   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
230   coi[0] = 0;
231 
232   /* set initial free space to be 3*pon*max( nnz(AP) per row) */
233   nnz           = 3*pon*ap->abnz_max + 1;
234   ierr          = PetscFreeSpaceGet(nnz,&free_space);
235   current_space = free_space;
236 
237   for (i=0; i<pon; i++) {
238     nnz  = 0;
239     pnz = poti[i+1] - poti[i];
240     j     = pnz;
241     ptJ   = potj + poti[i+1];
242     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
243       j--; ptJ--;
244       row  = *ptJ; /* row of AP == col of Pot */
245       apnz = api[row+1] - api[row];
246       Jptr   = apj + api[row];
247       /* add non-zero cols of AP into the sorted linked list lnk */
248       ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
249       nnz += nlnk;
250     }
251 
252     /* If free space is not available, double the total space in the list */
253     if (current_space->local_remaining<nnz) {
254       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
255     }
256 
257     /* Copy data into free space, and zero out denserows */
258     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
259     current_space->array           += nnz;
260     current_space->local_used      += nnz;
261     current_space->local_remaining -= nnz;
262     coi[i+1] = coi[i] + nnz;
263   }
264   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
265   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
266   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
267 
268   /* send j-array (coj) of Co to other processors */
269   /*----------------------------------------------*/
270   /* determine row ownership */
271   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
272   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
273   merge->rowmap->n = pn;
274   merge->rowmap->bs = 1;
275   ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
276   owners = merge->rowmap->range;
277 
278   /* determine the number of messages to send, their lengths */
279   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
280   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
281   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
282   len_s = merge->len_s;
283   merge->nsend = 0;
284 
285   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
286   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
287 
288   proc = 0;
289   for (i=0; i<pon; i++){
290     while (prmap[i] >= owners[proc+1]) proc++;
291     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
292     len_s[proc] += coi[i+1] - coi[i];
293   }
294 
295   len   = 0;  /* max length of buf_si[] */
296   owners_co[0] = 0;
297   for (proc=0; proc<size; proc++){
298     owners_co[proc+1] = owners_co[proc] + len_si[proc];
299     if (len_si[proc]){
300       merge->nsend++;
301       len_si[proc] = 2*(len_si[proc] + 1);
302       len += len_si[proc];
303     }
304   }
305 
306   /* determine the number and length of messages to receive for coi and coj  */
307   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
308   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
309 
310   /* post the Irecv and Isend of coj */
311   ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr);
312   ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
313 
314   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
315 
316   for (proc=0, k=0; proc<size; proc++){
317     if (!len_s[proc]) continue;
318     i = owners_co[proc];
319     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr);
320     k++;
321   }
322 
323   /* receives and sends of coj are complete */
324   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
325   i = merge->nrecv;
326   while (i--) {
327     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
328   }
329   ierr = PetscFree(rwaits);CHKERRQ(ierr);
330   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
331 
332   /* send and recv coi */
333   /*-------------------*/
334   ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
335 
336   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
337   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
338   for (proc=0,k=0; proc<size; proc++){
339     if (!len_s[proc]) continue;
340     /* form outgoing message for i-structure:
341          buf_si[0]:                 nrows to be sent
342                [1:nrows]:           row index (global)
343                [nrows+1:2*nrows+1]: i-structure index
344     */
345     /*-------------------------------------------*/
346     nrows = len_si[proc]/2 - 1;
347     buf_si_i    = buf_si + nrows+1;
348     buf_si[0]   = nrows;
349     buf_si_i[0] = 0;
350     nrows = 0;
351     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
352       nzi = coi[i+1] - coi[i];
353       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
354       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
355       nrows++;
356     }
357     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr);
358     k++;
359     buf_si += len_si[proc];
360   }
361   i = merge->nrecv;
362   while (i--) {
363     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
364   }
365   ierr = PetscFree(rwaits);CHKERRQ(ierr);
366   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
367   /*
368   ierr = PetscInfo2(A,"nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr);
369   for (i=0; i<merge->nrecv; i++){
370     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);
371   }
372   */
373   ierr = PetscFree(len_si);CHKERRQ(ierr);
374   ierr = PetscFree(len_ri);CHKERRQ(ierr);
375   ierr = PetscFree(swaits);CHKERRQ(ierr);
376   ierr = PetscFree(sstatus);CHKERRQ(ierr);
377   ierr = PetscFree(buf_s);CHKERRQ(ierr);
378 
379   /* compute the local portion of C (mpi mat) */
380   /*------------------------------------------*/
381   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
382 
383   /* allocate bi array and free space for accumulating nonzero column info */
384   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
385   bi[0] = 0;
386 
387   /* set initial free space to be 3*pn*max( nnz(AP) per row) */
388   nnz           = 3*pn*ap->abnz_max + 1;
389   ierr          = PetscFreeSpaceGet(nnz,&free_space);
390   current_space = free_space;
391 
392   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
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   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
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       ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
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         ierr = PetscLLAdd(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
421         nnz += nlnk;
422         nextrow[k]++; nextci[k]++;
423       }
424     }
425 
426     /* if free space is not available, make more free space */
427     if (current_space->local_remaining<nnz) {
428       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
429     }
430     /* copy data into free space, then initialize lnk */
431     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
432     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
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   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
439   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
440 
441   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
442   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
443   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
444 
445   /* create symbolic parallel matrix B_mpi */
446   /*---------------------------------------*/
447   ierr = MatCreate(comm,&B_mpi);CHKERRQ(ierr);
448   ierr = MatSetSizes(B_mpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
449   ierr = MatSetType(B_mpi,MATMPIAIJ);CHKERRQ(ierr);
450   ierr = MatMPIAIJSetPreallocation(B_mpi,0,dnz,0,onz);CHKERRQ(ierr);
451   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
452 
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;
462 
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;
467 
468   /* attach the supporting struct to B_mpi for reuse */
469   ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr);
470   ierr = PetscContainerSetPointer(container,merge);CHKERRQ(ierr);
471   ierr = PetscObjectCompose((PetscObject)B_mpi,"MatMergeSeqsToMPI",(PetscObject)container);CHKERRQ(ierr);
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     ierr = PetscInfo3(B_mpi,"Reallocs %D; Fill ratio: given %G needed %G when computing A*P.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
478     ierr = PetscInfo1(B_mpi,"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
479   } else {
480     ierr = PetscInfo(B_mpi,"Empty matrix product\n");CHKERRQ(ierr);
481   }
482 #endif
483   PetscFunctionReturn(0);
484 }
485 
486 #undef __FUNCT__
487 #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
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;
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=((PetscObject)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;
513 
514   PetscFunctionBegin;
515   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
516   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
517 
518   ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&cont_merge);CHKERRQ(ierr);
519   if (cont_merge) {
520     ierr  = PetscContainerGetPointer(cont_merge,(void **)&merge);CHKERRQ(ierr);
521   } else {
522     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses an object container");
523   }
524 
525   ierr = PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&cont_ptap);CHKERRQ(ierr);
526   if (cont_ptap) {
527     ierr  = PetscContainerGetPointer(cont_ptap,(void **)&ap);CHKERRQ(ierr);
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       ierr = MatGetBrowsOfAoCols(A,P,MAT_REUSE_MATRIX,&ap->startsj,&ap->startsj_r,&ap->bufa,&ap->B_oth);CHKERRQ(ierr);
532       ierr = MatGetLocalMat(P,MAT_REUSE_MATRIX,&ap->B_loc);CHKERRQ(ierr);
533     }
534   } else {
535     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container");
536   }
537 
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   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
546   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
547 
548   bi     = merge->bi; bj = merge->bj;
549   owners = merge->rowmap->range;
550   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
551   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
552 
553   /* get data from symbolic A*P */
554   ierr = PetscMalloc((ap->abnz_max+1)*sizeof(MatScalar),&apa);CHKERRQ(ierr);
555   ierr = PetscMemzero(apa,ap->abnz_max*sizeof(MatScalar));CHKERRQ(ierr);
556 
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       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
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       ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
593       aoa++;
594     }
595 
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       ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
612       pA++;
613     }
614 
615     /* zero the current row info for A*P */
616     ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
617   }
618   ierr = PetscFree(apa);CHKERRQ(ierr);
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   ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
626   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
627 
628   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&s_waits);CHKERRQ(ierr);
629   for (proc=0,k=0; proc<size; proc++){
630     if (!len_s[proc]) continue;
631     i = merge->owners_co[proc];
632     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
633     k++;
634   }
635   ierr = PetscMalloc(size*sizeof(MPI_Status),&status);CHKERRQ(ierr);
636   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
637   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
638   ierr = PetscFree(status);CHKERRQ(ierr);
639 
640   ierr = PetscFree(s_waits);CHKERRQ(ierr);
641   ierr = PetscFree(r_waits);CHKERRQ(ierr);
642   ierr = PetscFree(coa);CHKERRQ(ierr);
643 
644   /* insert local and received values into C */
645   /*-----------------------------------------*/
646   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
647 
648   for (k=0; k<merge->nrecv; k++){
649     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
650     nrows       = *(buf_ri_k[k]);
651     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
652     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
653   }
654 
655   for (i=0; i<cm; i++) {
656     row = owners[rank] + i; /* global row index of C_seq */
657     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
658     ba_i = ba + bi[i];
659     bnz  = bi[i+1] - bi[i];
660     /* add received vals into ba */
661     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
662       /* i-th row */
663       if (i == *nextrow[k]) {
664         cnz = *(nextci[k]+1) - *nextci[k];
665         cj  = buf_rj[k] + *(nextci[k]);
666         ca  = abuf_r[k] + *(nextci[k]);
667         nextcj = 0;
668         for (j=0; nextcj<cnz; j++){
669           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
670             ba_i[j] += ca[nextcj++];
671           }
672         }
673         nextrow[k]++; nextci[k]++;
674       }
675     }
676     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
677     ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
678   }
679   ierr = MatSetBlockSize(C,1);CHKERRQ(ierr);
680   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
681   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
682 
683   ierr = PetscFree(ba);CHKERRQ(ierr);
684   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
685   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
686   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
687   PetscFunctionReturn(0);
688 }
689