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