xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision 6ce94e8fe7f45068660bdbda61d05c45a2a71fa4)
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 /*
13 #define DEBUG_MATPTAP
14  */
15 
16 extern PetscErrorCode MatDestroy_MPIAIJ(Mat);
17 #undef __FUNCT__
18 #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP"
19 PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
20 {
21   PetscErrorCode       ierr;
22   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data;
23   Mat_PtAPMPI          *ptap=a->ptap;
24 
25   PetscFunctionBegin;
26   if (ptap){
27     Mat_Merge_SeqsToMPI  *merge=ptap->merge;
28     ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
29     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
30     ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
31     ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
32     ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
33     if (ptap->api){ierr = PetscFree(ptap->api);CHKERRQ(ierr);}
34     if (ptap->apj){ierr = PetscFree(ptap->apj);CHKERRQ(ierr);}
35     if (ptap->apa){ierr = PetscFree(ptap->apa);CHKERRQ(ierr);}
36     if (merge) {
37       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
38       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
39       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
40       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
41       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
42       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
43       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
44       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
45       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
46       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
47       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
48       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
49       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
50       ierr = merge->destroy(A);CHKERRQ(ierr);
51       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
52     }
53     ierr = PetscFree(ptap);CHKERRQ(ierr);
54   }
55 
56   PetscFunctionReturn(0);
57 }
58 
59 #undef __FUNCT__
60 #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP"
61 PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
62 {
63   PetscErrorCode       ierr;
64   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data;
65   Mat_PtAPMPI          *ptap = a->ptap;
66   Mat_Merge_SeqsToMPI  *merge = ptap->merge;
67 
68   PetscFunctionBegin;
69   ierr = (*merge->duplicate)(A,op,M);CHKERRQ(ierr);
70   (*M)->ops->destroy   = merge->destroy;
71   (*M)->ops->duplicate = merge->duplicate;
72   PetscFunctionReturn(0);
73 }
74 
75 #undef __FUNCT__
76 #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ"
77 PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
78 {
79   PetscErrorCode ierr;
80 
81   PetscFunctionBegin;
82   if (scall == MAT_INITIAL_MATRIX){
83     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
84     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
85     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
86   }
87   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
88   ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr);
89   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
90   PetscFunctionReturn(0);
91 }
92 
93 #undef __FUNCT__
94 #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
95 PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
96 {
97   PetscErrorCode       ierr;
98   Mat                  Cmpi;
99   Mat_PtAPMPI          *ptap;
100   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
101   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
102   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
103   Mat_SeqAIJ           *p_loc,*p_oth;
104   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
105   PetscInt             *adi=ad->i,*aj,*aoi=ao->i,nnz;
106   PetscInt             nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row;
107   PetscInt             am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
108   PetscBT              lnkbt;
109   MPI_Comm             comm=((PetscObject)A)->comm;
110   PetscMPIInt          size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
111   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
112   PetscInt             len,proc,*dnz,*onz,*owners;
113   PetscInt             nzi,*bi,*bj;
114   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
115   MPI_Request          *swaits,*rwaits;
116   MPI_Status           *sstatus,rstatus;
117   Mat_Merge_SeqsToMPI  *merge;
118   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
119   PetscReal            afill=1.0,afill_tmp;
120   PetscInt             rstart = P->cmap->rstart,rmax;
121   PetscScalar          *vals;
122 
123   PetscFunctionBegin;
124   /* check if matrix local sizes are compatible */
125   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend){
126     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Arow (%D, %D) != Prow (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
127   }
128   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend){
129     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Acol (%D, %D) != Prow (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend);
130   }
131 
132   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
133   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
134 
135   /* create struct Mat_PtAPMPI and attached it to C later */
136   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
137   ptap->reuse    = MAT_INITIAL_MATRIX;
138 
139 
140   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
141   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
142 #if defined(DEBUG_MATPTAP)
143    if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank);
144 #endif
145 
146   /* get P_loc by taking all local rows of P */
147   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
148 #if defined(DEBUG_MATPTAP)
149   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
150   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_loc is done...\n",rank);
151 #endif
152 
153   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
154   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
155   pi_loc = p_loc->i; pj_loc = p_loc->j;
156   pi_oth = p_oth->i; pj_oth = p_oth->j;
157 
158   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
159   /*-------------------------------------------------------------------*/
160   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&api);CHKERRQ(ierr);
161   api[0] = 0;
162 
163   /* create and initialize a linked list */
164   nlnk = pN+1;
165 #if defined(DEBUG_MATPTAP)
166   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
167   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call PetscLLCreate(), pN %D, P_rmax %D %D; Annz %D\n",rank,pN,p_loc->rmax,p_oth->rmax,adi[am]+aoi[am]);
168 #endif
169   ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
170 #if defined(DEBUG_MATPTAP)
171   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
172   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call PetscLLCreate() is done\n",rank);
173 #endif
174 
175   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
176   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&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     aj  = ad->j + adi[i];
184     for (j=0; j<nzi; j++){
185       row  = aj[j];
186       pnz  = pi_loc[row+1] - pi_loc[row];
187       Jptr = pj_loc + pi_loc[row];
188       /* add non-zero cols of P into the sorted linked list lnk */
189       ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
190       apnz += nlnk;
191     }
192     /* off-diagonal portion of A */
193     nzi = aoi[i+1] - aoi[i];
194     aj  = ao->j + aoi[i];
195     for (j=0; j<nzi; j++){
196       row = aj[j];
197       pnz = pi_oth[row+1] - pi_oth[row];
198       Jptr  = pj_oth + pi_oth[row];
199       ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
200       apnz += nlnk;
201     }
202 
203     api[i+1] = api[i] + apnz;
204     if (ap_rmax < apnz) ap_rmax = apnz;
205 
206     /* if free space is not available, double the total space in the list */
207     if (current_space->local_remaining<apnz) {
208       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
209       nspacedouble++;
210     }
211 
212     /* Copy data into free space, then initialize lnk */
213     ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
214     current_space->array           += apnz;
215     current_space->local_used      += apnz;
216     current_space->local_remaining -= apnz;
217   }
218   /* Allocate space for apj, initialize apj, and */
219   /* destroy list of free space and other temporary array(s) */
220   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&apj);CHKERRQ(ierr);
221   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
222   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
223   if (afill_tmp > afill) afill = afill_tmp;
224 
225   /* determine symbolic Co=(p->B)^T*AP - send to others */
226   /*----------------------------------------------------*/
227   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
228 
229   /* then, compute symbolic Co = (p->B)^T*AP */
230   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
231                          >= (num of nonzero rows of C_seq) - pn */
232   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
233   coi[0] = 0;
234 
235   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
236   nnz           = fill*(poti[pon] + api[am]);
237   ierr          = PetscFreeSpaceGet(nnz,&free_space);
238   current_space = free_space;
239 
240   for (i=0; i<pon; i++) {
241     nnz = 0;
242     pnz = poti[i+1] - poti[i];
243     ptJ = potj + poti[i];
244     for (j=0; j<pnz; j++){
245       row  = ptJ[j]; /* row of AP == col of Pot */
246       apnz = api[row+1] - api[row];
247       Jptr = apj + api[row];
248       /* add non-zero cols of AP into the sorted linked list lnk */
249       ierr = PetscLLAddSorted(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
250       nnz += nlnk;
251     }
252 
253     /* If free space is not available, double the total space in the list */
254     if (current_space->local_remaining<nnz) {
255       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
256       nspacedouble++;
257     }
258 
259     /* Copy data into free space, and zero out denserows */
260     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
261     current_space->array           += nnz;
262     current_space->local_used      += nnz;
263     current_space->local_remaining -= nnz;
264     coi[i+1] = coi[i] + nnz;
265   }
266   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
267   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
268   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
269   if (afill_tmp > afill) afill = afill_tmp;
270   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
271 
272 #if defined(DEBUG_MATPTAP)
273   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
274   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] local Co is done\n",rank);
275 #endif
276 
277   /* send j-array (coj) of Co to other processors */
278   /*----------------------------------------------*/
279   /* determine row ownership */
280   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
281   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
282   merge->rowmap->n = pn;
283   merge->rowmap->bs = 1;
284   ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
285   owners = merge->rowmap->range;
286 
287   /* determine the number of messages to send, their lengths */
288   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
289   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
290   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
291   len_s = merge->len_s;
292   merge->nsend = 0;
293 
294   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
295   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
296 
297   proc = 0;
298   for (i=0; i<pon; i++){
299     while (prmap[i] >= owners[proc+1]) proc++;
300     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
301     len_s[proc] += coi[i+1] - coi[i];
302   }
303 
304   len   = 0;  /* max length of buf_si[] */
305   owners_co[0] = 0;
306   for (proc=0; proc<size; proc++){
307     owners_co[proc+1] = owners_co[proc] + len_si[proc];
308     if (len_si[proc]){
309       merge->nsend++;
310       len_si[proc] = 2*(len_si[proc] + 1);
311       len += len_si[proc];
312     }
313   }
314 
315   /* determine the number and length of messages to receive for coi and coj  */
316   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
317   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
318 
319   /* post the Irecv and Isend of coj */
320   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
321   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
322   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
323   for (proc=0, k=0; proc<size; proc++){
324     if (!len_s[proc]) continue;
325     i = owners_co[proc];
326     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
327     k++;
328   }
329 
330   /* receives and sends of coj are complete */
331   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
332   for (i=0; i<merge->nrecv; i++){
333     PetscMPIInt icompleted;
334     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
335   }
336   ierr = PetscFree(rwaits);CHKERRQ(ierr);
337   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
338 
339   /* send and recv coi */
340   /*-------------------*/
341   ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
342   ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
343   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
344   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
345   for (proc=0,k=0; proc<size; proc++){
346     if (!len_s[proc]) continue;
347     /* form outgoing message for i-structure:
348          buf_si[0]:                 nrows to be sent
349                [1:nrows]:           row index (global)
350                [nrows+1:2*nrows+1]: i-structure index
351     */
352     /*-------------------------------------------*/
353     nrows = len_si[proc]/2 - 1;
354     buf_si_i    = buf_si + nrows+1;
355     buf_si[0]   = nrows;
356     buf_si_i[0] = 0;
357     nrows = 0;
358     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
359       nzi = coi[i+1] - coi[i];
360       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
361       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
362       nrows++;
363     }
364     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
365     k++;
366     buf_si += len_si[proc];
367   }
368   i = merge->nrecv;
369   while (i--) {
370     PetscMPIInt icompleted;
371     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
372   }
373   ierr = PetscFree(rwaits);CHKERRQ(ierr);
374   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
375   /*
376   ierr = PetscInfo2(A,"nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr);
377   for (i=0; i<merge->nrecv; i++){
378     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);
379   }
380   */
381   ierr = PetscFree(len_si);CHKERRQ(ierr);
382   ierr = PetscFree(len_ri);CHKERRQ(ierr);
383   ierr = PetscFree(swaits);CHKERRQ(ierr);
384   ierr = PetscFree(sstatus);CHKERRQ(ierr);
385   ierr = PetscFree(buf_s);CHKERRQ(ierr);
386 
387 #if defined(DEBUG_MATPTAP)
388   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
389   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Co is done\n",rank);
390 #endif
391 
392   /* compute the local portion of C (mpi mat) */
393   /*------------------------------------------*/
394   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
395 
396   /* allocate bi array and free space for accumulating nonzero column info */
397   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
398   bi[0] = 0;
399 
400   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
401   nnz           = fill*(pi_loc[pm] + api[am]);
402   ierr          = PetscFreeSpaceGet(nnz,&free_space);
403   current_space = free_space;
404 
405   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
406   for (k=0; k<merge->nrecv; k++){
407     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
408     nrows       = *buf_ri_k[k];
409     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
410     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
411   }
412   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
413   rmax = 0;
414   for (i=0; i<pn; i++) {
415     /* add pdt[i,:]*AP into lnk */
416     nnz = 0;
417     pnz = pdti[i+1] - pdti[i];
418     ptJ = pdtj + pdti[i];
419     for (j=0; j<pnz; j++){
420       row  = ptJ[j];  /* row of AP == col of Pt */
421       apnz = api[row+1] - api[row];
422       Jptr = apj + api[row];
423       /* add non-zero cols of AP into the sorted linked list lnk */
424       ierr = PetscLLAddSorted(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
425       nnz += nlnk;
426     }
427     /* add received col data into lnk */
428     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
429       if (i == *nextrow[k]) { /* i-th row */
430         nzi = *(nextci[k]+1) - *nextci[k];
431         Jptr  = buf_rj[k] + *nextci[k];
432         ierr = PetscLLAddSorted(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
433         nnz += nlnk;
434         nextrow[k]++; nextci[k]++;
435       }
436     }
437 
438     /* if free space is not available, make more free space */
439     if (current_space->local_remaining<nnz) {
440       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
441       nspacedouble++;
442     }
443     /* copy data into free space, then initialize lnk */
444     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
445     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
446     current_space->array           += nnz;
447     current_space->local_used      += nnz;
448     current_space->local_remaining -= nnz;
449     bi[i+1] = bi[i] + nnz;
450     if (nnz > rmax) rmax = nnz;
451   }
452   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
453   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
454 
455   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
456   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
457   afill_tmp = (PetscReal)bi[pn]/(pi_loc[pm] + api[am]+1);
458   if (afill_tmp > afill) afill = afill_tmp;
459   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
460 
461   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part - check runex55_SA ?  */
462   /*------------------------------------------------------------------------------------------------------*/
463   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
464   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
465 
466   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
467   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
468   ierr = MatSetBlockSizes(Cmpi,P->cmap->bs,P->cmap->bs);CHKERRQ(ierr);
469   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
470   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
471   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
472 
473   for (i=0; i<pn; i++){
474     row = i + rstart;
475     nnz = bi[i+1] - bi[i];
476     Jptr = bj + bi[i];
477     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
478   }
479   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
480   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
481   ierr = PetscFree(vals);CHKERRQ(ierr);
482 
483   merge->bi            = bi;
484   merge->bj            = bj;
485   merge->coi           = coi;
486   merge->coj           = coj;
487   merge->buf_ri        = buf_ri;
488   merge->buf_rj        = buf_rj;
489   merge->owners_co     = owners_co;
490   merge->destroy       = Cmpi->ops->destroy;
491   merge->duplicate     = Cmpi->ops->duplicate;
492 
493   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
494   /* Cmpi->assembled      = PETSC_FALSE; */
495   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
496   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
497 
498   /* attach the supporting struct to Cmpi for reuse */
499   c = (Mat_MPIAIJ*)Cmpi->data;
500   c->ptap        = ptap;
501   ptap->api      = api;
502   ptap->apj      = apj;
503   ptap->merge    = merge;
504   ptap->rmax     = ap_rmax;
505 
506   *C = Cmpi;
507 
508   /* flag 'scalable' determines which implementations to be used:
509        0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa;
510        1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */
511   /* set default scalable */
512   ptap->scalable = PETSC_TRUE;
513   ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,PETSC_NULL);CHKERRQ(ierr);
514   if (!ptap->scalable){  /* Do dense axpy */
515     ierr = PetscMalloc(pN*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr);
516   } else {
517     ierr = PetscMalloc((ap_rmax+1)*sizeof(PetscScalar),&ptap->apa);CHKERRQ(ierr);
518   }
519 
520 #if defined(PETSC_USE_INFO)
521   if (bi[pn] != 0) {
522     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
523     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
524   } else {
525     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
526   }
527 #endif
528   PetscFunctionReturn(0);
529 }
530 
531 #undef __FUNCT__
532 #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
533 PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
534 {
535   PetscErrorCode       ierr;
536   Mat_Merge_SeqsToMPI  *merge;
537   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
538   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
539   Mat_SeqAIJ           *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
540   Mat_SeqAIJ           *p_loc,*p_oth;
541   Mat_PtAPMPI          *ptap;
542   PetscInt             *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
543   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
544   PetscInt             i,j,k,anz,pnz,apnz,nextap,row,*cj;
545   MatScalar            *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
546   PetscInt             am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
547   MPI_Comm             comm=((PetscObject)C)->comm;
548   PetscMPIInt          size,rank,taga,*len_s;
549   PetscInt             *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
550   PetscInt             **buf_ri,**buf_rj;
551   PetscInt             cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
552   MPI_Request          *s_waits,*r_waits;
553   MPI_Status           *status;
554   MatScalar            **abuf_r,*ba_i,*pA,*coa,*ba;
555   PetscInt             *api,*apj,*coi,*coj;
556   PetscInt             *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
557   PetscBool            scalable;
558 
559   PetscFunctionBegin;
560   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
561   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
562 
563   ptap = c->ptap;
564   if (!ptap) SETERRQ(((PetscObject)C)->comm,PETSC_ERR_ARG_INCOMP,"MatPtAP() has not been called to create matrix C yet, cannot use MAT_REUSE_MATRIX");
565   merge    = ptap->merge;
566   apa      = ptap->apa;
567   scalable = ptap->scalable;
568 
569   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
570   /*--------------------------------------------------*/
571   if (ptap->reuse == MAT_INITIAL_MATRIX){
572     ptap->reuse = MAT_REUSE_MATRIX;
573   } else { /* update numerical values of P_oth and P_loc */
574     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
575     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
576   }
577 
578   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
579   /*--------------------------------------------------------------*/
580   /* get data from symbolic products */
581   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
582   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
583   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a;
584   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
585 
586   coi = merge->coi; coj = merge->coj;
587   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
588   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
589 
590   bi     = merge->bi; bj = merge->bj;
591   owners = merge->rowmap->range;
592   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
593   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
594 
595   api = ptap->api; apj = ptap->apj;
596 
597   if (!scalable){  /* Do dense axpy */
598     /*--------------------------------------------------*/
599     /* apa (length of pN) stores dense row A[i,:]*P - nonscalable! */
600     ierr = PetscMemzero(apa,P->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
601 
602     for (i=0; i<am; i++) {
603       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
604       /*------------------------------------------------------------*/
605       apJ = apj + api[i];
606 
607       /* diagonal portion of A */
608       anz = adi[i+1] - adi[i];
609       adj = ad->j + adi[i];
610       ada = ad->a + adi[i];
611       for (j=0; j<anz; j++) {
612         row = adj[j];
613         pnz = pi_loc[row+1] - pi_loc[row];
614         pj  = pj_loc + pi_loc[row];
615         pa  = pa_loc + pi_loc[row];
616 
617         /* perform dense axpy */
618         valtmp = ada[j];
619         for (k=0; k<pnz; k++){
620           apa[pj[k]] += valtmp*pa[k];
621         }
622         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
623       }
624 
625       /* off-diagonal portion of A */
626       anz = aoi[i+1] - aoi[i];
627       aoj = ao->j + aoi[i];
628       aoa = ao->a + aoi[i];
629       for (j=0; j<anz; j++) {
630         row = aoj[j];
631         pnz = pi_oth[row+1] - pi_oth[row];
632         pj  = pj_oth + pi_oth[row];
633         pa  = pa_oth + pi_oth[row];
634 
635         /* perform dense axpy */
636         valtmp = aoa[j];
637         for (k=0; k<pnz; k++){
638           apa[pj[k]] += valtmp*pa[k];
639         }
640         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
641       }
642 
643       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
644       /*--------------------------------------------------------------*/
645       apnz = api[i+1] - api[i];
646       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
647       pnz = po->i[i+1] - po->i[i];
648       poJ = po->j + po->i[i];
649       pA  = po->a + po->i[i];
650       for (j=0; j<pnz; j++){
651         row = poJ[j];
652         cnz = coi[row+1] - coi[row];
653         cj  = coj + coi[row];
654         ca  = coa + coi[row];
655         /* perform dense axpy */
656         valtmp = pA[j];
657         for (k=0; k<cnz; k++) {
658           ca[k] += valtmp*apa[cj[k]];
659         }
660         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
661       }
662 
663       /* put the value into Cd (diagonal part) */
664       pnz = pd->i[i+1] - pd->i[i];
665       pdJ = pd->j + pd->i[i];
666       pA  = pd->a + pd->i[i];
667       for (j=0; j<pnz; j++){
668         row = pdJ[j];
669         cnz = bi[row+1] - bi[row];
670         cj  = bj + bi[row];
671         ca  = ba + bi[row];
672         /* perform dense axpy */
673         valtmp = pA[j];
674         for (k=0; k<cnz; k++) {
675           ca[k] += valtmp*apa[cj[k]];
676         }
677         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
678       }
679 
680       /* zero the current row of A*P */
681       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
682     }
683   } else {/* Perform sparse axpy */
684     /*----------------------------------------------------*/
685     /* apa (length ap_rmax) stores sparse row A[i,:]*P */
686     ierr = PetscMemzero(apa,ptap->rmax*sizeof(MatScalar));CHKERRQ(ierr);
687 
688     pA=pa_loc;
689     for (i=0; i<am; i++) {
690       /* form i-th sparse row of A*P */
691       apnz = api[i+1] - api[i];
692       apJ  = apj + api[i];
693       /* diagonal portion of A */
694       anz = adi[i+1] - adi[i];
695       adj = ad->j + adi[i];
696       ada = ad->a + adi[i];
697       for (j=0; j<anz; j++) {
698         row = adj[j];
699         pnz = pi_loc[row+1] - pi_loc[row];
700         pj  = pj_loc + pi_loc[row];
701         pa  = pa_loc + pi_loc[row];
702         valtmp = ada[j];
703         nextp  = 0;
704         for (k=0; nextp<pnz; k++) {
705           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
706             apa[k] += valtmp*pa[nextp++];
707           }
708         }
709         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
710       }
711       /* off-diagonal portion of A */
712       anz = aoi[i+1] - aoi[i];
713       aoj = ao->j + aoi[i];
714       aoa = ao->a + aoi[i];
715       for (j=0; j<anz; j++) {
716         row = aoj[j];
717         pnz = pi_oth[row+1] - pi_oth[row];
718         pj  = pj_oth + pi_oth[row];
719         pa  = pa_oth + pi_oth[row];
720         valtmp = aoa[j];
721         nextp  = 0;
722         for (k=0; nextp<pnz; k++) {
723           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
724             apa[k] += valtmp*pa[nextp++];
725           }
726         }
727         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
728       }
729 
730       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
731       /*--------------------------------------------------------------*/
732       pnz = pi_loc[i+1] - pi_loc[i];
733       pJ  = pj_loc + pi_loc[i];
734       for (j=0; j<pnz; j++) {
735         nextap = 0;
736         row    = pJ[j]; /* global index */
737         if (row < pcstart || row >=pcend) { /* put the value into Co */
738           row = *poJ;
739           cj  = coj + coi[row];
740           ca  = coa + coi[row]; poJ++;
741         } else {                            /* put the value into Cd */
742           row  = *pdJ;
743           cj   = bj + bi[row];
744           ca   = ba + bi[row]; pdJ++;
745         }
746         valtmp = pA[j];
747         for (k=0; nextap<apnz; k++) {
748           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
749         }
750         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
751       }
752       pA += pnz;
753       /* zero the current row info for A*P */
754       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
755     }
756   }
757 
758   /* 3) send and recv matrix values coa */
759   /*------------------------------------*/
760   buf_ri = merge->buf_ri;
761   buf_rj = merge->buf_rj;
762   len_s  = merge->len_s;
763   ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
764   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
765 
766   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
767   for (proc=0,k=0; proc<size; proc++){
768     if (!len_s[proc]) continue;
769     i = merge->owners_co[proc];
770     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
771     k++;
772   }
773   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
774   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
775 
776   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
777   ierr = PetscFree(r_waits);CHKERRQ(ierr);
778   ierr = PetscFree(coa);CHKERRQ(ierr);
779 
780   /* 4) insert local Cseq and received values into Cmpi */
781   /*------------------------------------------------------*/
782   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
783   for (k=0; k<merge->nrecv; k++){
784     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
785     nrows       = *(buf_ri_k[k]);
786     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
787     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
788   }
789 
790   for (i=0; i<cm; i++) {
791     row = owners[rank] + i; /* global row index of C_seq */
792     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
793     ba_i = ba + bi[i];
794     bnz  = bi[i+1] - bi[i];
795     /* add received vals into ba */
796     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
797       /* i-th row */
798       if (i == *nextrow[k]) {
799         cnz = *(nextci[k]+1) - *nextci[k];
800         cj  = buf_rj[k] + *(nextci[k]);
801         ca  = abuf_r[k] + *(nextci[k]);
802         nextcj = 0;
803         for (j=0; nextcj<cnz; j++){
804           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
805             ba_i[j] += ca[nextcj++];
806           }
807         }
808         nextrow[k]++; nextci[k]++;
809         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
810       }
811     }
812     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
813   }
814   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
815   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
816 
817   ierr = PetscFree(ba);CHKERRQ(ierr);
818   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
819   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
820   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
821 #if defined(PETSC_USE_INFO)
822   if (scalable){
823     ierr = PetscInfo(C,"Use scalable sparse axpy\n");CHKERRQ(ierr);
824   } else {
825     ierr = PetscInfo(C,"Use non-scalable dense axpy\n");CHKERRQ(ierr);
826   }
827 #endif
828   PetscFunctionReturn(0);
829 }
830