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