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