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