xref: /petsc/src/mat/impls/aij/mpi/mpiptap.c (revision cf3ca8ceebd227328df08578a8d35a5b419c4270)
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 
12fe615159SHong Zhang /*
13fe615159SHong Zhang #define DEBUG_MATPTAP
14fe615159SHong Zhang  */
15fe615159SHong Zhang 
1609573ac7SBarry Smith extern PetscErrorCode MatDestroy_MPIAIJ(Mat);
17cc6cb787SHong Zhang #undef __FUNCT__
18f8487c73SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP"
19f8487c73SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A)
20cc6cb787SHong Zhang {
21cc6cb787SHong Zhang   PetscErrorCode       ierr;
22f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data;
23f8487c73SHong Zhang   Mat_PtAPMPI          *ptap=a->ptap;
24cc6cb787SHong Zhang 
25cc6cb787SHong Zhang   PetscFunctionBegin;
26f8487c73SHong Zhang   if (ptap){
27c8b0795cSMark F. Adams     Mat_Merge_SeqsToMPI  *merge=ptap->merge;
28b7f45c76SHong Zhang     ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr);
29f8487c73SHong Zhang     ierr = PetscFree(ptap->bufa);CHKERRQ(ierr);
30a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr);
31a1a4e84aSHong Zhang     ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr);
32c5af039cSHong Zhang     ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */
33c5af039cSHong Zhang     if (ptap->api){ierr = PetscFree(ptap->api);CHKERRQ(ierr);}
34c5af039cSHong Zhang     if (ptap->apj){ierr = PetscFree(ptap->apj);CHKERRQ(ierr);}
35c8b0795cSMark F. Adams     if (merge) {
36cc6cb787SHong Zhang       ierr = PetscFree(merge->id_r);CHKERRQ(ierr);
37cc6cb787SHong Zhang       ierr = PetscFree(merge->len_s);CHKERRQ(ierr);
38cc6cb787SHong Zhang       ierr = PetscFree(merge->len_r);CHKERRQ(ierr);
39cc6cb787SHong Zhang       ierr = PetscFree(merge->bi);CHKERRQ(ierr);
40cc6cb787SHong Zhang       ierr = PetscFree(merge->bj);CHKERRQ(ierr);
41c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr);
42cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr);
43c05d87d6SBarry Smith       ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr);
44cc6cb787SHong Zhang       ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr);
4505b42c5fSBarry Smith       ierr = PetscFree(merge->coi);CHKERRQ(ierr);
4605b42c5fSBarry Smith       ierr = PetscFree(merge->coj);CHKERRQ(ierr);
4705b42c5fSBarry Smith       ierr = PetscFree(merge->owners_co);CHKERRQ(ierr);
486bf464f9SBarry Smith       ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr);
49dce485f0SBarry Smith       ierr = merge->destroy(A);CHKERRQ(ierr);
50f8487c73SHong Zhang       ierr = PetscFree(ptap->merge);CHKERRQ(ierr);
51bf0cc555SLisandro Dalcin     }
52f8487c73SHong Zhang     ierr = PetscFree(ptap);CHKERRQ(ierr);
53c8b0795cSMark F. Adams   }
54c8b0795cSMark F. Adams 
55cc6cb787SHong Zhang   PetscFunctionReturn(0);
56cc6cb787SHong Zhang }
57cc6cb787SHong Zhang 
58cc6cb787SHong Zhang #undef __FUNCT__
59cc6cb787SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP"
60f0c0a3a6SBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M)
61f0c0a3a6SBarry Smith {
62cc6cb787SHong Zhang   PetscErrorCode       ierr;
6311a6856fSHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data;
6411a6856fSHong Zhang   Mat_PtAPMPI          *ptap = a->ptap;
6511a6856fSHong Zhang   Mat_Merge_SeqsToMPI  *merge = ptap->merge;
66cc6cb787SHong Zhang 
67cc6cb787SHong Zhang   PetscFunctionBegin;
68dce485f0SBarry Smith   ierr = (*merge->duplicate)(A,op,M);CHKERRQ(ierr);
6911a6856fSHong Zhang   (*M)->ops->destroy   = merge->destroy;
7011a6856fSHong Zhang   (*M)->ops->duplicate = merge->duplicate;
71cc6cb787SHong Zhang   PetscFunctionReturn(0);
72cc6cb787SHong Zhang }
73cc6cb787SHong Zhang 
7442c57489SHong Zhang #undef __FUNCT__
75*cf3ca8ceSHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ"
76*cf3ca8ceSHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C)
7742c57489SHong Zhang {
7842c57489SHong Zhang   PetscErrorCode ierr;
7942c57489SHong Zhang 
8042c57489SHong Zhang   PetscFunctionBegin;
81*cf3ca8ceSHong Zhang   if (scall == MAT_INITIAL_MATRIX){
82*cf3ca8ceSHong Zhang     ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
83*cf3ca8ceSHong Zhang     ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);
84*cf3ca8ceSHong Zhang     ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr);
857a7894deSKris Buschelman   }
86*cf3ca8ceSHong Zhang   ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
87*cf3ca8ceSHong Zhang   ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr);
88*cf3ca8ceSHong Zhang   ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr);
8942c57489SHong Zhang   PetscFunctionReturn(0);
9042c57489SHong Zhang }
9142c57489SHong Zhang 
9242c57489SHong Zhang #undef __FUNCT__
9342c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ"
9442c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C)
9542c57489SHong Zhang {
9642c57489SHong Zhang   PetscErrorCode       ierr;
9777584fe6SHong Zhang   Mat                  Cmpi;
98f8487c73SHong Zhang   Mat_PtAPMPI          *ptap;
99a1a86e44SBarry Smith   PetscFreeSpaceList   free_space=PETSC_NULL,current_space=PETSC_NULL;
100f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c;
10142c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
10242c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
103ded4f561SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ;
10477584fe6SHong Zhang   PetscInt             *adi=ad->i,*aj,*aoi=ao->i,nnz;
10513f74950SBarry Smith   PetscInt             nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row;
106d16ca5f0SHong Zhang   PetscInt             am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n;
10742c57489SHong Zhang   PetscBT              lnkbt;
1087adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)A)->comm;
109529bc5dcSHong Zhang   PetscMPIInt          size,rank,tagi,tagj,*len_si,*len_s,*len_ri;
11042c57489SHong Zhang   PetscInt             **buf_rj,**buf_ri,**buf_ri_k;
11142c57489SHong Zhang   PetscInt             len,proc,*dnz,*onz,*owners;
112ded4f561SHong Zhang   PetscInt             nzi,*bi,*bj;
11342c57489SHong Zhang   PetscInt             nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci;
114ded4f561SHong Zhang   MPI_Request          *swaits,*rwaits;
115ded4f561SHong Zhang   MPI_Status           *sstatus,rstatus;
11642c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
11777584fe6SHong Zhang   PetscInt             *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0;
118d16ca5f0SHong Zhang   PetscReal            afill=1.0,afill_tmp;
11908cb4532SHong Zhang   PetscInt             rstart = P->cmap->rstart,rmax;
12008cb4532SHong Zhang   PetscScalar          *vals;
12142c57489SHong Zhang 
12242c57489SHong Zhang   PetscFunctionBegin;
123c5af039cSHong Zhang   /* check if matrix local sizes are compatible */
124c5af039cSHong Zhang   if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend){
125c5af039cSHong Zhang     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Arow (%D, %D) != Prow (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend);
126c5af039cSHong Zhang   }
127c5af039cSHong Zhang   if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend){
128c5af039cSHong Zhang     SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Acol (%D, %D) != Prow (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend);
129c5af039cSHong Zhang   }
130c5af039cSHong Zhang 
13183445d95SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
13283445d95SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
13383445d95SHong Zhang 
134f8487c73SHong Zhang   /* create struct Mat_PtAPMPI and attached it to C later */
135f8487c73SHong Zhang   ierr = PetscNew(Mat_PtAPMPI,&ptap);CHKERRQ(ierr);
136f8487c73SHong Zhang   ptap->reuse    = MAT_INITIAL_MATRIX;
1376c6e00e0SHong Zhang 
13877584fe6SHong Zhang 
1396c6e00e0SHong Zhang   /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */
140b7f45c76SHong Zhang   ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
141fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
142fe615159SHong Zhang    if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_oth is done...\n",rank);
143fe615159SHong Zhang #endif
144fe615159SHong Zhang 
1456c6e00e0SHong Zhang   /* get P_loc by taking all local rows of P */
146a1a4e84aSHong Zhang   ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
147fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
148fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
149fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] P_loc is done...\n",rank);
150fe615159SHong Zhang #endif
1516c6e00e0SHong Zhang 
152a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
153a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
1546c6e00e0SHong Zhang   pi_loc = p_loc->i; pj_loc = p_loc->j;
1556c6e00e0SHong Zhang   pi_oth = p_oth->i; pj_oth = p_oth->j;
15642c57489SHong Zhang 
157fd0ff01cSHong Zhang   /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */
158fd0ff01cSHong Zhang   /*-------------------------------------------------------------------*/
15977584fe6SHong Zhang   ierr = PetscMalloc((am+1)*sizeof(PetscInt),&api);CHKERRQ(ierr);
16082412ba7SHong Zhang   api[0] = 0;
16183445d95SHong Zhang 
16283445d95SHong Zhang   /* create and initialize a linked list */
16383445d95SHong Zhang   nlnk = pN+1;
164fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
165fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
166fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call PetscLLCreate(), pN %D, P_rmax %D %D; Annz %D\n",rank,pN,p_loc->rmax,p_oth->rmax,adi[am]+aoi[am]);
167fe615159SHong Zhang #endif
16883445d95SHong Zhang   ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
169fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
170fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
171fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] call PetscLLCreate() is done\n",rank);
172fe615159SHong Zhang #endif
17383445d95SHong Zhang 
174d16ca5f0SHong Zhang   /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */
175d16ca5f0SHong Zhang   ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr);
176f4a743e1SHong Zhang   current_space = free_space;
177f4a743e1SHong Zhang 
178f4a743e1SHong Zhang   for (i=0; i<am; i++) {
17982412ba7SHong Zhang     apnz = 0;
180f4a743e1SHong Zhang     /* diagonal portion of A */
181ded4f561SHong Zhang     nzi = adi[i+1] - adi[i];
18277584fe6SHong Zhang     aj  = ad->j + adi[i];
183ded4f561SHong Zhang     for (j=0; j<nzi; j++){
18477584fe6SHong Zhang       row  = aj[j];
18582412ba7SHong Zhang       pnz  = pi_loc[row+1] - pi_loc[row];
186ded4f561SHong Zhang       Jptr = pj_loc + pi_loc[row];
18783445d95SHong Zhang       /* add non-zero cols of P into the sorted linked list lnk */
188dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
18982412ba7SHong Zhang       apnz += nlnk;
190f4a743e1SHong Zhang     }
191f4a743e1SHong Zhang     /* off-diagonal portion of A */
192ded4f561SHong Zhang     nzi = aoi[i+1] - aoi[i];
19377584fe6SHong Zhang     aj  = ao->j + aoi[i];
194ded4f561SHong Zhang     for (j=0; j<nzi; j++){
19577584fe6SHong Zhang       row = aj[j];
19682412ba7SHong Zhang       pnz = pi_oth[row+1] - pi_oth[row];
197ded4f561SHong Zhang       Jptr  = pj_oth + pi_oth[row];
198dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
19982412ba7SHong Zhang       apnz += nlnk;
200f4a743e1SHong Zhang     }
201f4a743e1SHong Zhang 
20282412ba7SHong Zhang     api[i+1] = api[i] + apnz;
20377584fe6SHong Zhang     if (ap_rmax < apnz) ap_rmax = apnz;
204f4a743e1SHong Zhang 
20583445d95SHong Zhang     /* if free space is not available, double the total space in the list */
20682412ba7SHong Zhang     if (current_space->local_remaining<apnz) {
2072ba1acd5SHong Zhang       ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
208f2b054eeSHong Zhang       nspacedouble++;
209f4a743e1SHong Zhang     }
210f4a743e1SHong Zhang 
211f4a743e1SHong Zhang     /* Copy data into free space, then initialize lnk */
21282412ba7SHong Zhang     ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
21382412ba7SHong Zhang     current_space->array           += apnz;
21482412ba7SHong Zhang     current_space->local_used      += apnz;
21582412ba7SHong Zhang     current_space->local_remaining -= apnz;
216f4a743e1SHong Zhang   }
21782412ba7SHong Zhang   /* Allocate space for apj, initialize apj, and */
218f4a743e1SHong Zhang   /* destroy list of free space and other temporary array(s) */
21977584fe6SHong Zhang   ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&apj);CHKERRQ(ierr);
22077584fe6SHong Zhang   ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr);
221118727c9SMark F. Adams   afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1);
222d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
223f4a743e1SHong Zhang 
224ded4f561SHong Zhang   /* determine symbolic Co=(p->B)^T*AP - send to others */
225ded4f561SHong Zhang   /*----------------------------------------------------*/
226de0260b3SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
22742c57489SHong Zhang 
228ded4f561SHong Zhang   /* then, compute symbolic Co = (p->B)^T*AP */
229d0f46423SBarry Smith   pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors
23083445d95SHong Zhang                          >= (num of nonzero rows of C_seq) - pn */
231de0260b3SHong Zhang   ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr);
232de0260b3SHong Zhang   coi[0] = 0;
23342c57489SHong Zhang 
234d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */
235d16ca5f0SHong Zhang   nnz           = fill*(poti[pon] + api[am]);
236a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
23742c57489SHong Zhang   current_space = free_space;
23842c57489SHong Zhang 
239de0260b3SHong Zhang   for (i=0; i<pon; i++) {
240ded4f561SHong Zhang     nnz = 0;
24182412ba7SHong Zhang     pnz = poti[i+1] - poti[i];
24277584fe6SHong Zhang     ptJ = potj + poti[i];
24377584fe6SHong Zhang     for (j=0; j<pnz; j++){
24477584fe6SHong Zhang       row  = ptJ[j]; /* row of AP == col of Pot */
24582412ba7SHong Zhang       apnz = api[row+1] - api[row];
246ded4f561SHong Zhang       Jptr = apj + api[row];
24783445d95SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
248dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
249ded4f561SHong Zhang       nnz += nlnk;
25042c57489SHong Zhang     }
25142c57489SHong Zhang 
25283445d95SHong Zhang     /* If free space is not available, double the total space in the list */
253ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
2544238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
255d16ca5f0SHong Zhang       nspacedouble++;
25642c57489SHong Zhang     }
25742c57489SHong Zhang 
25842c57489SHong Zhang     /* Copy data into free space, and zero out denserows */
259ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
260ded4f561SHong Zhang     current_space->array           += nnz;
261ded4f561SHong Zhang     current_space->local_used      += nnz;
262ded4f561SHong Zhang     current_space->local_remaining -= nnz;
263ded4f561SHong Zhang     coi[i+1] = coi[i] + nnz;
26442c57489SHong Zhang   }
265de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr);
266a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr);
267118727c9SMark F. Adams   afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1);
268d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
269de0260b3SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr);
27042c57489SHong Zhang 
271fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
272fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
273fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] local Co is done\n",rank);
274fe615159SHong Zhang #endif
275fe615159SHong Zhang 
276ded4f561SHong Zhang   /* send j-array (coj) of Co to other processors */
277ded4f561SHong Zhang   /*----------------------------------------------*/
27842c57489SHong Zhang   /* determine row ownership */
27983445d95SHong Zhang   ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr);
28026283091SBarry Smith   ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr);
2817a2fc3feSBarry Smith   merge->rowmap->n = pn;
2827a2fc3feSBarry Smith   merge->rowmap->bs = 1;
28326283091SBarry Smith   ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr);
2847a2fc3feSBarry Smith   owners = merge->rowmap->range;
28542c57489SHong Zhang 
28642c57489SHong Zhang   /* determine the number of messages to send, their lengths */
28742c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr);
28883445d95SHong Zhang   ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
28942c57489SHong Zhang   ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr);
29042c57489SHong Zhang   len_s = merge->len_s;
29142c57489SHong Zhang   merge->nsend = 0;
29283445d95SHong Zhang 
293de0260b3SHong Zhang   ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr);
294de0260b3SHong Zhang   ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr);
295de0260b3SHong Zhang 
29683445d95SHong Zhang   proc = 0;
297de0260b3SHong Zhang   for (i=0; i<pon; i++){
298de0260b3SHong Zhang     while (prmap[i] >= owners[proc+1]) proc++;
299de0260b3SHong Zhang     len_si[proc]++;  /* num of rows in Co to be sent to [proc] */
300de0260b3SHong Zhang     len_s[proc] += coi[i+1] - coi[i];
30183445d95SHong Zhang   }
302de0260b3SHong Zhang 
303de0260b3SHong Zhang   len   = 0;  /* max length of buf_si[] */
304de0260b3SHong Zhang   owners_co[0] = 0;
30542c57489SHong Zhang   for (proc=0; proc<size; proc++){
306de0260b3SHong Zhang     owners_co[proc+1] = owners_co[proc] + len_si[proc];
30783445d95SHong Zhang     if (len_si[proc]){
30842c57489SHong Zhang       merge->nsend++;
30983445d95SHong Zhang       len_si[proc] = 2*(len_si[proc] + 1);
31042c57489SHong Zhang       len += len_si[proc];
31142c57489SHong Zhang     }
31242c57489SHong Zhang   }
31342c57489SHong Zhang 
314ded4f561SHong Zhang   /* determine the number and length of messages to receive for coi and coj  */
31542c57489SHong Zhang   ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr);
31642c57489SHong Zhang   ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr);
31742c57489SHong Zhang 
318ded4f561SHong Zhang   /* post the Irecv and Isend of coj */
319529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr);
320529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr);
321ded4f561SHong Zhang   ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr);
32242c57489SHong Zhang   for (proc=0, k=0; proc<size; proc++){
32342c57489SHong Zhang     if (!len_s[proc]) continue;
324de0260b3SHong Zhang     i = owners_co[proc];
325529bc5dcSHong Zhang     ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr);
32642c57489SHong Zhang     k++;
32742c57489SHong Zhang   }
32842c57489SHong Zhang 
329ded4f561SHong Zhang   /* receives and sends of coj are complete */
330ded4f561SHong Zhang   ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr);
33177584fe6SHong Zhang   for (i=0; i<merge->nrecv; i++){
332c280ed6aSJed Brown     PetscMPIInt icompleted;
333c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
334ded4f561SHong Zhang   }
335ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3360c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
33782412ba7SHong Zhang 
338ded4f561SHong Zhang   /* send and recv coi */
339ded4f561SHong Zhang   /*-------------------*/
340529bc5dcSHong Zhang   ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr);
341529bc5dcSHong Zhang   ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr);
34242c57489SHong Zhang   ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr);
34342c57489SHong Zhang   buf_si = buf_s;  /* points to the beginning of k-th msg to be sent */
34442c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
34542c57489SHong Zhang     if (!len_s[proc]) continue;
34642c57489SHong Zhang     /* form outgoing message for i-structure:
34742c57489SHong Zhang          buf_si[0]:                 nrows to be sent
34842c57489SHong Zhang                [1:nrows]:           row index (global)
34942c57489SHong Zhang                [nrows+1:2*nrows+1]: i-structure index
35042c57489SHong Zhang     */
35142c57489SHong Zhang     /*-------------------------------------------*/
35242c57489SHong Zhang     nrows = len_si[proc]/2 - 1;
35342c57489SHong Zhang     buf_si_i    = buf_si + nrows+1;
35442c57489SHong Zhang     buf_si[0]   = nrows;
35542c57489SHong Zhang     buf_si_i[0] = 0;
35642c57489SHong Zhang     nrows = 0;
357de0260b3SHong Zhang     for (i=owners_co[proc]; i<owners_co[proc+1]; i++){
358ded4f561SHong Zhang       nzi = coi[i+1] - coi[i];
359ded4f561SHong Zhang       buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */
360de0260b3SHong Zhang       buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */
36142c57489SHong Zhang       nrows++;
36242c57489SHong Zhang     }
363529bc5dcSHong Zhang     ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr);
36442c57489SHong Zhang     k++;
36542c57489SHong Zhang     buf_si += len_si[proc];
36642c57489SHong Zhang   }
367ded4f561SHong Zhang   i = merge->nrecv;
368ded4f561SHong Zhang   while (i--) {
369c280ed6aSJed Brown     PetscMPIInt icompleted;
370c280ed6aSJed Brown     ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr);
371ded4f561SHong Zhang   }
372ded4f561SHong Zhang   ierr = PetscFree(rwaits);CHKERRQ(ierr);
3730c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);}
374f2b054eeSHong Zhang   /*
375ae15b995SBarry Smith   ierr = PetscInfo2(A,"nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr);
37642c57489SHong Zhang   for (i=0; i<merge->nrecv; i++){
377ae15b995SBarry 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);
37842c57489SHong Zhang   }
379f2b054eeSHong Zhang   */
38042c57489SHong Zhang   ierr = PetscFree(len_si);CHKERRQ(ierr);
38142c57489SHong Zhang   ierr = PetscFree(len_ri);CHKERRQ(ierr);
382ded4f561SHong Zhang   ierr = PetscFree(swaits);CHKERRQ(ierr);
383ded4f561SHong Zhang   ierr = PetscFree(sstatus);CHKERRQ(ierr);
38442c57489SHong Zhang   ierr = PetscFree(buf_s);CHKERRQ(ierr);
38542c57489SHong Zhang 
386fe615159SHong Zhang #if defined(DEBUG_MATPTAP)
387fe615159SHong Zhang   ierr = MPI_Barrier(comm);CHKERRQ(ierr);
388fe615159SHong Zhang   if (!rank) ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] Co is done\n",rank);
389fe615159SHong Zhang #endif
390fe615159SHong Zhang 
391de0260b3SHong Zhang   /* compute the local portion of C (mpi mat) */
392de0260b3SHong Zhang   /*------------------------------------------*/
393ded4f561SHong Zhang   ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
394ded4f561SHong Zhang 
39542c57489SHong Zhang   /* allocate bi array and free space for accumulating nonzero column info */
39642c57489SHong Zhang   ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr);
39742c57489SHong Zhang   bi[0] = 0;
39842c57489SHong Zhang 
399d16ca5f0SHong Zhang   /* set initial free space to be fill*(nnz(P) + nnz(AP)) */
400d16ca5f0SHong Zhang   nnz           = fill*(pi_loc[pm] + api[am]);
401a1a86e44SBarry Smith   ierr          = PetscFreeSpaceGet(nnz,&free_space);
40242c57489SHong Zhang   current_space = free_space;
40342c57489SHong Zhang 
404687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
40542c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
40642c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
40742c57489SHong Zhang     nrows       = *buf_ri_k[k];
40842c57489SHong Zhang     nextrow[k]  = buf_ri_k[k] + 1;  /* next row number of k-th recved i-structure */
40942c57489SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
41042c57489SHong Zhang   }
41142c57489SHong Zhang   ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr);
41208cb4532SHong Zhang   rmax = 0;
41342c57489SHong Zhang   for (i=0; i<pn; i++) {
414ded4f561SHong Zhang     /* add pdt[i,:]*AP into lnk */
415ded4f561SHong Zhang     nnz = 0;
416ded4f561SHong Zhang     pnz = pdti[i+1] - pdti[i];
41777584fe6SHong Zhang     ptJ = pdtj + pdti[i];
41877584fe6SHong Zhang     for (j=0; j<pnz; j++){
41977584fe6SHong Zhang       row  = ptJ[j];  /* row of AP == col of Pt */
420ded4f561SHong Zhang       apnz = api[row+1] - api[row];
421ded4f561SHong Zhang       Jptr = apj + api[row];
422ded4f561SHong Zhang       /* add non-zero cols of AP into the sorted linked list lnk */
423dadf0e6bSHong Zhang       ierr = PetscLLAddSorted(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
424ded4f561SHong Zhang       nnz += nlnk;
425ded4f561SHong Zhang     }
42642c57489SHong Zhang     /* add received col data into lnk */
42742c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
42842c57489SHong Zhang       if (i == *nextrow[k]) { /* i-th row */
429ded4f561SHong Zhang         nzi = *(nextci[k]+1) - *nextci[k];
430ded4f561SHong Zhang         Jptr  = buf_rj[k] + *nextci[k];
431dadf0e6bSHong Zhang         ierr = PetscLLAddSorted(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr);
432ded4f561SHong Zhang         nnz += nlnk;
43342c57489SHong Zhang         nextrow[k]++; nextci[k]++;
43442c57489SHong Zhang       }
43542c57489SHong Zhang     }
43642c57489SHong Zhang 
43742c57489SHong Zhang     /* if free space is not available, make more free space */
438ded4f561SHong Zhang     if (current_space->local_remaining<nnz) {
4394238b7adSHong Zhang       ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,&current_space);CHKERRQ(ierr);
440d16ca5f0SHong Zhang       nspacedouble++;
44142c57489SHong Zhang     }
44242c57489SHong Zhang     /* copy data into free space, then initialize lnk */
443ded4f561SHong Zhang     ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr);
444ded4f561SHong Zhang     ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr);
445ded4f561SHong Zhang     current_space->array           += nnz;
446ded4f561SHong Zhang     current_space->local_used      += nnz;
447ded4f561SHong Zhang     current_space->local_remaining -= nnz;
448ded4f561SHong Zhang     bi[i+1] = bi[i] + nnz;
44908cb4532SHong Zhang     if (nnz > rmax) rmax = nnz;
45042c57489SHong Zhang   }
451ded4f561SHong Zhang   ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr);
4520572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
45342c57489SHong Zhang 
45442c57489SHong Zhang   ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr);
455a1a86e44SBarry Smith   ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr);
456118727c9SMark F. Adams   afill_tmp = (PetscReal)bi[pn]/(pi_loc[pm] + api[am]+1);
457d16ca5f0SHong Zhang   if (afill_tmp > afill) afill = afill_tmp;
45842c57489SHong Zhang   ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr);
45942c57489SHong Zhang 
46008cb4532SHong Zhang   /* create symbolic parallel matrix Cmpi - why cannot be assembled in Numeric part - check runex55_SA ?  */
46108cb4532SHong Zhang   /*------------------------------------------------------------------------------------------------------*/
46208cb4532SHong Zhang   ierr = PetscMalloc((rmax+1)*sizeof(PetscScalar),&vals);CHKERRQ(ierr);
46308cb4532SHong Zhang   ierr = PetscMemzero(vals,rmax*sizeof(PetscScalar));CHKERRQ(ierr);
46408cb4532SHong Zhang 
46577584fe6SHong Zhang   ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr);
46677584fe6SHong Zhang   ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
467a2f3521dSMark F. Adams   ierr = MatSetBlockSizes(Cmpi,P->cmap->bs,P->cmap->bs);CHKERRQ(ierr);
46877584fe6SHong Zhang   ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr);
46977584fe6SHong Zhang   ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr);
47042c57489SHong Zhang   ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr);
471a2f3521dSMark F. Adams 
47208cb4532SHong Zhang   for (i=0; i<pn; i++){
47308cb4532SHong Zhang     row = i + rstart;
47408cb4532SHong Zhang     nnz = bi[i+1] - bi[i];
47508cb4532SHong Zhang     Jptr = bj + bi[i];
47608cb4532SHong Zhang     ierr = MatSetValues(Cmpi,1,&row,nnz,Jptr,vals,INSERT_VALUES);CHKERRQ(ierr);
47708cb4532SHong Zhang   }
47808cb4532SHong Zhang   ierr = MatAssemblyBegin(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
47908cb4532SHong Zhang   ierr = MatAssemblyEnd(Cmpi,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
48008cb4532SHong Zhang   ierr = PetscFree(vals);CHKERRQ(ierr);
48142c57489SHong Zhang 
48242c57489SHong Zhang   merge->bi            = bi;
48342c57489SHong Zhang   merge->bj            = bj;
484de0260b3SHong Zhang   merge->coi           = coi;
485de0260b3SHong Zhang   merge->coj           = coj;
48642c57489SHong Zhang   merge->buf_ri        = buf_ri;
48742c57489SHong Zhang   merge->buf_rj        = buf_rj;
488de0260b3SHong Zhang   merge->owners_co     = owners_co;
48977584fe6SHong Zhang   merge->destroy       = Cmpi->ops->destroy;
49077584fe6SHong Zhang   merge->duplicate     = Cmpi->ops->duplicate;
491cc6cb787SHong Zhang 
49277584fe6SHong Zhang   /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */
49308cb4532SHong Zhang   /* Cmpi->assembled      = PETSC_FALSE; */
49477584fe6SHong Zhang   Cmpi->ops->destroy   = MatDestroy_MPIAIJ_PtAP;
49577584fe6SHong Zhang   Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP;
49642c57489SHong Zhang 
49777584fe6SHong Zhang   /* attach the supporting struct to Cmpi for reuse */
49877584fe6SHong Zhang   c = (Mat_MPIAIJ*)Cmpi->data;
499f8487c73SHong Zhang   c->ptap        = ptap;
50077584fe6SHong Zhang   ptap->api      = api;
50177584fe6SHong Zhang   ptap->apj      = apj;
502f8487c73SHong Zhang   ptap->merge    = merge;
503d6ab1dc0SHong Zhang   ptap->rmax     = ap_rmax;
504f8487c73SHong Zhang 
50577584fe6SHong Zhang   *C = Cmpi;
506f2b054eeSHong Zhang #if defined(PETSC_USE_INFO)
507f2b054eeSHong Zhang   if (bi[pn] != 0) {
50877584fe6SHong Zhang     ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %G needed %G.\n",nspacedouble,fill,afill);CHKERRQ(ierr);
50977584fe6SHong Zhang     ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr);
510f2b054eeSHong Zhang   } else {
51177584fe6SHong Zhang     ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr);
512f2b054eeSHong Zhang   }
513f2b054eeSHong Zhang #endif
51442c57489SHong Zhang   PetscFunctionReturn(0);
51542c57489SHong Zhang }
51642c57489SHong Zhang 
51742c57489SHong Zhang #undef __FUNCT__
51842c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ"
51942c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C)
52042c57489SHong Zhang {
52142c57489SHong Zhang   PetscErrorCode       ierr;
52242c57489SHong Zhang   Mat_Merge_SeqsToMPI  *merge;
523f8487c73SHong Zhang   Mat_MPIAIJ           *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data;
52442c57489SHong Zhang   Mat_SeqAIJ           *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data;
525de0260b3SHong Zhang   Mat_SeqAIJ           *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data;
52642c57489SHong Zhang   Mat_SeqAIJ           *p_loc,*p_oth;
527f8487c73SHong Zhang   Mat_PtAPMPI          *ptap;
5281d633602SHong Zhang   PetscInt             *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp;
52982412ba7SHong Zhang   PetscInt             *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj;
53082412ba7SHong Zhang   PetscInt             i,j,k,anz,pnz,apnz,nextap,row,*cj;
531e900a757SHong Zhang   MatScalar            *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp;
532d0f46423SBarry Smith   PetscInt             am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n;
5337adad957SLisandro Dalcin   MPI_Comm             comm=((PetscObject)C)->comm;
53442c57489SHong Zhang   PetscMPIInt          size,rank,taga,*len_s;
53582412ba7SHong Zhang   PetscInt             *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci;
53642c57489SHong Zhang   PetscInt             **buf_ri,**buf_rj;
53750f5bf86SHong Zhang   PetscInt             cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */
53842c57489SHong Zhang   MPI_Request          *s_waits,*r_waits;
53942c57489SHong Zhang   MPI_Status           *status;
54082412ba7SHong Zhang   MatScalar            **abuf_r,*ba_i,*pA,*coa,*ba;
54182412ba7SHong Zhang   PetscInt             *api,*apj,*coi,*coj;
542d0f46423SBarry Smith   PetscInt             *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend;
543a9d06231SHong Zhang   PetscInt             sparse_axpy;
54442c57489SHong Zhang 
54542c57489SHong Zhang   PetscFunctionBegin;
54642c57489SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
54742c57489SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
54842c57489SHong Zhang 
549f8487c73SHong Zhang   ptap  = c->ptap;
5501c4805f8SJed Brown   if (!ptap) SETERRQ(((PetscObject)C)->comm,PETSC_ERR_ARG_INCOMP,"MatPtAP() has not been called to create matrix C yet, cannot use MAT_REUSE_MATRIX");
551f8487c73SHong Zhang   merge = ptap->merge;
5526c6e00e0SHong Zhang 
553a1a4e84aSHong Zhang   /* 1) get P_oth = ptap->P_oth  and P_loc = ptap->P_loc */
554f9473cf7SHong Zhang   /*--------------------------------------------------*/
555f8487c73SHong Zhang   if (ptap->reuse == MAT_INITIAL_MATRIX){
556f8487c73SHong Zhang     ptap->reuse = MAT_REUSE_MATRIX;
5576c6e00e0SHong Zhang   } else { /* update numerical values of P_oth and P_loc */
558b7f45c76SHong Zhang     ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr);
559a1a4e84aSHong Zhang     ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr);
5606c6e00e0SHong Zhang   }
561f8487c73SHong Zhang 
562f9473cf7SHong Zhang   /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */
563f9473cf7SHong Zhang   /*--------------------------------------------------------------*/
56442c57489SHong Zhang   /* get data from symbolic products */
565a1a4e84aSHong Zhang   p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data;
566a1a4e84aSHong Zhang   p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data;
567a9d06231SHong Zhang   pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a;
56842c57489SHong Zhang   pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a;
56942c57489SHong Zhang 
570de0260b3SHong Zhang   coi = merge->coi; coj = merge->coj;
571de0260b3SHong Zhang   ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr);
572de0260b3SHong Zhang   ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr);
573de0260b3SHong Zhang 
574de0260b3SHong Zhang   bi     = merge->bi; bj = merge->bj;
5757a2fc3feSBarry Smith   owners = merge->rowmap->range;
576de0260b3SHong Zhang   ierr   = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr);
577de0260b3SHong Zhang   ierr   = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr);
57842c57489SHong Zhang 
579a1a4e84aSHong Zhang   api = ptap->api; apj = ptap->apj;
5800496f32cSHong Zhang   /* flag 'sparse_axpy' determines which implementations to be used:
5810496f32cSHong Zhang        0: do dense axpy in MatPtAPNumeric() - fastest, but requires storage of a dense array apa; (default)
5820496f32cSHong Zhang        1: do one sparse axpy - uses same memory as sparse_axpy=0 and might execute less flops
5830496f32cSHong Zhang           (apnz vs. cnz in the outerproduct), slower than case '0' when cnz is not too large than apnz;
5840496f32cSHong Zhang        2: do two sparse axpy in MatPtAPNumeric() - slowest, uses a sparse array apa */
585a9d06231SHong Zhang   /* set default sparse_axpy */
586e900a757SHong Zhang   sparse_axpy = 0;
5870496f32cSHong Zhang   ierr = PetscOptionsGetInt(PETSC_NULL,"-matptap_sparseaxpy",&sparse_axpy,PETSC_NULL);CHKERRQ(ierr);
588a9d06231SHong Zhang   if (sparse_axpy == 0){  /* Do not perform sparse axpy */
589a9d06231SHong Zhang     /*--------------------------------------------------*/
590a9d06231SHong Zhang     /* malloc apa to store dense row A[i,:]*P */
591a9d06231SHong Zhang     ierr = PetscMalloc((P->cmap->N)*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
592a9d06231SHong Zhang     ierr = PetscMemzero(apa,P->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
593a9d06231SHong Zhang 
594a9d06231SHong Zhang     for (i=0; i<am; i++) {
595a9d06231SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
596a9d06231SHong Zhang       /*------------------------------------------------------------*/
597a9d06231SHong Zhang       apJ = apj + api[i];
598a9d06231SHong Zhang 
599a9d06231SHong Zhang       /* diagonal portion of A */
600a9d06231SHong Zhang       anz = adi[i+1] - adi[i];
601a9d06231SHong Zhang       adj = ad->j + adi[i];
602a9d06231SHong Zhang       ada = ad->a + adi[i];
603a9d06231SHong Zhang       for (j=0; j<anz; j++) {
604a9d06231SHong Zhang         row = adj[j];
605a9d06231SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
606a9d06231SHong Zhang         pj  = pj_loc + pi_loc[row];
607a9d06231SHong Zhang         pa  = pa_loc + pi_loc[row];
608a9d06231SHong Zhang 
609a9d06231SHong Zhang         /* perform dense axpy */
610e900a757SHong Zhang         valtmp = ada[j];
611a9d06231SHong Zhang         for (k=0; k<pnz; k++){
612e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
613a9d06231SHong Zhang         }
614a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
615a9d06231SHong Zhang       }
616a9d06231SHong Zhang 
617a9d06231SHong Zhang       /* off-diagonal portion of A */
618a9d06231SHong Zhang       anz = aoi[i+1] - aoi[i];
619a9d06231SHong Zhang       aoj = ao->j + aoi[i];
620a9d06231SHong Zhang       aoa = ao->a + aoi[i];
621a9d06231SHong Zhang       for (j=0; j<anz; j++) {
622a9d06231SHong Zhang         row = aoj[j];
623a9d06231SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
624a9d06231SHong Zhang         pj  = pj_oth + pi_oth[row];
625a9d06231SHong Zhang         pa  = pa_oth + pi_oth[row];
626a9d06231SHong Zhang 
627a9d06231SHong Zhang         /* perform dense axpy */
628e900a757SHong Zhang         valtmp = aoa[j];
629a9d06231SHong Zhang         for (k=0; k<pnz; k++){
630e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
631a9d06231SHong Zhang         }
632a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
633a9d06231SHong Zhang       }
634a9d06231SHong Zhang 
635a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
636a9d06231SHong Zhang       /*--------------------------------------------------------------*/
637a9d06231SHong Zhang       apnz = api[i+1] - api[i];
638a9d06231SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
639a9d06231SHong Zhang       pnz = po->i[i+1] - po->i[i];
640e900a757SHong Zhang       poJ = po->j + po->i[i];
641a9d06231SHong Zhang       pA  = po->a + po->i[i];
642a9d06231SHong Zhang       for (j=0; j<pnz; j++){
643e900a757SHong Zhang         row = poJ[j];
644e900a757SHong Zhang         cnz = coi[row+1] - coi[row];
645e900a757SHong Zhang         cj  = coj + coi[row];
646e900a757SHong Zhang         ca  = coa + coi[row];
647a9d06231SHong Zhang         /* perform dense axpy */
648e900a757SHong Zhang         valtmp = pA[j];
649a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
650e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
651a9d06231SHong Zhang         }
652a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
653a9d06231SHong Zhang       }
654a9d06231SHong Zhang 
655a9d06231SHong Zhang       /* put the value into Cd (diagonal part) */
656a9d06231SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
657e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
658a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
659a9d06231SHong Zhang       for (j=0; j<pnz; j++){
660e900a757SHong Zhang         row = pdJ[j];
661e900a757SHong Zhang         cnz = bi[row+1] - bi[row];
662e900a757SHong Zhang         cj  = bj + bi[row];
663e900a757SHong Zhang         ca  = ba + bi[row];
664a9d06231SHong Zhang         /* perform dense axpy */
665e900a757SHong Zhang         valtmp = pA[j];
666a9d06231SHong Zhang         for (k=0; k<cnz; k++) {
667e900a757SHong Zhang           ca[k] += valtmp*apa[cj[k]];
668a9d06231SHong Zhang         }
669a9d06231SHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
670a9d06231SHong Zhang       }
671a9d06231SHong Zhang 
672a9d06231SHong Zhang       /* zero the current row of A*P */
673a9d06231SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
674a9d06231SHong Zhang     }
675a9d06231SHong Zhang   } else if (sparse_axpy == 1){  /* Perform one sparse axpy */
676a9d06231SHong Zhang     /*------------------------------------------------------*/
6771d633602SHong Zhang     /* malloc apa to store dense row A[i,:]*P */
6781d633602SHong Zhang     ierr = PetscMalloc((P->cmap->N)*sizeof(PetscScalar),&apa);CHKERRQ(ierr);
6791d633602SHong Zhang     ierr = PetscMemzero(apa,P->cmap->N*sizeof(PetscScalar));CHKERRQ(ierr);
6801d633602SHong Zhang 
6811d633602SHong Zhang     for (i=0; i<am; i++) {
682f9473cf7SHong Zhang       /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */
683f9473cf7SHong Zhang       /*------------------------------------------------------------*/
6841d633602SHong Zhang       apJ = apj + api[i];
6851d633602SHong Zhang 
6861d633602SHong Zhang       /* diagonal portion of A */
6871d633602SHong Zhang       anz = adi[i+1] - adi[i];
6881d633602SHong Zhang       adj = ad->j + adi[i];
6891d633602SHong Zhang       ada = ad->a + adi[i];
6901d633602SHong Zhang       for (j=0; j<anz; j++) {
6911d633602SHong Zhang         row = adj[j];
6921d633602SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
6931d633602SHong Zhang         pj  = pj_loc + pi_loc[row];
6941d633602SHong Zhang         pa  = pa_loc + pi_loc[row];
6951d633602SHong Zhang 
6961d633602SHong Zhang         /* perform dense axpy */
697e900a757SHong Zhang         valtmp = ada[j];
6981d633602SHong Zhang         for (k=0; k<pnz; k++){
699e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
7001d633602SHong Zhang         }
7011d633602SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
7021d633602SHong Zhang       }
7031d633602SHong Zhang 
7041d633602SHong Zhang       /* off-diagonal portion of A */
7051d633602SHong Zhang       anz = aoi[i+1] - aoi[i];
7061d633602SHong Zhang       aoj = ao->j + aoi[i];
7071d633602SHong Zhang       aoa = ao->a + aoi[i];
7081d633602SHong Zhang       for (j=0; j<anz; j++) {
7091d633602SHong Zhang         row = aoj[j];
7101d633602SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
7111d633602SHong Zhang         pj  = pj_oth + pi_oth[row];
7121d633602SHong Zhang         pa  = pa_oth + pi_oth[row];
7131d633602SHong Zhang 
7141d633602SHong Zhang         /* perform dense axpy */
715e900a757SHong Zhang         valtmp = aoa[j];
7161d633602SHong Zhang         for (k=0; k<pnz; k++){
717e900a757SHong Zhang           apa[pj[k]] += valtmp*pa[k];
7181d633602SHong Zhang         }
7191d633602SHong Zhang         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
7201d633602SHong Zhang       }
7211d633602SHong Zhang 
722f9473cf7SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
723f9473cf7SHong Zhang       /*--------------------------------------------------------------*/
7241d633602SHong Zhang       apnz = api[i+1] - api[i];
725f9473cf7SHong Zhang       /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */
726f9473cf7SHong Zhang       pnz = po->i[i+1] - po->i[i];
727e900a757SHong Zhang       poJ = po->j + po->i[i];
728a9d06231SHong Zhang       pA  = po->a + po->i[i];
7291d633602SHong Zhang       for (j=0; j<pnz; j++){
730e900a757SHong Zhang         row    = poJ[j];
731e900a757SHong Zhang         cj     = coj + coi[row];
732e900a757SHong Zhang         ca     = coa + coi[row];
733e900a757SHong Zhang         valtmp = pA[j];
7341d633602SHong Zhang         /* perform sparse axpy */
7351d633602SHong Zhang         nextap = 0;
7361d633602SHong Zhang         for (k=0; nextap<apnz; k++) {
7371d633602SHong Zhang           if (cj[k]==apJ[nextap]) { /* global column index */
738e900a757SHong Zhang             ca[k] += valtmp*apa[cj[k]]; nextap++;
7391d633602SHong Zhang           }
7401d633602SHong Zhang         }
7411d633602SHong Zhang         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
7421d633602SHong Zhang       }
743f9473cf7SHong Zhang 
744f9473cf7SHong Zhang       /* put the value into Cd (diagonal part) */
745f9473cf7SHong Zhang       pnz = pd->i[i+1] - pd->i[i];
746e900a757SHong Zhang       pdJ = pd->j + pd->i[i];
747a9d06231SHong Zhang       pA  = pd->a + pd->i[i];
748f9473cf7SHong Zhang       for (j=0; j<pnz; j++){
749e900a757SHong Zhang         row    = pdJ[j];
750e900a757SHong Zhang         cj     = bj + bi[row];
751e900a757SHong Zhang         ca     = ba + bi[row];
752e900a757SHong Zhang         valtmp = pA[j];
753f9473cf7SHong Zhang         /* perform sparse axpy */
754f9473cf7SHong Zhang         nextap = 0;
755f9473cf7SHong Zhang         for (k=0; nextap<apnz; k++) {
756f9473cf7SHong Zhang           if (cj[k]==apJ[nextap]) { /* global column index */
757e900a757SHong Zhang             ca[k] += valtmp*apa[cj[k]];
758a9d06231SHong Zhang             nextap++;
759f9473cf7SHong Zhang           }
760f9473cf7SHong Zhang         }
761f9473cf7SHong Zhang         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
762f9473cf7SHong Zhang       }
763f9473cf7SHong Zhang 
764f9473cf7SHong Zhang       /* zero the current row of A*P */
7651d633602SHong Zhang       for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0;
7661d633602SHong Zhang     }
7670496f32cSHong Zhang   } else if (sparse_axpy == 2){/* Perform two sparse axpy */
768a9d06231SHong Zhang     /*----------------------------------------------------*/
7691d633602SHong Zhang     /* malloc apa to store sparse row A[i,:]*P */
770d6ab1dc0SHong Zhang     ierr = PetscMalloc((ptap->rmax+1)*sizeof(MatScalar),&apa);CHKERRQ(ierr);
771d6ab1dc0SHong Zhang     ierr = PetscMemzero(apa,ptap->rmax*sizeof(MatScalar));CHKERRQ(ierr);
77242c57489SHong Zhang 
773a9d06231SHong Zhang     pA=pa_loc;
77442c57489SHong Zhang     for (i=0; i<am; i++) {
775f4a743e1SHong Zhang       /* form i-th sparse row of A*P */
77682412ba7SHong Zhang       apnz = api[i+1] - api[i];
77782412ba7SHong Zhang       apJ  = apj + api[i];
77842c57489SHong Zhang       /* diagonal portion of A */
77982412ba7SHong Zhang       anz = adi[i+1] - adi[i];
7801d633602SHong Zhang       adj = ad->j + adi[i];
7811d633602SHong Zhang       ada = ad->a + adi[i];
78282412ba7SHong Zhang       for (j=0; j<anz; j++) {
7831d633602SHong Zhang         row = adj[j];
78482412ba7SHong Zhang         pnz = pi_loc[row+1] - pi_loc[row];
78582412ba7SHong Zhang         pj  = pj_loc + pi_loc[row];
78682412ba7SHong Zhang         pa  = pa_loc + pi_loc[row];
787e900a757SHong Zhang         valtmp = ada[j];
788f4a743e1SHong Zhang         nextp  = 0;
78982412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
79082412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
791e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
79242c57489SHong Zhang           }
79342c57489SHong Zhang         }
794dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
79542c57489SHong Zhang       }
79642c57489SHong Zhang       /* off-diagonal portion of A */
79782412ba7SHong Zhang       anz = aoi[i+1] - aoi[i];
7981d633602SHong Zhang       aoj = ao->j + aoi[i];
7991d633602SHong Zhang       aoa = ao->a + aoi[i];
80082412ba7SHong Zhang       for (j=0; j<anz; j++) {
8011d633602SHong Zhang         row = aoj[j];
80282412ba7SHong Zhang         pnz = pi_oth[row+1] - pi_oth[row];
80382412ba7SHong Zhang         pj  = pj_oth + pi_oth[row];
80482412ba7SHong Zhang         pa  = pa_oth + pi_oth[row];
805e900a757SHong Zhang         valtmp = aoa[j];
806f4a743e1SHong Zhang         nextp  = 0;
80782412ba7SHong Zhang         for (k=0; nextp<pnz; k++) {
80882412ba7SHong Zhang           if (apJ[k] == pj[nextp]) { /* col of AP == col of P */
809e900a757SHong Zhang             apa[k] += valtmp*pa[nextp++];
81042c57489SHong Zhang           }
81142c57489SHong Zhang         }
812dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr);
81342c57489SHong Zhang       }
81442c57489SHong Zhang 
815a9d06231SHong Zhang       /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */
816a9d06231SHong Zhang       /*--------------------------------------------------------------*/
81782412ba7SHong Zhang       pnz = pi_loc[i+1] - pi_loc[i];
818f9473cf7SHong Zhang       pJ  = pj_loc + pi_loc[i];
81982412ba7SHong Zhang       for (j=0; j<pnz; j++) {
82042c57489SHong Zhang         nextap = 0;
821f9473cf7SHong Zhang         row    = pJ[j]; /* global index */
82282412ba7SHong Zhang         if (row < pcstart || row >=pcend) { /* put the value into Co */
823e900a757SHong Zhang           row = *poJ;
824e900a757SHong Zhang           cj  = coj + coi[row];
825e900a757SHong Zhang           ca  = coa + coi[row]; poJ++;
82682412ba7SHong Zhang         } else {                            /* put the value into Cd */
827e900a757SHong Zhang           row  = *pdJ;
828e900a757SHong Zhang           cj   = bj + bi[row];
829e900a757SHong Zhang           ca   = ba + bi[row]; pdJ++;
83042c57489SHong Zhang         }
831e900a757SHong Zhang         valtmp = pA[j];
83282412ba7SHong Zhang         for (k=0; nextap<apnz; k++) {
833e900a757SHong Zhang           if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++];
83442c57489SHong Zhang         }
835dc0b31edSSatish Balay         ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr);
836de0260b3SHong Zhang       }
8370496f32cSHong Zhang       pA += pnz;
83842c57489SHong Zhang       /* zero the current row info for A*P */
83982412ba7SHong Zhang       ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr);
84042c57489SHong Zhang     }
841a9d06231SHong Zhang   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"sparse_axpy only takes values 0, 1 and 2");
84242c57489SHong Zhang   ierr = PetscFree(apa);CHKERRQ(ierr);
84342c57489SHong Zhang 
844a9d06231SHong Zhang   /* 3) send and recv matrix values coa */
845a9d06231SHong Zhang   /*------------------------------------*/
84642c57489SHong Zhang   buf_ri = merge->buf_ri;
84742c57489SHong Zhang   buf_rj = merge->buf_rj;
84842c57489SHong Zhang   len_s  = merge->len_s;
849fc42d0c8SSatish Balay   ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr);
85042c57489SHong Zhang   ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr);
85142c57489SHong Zhang 
852a9d06231SHong Zhang   ierr = PetscMalloc2(merge->nsend+1,MPI_Request,&s_waits,size,MPI_Status,&status);CHKERRQ(ierr);
85342c57489SHong Zhang   for (proc=0,k=0; proc<size; proc++){
85442c57489SHong Zhang     if (!len_s[proc]) continue;
855de0260b3SHong Zhang     i = merge->owners_co[proc];
856de0260b3SHong Zhang     ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr);
85742c57489SHong Zhang     k++;
85842c57489SHong Zhang   }
8590c468ba9SBarry Smith   if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);}
8600c468ba9SBarry Smith   if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);}
86142c57489SHong Zhang 
862a9d06231SHong Zhang   ierr = PetscFree2(s_waits,status);CHKERRQ(ierr);
86342c57489SHong Zhang   ierr = PetscFree(r_waits);CHKERRQ(ierr);
86482412ba7SHong Zhang   ierr = PetscFree(coa);CHKERRQ(ierr);
86542c57489SHong Zhang 
866a9d06231SHong Zhang   /* 4) insert local Cseq and received values into Cmpi */
867a9d06231SHong Zhang   /*------------------------------------------------------*/
868687f1162SBarry Smith   ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr);
86942c57489SHong Zhang   for (k=0; k<merge->nrecv; k++){
87042c57489SHong Zhang     buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */
87142c57489SHong Zhang     nrows       = *(buf_ri_k[k]);
87242c57489SHong Zhang     nextrow[k]  = buf_ri_k[k]+1;  /* next row number of k-th recved i-structure */
87382412ba7SHong Zhang     nextci[k]   = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure  */
87442c57489SHong Zhang   }
87542c57489SHong Zhang 
876de0260b3SHong Zhang   for (i=0; i<cm; i++) {
87782412ba7SHong Zhang     row = owners[rank] + i; /* global row index of C_seq */
87842c57489SHong Zhang     bj_i = bj + bi[i];  /* col indices of the i-th row of C */
879de0260b3SHong Zhang     ba_i = ba + bi[i];
88082412ba7SHong Zhang     bnz  = bi[i+1] - bi[i];
88142c57489SHong Zhang     /* add received vals into ba */
88242c57489SHong Zhang     for (k=0; k<merge->nrecv; k++){ /* k-th received message */
88342c57489SHong Zhang       /* i-th row */
88442c57489SHong Zhang       if (i == *nextrow[k]) {
88582412ba7SHong Zhang         cnz = *(nextci[k]+1) - *nextci[k];
88682412ba7SHong Zhang         cj  = buf_rj[k] + *(nextci[k]);
88782412ba7SHong Zhang         ca  = abuf_r[k] + *(nextci[k]);
88882412ba7SHong Zhang         nextcj = 0;
88982412ba7SHong Zhang         for (j=0; nextcj<cnz; j++){
89082412ba7SHong Zhang           if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */
89182412ba7SHong Zhang             ba_i[j] += ca[nextcj++];
89242c57489SHong Zhang           }
89342c57489SHong Zhang         }
89482412ba7SHong Zhang         nextrow[k]++; nextci[k]++;
895c5af039cSHong Zhang         ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr);
89642c57489SHong Zhang       }
89742c57489SHong Zhang     }
89882412ba7SHong Zhang     ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr);
89942c57489SHong Zhang   }
90042c57489SHong Zhang   ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90142c57489SHong Zhang   ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90242c57489SHong Zhang 
903de0260b3SHong Zhang   ierr = PetscFree(ba);CHKERRQ(ierr);
904c05d87d6SBarry Smith   ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr);
90542c57489SHong Zhang   ierr = PetscFree(abuf_r);CHKERRQ(ierr);
9060572522cSBarry Smith   ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr);
90742c57489SHong Zhang   PetscFunctionReturn(0);
90842c57489SHong Zhang }
909