xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision d16ca5f04817f8f31fe12367dc76ec08c930022f)
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__
14f8487c73SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP"
15f8487c73SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
16cc6cb787SHong Zhang {
17cc6cb787SHong Zhang   PetscErrorCode       ierr;
18f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data;
19f8487c73SHong Zhang   Mat_PtAPMPI          *ptap=a->ptap;
20cc6cb787SHong Zhang 
21cc6cb787SHong Zhang   PetscFunctionBegin;
22f8487c73SHong Zhang   if (ptap){
23f8487c73SHong Zhang     ierr = PetscFree2(ptap->startsj,ptap->startsj_r);CHKERRQ(ierr);
24f8487c73SHong Zhang     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
25a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
26a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
27a1a4e84aSHong Zhang     ierr = PetscFree(ptap->api);CHKERRQ(ierr);
28a1a4e84aSHong Zhang     ierr = PetscFree(ptap->apj);CHKERRQ(ierr);
29f8487c73SHong Zhang   }
30f8487c73SHong Zhang   if (ptap->merge) {
31f8487c73SHong 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);
46f8487c73SHong Zhang     ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
47bf0cc555SLisandro Dalcin   }
48f8487c73SHong 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;
101f8487c73SHong Zhang   Mat_PtAPMPI          *ptap;
102a1a86e44SBarry Smith   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
103f8487c73SHong 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;
109*d16ca5f0SHong Zhang   PetscInt             am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->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;
120598bc09dSHong Zhang   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j;
121*d16ca5f0SHong Zhang   PetscReal            afill=1.0,afill_tmp;
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 
127f8487c73SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
128f8487c73SHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
129a1a4e84aSHong Zhang   ptap->api=PETSC_NULL; ptap->apj=PETSC_NULL;
130f8487c73SHong Zhang   ptap->abnz_max = 0;
131f8487c73SHong 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 */
134a1a4e84aSHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
1356c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
136a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
1376c6e00e0SHong Zhang 
138a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
139a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_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);
146a1a4e84aSHong Zhang   ptap->api = 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 
153*d16ca5f0SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
154*d16ca5f0SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&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;
180f8487c73SHong 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) */
196a1a4e84aSHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ptap->apj);CHKERRQ(ierr);
197a1a4e84aSHong Zhang   apj = ptap->apj;
198a1a4e84aSHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,ptap->apj);CHKERRQ(ierr);
199*d16ca5f0SHong Zhang   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]);
200*d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
201f4a743e1SHong Zhang 
202ded4f561SHong Zhang   /* determine symbolic Co=(p->B)^T*AP - send to others */
203ded4f561SHong Zhang   /*----------------------------------------------------*/
204de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
20542c57489SHong Zhang 
206ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
207d0f46423SBarry Smith   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
20883445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
209de0260b3SHong Zhang   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
210de0260b3SHong Zhang   coi[0] = 0;
21142c57489SHong Zhang 
212*d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
213*d16ca5f0SHong Zhang   nnz           = fill*(poti[pon] + api[am]);
214a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
21542c57489SHong Zhang   current_space = free_space;
21642c57489SHong Zhang 
217de0260b3SHong Zhang   for (i=0; i<pon; i++) {
218ded4f561SHong Zhang     nnz  = 0;
21982412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
22082412ba7SHong Zhang     j     = pnz;
221de0260b3SHong Zhang     ptJ   = potj + poti[i+1];
22283445d95SHong Zhang     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
22383445d95SHong Zhang       j--; ptJ--;
22482412ba7SHong Zhang       row  = *ptJ; /* row of AP == col of Pot */
22582412ba7SHong Zhang       apnz = api[row+1] - api[row];
226ded4f561SHong Zhang       Jptr   = apj + api[row];
22783445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
228ded4f561SHong Zhang       ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
229ded4f561SHong Zhang       nnz += nlnk;
23042c57489SHong Zhang     }
23142c57489SHong Zhang 
23283445d95SHong Zhang     /* If free space is not available, double the total space in the list */
233ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
2344238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
235*d16ca5f0SHong Zhang       nspacedouble++;
23642c57489SHong Zhang     }
23742c57489SHong Zhang 
23842c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
239ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
240ded4f561SHong Zhang     current_space->array           += nnz;
241ded4f561SHong Zhang     current_space->local_used      += nnz;
242ded4f561SHong Zhang     current_space->local_remaining -= nnz;
243ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
24442c57489SHong Zhang   }
245de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
246a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
247*d16ca5f0SHong Zhang   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]);
248*d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
249de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
25042c57489SHong Zhang 
251ded4f561SHong Zhang   /* send j-array (coj) of Co to other processors */
252ded4f561SHong Zhang   /*----------------------------------------------*/
25342c57489SHong Zhang   /* determine row ownership */
25483445d95SHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
25526283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
2567a2fc3feSBarry Smith   merge->rowmap->n = pn;
2577a2fc3feSBarry Smith   merge->rowmap->bs = 1;
25826283091SBarry Smith   ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
2597a2fc3feSBarry Smith   owners = merge->rowmap->range;
26042c57489SHong Zhang 
26142c57489SHong Zhang   /* determine the number of messages to send, their lengths */
26242c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
26383445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
26442c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
26542c57489SHong Zhang   len_s = merge->len_s;
26642c57489SHong Zhang   merge->nsend = 0;
26783445d95SHong Zhang 
268de0260b3SHong Zhang   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
269de0260b3SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
270de0260b3SHong Zhang 
27183445d95SHong Zhang   proc = 0;
272de0260b3SHong Zhang   for (i=0; i<pon; i++){
273de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
274de0260b3SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
275de0260b3SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
27683445d95SHong Zhang   }
277de0260b3SHong Zhang 
278de0260b3SHong Zhang   len   = 0;  /* max length of buf_si[] */
279de0260b3SHong Zhang   owners_co[0] = 0;
28042c57489SHong Zhang   for (proc=0; proc<size; proc++){
281de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
28283445d95SHong Zhang     if (len_si[proc]){
28342c57489SHong Zhang       merge->nsend++;
28483445d95SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
28542c57489SHong Zhang       len += len_si[proc];
28642c57489SHong Zhang     }
28742c57489SHong Zhang   }
28842c57489SHong Zhang 
289ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
29042c57489SHong Zhang   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
29142c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
29242c57489SHong Zhang 
293ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
294fc42d0c8SSatish Balay   ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr);
295ded4f561SHong Zhang   ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
29642c57489SHong Zhang 
297ded4f561SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
29842c57489SHong Zhang 
29942c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++){
30042c57489SHong Zhang     if (!len_s[proc]) continue;
301de0260b3SHong Zhang     i = owners_co[proc];
302ded4f561SHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr);
30342c57489SHong Zhang     k++;
30442c57489SHong Zhang   }
30542c57489SHong Zhang 
306ded4f561SHong Zhang   /* receives and sends of coj are complete */
307ded4f561SHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
308ded4f561SHong Zhang   i = merge->nrecv;
309ded4f561SHong Zhang   while (i--) {
310ded4f561SHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
311ded4f561SHong Zhang   }
312ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3130c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
31482412ba7SHong Zhang 
315ded4f561SHong Zhang   /* send and recv coi */
316ded4f561SHong Zhang   /*-------------------*/
317ded4f561SHong Zhang   ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
31842c57489SHong Zhang 
31942c57489SHong Zhang   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
32042c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
32142c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
32242c57489SHong Zhang     if (!len_s[proc]) continue;
32342c57489SHong Zhang     /* form outgoing message for i-structure:
32442c57489SHong Zhang          buf_si[0]:                 nrows to be sent
32542c57489SHong Zhang                [1:nrows]:           row index (global)
32642c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
32742c57489SHong Zhang     */
32842c57489SHong Zhang     /*-------------------------------------------*/
32942c57489SHong Zhang     nrows = len_si[proc]/2 - 1;
33042c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
33142c57489SHong Zhang     buf_si[0]   = nrows;
33242c57489SHong Zhang     buf_si_i[0] = 0;
33342c57489SHong Zhang     nrows = 0;
334de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
335ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
336ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
337de0260b3SHong Zhang       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
33842c57489SHong Zhang       nrows++;
33942c57489SHong Zhang     }
340ded4f561SHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr);
34142c57489SHong Zhang     k++;
34242c57489SHong Zhang     buf_si += len_si[proc];
34342c57489SHong Zhang   }
344ded4f561SHong Zhang   i = merge->nrecv;
345ded4f561SHong Zhang   while (i--) {
346ded4f561SHong Zhang     ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr);
347ded4f561SHong Zhang   }
348ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3490c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
350f2b054eeSHong Zhang   /*
351ae15b995SBarry Smith   ierr = PetscInfo2(A,"nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr);
35242c57489SHong Zhang   for (i=0; i<merge->nrecv; i++){
353ae15b995SBarry 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);
35442c57489SHong Zhang   }
355f2b054eeSHong Zhang   */
35642c57489SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
35742c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
358ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
359ded4f561SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
36042c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
36142c57489SHong Zhang 
362de0260b3SHong Zhang   /* compute the local portion of C (mpi mat) */
363de0260b3SHong Zhang   /*------------------------------------------*/
364ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
365ded4f561SHong Zhang 
36642c57489SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
36742c57489SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
36842c57489SHong Zhang   bi[0] = 0;
36942c57489SHong Zhang 
370*d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
371*d16ca5f0SHong Zhang   nnz           = fill*(pi_loc[pm] + api[am]);
372a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
37342c57489SHong Zhang   current_space = free_space;
37442c57489SHong Zhang 
375687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
37642c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
37742c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
37842c57489SHong Zhang     nrows       = *buf_ri_k[k];
37942c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
38042c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
38142c57489SHong Zhang   }
38242c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
38342c57489SHong Zhang   for (i=0; i<pn; i++) {
384ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
385ded4f561SHong Zhang     nnz = 0;
386ded4f561SHong Zhang     pnz  = pdti[i+1] - pdti[i];
387ded4f561SHong Zhang     j    = pnz;
388ded4f561SHong Zhang     ptJ  = pdtj + pdti[i+1];
389ded4f561SHong Zhang     while (j){/* assume cols are almost in increasing order, starting from its end saves computation */
390ded4f561SHong Zhang       j--; ptJ--;
391ded4f561SHong Zhang       row  = *ptJ; /* row of AP == col of Pt */
392ded4f561SHong Zhang       apnz = api[row+1] - api[row];
393ded4f561SHong Zhang       Jptr   = apj + api[row];
394ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
395ded4f561SHong Zhang       ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
396ded4f561SHong Zhang       nnz += nlnk;
397ded4f561SHong Zhang     }
39842c57489SHong Zhang     /* add received col data into lnk */
39942c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
40042c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
401ded4f561SHong Zhang         nzi = *(nextci[k]+1) - *nextci[k];
402ded4f561SHong Zhang         Jptr  = buf_rj[k] + *nextci[k];
403ded4f561SHong Zhang         ierr = PetscLLAdd(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
404ded4f561SHong Zhang         nnz += nlnk;
40542c57489SHong Zhang         nextrow[k]++; nextci[k]++;
40642c57489SHong Zhang       }
40742c57489SHong Zhang     }
40842c57489SHong Zhang 
40942c57489SHong Zhang     /* if free space is not available, make more free space */
410ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
4114238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
412*d16ca5f0SHong Zhang       nspacedouble++;
41342c57489SHong Zhang     }
41442c57489SHong Zhang     /* copy data into free space, then initialize lnk */
415ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
416ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
417ded4f561SHong Zhang     current_space->array           += nnz;
418ded4f561SHong Zhang     current_space->local_used      += nnz;
419ded4f561SHong Zhang     current_space->local_remaining -= nnz;
420ded4f561SHong Zhang     bi[i+1] = bi[i] + nnz;
42142c57489SHong Zhang   }
422ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
4230572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
42442c57489SHong Zhang 
42542c57489SHong Zhang   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
426a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
427*d16ca5f0SHong Zhang   afill_tmp = (PetscReal)bi[pn]/(pi_loc[pm] + api[am]);
428*d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
42942c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
43042c57489SHong Zhang 
43142c57489SHong Zhang   /* create symbolic parallel matrix B_mpi */
43242c57489SHong Zhang   /*---------------------------------------*/
433f69a0ea3SMatthew Knepley   ierr = MatCreate(comm,&B_mpi);CHKERRQ(ierr);
434f69a0ea3SMatthew Knepley   ierr = MatSetSizes(B_mpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
43542c57489SHong Zhang   ierr = MatSetType(B_mpi,MATMPIAIJ);CHKERRQ(ierr);
43642c57489SHong Zhang   ierr = MatMPIAIJSetPreallocation(B_mpi,0,dnz,0,onz);CHKERRQ(ierr);
43742c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
43842c57489SHong Zhang 
43942c57489SHong Zhang   merge->bi            = bi;
44042c57489SHong Zhang   merge->bj            = bj;
441de0260b3SHong Zhang   merge->coi           = coi;
442de0260b3SHong Zhang   merge->coj           = coj;
44342c57489SHong Zhang   merge->buf_ri        = buf_ri;
44442c57489SHong Zhang   merge->buf_rj        = buf_rj;
445de0260b3SHong Zhang   merge->owners_co     = owners_co;
446dce485f0SBarry Smith   merge->destroy       = B_mpi->ops->destroy;
447dce485f0SBarry Smith   merge->duplicate     = B_mpi->ops->duplicate;
448cc6cb787SHong Zhang 
449cc6cb787SHong Zhang   /* B_mpi is not ready for use - assembly will be done by MatPtAPNumeric() */
450cc6cb787SHong Zhang   B_mpi->assembled      = PETSC_FALSE;
451f8487c73SHong Zhang   B_mpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
452cc6cb787SHong Zhang   B_mpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
453f8487c73SHong Zhang   ierr = MatSetBlockSize(B_mpi,1);CHKERRQ(ierr);
45442c57489SHong Zhang 
45542c57489SHong Zhang   /* attach the supporting struct to B_mpi for reuse */
456f8487c73SHong Zhang   c = (Mat_MPIAIJ*)B_mpi->data;
457f8487c73SHong Zhang   c->ptap     = ptap;
458f8487c73SHong Zhang   ptap->merge = merge;
459f8487c73SHong Zhang 
46042c57489SHong Zhang   *C = B_mpi;
461f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
462f2b054eeSHong Zhang   if (bi[pn] != 0) {
463*d16ca5f0SHong Zhang     ierr = PetscInfo3(B_mpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
464f2b054eeSHong Zhang     ierr = PetscInfo1(B_mpi,"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
465f2b054eeSHong Zhang   } else {
466f2b054eeSHong Zhang     ierr = PetscInfo(B_mpi,"Empty matrix product\n");CHKERRQ(ierr);
467f2b054eeSHong Zhang   }
468f2b054eeSHong Zhang #endif
46942c57489SHong Zhang   PetscFunctionReturn(0);
47042c57489SHong Zhang }
47142c57489SHong Zhang 
47242c57489SHong Zhang #undef __FUNCT__
47342c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
47442c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
47542c57489SHong Zhang {
47642c57489SHong Zhang   PetscErrorCode       ierr;
47742c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
478f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
47942c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
480de0260b3SHong Zhang   Mat_SeqAIJ           *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
48142c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
482f8487c73SHong Zhang   Mat_PtAPMPI          *ptap;
4831d633602SHong Zhang   PetscInt             *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
48482412ba7SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
48582412ba7SHong Zhang   PetscInt             i,j,k,anz,pnz,apnz,nextap,row,*cj;
486e900a757SHong Zhang   MatScalar            *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
487d0f46423SBarry Smith   PetscInt             am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
4887adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)C)->comm;
48942c57489SHong Zhang   PetscMPIInt          size,rank,taga,*len_s;
49082412ba7SHong Zhang   PetscInt             *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
49142c57489SHong Zhang   PetscInt             **buf_ri,**buf_rj;
49250f5bf86SHong Zhang   PetscInt             cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
49342c57489SHong Zhang   MPI_Request          *s_waits,*r_waits;
49442c57489SHong Zhang   MPI_Status           *status;
49582412ba7SHong Zhang   MatScalar            **abuf_r,*ba_i,*pA,*coa,*ba;
49682412ba7SHong Zhang   PetscInt             *api,*apj,*coi,*coj;
497d0f46423SBarry Smith   PetscInt             *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
498a9d06231SHong Zhang   PetscInt             sparse_axpy;
499f9473cf7SHong 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;
50042c57489SHong Zhang 
50142c57489SHong Zhang   PetscFunctionBegin;
502f8487c73SHong Zhang   /* ierr = MPI_Barrier(comm);CHKERRQ(ierr); */
503f9473cf7SHong Zhang   ierr = PetscGetTime(&t00);CHKERRQ(ierr);
50442c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
50542c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
50642c57489SHong Zhang 
507f9473cf7SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
508f8487c73SHong Zhang   ptap  = c->ptap;
509f8487c73SHong Zhang   merge = ptap->merge;
5106c6e00e0SHong Zhang 
511a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
512f9473cf7SHong Zhang   /*--------------------------------------------------*/
513f8487c73SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX){
514f8487c73SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
5156c6e00e0SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
516a1a4e84aSHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
517a1a4e84aSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
5186c6e00e0SHong Zhang   }
519f8487c73SHong Zhang 
520f9473cf7SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
521f9473cf7SHong Zhang   time_matupdate += tf-t0;
52242c57489SHong Zhang 
523f9473cf7SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
524f9473cf7SHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
525f9473cf7SHong Zhang   /*--------------------------------------------------------------*/
52642c57489SHong Zhang   /* get data from symbolic products */
527a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
528a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
529a9d06231SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a;
53042c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
53142c57489SHong Zhang 
532de0260b3SHong Zhang   coi = merge->coi; coj = merge->coj;
533de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
534de0260b3SHong Zhang   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
535de0260b3SHong Zhang 
536de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
5377a2fc3feSBarry Smith   owners = merge->rowmap->range;
538de0260b3SHong Zhang   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
539de0260b3SHong Zhang   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
540f9473cf7SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
541f9473cf7SHong Zhang   time_malloc += tf-t0;
54242c57489SHong Zhang 
543a1a4e84aSHong Zhang   api = ptap->api; apj = ptap->apj;
5440496f32cSHong Zhang   /* flag 'sparse_axpy' determines which implementations to be used:
5450496f32cSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fastest, but requires storage of a dense array apa; (default)
5460496f32cSHong Zhang        1: do one sparse axpy - uses same memory as sparse_axpy=0 and might execute less flops
5470496f32cSHong Zhang           (apnz vs. cnz in the outerproduct), slower than case '0' when cnz is not too large than apnz;
5480496f32cSHong Zhang        2: do two sparse axpy in MatPtAPNumeric() - slowest, uses a sparse array apa */
549a9d06231SHong Zhang   /* set default sparse_axpy */
550e900a757SHong Zhang   sparse_axpy = 0;
5510496f32cSHong Zhang   ierr = PetscOptionsGetInt(PETSC_NULL,"-matptap_sparseaxpy",&sparse_axpy,PETSC_NULL);CHKERRQ(ierr);
552a9d06231SHong Zhang   if (sparse_axpy == 0){  /* Do not perform sparse axpy */
553a9d06231SHong Zhang     /*--------------------------------------------------*/
554a9d06231SHong Zhang     /* malloc apa to store dense row A[i,:]*P */
555a9d06231SHong Zhang     ierr = PetscMalloc((P->cmap->N)*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
556a9d06231SHong Zhang     ierr = PetscMemzero(apa,P->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
557a9d06231SHong Zhang 
558a9d06231SHong Zhang     for (i=0; i<am; i++) {
559a9d06231SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
560a9d06231SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
561a9d06231SHong Zhang       /*------------------------------------------------------------*/
562a9d06231SHong Zhang       apJ = apj + api[i];
563a9d06231SHong Zhang 
564a9d06231SHong Zhang       /* diagonal portion of A */
565a9d06231SHong Zhang       anz = adi[i+1] - adi[i];
566a9d06231SHong Zhang       adj = ad->j + adi[i];
567a9d06231SHong Zhang       ada = ad->a + adi[i];
568a9d06231SHong Zhang       for (j=0; j<anz; j++) {
569a9d06231SHong Zhang         row = adj[j];
570a9d06231SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
571a9d06231SHong Zhang         pj  = pj_loc + pi_loc[row];
572a9d06231SHong Zhang         pa  = pa_loc + pi_loc[row];
573a9d06231SHong Zhang 
574a9d06231SHong Zhang         /* perform dense axpy */
575e900a757SHong Zhang         valtmp = ada[j];
576a9d06231SHong Zhang         for (k=0; k<pnz; k++){
577e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
578a9d06231SHong Zhang         }
579a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
580a9d06231SHong Zhang       }
581a9d06231SHong Zhang 
582a9d06231SHong Zhang       /* off-diagonal portion of A */
583a9d06231SHong Zhang       anz = aoi[i+1] - aoi[i];
584a9d06231SHong Zhang       aoj = ao->j + aoi[i];
585a9d06231SHong Zhang       aoa = ao->a + aoi[i];
586a9d06231SHong Zhang       for (j=0; j<anz; j++) {
587a9d06231SHong Zhang         row = aoj[j];
588a9d06231SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
589a9d06231SHong Zhang         pj  = pj_oth + pi_oth[row];
590a9d06231SHong Zhang         pa  = pa_oth + pi_oth[row];
591a9d06231SHong Zhang 
592a9d06231SHong Zhang         /* perform dense axpy */
593e900a757SHong Zhang         valtmp = aoa[j];
594a9d06231SHong Zhang         for (k=0; k<pnz; k++){
595e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
596a9d06231SHong Zhang         }
597a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
598a9d06231SHong Zhang       }
599a9d06231SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
600a9d06231SHong Zhang       time_Cseq0 += tf - t0;
601a9d06231SHong Zhang 
602a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
603a9d06231SHong Zhang       /*--------------------------------------------------------------*/
604a9d06231SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
605a9d06231SHong Zhang       apnz = api[i+1] - api[i];
606a9d06231SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
607a9d06231SHong Zhang       pnz = po->i[i+1] - po->i[i];
608e900a757SHong Zhang       poJ = po->j + po->i[i];
609a9d06231SHong Zhang       pA  = po->a + po->i[i];
610a9d06231SHong Zhang       for (j=0; j<pnz; j++){
611e900a757SHong Zhang         row = poJ[j];
612e900a757SHong Zhang         cnz = coi[row+1] - coi[row];
613e900a757SHong Zhang         cj  = coj + coi[row];
614e900a757SHong Zhang         ca  = coa + coi[row];
615a9d06231SHong Zhang         /* perform dense axpy */
616e900a757SHong Zhang         valtmp = pA[j];
617a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
618e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
619a9d06231SHong Zhang         }
620a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
621a9d06231SHong Zhang       }
622a9d06231SHong Zhang 
623a9d06231SHong Zhang       /* put the value into Cd (diagonal part) */
624a9d06231SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
625e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
626a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
627a9d06231SHong Zhang       for (j=0; j<pnz; j++){
628e900a757SHong Zhang         row = pdJ[j];
629e900a757SHong Zhang         cnz = bi[row+1] - bi[row];
630e900a757SHong Zhang         cj  = bj + bi[row];
631e900a757SHong Zhang         ca  = ba + bi[row];
632a9d06231SHong Zhang         /* perform dense axpy */
633e900a757SHong Zhang         valtmp = pA[j];
634a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
635e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
636a9d06231SHong Zhang         }
637a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
638a9d06231SHong Zhang       }
639a9d06231SHong Zhang 
640a9d06231SHong Zhang       /* zero the current row of A*P */
641a9d06231SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
642a9d06231SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
643a9d06231SHong Zhang       time_Cseq1 += tf - t0;
644a9d06231SHong Zhang     }
645a9d06231SHong Zhang   } else if (sparse_axpy == 1){  /* Perform one sparse axpy */
646a9d06231SHong Zhang     /*------------------------------------------------------*/
6471d633602SHong Zhang     /* malloc apa to store dense row A[i,:]*P */
6481d633602SHong Zhang     ierr = PetscMalloc((P->cmap->N)*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
6491d633602SHong Zhang     ierr = PetscMemzero(apa,P->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
6501d633602SHong Zhang 
6511d633602SHong Zhang     for (i=0; i<am; i++) {
652f9473cf7SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
653f9473cf7SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
654f9473cf7SHong Zhang       /*------------------------------------------------------------*/
6551d633602SHong Zhang       apJ = apj + api[i];
6561d633602SHong Zhang 
6571d633602SHong Zhang       /* diagonal portion of A */
6581d633602SHong Zhang       anz = adi[i+1] - adi[i];
6591d633602SHong Zhang       adj = ad->j + adi[i];
6601d633602SHong Zhang       ada = ad->a + adi[i];
6611d633602SHong Zhang       for (j=0; j<anz; j++) {
6621d633602SHong Zhang         row = adj[j];
6631d633602SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
6641d633602SHong Zhang         pj  = pj_loc + pi_loc[row];
6651d633602SHong Zhang         pa  = pa_loc + pi_loc[row];
6661d633602SHong Zhang 
6671d633602SHong Zhang         /* perform dense axpy */
668e900a757SHong Zhang         valtmp = ada[j];
6691d633602SHong Zhang         for (k=0; k<pnz; k++){
670e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
6711d633602SHong Zhang         }
6721d633602SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
6731d633602SHong Zhang       }
6741d633602SHong Zhang 
6751d633602SHong Zhang       /* off-diagonal portion of A */
6761d633602SHong Zhang       anz = aoi[i+1] - aoi[i];
6771d633602SHong Zhang       aoj = ao->j + aoi[i];
6781d633602SHong Zhang       aoa = ao->a + aoi[i];
6791d633602SHong Zhang       for (j=0; j<anz; j++) {
6801d633602SHong Zhang         row = aoj[j];
6811d633602SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
6821d633602SHong Zhang         pj  = pj_oth + pi_oth[row];
6831d633602SHong Zhang         pa  = pa_oth + pi_oth[row];
6841d633602SHong Zhang 
6851d633602SHong Zhang         /* perform dense axpy */
686e900a757SHong Zhang         valtmp = aoa[j];
6871d633602SHong Zhang         for (k=0; k<pnz; k++){
688e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
6891d633602SHong Zhang         }
6901d633602SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
6911d633602SHong Zhang       }
692f9473cf7SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
693f9473cf7SHong Zhang       time_Cseq0 += tf - t0;
6941d633602SHong Zhang 
695f9473cf7SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
696f9473cf7SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
697f9473cf7SHong Zhang       /*--------------------------------------------------------------*/
6981d633602SHong Zhang       apnz = api[i+1] - api[i];
699f9473cf7SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
700f9473cf7SHong Zhang       pnz = po->i[i+1] - po->i[i];
701e900a757SHong Zhang       poJ = po->j + po->i[i];
702a9d06231SHong Zhang       pA  = po->a + po->i[i];
7031d633602SHong Zhang       for (j=0; j<pnz; j++){
704e900a757SHong Zhang         row    = poJ[j];
705e900a757SHong Zhang         cj     = coj + coi[row];
706e900a757SHong Zhang         ca     = coa + coi[row];
707e900a757SHong Zhang         valtmp = pA[j];
7081d633602SHong Zhang         /* perform sparse axpy */
7091d633602SHong Zhang         nextap = 0;
7101d633602SHong Zhang         for (k=0; nextap<apnz; k++) {
7111d633602SHong Zhang           if (cj[k]==apJ[nextap]) { /* global column index */
712e900a757SHong Zhang             ca[k] += valtmp*apa[cj[k]]; nextap++;
7131d633602SHong Zhang           }
7141d633602SHong Zhang         }
7151d633602SHong Zhang         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
7161d633602SHong Zhang       }
717f9473cf7SHong Zhang 
718f9473cf7SHong Zhang       /* put the value into Cd (diagonal part) */
719f9473cf7SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
720e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
721a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
722f9473cf7SHong Zhang       for (j=0; j<pnz; j++){
723e900a757SHong Zhang         row    = pdJ[j];
724e900a757SHong Zhang         cj     = bj + bi[row];
725e900a757SHong Zhang         ca     = ba + bi[row];
726e900a757SHong Zhang         valtmp = pA[j];
727f9473cf7SHong Zhang         /* perform sparse axpy */
728f9473cf7SHong Zhang         nextap = 0;
729f9473cf7SHong Zhang         for (k=0; nextap<apnz; k++) {
730f9473cf7SHong Zhang           if (cj[k]==apJ[nextap]) { /* global column index */
731e900a757SHong Zhang             ca[k] += valtmp*apa[cj[k]];
732a9d06231SHong Zhang             nextap++;
733f9473cf7SHong Zhang           }
734f9473cf7SHong Zhang         }
735f9473cf7SHong Zhang         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
736f9473cf7SHong Zhang       }
737f9473cf7SHong Zhang 
738f9473cf7SHong Zhang       /* zero the current row of A*P */
7391d633602SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
740f9473cf7SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
741f9473cf7SHong Zhang       time_Cseq1 += tf - t0;
7421d633602SHong Zhang     }
7430496f32cSHong Zhang   } else if (sparse_axpy == 2){/* Perform two sparse axpy */
744a9d06231SHong Zhang     /*----------------------------------------------------*/
7451d633602SHong Zhang     /* malloc apa to store sparse row A[i,:]*P */
746f8487c73SHong Zhang     ierr = PetscMalloc((ptap->abnz_max+1)*sizeof(MatScalar),&apa);CHKERRQ(ierr);
747f8487c73SHong Zhang     ierr = PetscMemzero(apa,ptap->abnz_max*sizeof(MatScalar));CHKERRQ(ierr);
74842c57489SHong Zhang 
749a9d06231SHong Zhang     pA=pa_loc;
75042c57489SHong Zhang     for (i=0; i<am; i++) {
751a9d06231SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
752f4a743e1SHong Zhang       /* form i-th sparse row of A*P */
75382412ba7SHong Zhang       apnz = api[i+1] - api[i];
75482412ba7SHong Zhang       apJ  = apj + api[i];
75542c57489SHong Zhang       /* diagonal portion of A */
75682412ba7SHong Zhang       anz = adi[i+1] - adi[i];
7571d633602SHong Zhang       adj = ad->j + adi[i];
7581d633602SHong Zhang       ada = ad->a + adi[i];
75982412ba7SHong Zhang       for (j=0; j<anz; j++) {
7601d633602SHong Zhang         row = adj[j];
76182412ba7SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
76282412ba7SHong Zhang         pj  = pj_loc + pi_loc[row];
76382412ba7SHong Zhang         pa  = pa_loc + pi_loc[row];
764e900a757SHong Zhang         valtmp = ada[j];
765f4a743e1SHong Zhang         nextp  = 0;
76682412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
76782412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
768e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
76942c57489SHong Zhang           }
77042c57489SHong Zhang         }
771dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
77242c57489SHong Zhang       }
77342c57489SHong Zhang       /* off-diagonal portion of A */
77482412ba7SHong Zhang       anz = aoi[i+1] - aoi[i];
7751d633602SHong Zhang       aoj = ao->j + aoi[i];
7761d633602SHong Zhang       aoa = ao->a + aoi[i];
77782412ba7SHong Zhang       for (j=0; j<anz; j++) {
7781d633602SHong Zhang         row = aoj[j];
77982412ba7SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
78082412ba7SHong Zhang         pj  = pj_oth + pi_oth[row];
78182412ba7SHong Zhang         pa  = pa_oth + pi_oth[row];
782e900a757SHong Zhang         valtmp = aoa[j];
783f4a743e1SHong Zhang         nextp  = 0;
78482412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
78582412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
786e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
78742c57489SHong Zhang           }
78842c57489SHong Zhang         }
789dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
79042c57489SHong Zhang       }
791a9d06231SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
792a9d06231SHong Zhang       time_Cseq0 += tf - t0;
79342c57489SHong Zhang 
794a9d06231SHong Zhang       ierr = PetscGetTime(&t0);CHKERRQ(ierr);
795a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
796a9d06231SHong Zhang       /*--------------------------------------------------------------*/
79782412ba7SHong Zhang       pnz = pi_loc[i+1] - pi_loc[i];
798f9473cf7SHong Zhang       pJ  = pj_loc + pi_loc[i];
79982412ba7SHong Zhang       for (j=0; j<pnz; j++) {
80042c57489SHong Zhang         nextap = 0;
801f9473cf7SHong Zhang         row    = pJ[j]; /* global index */
80282412ba7SHong Zhang         if (row < pcstart || row >=pcend) { /* put the value into Co */
803e900a757SHong Zhang           row = *poJ;
804e900a757SHong Zhang           cj  = coj + coi[row];
805e900a757SHong Zhang           ca  = coa + coi[row]; poJ++;
80682412ba7SHong Zhang         } else {                            /* put the value into Cd */
807e900a757SHong Zhang           row  = *pdJ;
808e900a757SHong Zhang           cj   = bj + bi[row];
809e900a757SHong Zhang           ca   = ba + bi[row]; pdJ++;
81042c57489SHong Zhang         }
811e900a757SHong Zhang         valtmp = pA[j];
81282412ba7SHong Zhang         for (k=0; nextap<apnz; k++) {
813e900a757SHong Zhang           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
81442c57489SHong Zhang         }
815dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
816de0260b3SHong Zhang       }
8170496f32cSHong Zhang       pA += pnz;
81842c57489SHong Zhang       /* zero the current row info for A*P */
81982412ba7SHong Zhang       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
820a9d06231SHong Zhang       ierr = PetscGetTime(&tf);CHKERRQ(ierr);
821a9d06231SHong Zhang       time_Cseq1 += tf - t0;
82242c57489SHong Zhang     }
823a9d06231SHong Zhang   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"sparse_axpy only takes values 0, 1 and 2");
82442c57489SHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
82542c57489SHong Zhang 
826a9d06231SHong Zhang   /* 3) send and recv matrix values coa */
827a9d06231SHong Zhang   /*------------------------------------*/
82842c57489SHong Zhang   buf_ri = merge->buf_ri;
82942c57489SHong Zhang   buf_rj = merge->buf_rj;
83042c57489SHong Zhang   len_s  = merge->len_s;
831fc42d0c8SSatish Balay   ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
83242c57489SHong Zhang   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
83342c57489SHong Zhang 
834a9d06231SHong Zhang   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
83542c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
83642c57489SHong Zhang     if (!len_s[proc]) continue;
837de0260b3SHong Zhang     i = merge->owners_co[proc];
838de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
83942c57489SHong Zhang     k++;
84042c57489SHong Zhang   }
8410c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
8420c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
84342c57489SHong Zhang 
844a9d06231SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
84542c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
84682412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
84742c57489SHong Zhang 
848a9d06231SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
849a9d06231SHong Zhang   /*------------------------------------------------------*/
850f9473cf7SHong Zhang   ierr = PetscGetTime(&t0);CHKERRQ(ierr);
851687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
85242c57489SHong Zhang 
85342c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
85442c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
85542c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
85642c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
85782412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
85842c57489SHong Zhang   }
85942c57489SHong Zhang 
860de0260b3SHong Zhang   for (i=0; i<cm; i++) {
86182412ba7SHong Zhang     row = owners[rank] + i; /* global row index of C_seq */
86242c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
863de0260b3SHong Zhang     ba_i = ba + bi[i];
86482412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
86542c57489SHong Zhang     /* add received vals into ba */
86642c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
86742c57489SHong Zhang       /* i-th row */
86842c57489SHong Zhang       if (i == *nextrow[k]) {
86982412ba7SHong Zhang         cnz = *(nextci[k]+1) - *nextci[k];
87082412ba7SHong Zhang         cj  = buf_rj[k] + *(nextci[k]);
87182412ba7SHong Zhang         ca  = abuf_r[k] + *(nextci[k]);
87282412ba7SHong Zhang         nextcj = 0;
87382412ba7SHong Zhang         for (j=0; nextcj<cnz; j++){
87482412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
87582412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
87642c57489SHong Zhang           }
87742c57489SHong Zhang         }
87882412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
87942c57489SHong Zhang       }
88042c57489SHong Zhang     }
88182412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
882dc0b31edSSatish Balay     ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
88342c57489SHong Zhang   }
88442c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
88542c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
886f9473cf7SHong Zhang   ierr = PetscGetTime(&tf);CHKERRQ(ierr);
887f9473cf7SHong Zhang   time_setvals += tf-t0;
88842c57489SHong Zhang 
889de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
890c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
89142c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
8920572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
893f9473cf7SHong Zhang 
894f9473cf7SHong Zhang   ierr = PetscGetTime(&tff);CHKERRQ(ierr);
895f9473cf7SHong Zhang   etime += tff - t00;
896f9473cf7SHong Zhang   /*
897a9d06231SHong Zhang   PetscInt prid=0;
898a9d06231SHong Zhang   if (rank == prid){
899f9473cf7SHong 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);
900a9d06231SHong Zhang   }
901f9473cf7SHong Zhang    */
90242c57489SHong Zhang   PetscFunctionReturn(0);
90342c57489SHong Zhang }
904