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