xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision f8487c7373c1533f3bf25724e12928e7f7a4f17b)
182412ba7SHong Zhang 
242c57489SHong Zhang /*
342c57489SHong Zhang   Defines projective product routines where A is a MPIAIJ matrix
442c57489SHong Zhang           C = P^T * A * P
542c57489SHong Zhang */
642c57489SHong Zhang 
7c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h>   /*I "petscmat.h" I*/
8c6db04a5SJed Brown #include <../src/mat/utils/freespace.h>
9c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h>
10c6db04a5SJed Brown #include <petscbt.h>
1142c57489SHong Zhang 
1209573ac7SBarry Smith extern PetscErrorCode MatDestroy_MPIAIJ(Mat);
13cc6cb787SHong Zhang #undef __FUNCT__
14*f8487c73SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP"
15*f8487c73SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
16cc6cb787SHong Zhang {
17cc6cb787SHong Zhang   PetscErrorCode       ierr;
18*f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data;
19*f8487c73SHong Zhang   Mat_PtAPMPI          *ptap=a->ptap;
20cc6cb787SHong Zhang 
21cc6cb787SHong Zhang   PetscFunctionBegin;
22*f8487c73SHong Zhang   if (ptap){
23*f8487c73SHong Zhang     ierr = PetscFree2(ptap->startsj,ptap->startsj_r);CHKERRQ(ierr);
24*f8487c73SHong Zhang     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
25*f8487c73SHong Zhang     ierr = MatDestroy(&ptap->B_loc);CHKERRQ(ierr);
26*f8487c73SHong Zhang     ierr = MatDestroy(&ptap->B_oth);CHKERRQ(ierr);
27*f8487c73SHong Zhang     ierr = PetscFree(ptap->abi);CHKERRQ(ierr);
28*f8487c73SHong Zhang     ierr = PetscFree(ptap->abj);CHKERRQ(ierr);
29*f8487c73SHong Zhang   }
30*f8487c73SHong Zhang   if (ptap->merge) {
31*f8487c73SHong Zhang     Mat_Merge_SeqsToMPI  *merge=ptap->merge;
32cc6cb787SHong Zhang     ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
33cc6cb787SHong Zhang     ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
34cc6cb787SHong Zhang     ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
35cc6cb787SHong Zhang     ierr = PetscFree(merge->bi);CHKERRQ(ierr);
36cc6cb787SHong Zhang     ierr = PetscFree(merge->bj);CHKERRQ(ierr);
37c05d87d6SBarry Smith     ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
38cc6cb787SHong Zhang     ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
39c05d87d6SBarry Smith     ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
40cc6cb787SHong Zhang     ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
4105b42c5fSBarry Smith     ierr = PetscFree(merge->coi);CHKERRQ(ierr);
4205b42c5fSBarry Smith     ierr = PetscFree(merge->coj);CHKERRQ(ierr);
4305b42c5fSBarry Smith     ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
446bf464f9SBarry Smith     ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
45dce485f0SBarry Smith     ierr = merge->destroy(A);CHKERRQ(ierr);
46*f8487c73SHong Zhang     ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
47bf0cc555SLisandro Dalcin   }
48*f8487c73SHong Zhang   ierr = PetscFree(ptap);CHKERRQ(ierr);
49cc6cb787SHong Zhang   PetscFunctionReturn(0);
50cc6cb787SHong Zhang }
51cc6cb787SHong Zhang 
52cc6cb787SHong Zhang #undef __FUNCT__
53cc6cb787SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP"
54f0c0a3a6SBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
55f0c0a3a6SBarry Smith {
56cc6cb787SHong Zhang   PetscErrorCode       ierr;
57cc6cb787SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
58776b82aeSLisandro Dalcin   PetscContainer       container;
59cc6cb787SHong Zhang 
60cc6cb787SHong Zhang   PetscFunctionBegin;
61cc6cb787SHong Zhang   ierr = PetscObjectQuery((PetscObject)A,"MatMergeSeqsToMPI",(PetscObject *)&container);CHKERRQ(ierr);
62cc6cb787SHong Zhang   if (container) {
63776b82aeSLisandro Dalcin     ierr  = PetscContainerGetPointer(container,(void **)&merge);CHKERRQ(ierr);
64e7e72b3dSBarry Smith   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit");
65dce485f0SBarry Smith   ierr = (*merge->duplicate)(A,op,M);CHKERRQ(ierr);
66dce485f0SBarry Smith   (*M)->ops->destroy   = merge->destroy;   /* =MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */
67dce485f0SBarry Smith   (*M)->ops->duplicate = merge->duplicate; /* =MatDuplicate_ MPIAIJ */
68cc6cb787SHong Zhang   PetscFunctionReturn(0);
69cc6cb787SHong Zhang }
70cc6cb787SHong Zhang 
7142c57489SHong Zhang #undef __FUNCT__
727a7894deSKris Buschelman #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ"
737a7894deSKris Buschelman PetscErrorCode MatPtAPSymbolic_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
7442c57489SHong Zhang {
7542c57489SHong Zhang   PetscErrorCode ierr;
7642c57489SHong Zhang 
7742c57489SHong Zhang   PetscFunctionBegin;
7865e19b50SBarry Smith   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);
797a7894deSKris Buschelman   ierr = (*P->ops->ptapsymbolic_mpiaij)(A,P,fill,C);CHKERRQ(ierr);
807a7894deSKris Buschelman   PetscFunctionReturn(0);
817a7894deSKris Buschelman }
827a7894deSKris Buschelman 
837a7894deSKris Buschelman #undef __FUNCT__
847a7894deSKris Buschelman #define __FUNCT__ "MatPtAPNumeric_MPIAIJ"
857a7894deSKris Buschelman PetscErrorCode MatPtAPNumeric_MPIAIJ(Mat A,Mat P,Mat C)
867a7894deSKris Buschelman {
877a7894deSKris Buschelman   PetscErrorCode ierr;
887a7894deSKris Buschelman 
897a7894deSKris Buschelman   PetscFunctionBegin;
9065e19b50SBarry Smith   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);
917a7894deSKris Buschelman   ierr = (*P->ops->ptapnumeric_mpiaij)(A,P,C);CHKERRQ(ierr);
9242c57489SHong Zhang   PetscFunctionReturn(0);
9342c57489SHong Zhang }
9442c57489SHong Zhang 
9542c57489SHong Zhang #undef __FUNCT__
9642c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
9742c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
9842c57489SHong Zhang {
9942c57489SHong Zhang   PetscErrorCode       ierr;
1006c6e00e0SHong Zhang   Mat                  B_mpi;
101*f8487c73SHong Zhang   Mat_PtAPMPI          *ptap;
102a1a86e44SBarry Smith   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
103*f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
10442c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
10542c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
106ded4f561SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
10783445d95SHong Zhang   PetscInt             *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz;
10813f74950SBarry Smith   PetscInt             nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row;
109d0f46423SBarry Smith   PetscInt             am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n;
11042c57489SHong Zhang   PetscBT              lnkbt;
1117adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)A)->comm;
112ded4f561SHong Zhang   PetscMPIInt          size,rank,tag,*len_si,*len_s,*len_ri;
11342c57489SHong Zhang   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
11442c57489SHong Zhang   PetscInt             len,proc,*dnz,*onz,*owners;
115ded4f561SHong Zhang   PetscInt             nzi,*bi,*bj;
11642c57489SHong Zhang   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
117ded4f561SHong Zhang   MPI_Request          *swaits,*rwaits;
118ded4f561SHong Zhang   MPI_Status           *sstatus,rstatus;
11942c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
120f2b054eeSHong Zhang   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0;
12113f74950SBarry Smith   PetscMPIInt          j;
12242c57489SHong Zhang 
12342c57489SHong Zhang   PetscFunctionBegin;
12483445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
12583445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
12683445d95SHong Zhang 
127*f8487c73SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
128*f8487c73SHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
129*f8487c73SHong Zhang   ptap->abi=PETSC_NULL; ptap->abj=PETSC_NULL;
130*f8487c73SHong Zhang   ptap->abnz_max = 0;
131*f8487c73SHong Zhang   ptap->reuse     = MAT_INITIAL_MATRIX;
1326c6e00e0SHong Zhang 
1336c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
134*f8487c73SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->B_oth);CHKERRQ(ierr);
1356c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
136*f8487c73SHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->B_loc);CHKERRQ(ierr);
1376c6e00e0SHong Zhang 
138*f8487c73SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->B_loc)->data;
139*f8487c73SHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->B_oth)->data;
1406c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
1416c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
14242c57489SHong Zhang 
143fd0ff01cSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
144fd0ff01cSHong Zhang   /*-------------------------------------------------------------------*/
14582412ba7SHong Zhang   ierr  = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr);
146*f8487c73SHong Zhang   ptap->abi = api;
14782412ba7SHong Zhang   api[0] = 0;
14883445d95SHong Zhang 
14983445d95SHong Zhang   /* create and initialize a linked list */
15083445d95SHong Zhang   nlnk = pN+1;
15183445d95SHong Zhang   ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
15283445d95SHong Zhang 
15382412ba7SHong Zhang   /* Initial FreeSpace size is fill*nnz(A) */
154a1a86e44SBarry Smith   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am])),&free_space);CHKERRQ(ierr);
155f4a743e1SHong Zhang   current_space = free_space;
156f4a743e1SHong Zhang 
157f4a743e1SHong Zhang   for (i=0;i<am;i++) {
15882412ba7SHong Zhang     apnz = 0;
159f4a743e1SHong Zhang     /* diagonal portion of A */
160ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
161ded4f561SHong Zhang     for (j=0; j<nzi; j++){
16282412ba7SHong Zhang       row = *adj++;
16382412ba7SHong Zhang       pnz = pi_loc[row+1] - pi_loc[row];
164ded4f561SHong Zhang       Jptr  = pj_loc + pi_loc[row];
16583445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
166ded4f561SHong Zhang       ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
16782412ba7SHong Zhang       apnz += nlnk;
168f4a743e1SHong Zhang     }
169f4a743e1SHong Zhang     /* off-diagonal portion of A */
170ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
171ded4f561SHong Zhang     for (j=0; j<nzi; j++){
17282412ba7SHong Zhang       row = *aoj++;
17382412ba7SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
174ded4f561SHong Zhang       Jptr  = pj_oth + pi_oth[row];
175ded4f561SHong Zhang       ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
17682412ba7SHong Zhang       apnz += nlnk;
177f4a743e1SHong Zhang     }
178f4a743e1SHong Zhang 
17982412ba7SHong Zhang     api[i+1] = api[i] + apnz;
180*f8487c73SHong Zhang     if (ptap->abnz_max < apnz) ptap->abnz_max = apnz;
181f4a743e1SHong Zhang 
18283445d95SHong Zhang     /* if free space is not available, double the total space in the list */
18382412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
1842ba1acd5SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
185f2b054eeSHong Zhang       nspacedouble++;
186f4a743e1SHong Zhang     }
187f4a743e1SHong Zhang 
188f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
18982412ba7SHong Zhang     ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
19082412ba7SHong Zhang     current_space->array           += apnz;
19182412ba7SHong Zhang     current_space->local_used      += apnz;
19282412ba7SHong Zhang     current_space->local_remaining -= apnz;
193f4a743e1SHong Zhang   }
19482412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
195f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
196*f8487c73SHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->abj);CHKERRQ(ierr);
197*f8487c73SHong Zhang   apj = ptap->abj;
198*f8487c73SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->abj);CHKERRQ(ierr);
199f4a743e1SHong Zhang 
200ded4f561SHong Zhang   /* determine symbolic Co=(p->B)^T*AP - send to others */
201ded4f561SHong Zhang   /*----------------------------------------------------*/
202de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
20342c57489SHong Zhang 
204ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
205d0f46423SBarry Smith   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
20683445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
207de0260b3SHong Zhang   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
208de0260b3SHong Zhang   coi[0] = 0;
20942c57489SHong Zhang 
21082412ba7SHong Zhang   /* set initial free space to be 3*pon*max( nnz(AP) per row) */
211*f8487c73SHong Zhang   nnz           = 3*pon*ptap->abnz_max + 1;
212a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
21342c57489SHong Zhang   current_space = free_space;
21442c57489SHong Zhang 
215de0260b3SHong Zhang   for (i=0; i<pon; i++) {
216ded4f561SHong Zhang     nnz  = 0;
21782412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
21882412ba7SHong Zhang     j     = pnz;
219de0260b3SHong Zhang     ptJ   = potj + poti[i+1];
22083445d95SHong Zhang     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
22183445d95SHong Zhang       j--; ptJ--;
22282412ba7SHong Zhang       row  = *ptJ; /* row of AP == col of Pot */
22382412ba7SHong Zhang       apnz = api[row+1] - api[row];
224ded4f561SHong Zhang       Jptr   = apj + api[row];
22583445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
226ded4f561SHong Zhang       ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
227ded4f561SHong Zhang       nnz += nlnk;
22842c57489SHong Zhang     }
22942c57489SHong Zhang 
23083445d95SHong Zhang     /* If free space is not available, double the total space in the list */
231ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
2324238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
23342c57489SHong Zhang     }
23442c57489SHong Zhang 
23542c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
236ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
237ded4f561SHong Zhang     current_space->array           += nnz;
238ded4f561SHong Zhang     current_space->local_used      += nnz;
239ded4f561SHong Zhang     current_space->local_remaining -= nnz;
240ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
24142c57489SHong Zhang   }
242de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
243a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
244de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
24542c57489SHong Zhang 
246ded4f561SHong Zhang   /* send j-array (coj) of Co to other processors */
247ded4f561SHong Zhang   /*----------------------------------------------*/
24842c57489SHong Zhang   /* determine row ownership */
24983445d95SHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
25026283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
2517a2fc3feSBarry Smith   merge->rowmap->n = pn;
2527a2fc3feSBarry Smith   merge->rowmap->bs = 1;
25326283091SBarry Smith   ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
2547a2fc3feSBarry Smith   owners = merge->rowmap->range;
25542c57489SHong Zhang 
25642c57489SHong Zhang   /* determine the number of messages to send, their lengths */
25742c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
25883445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
25942c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
26042c57489SHong Zhang   len_s = merge->len_s;
26142c57489SHong Zhang   merge->nsend = 0;
26283445d95SHong Zhang 
263de0260b3SHong Zhang   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
264de0260b3SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
265de0260b3SHong Zhang 
26683445d95SHong Zhang   proc = 0;
267de0260b3SHong Zhang   for (i=0; i<pon; i++){
268de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
269de0260b3SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
270de0260b3SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
27183445d95SHong Zhang   }
272de0260b3SHong Zhang 
273de0260b3SHong Zhang   len   = 0;  /* max length of buf_si[] */
274de0260b3SHong Zhang   owners_co[0] = 0;
27542c57489SHong Zhang   for (proc=0; proc<size; proc++){
276de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
27783445d95SHong Zhang     if (len_si[proc]){
27842c57489SHong Zhang       merge->nsend++;
27983445d95SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
28042c57489SHong Zhang       len += len_si[proc];
28142c57489SHong Zhang     }
28242c57489SHong Zhang   }
28342c57489SHong Zhang 
284ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
28542c57489SHong Zhang   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
28642c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
28742c57489SHong Zhang 
288ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
289fc42d0c8SSatish Balay   ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr);
290ded4f561SHong Zhang   ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
29142c57489SHong Zhang 
292ded4f561SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
29342c57489SHong Zhang 
29442c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++){
29542c57489SHong Zhang     if (!len_s[proc]) continue;
296de0260b3SHong Zhang     i = owners_co[proc];
297ded4f561SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr);
29842c57489SHong Zhang     k++;
29942c57489SHong Zhang   }
30042c57489SHong Zhang 
301ded4f561SHong Zhang   /* receives and sends of coj are complete */
302ded4f561SHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
303ded4f561SHong Zhang   i = merge->nrecv;
304ded4f561SHong Zhang   while (i--) {
305ded4f561SHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
306ded4f561SHong Zhang   }
307ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3080c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
30982412ba7SHong Zhang 
310ded4f561SHong Zhang   /* send and recv coi */
311ded4f561SHong Zhang   /*-------------------*/
312ded4f561SHong Zhang   ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
31342c57489SHong Zhang 
31442c57489SHong Zhang   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
31542c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
31642c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
31742c57489SHong Zhang     if (!len_s[proc]) continue;
31842c57489SHong Zhang     /* form outgoing message for i-structure:
31942c57489SHong Zhang          buf_si[0]:                 nrows to be sent
32042c57489SHong Zhang                [1:nrows]:           row index (global)
32142c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
32242c57489SHong Zhang     */
32342c57489SHong Zhang     /*-------------------------------------------*/
32442c57489SHong Zhang     nrows = len_si[proc]/2 - 1;
32542c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
32642c57489SHong Zhang     buf_si[0]   = nrows;
32742c57489SHong Zhang     buf_si_i[0] = 0;
32842c57489SHong Zhang     nrows = 0;
329de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
330ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
331ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
332de0260b3SHong Zhang       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
33342c57489SHong Zhang       nrows++;
33442c57489SHong Zhang     }
335ded4f561SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr);
33642c57489SHong Zhang     k++;
33742c57489SHong Zhang     buf_si += len_si[proc];
33842c57489SHong Zhang   }
339ded4f561SHong Zhang   i = merge->nrecv;
340ded4f561SHong Zhang   while (i--) {
341ded4f561SHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
342ded4f561SHong Zhang   }
343ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3440c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
345f2b054eeSHong Zhang   /*
346ae15b995SBarry Smith   ierr = PetscInfo2(A,"nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr);
34742c57489SHong Zhang   for (i=0; i<merge->nrecv; i++){
348ae15b995SBarry Smith     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);
34942c57489SHong Zhang   }
350f2b054eeSHong Zhang   */
35142c57489SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
35242c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
353ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
354ded4f561SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
35542c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
35642c57489SHong Zhang 
357de0260b3SHong Zhang   /* compute the local portion of C (mpi mat) */
358de0260b3SHong Zhang   /*------------------------------------------*/
359ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
360ded4f561SHong Zhang 
36142c57489SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
36242c57489SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
36342c57489SHong Zhang   bi[0] = 0;
36442c57489SHong Zhang 
365ded4f561SHong Zhang   /* set initial free space to be 3*pn*max( nnz(AP) per row) */
366*f8487c73SHong Zhang   nnz           = 3*pn*ptap->abnz_max + 1;
367a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
36842c57489SHong Zhang   current_space = free_space;
36942c57489SHong Zhang 
370687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
37142c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
37242c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
37342c57489SHong Zhang     nrows       = *buf_ri_k[k];
37442c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
37542c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
37642c57489SHong Zhang   }
37742c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
37842c57489SHong Zhang   for (i=0; i<pn; i++) {
379ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
380ded4f561SHong Zhang     nnz = 0;
381ded4f561SHong Zhang     pnz  = pdti[i+1] - pdti[i];
382ded4f561SHong Zhang     j    = pnz;
383ded4f561SHong Zhang     ptJ  = pdtj + pdti[i+1];
384ded4f561SHong Zhang     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
385ded4f561SHong Zhang       j--; ptJ--;
386ded4f561SHong Zhang       row  = *ptJ; /* row of AP == col of Pt */
387ded4f561SHong Zhang       apnz = api[row+1] - api[row];
388ded4f561SHong Zhang       Jptr   = apj + api[row];
389ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
390ded4f561SHong Zhang       ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
391ded4f561SHong Zhang       nnz += nlnk;
392ded4f561SHong Zhang     }
39342c57489SHong Zhang     /* add received col data into lnk */
39442c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
39542c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
396ded4f561SHong Zhang         nzi = *(nextci[k]+1) - *nextci[k];
397ded4f561SHong Zhang         Jptr  = buf_rj[k] + *nextci[k];
398ded4f561SHong Zhang         ierr = PetscLLAdd(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
399ded4f561SHong Zhang         nnz += nlnk;
40042c57489SHong Zhang         nextrow[k]++; nextci[k]++;
40142c57489SHong Zhang       }
40242c57489SHong Zhang     }
40342c57489SHong Zhang 
40442c57489SHong Zhang     /* if free space is not available, make more free space */
405ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
4064238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
40742c57489SHong Zhang     }
40842c57489SHong Zhang     /* copy data into free space, then initialize lnk */
409ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
410ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
411ded4f561SHong Zhang     current_space->array           += nnz;
412ded4f561SHong Zhang     current_space->local_used      += nnz;
413ded4f561SHong Zhang     current_space->local_remaining -= nnz;
414ded4f561SHong Zhang     bi[i+1] = bi[i] + nnz;
41542c57489SHong Zhang   }
416ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
4170572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
41842c57489SHong Zhang 
41942c57489SHong Zhang   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
420a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
42142c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
42242c57489SHong Zhang 
42342c57489SHong Zhang   /* create symbolic parallel matrix B_mpi */
42442c57489SHong Zhang   /*---------------------------------------*/
425f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,&B_mpi);CHKERRQ(ierr);
426f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B_mpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
42742c57489SHong Zhang   ierr = MatSetType(B_mpi,MATMPIAIJ);CHKERRQ(ierr);
42842c57489SHong Zhang   ierr = MatMPIAIJSetPreallocation(B_mpi,0,dnz,0,onz);CHKERRQ(ierr);
42942c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
43042c57489SHong Zhang 
43142c57489SHong Zhang   merge->bi            = bi;
43242c57489SHong Zhang   merge->bj            = bj;
433de0260b3SHong Zhang   merge->coi           = coi;
434de0260b3SHong Zhang   merge->coj           = coj;
43542c57489SHong Zhang   merge->buf_ri        = buf_ri;
43642c57489SHong Zhang   merge->buf_rj        = buf_rj;
437de0260b3SHong Zhang   merge->owners_co     = owners_co;
438dce485f0SBarry Smith   merge->destroy       = B_mpi->ops->destroy;
439dce485f0SBarry Smith   merge->duplicate     = B_mpi->ops->duplicate;
440cc6cb787SHong Zhang 
441cc6cb787SHong Zhang   /* B_mpi is not ready for use - assembly will be done by MatPtAPNumeric() */
442cc6cb787SHong Zhang   B_mpi->assembled      = PETSC_FALSE;
443*f8487c73SHong Zhang   B_mpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
444cc6cb787SHong Zhang   B_mpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
445*f8487c73SHong Zhang   ierr = MatSetBlockSize(B_mpi,1);CHKERRQ(ierr);
44642c57489SHong Zhang 
44742c57489SHong Zhang   /* attach the supporting struct to B_mpi for reuse */
448*f8487c73SHong Zhang   c = (Mat_MPIAIJ*)B_mpi->data;
449*f8487c73SHong Zhang   c->ptap     = ptap;
450*f8487c73SHong Zhang   ptap->merge = merge;
451*f8487c73SHong Zhang 
45242c57489SHong Zhang   *C = B_mpi;
453f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
454f2b054eeSHong Zhang   if (bi[pn] != 0) {
455f2b054eeSHong Zhang     PetscReal afill = ((PetscReal)bi[pn])/(adi[am]+aoi[am]);
456f2b054eeSHong Zhang     if (afill < 1.0) afill = 1.0;
457f2b054eeSHong Zhang     ierr = PetscInfo3(B_mpi,"Reallocs %D; Fill ratio: given %G needed %G when computing A*P.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
458f2b054eeSHong Zhang     ierr = PetscInfo1(B_mpi,"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
459f2b054eeSHong Zhang   } else {
460f2b054eeSHong Zhang     ierr = PetscInfo(B_mpi,"Empty matrix product\n");CHKERRQ(ierr);
461f2b054eeSHong Zhang   }
462f2b054eeSHong Zhang #endif
46342c57489SHong Zhang   PetscFunctionReturn(0);
46442c57489SHong Zhang }
46542c57489SHong Zhang 
46642c57489SHong Zhang #undef __FUNCT__
46742c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
46842c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
46942c57489SHong Zhang {
47042c57489SHong Zhang   PetscErrorCode       ierr;
47142c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
472*f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
47342c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
474de0260b3SHong Zhang   Mat_SeqAIJ           *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
47542c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
476*f8487c73SHong Zhang   Mat_PtAPMPI          *ptap;
4771d633602SHong Zhang   PetscInt             *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
47882412ba7SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
47982412ba7SHong Zhang   PetscInt             i,j,k,anz,pnz,apnz,nextap,row,*cj;
4801d633602SHong Zhang   MatScalar            *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth;
481d0f46423SBarry Smith   PetscInt             am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
4827adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)C)->comm;
48342c57489SHong Zhang   PetscMPIInt          size,rank,taga,*len_s;
48482412ba7SHong Zhang   PetscInt             *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
48542c57489SHong Zhang   PetscInt             **buf_ri,**buf_rj;
48650f5bf86SHong Zhang   PetscInt             cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
48742c57489SHong Zhang   MPI_Request          *s_waits,*r_waits;
48842c57489SHong Zhang   MPI_Status           *status;
48982412ba7SHong Zhang   MatScalar            **abuf_r,*ba_i,*pA,*coa,*ba;
49082412ba7SHong Zhang   PetscInt             *api,*apj,*coi,*coj;
491d0f46423SBarry Smith   PetscInt             *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
492a9d06231SHong Zhang   PetscInt             sparse_axpy;
493f9473cf7SHong Zhang   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;
49442c57489SHong Zhang 
49542c57489SHong Zhang   PetscFunctionBegin;
496*f8487c73SHong Zhang   /* ierr = MPI_Barrier(comm);CHKERRQ(ierr); */
497f9473cf7SHong Zhang   ierr = PetscGetTime(&t00);CHKERRQ(ierr);
49842c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
49942c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
50042c57489SHong Zhang 
501f9473cf7SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
502*f8487c73SHong Zhang   ptap  = c->ptap;
503*f8487c73SHong Zhang   merge = ptap->merge;
5046c6e00e0SHong Zhang 
505*f8487c73SHong Zhang   /* 1) get P_oth = ptap->B_oth  and P_loc = ptap->B_loc */
506f9473cf7SHong Zhang   /*--------------------------------------------------*/
507*f8487c73SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX){
508*f8487c73SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
5096c6e00e0SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
510*f8487c73SHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->B_oth);CHKERRQ(ierr);
511*f8487c73SHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->B_loc);CHKERRQ(ierr);
5126c6e00e0SHong Zhang   }
513*f8487c73SHong Zhang 
514f9473cf7SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
515f9473cf7SHong Zhang   time_matupdate += tf-t0;
51642c57489SHong Zhang 
517f9473cf7SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
518f9473cf7SHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
519f9473cf7SHong Zhang   /*--------------------------------------------------------------*/
52042c57489SHong Zhang   /* get data from symbolic products */
521*f8487c73SHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->B_loc)->data;
522*f8487c73SHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->B_oth)->data;
523a9d06231SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a;
52442c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
52542c57489SHong Zhang 
526de0260b3SHong Zhang   coi = merge->coi; coj = merge->coj;
527de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
528de0260b3SHong Zhang   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
529de0260b3SHong Zhang 
530de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
5317a2fc3feSBarry Smith   owners = merge->rowmap->range;
532de0260b3SHong Zhang   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
533de0260b3SHong Zhang   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
534f9473cf7SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
535f9473cf7SHong Zhang   time_malloc += tf-t0;
53642c57489SHong Zhang 
537*f8487c73SHong Zhang   api = ptap->abi; apj = ptap->abj;
5380496f32cSHong Zhang   /* flag 'sparse_axpy' determines which implementations to be used:
5390496f32cSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fastest, but requires storage of a dense array apa; (default)
5400496f32cSHong Zhang        1: do one sparse axpy - uses same memory as sparse_axpy=0 and might execute less flops
5410496f32cSHong Zhang           (apnz vs. cnz in the outerproduct), slower than case '0' when cnz is not too large than apnz;
5420496f32cSHong Zhang        2: do two sparse axpy in MatPtAPNumeric() - slowest, uses a sparse array apa */
543a9d06231SHong Zhang   /* set default sparse_axpy */
544a9d06231SHong Zhang   sparse_axpy = 1;
5450496f32cSHong Zhang   ierr = PetscOptionsGetInt(PETSC_NULL,"-matptap_sparseaxpy",&sparse_axpy,PETSC_NULL);CHKERRQ(ierr);
546a9d06231SHong Zhang   if (sparse_axpy == 0){  /* Do not perform sparse axpy */
547a9d06231SHong Zhang     /*--------------------------------------------------*/
548a9d06231SHong Zhang     /* malloc apa to store dense row A[i,:]*P */
549a9d06231SHong Zhang     ierr = PetscMalloc((P->cmap->N)*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
550a9d06231SHong Zhang     ierr = PetscMemzero(apa,P->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
551a9d06231SHong Zhang 
552a9d06231SHong Zhang     for (i=0; i<am; i++) {
553a9d06231SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
554a9d06231SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
555a9d06231SHong Zhang       /*------------------------------------------------------------*/
556a9d06231SHong Zhang       apJ = apj + api[i];
557a9d06231SHong Zhang 
558a9d06231SHong Zhang       /* diagonal portion of A */
559a9d06231SHong Zhang       anz = adi[i+1] - adi[i];
560a9d06231SHong Zhang       adj = ad->j + adi[i];
561a9d06231SHong Zhang       ada = ad->a + adi[i];
562a9d06231SHong Zhang       for (j=0; j<anz; j++) {
563a9d06231SHong Zhang         row = adj[j];
564a9d06231SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
565a9d06231SHong Zhang         pj  = pj_loc + pi_loc[row];
566a9d06231SHong Zhang         pa  = pa_loc + pi_loc[row];
567a9d06231SHong Zhang 
568a9d06231SHong Zhang         /* perform dense axpy */
569a9d06231SHong Zhang         for (k=0; k<pnz; k++){
570a9d06231SHong Zhang           apa[pj[k]] += ada[j]*pa[k];
571a9d06231SHong Zhang         }
572a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
573a9d06231SHong Zhang       }
574a9d06231SHong Zhang 
575a9d06231SHong Zhang       /* off-diagonal portion of A */
576a9d06231SHong Zhang       anz = aoi[i+1] - aoi[i];
577a9d06231SHong Zhang       aoj = ao->j + aoi[i];
578a9d06231SHong Zhang       aoa = ao->a + aoi[i];
579a9d06231SHong Zhang       for (j=0; j<anz; j++) {
580a9d06231SHong Zhang         row = aoj[j];
581a9d06231SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
582a9d06231SHong Zhang         pj  = pj_oth + pi_oth[row];
583a9d06231SHong Zhang         pa  = pa_oth + pi_oth[row];
584a9d06231SHong Zhang 
585a9d06231SHong Zhang         /* perform dense axpy */
586a9d06231SHong Zhang         for (k=0; k<pnz; k++){
587a9d06231SHong Zhang           apa[pj[k]] += aoa[j]*pa[k];
588a9d06231SHong Zhang         }
589a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
590a9d06231SHong Zhang       }
591a9d06231SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
592a9d06231SHong Zhang       time_Cseq0 += tf - t0;
593a9d06231SHong Zhang 
594a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
595a9d06231SHong Zhang       /*--------------------------------------------------------------*/
596a9d06231SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
597a9d06231SHong Zhang       apnz = api[i+1] - api[i];
598a9d06231SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
599a9d06231SHong Zhang       pnz = po->i[i+1] - po->i[i];
600a9d06231SHong Zhang       pA  = po->a + po->i[i];
601a9d06231SHong Zhang       for (j=0; j<pnz; j++){
602a9d06231SHong Zhang         cnz = coi[*poJ+1] - coi[*poJ];
603a9d06231SHong Zhang         cj  = coj + coi[*poJ];
604a9d06231SHong Zhang         ca  = coa + coi[*poJ++];
605a9d06231SHong Zhang         /* perform dense axpy */
606a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
607a9d06231SHong Zhang           ca[k] += pA[j]*apa[cj[k]];
608a9d06231SHong Zhang         }
609a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
610a9d06231SHong Zhang       }
611a9d06231SHong Zhang 
612a9d06231SHong Zhang       /* put the value into Cd (diagonal part) */
613a9d06231SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
614a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
615a9d06231SHong Zhang       for (j=0; j<pnz; j++){
616a9d06231SHong Zhang         cnz = bi[*pdJ+1] - bi[*pdJ];
617a9d06231SHong Zhang         cj  = bj + bi[*pdJ];
618a9d06231SHong Zhang         ca  = ba + bi[*pdJ++];
619a9d06231SHong Zhang         /* perform dense axpy */
620a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
621a9d06231SHong Zhang           ca[k] += pA[j]*apa[cj[k]];
622a9d06231SHong Zhang         }
623a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
624a9d06231SHong Zhang       }
625a9d06231SHong Zhang 
626a9d06231SHong Zhang       /* zero the current row of A*P */
627a9d06231SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
628a9d06231SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
629a9d06231SHong Zhang       time_Cseq1 += tf - t0;
630a9d06231SHong Zhang     }
631a9d06231SHong Zhang   } else if (sparse_axpy == 1){  /* Perform one sparse axpy */
632a9d06231SHong Zhang     /*------------------------------------------------------*/
6331d633602SHong Zhang     /* malloc apa to store dense row A[i,:]*P */
6341d633602SHong Zhang     ierr = PetscMalloc((P->cmap->N)*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
6351d633602SHong Zhang     ierr = PetscMemzero(apa,P->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
6361d633602SHong Zhang 
6371d633602SHong Zhang     for (i=0; i<am; i++) {
638f9473cf7SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
639f9473cf7SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
640f9473cf7SHong Zhang       /*------------------------------------------------------------*/
6411d633602SHong Zhang       apJ = apj + api[i];
6421d633602SHong Zhang 
6431d633602SHong Zhang       /* diagonal portion of A */
6441d633602SHong Zhang       anz = adi[i+1] - adi[i];
6451d633602SHong Zhang       adj = ad->j + adi[i];
6461d633602SHong Zhang       ada = ad->a + adi[i];
6471d633602SHong Zhang       for (j=0; j<anz; j++) {
6481d633602SHong Zhang         row = adj[j];
6491d633602SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
6501d633602SHong Zhang         pj  = pj_loc + pi_loc[row];
6511d633602SHong Zhang         pa  = pa_loc + pi_loc[row];
6521d633602SHong Zhang 
6531d633602SHong Zhang         /* perform dense axpy */
6541d633602SHong Zhang         for (k=0; k<pnz; k++){
6551d633602SHong Zhang           apa[pj[k]] += ada[j]*pa[k];
6561d633602SHong Zhang         }
6571d633602SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
6581d633602SHong Zhang       }
6591d633602SHong Zhang 
6601d633602SHong Zhang       /* off-diagonal portion of A */
6611d633602SHong Zhang       anz = aoi[i+1] - aoi[i];
6621d633602SHong Zhang       aoj = ao->j + aoi[i];
6631d633602SHong Zhang       aoa = ao->a + aoi[i];
6641d633602SHong Zhang       for (j=0; j<anz; j++) {
6651d633602SHong Zhang         row = aoj[j];
6661d633602SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
6671d633602SHong Zhang         pj  = pj_oth + pi_oth[row];
6681d633602SHong Zhang         pa  = pa_oth + pi_oth[row];
6691d633602SHong Zhang 
6701d633602SHong Zhang         /* perform dense axpy */
6711d633602SHong Zhang         for (k=0; k<pnz; k++){
6721d633602SHong Zhang           apa[pj[k]] += aoa[j]*pa[k];
6731d633602SHong Zhang         }
6741d633602SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
6751d633602SHong Zhang       }
676f9473cf7SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
677f9473cf7SHong Zhang       time_Cseq0 += tf - t0;
6781d633602SHong Zhang 
679f9473cf7SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
680f9473cf7SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
681f9473cf7SHong Zhang       /*--------------------------------------------------------------*/
6821d633602SHong Zhang       apnz = api[i+1] - api[i];
683f9473cf7SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
684f9473cf7SHong Zhang       pnz = po->i[i+1] - po->i[i];
685a9d06231SHong Zhang       pA  = po->a + po->i[i];
6861d633602SHong Zhang       for (j=0; j<pnz; j++){
6871d633602SHong Zhang         cj  = coj + coi[*poJ];
6881d633602SHong Zhang         ca  = coa + coi[*poJ++];
6891d633602SHong Zhang         /* perform sparse axpy */
6901d633602SHong Zhang         nextap = 0;
6911d633602SHong Zhang         for (k=0; nextap<apnz; k++) {
6921d633602SHong Zhang           if (cj[k]==apJ[nextap]) { /* global column index */
693a9d06231SHong Zhang             ca[k] += pA[j]*apa[cj[k]]; nextap++;
6941d633602SHong Zhang           }
6951d633602SHong Zhang         }
6961d633602SHong Zhang         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
6971d633602SHong Zhang       }
698f9473cf7SHong Zhang 
699f9473cf7SHong Zhang       /* put the value into Cd (diagonal part) */
700f9473cf7SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
701a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
702f9473cf7SHong Zhang       for (j=0; j<pnz; j++){
703f9473cf7SHong Zhang         cj  = bj + bi[*pdJ];
704f9473cf7SHong Zhang         ca  = ba + bi[*pdJ++];
705f9473cf7SHong Zhang         /* perform sparse axpy */
706f9473cf7SHong Zhang         nextap = 0;
707f9473cf7SHong Zhang         for (k=0; nextap<apnz; k++) {
708f9473cf7SHong Zhang           if (cj[k]==apJ[nextap]) { /* global column index */
709a9d06231SHong Zhang             ca[k] += pA[j]*apa[cj[k]];
710a9d06231SHong Zhang             nextap++;
711f9473cf7SHong Zhang           }
712f9473cf7SHong Zhang         }
713f9473cf7SHong Zhang         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
714f9473cf7SHong Zhang       }
715f9473cf7SHong Zhang 
716f9473cf7SHong Zhang       /* zero the current row of A*P */
7171d633602SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
718f9473cf7SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
719f9473cf7SHong Zhang       time_Cseq1 += tf - t0;
7201d633602SHong Zhang     }
7210496f32cSHong Zhang   } else if (sparse_axpy == 2){/* Perform two sparse axpy */
722a9d06231SHong Zhang     /*----------------------------------------------------*/
7231d633602SHong Zhang     /* malloc apa to store sparse row A[i,:]*P */
724*f8487c73SHong Zhang     ierr = PetscMalloc((ptap->abnz_max+1)*sizeof(MatScalar),&apa);CHKERRQ(ierr);
725*f8487c73SHong Zhang     ierr = PetscMemzero(apa,ptap->abnz_max*sizeof(MatScalar));CHKERRQ(ierr);
72642c57489SHong Zhang 
727a9d06231SHong Zhang     pA=pa_loc;
72842c57489SHong Zhang     for (i=0; i<am; i++) {
729a9d06231SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
730f4a743e1SHong Zhang       /* form i-th sparse row of A*P */
73182412ba7SHong Zhang       apnz = api[i+1] - api[i];
73282412ba7SHong Zhang       apJ  = apj + api[i];
73342c57489SHong Zhang       /* diagonal portion of A */
73482412ba7SHong Zhang       anz = adi[i+1] - adi[i];
7351d633602SHong Zhang       adj = ad->j + adi[i];
7361d633602SHong Zhang       ada = ad->a + adi[i];
73782412ba7SHong Zhang       for (j=0; j<anz; j++) {
7381d633602SHong Zhang         row = adj[j];
73982412ba7SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
74082412ba7SHong Zhang         pj  = pj_loc + pi_loc[row];
74182412ba7SHong Zhang         pa  = pa_loc + pi_loc[row];
742f4a743e1SHong Zhang         nextp = 0;
74382412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
74482412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
7451d633602SHong Zhang             apa[k] += ada[j]*pa[nextp++];
74642c57489SHong Zhang           }
74742c57489SHong Zhang         }
748dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
74942c57489SHong Zhang       }
75042c57489SHong Zhang       /* off-diagonal portion of A */
75182412ba7SHong Zhang       anz = aoi[i+1] - aoi[i];
7521d633602SHong Zhang       aoj = ao->j + aoi[i];
7531d633602SHong Zhang       aoa = ao->a + aoi[i];
75482412ba7SHong Zhang       for (j=0; j<anz; j++) {
7551d633602SHong Zhang         row = aoj[j];
75682412ba7SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
75782412ba7SHong Zhang         pj  = pj_oth + pi_oth[row];
75882412ba7SHong Zhang         pa  = pa_oth + pi_oth[row];
759f4a743e1SHong Zhang         nextp = 0;
76082412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
76182412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
7621d633602SHong Zhang             apa[k] += aoa[j]*pa[nextp++];
76342c57489SHong Zhang           }
76442c57489SHong Zhang         }
765dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
76642c57489SHong Zhang       }
767a9d06231SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
768a9d06231SHong Zhang       time_Cseq0 += tf - t0;
76942c57489SHong Zhang 
770a9d06231SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
771a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
772a9d06231SHong Zhang       /*--------------------------------------------------------------*/
77382412ba7SHong Zhang       pnz = pi_loc[i+1] - pi_loc[i];
774f9473cf7SHong Zhang       pJ  = pj_loc + pi_loc[i];
77582412ba7SHong Zhang       for (j=0; j<pnz; j++) {
77642c57489SHong Zhang         nextap = 0;
777f9473cf7SHong Zhang         row    = pJ[j]; /* global index */
77882412ba7SHong Zhang         if (row < pcstart || row >=pcend) { /* put the value into Co */
77982412ba7SHong Zhang           cj  = coj + coi[*poJ];
78082412ba7SHong Zhang           ca  = coa + coi[*poJ++];
78182412ba7SHong Zhang         } else {                            /* put the value into Cd */
78282412ba7SHong Zhang           cj   = bj + bi[*pdJ];
78382412ba7SHong Zhang           ca   = ba + bi[*pdJ++];
78442c57489SHong Zhang         }
78582412ba7SHong Zhang         for (k=0; nextap<apnz; k++) {
7860496f32cSHong Zhang           if (cj[k]==apJ[nextap]) ca[k] += pA[j]*apa[nextap++];
78742c57489SHong Zhang         }
788dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
789de0260b3SHong Zhang       }
7900496f32cSHong Zhang       pA += pnz;
79142c57489SHong Zhang       /* zero the current row info for A*P */
79282412ba7SHong Zhang       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
793a9d06231SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
794a9d06231SHong Zhang       time_Cseq1 += tf - t0;
79542c57489SHong Zhang     }
796a9d06231SHong Zhang   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"sparse_axpy only takes values 0, 1 and 2");
79742c57489SHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
79842c57489SHong Zhang 
799a9d06231SHong Zhang   /* 3) send and recv matrix values coa */
800a9d06231SHong Zhang   /*------------------------------------*/
80142c57489SHong Zhang   buf_ri = merge->buf_ri;
80242c57489SHong Zhang   buf_rj = merge->buf_rj;
80342c57489SHong Zhang   len_s  = merge->len_s;
804fc42d0c8SSatish Balay   ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
80542c57489SHong Zhang   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
80642c57489SHong Zhang 
807a9d06231SHong Zhang   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
80842c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
80942c57489SHong Zhang     if (!len_s[proc]) continue;
810de0260b3SHong Zhang     i = merge->owners_co[proc];
811de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
81242c57489SHong Zhang     k++;
81342c57489SHong Zhang   }
8140c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
8150c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
81642c57489SHong Zhang 
817a9d06231SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
81842c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
81982412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
82042c57489SHong Zhang 
821a9d06231SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
822a9d06231SHong Zhang   /*------------------------------------------------------*/
823f9473cf7SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
824687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
82542c57489SHong Zhang 
82642c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
82742c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
82842c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
82942c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
83082412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
83142c57489SHong Zhang   }
83242c57489SHong Zhang 
833de0260b3SHong Zhang   for (i=0; i<cm; i++) {
83482412ba7SHong Zhang     row = owners[rank] + i; /* global row index of C_seq */
83542c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
836de0260b3SHong Zhang     ba_i = ba + bi[i];
83782412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
83842c57489SHong Zhang     /* add received vals into ba */
83942c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
84042c57489SHong Zhang       /* i-th row */
84142c57489SHong Zhang       if (i == *nextrow[k]) {
84282412ba7SHong Zhang         cnz = *(nextci[k]+1) - *nextci[k];
84382412ba7SHong Zhang         cj  = buf_rj[k] + *(nextci[k]);
84482412ba7SHong Zhang         ca  = abuf_r[k] + *(nextci[k]);
84582412ba7SHong Zhang         nextcj = 0;
84682412ba7SHong Zhang         for (j=0; nextcj<cnz; j++){
84782412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
84882412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
84942c57489SHong Zhang           }
85042c57489SHong Zhang         }
85182412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
85242c57489SHong Zhang       }
85342c57489SHong Zhang     }
85482412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
855dc0b31edSSatish Balay     ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
85642c57489SHong Zhang   }
85742c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
85842c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
859f9473cf7SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
860f9473cf7SHong Zhang   time_setvals += tf-t0;
86142c57489SHong Zhang 
862de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
863c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
86442c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
8650572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
866f9473cf7SHong Zhang 
867f9473cf7SHong Zhang   ierr = PetscGetTime(&tff);CHKERRQ(ierr);
868f9473cf7SHong Zhang   etime += tff - t00;
869f9473cf7SHong Zhang   /*
870a9d06231SHong Zhang   PetscInt prid=0;
871a9d06231SHong Zhang   if (rank == prid){
872f9473cf7SHong Zhang    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);
873a9d06231SHong Zhang   }
874f9473cf7SHong Zhang    */
87542c57489SHong Zhang   PetscFunctionReturn(0);
87642c57489SHong Zhang }
877