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 742c57489SHong Zhang #include "src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/ 842c57489SHong Zhang #include "src/mat/utils/freespace.h" 942c57489SHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h" 1042c57489SHong Zhang #include "petscbt.h" 1142c57489SHong Zhang 1242c57489SHong Zhang #undef __FUNCT__ 13*7a7894deSKris Buschelman #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ" 14*7a7894deSKris Buschelman PetscErrorCode MatPtAPSymbolic_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 1542c57489SHong Zhang { 1642c57489SHong Zhang PetscErrorCode ierr; 1742c57489SHong Zhang 1842c57489SHong Zhang PetscFunctionBegin; 19*7a7894deSKris Buschelman if (!P->ops->ptapsymbolic_mpiaij) { 20*7a7894deSKris Buschelman SETERRQ2(PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",A->type_name,P->type_name); 2142c57489SHong Zhang } 22*7a7894deSKris Buschelman ierr = (*P->ops->ptapsymbolic_mpiaij)(A,P,fill,C);CHKERRQ(ierr); 23*7a7894deSKris Buschelman PetscFunctionReturn(0); 24*7a7894deSKris Buschelman } 25*7a7894deSKris Buschelman 26*7a7894deSKris Buschelman #undef __FUNCT__ 27*7a7894deSKris Buschelman #define __FUNCT__ "MatPtAPNumeric_MPIAIJ" 28*7a7894deSKris Buschelman PetscErrorCode MatPtAPNumeric_MPIAIJ(Mat A,Mat P,Mat C) 29*7a7894deSKris Buschelman { 30*7a7894deSKris Buschelman PetscErrorCode ierr; 31*7a7894deSKris Buschelman 32*7a7894deSKris Buschelman PetscFunctionBegin; 33*7a7894deSKris Buschelman if (!P->ops->ptapnumeric_mpiaij) { 34*7a7894deSKris Buschelman SETERRQ2(PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",A->type_name,P->type_name); 35*7a7894deSKris Buschelman } 36*7a7894deSKris Buschelman ierr = (*P->ops->ptapnumeric_mpiaij)(A,P,C);CHKERRQ(ierr); 3742c57489SHong Zhang PetscFunctionReturn(0); 3842c57489SHong Zhang } 3942c57489SHong Zhang 4042c57489SHong Zhang #undef __FUNCT__ 4142c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ" 4242c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 4342c57489SHong Zhang { 4442c57489SHong Zhang PetscErrorCode ierr; 456c6e00e0SHong Zhang Mat B_mpi; 46de0260b3SHong Zhang Mat_MatMatMultMPI *ap; 4742c57489SHong Zhang PetscObjectContainer container; 4842c57489SHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 4983445d95SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 5042c57489SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 5142c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 52ded4f561SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ; 5383445d95SHong Zhang PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz; 5413f74950SBarry Smith PetscInt nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row; 55a6b2eed2SHong Zhang PetscInt am=A->m,pN=P->N,pn=P->n; 5642c57489SHong Zhang PetscBT lnkbt; 5742c57489SHong Zhang MPI_Comm comm=A->comm; 58ded4f561SHong Zhang PetscMPIInt size,rank,tag,*len_si,*len_s,*len_ri; 5942c57489SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 6042c57489SHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 61ded4f561SHong Zhang PetscInt nzi,*bi,*bj; 6242c57489SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 63ded4f561SHong Zhang MPI_Request *swaits,*rwaits; 64ded4f561SHong Zhang MPI_Status *sstatus,rstatus; 6542c57489SHong Zhang Mat_Merge_SeqsToMPI *merge; 66ded4f561SHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,pon; 6713f74950SBarry Smith PetscMPIInt j; 6842c57489SHong Zhang 6942c57489SHong Zhang PetscFunctionBegin; 7083445d95SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 7183445d95SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 7283445d95SHong Zhang 736c6e00e0SHong Zhang /* destroy the container 'Mat_MatMatMultMPI' in case that P is attached to it */ 746c6e00e0SHong Zhang ierr = PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 756c6e00e0SHong Zhang if (container) { 766c6e00e0SHong Zhang ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr); 776c6e00e0SHong Zhang ierr = PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",0);CHKERRQ(ierr); 786c6e00e0SHong Zhang } 796c6e00e0SHong Zhang 806c6e00e0SHong Zhang /* create the container 'Mat_MatMatMultMPI' and attach it to P */ 816c6e00e0SHong Zhang ierr = PetscNew(Mat_MatMatMultMPI,&ap);CHKERRQ(ierr); 826c6e00e0SHong Zhang ap->abi=PETSC_NULL; ap->abj=PETSC_NULL; 836c6e00e0SHong Zhang ap->abnz_max = 0; 846c6e00e0SHong Zhang 856c6e00e0SHong Zhang ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 866c6e00e0SHong Zhang ierr = PetscObjectContainerSetPointer(container,ap);CHKERRQ(ierr); 876c6e00e0SHong Zhang ierr = PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr); 886c6e00e0SHong Zhang ap->MatDestroy = P->ops->destroy; 896c6e00e0SHong Zhang P->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 906c6e00e0SHong Zhang ap->reuse = MAT_INITIAL_MATRIX; 916c6e00e0SHong Zhang ierr = PetscObjectContainerSetUserDestroy(container,PetscObjectContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr); 926c6e00e0SHong Zhang 936c6e00e0SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 946c6e00e0SHong Zhang ierr = MatGetBrowsOfAoCols(A,P,MAT_INITIAL_MATRIX,&ap->startsj,&ap->bufa,&ap->B_oth);CHKERRQ(ierr); 956c6e00e0SHong Zhang /* get P_loc by taking all local rows of P */ 966c6e00e0SHong Zhang ierr = MatGetLocalMat(P,MAT_INITIAL_MATRIX,&ap->B_loc);CHKERRQ(ierr); 976c6e00e0SHong Zhang 986c6e00e0SHong Zhang p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data; 996c6e00e0SHong Zhang p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data; 1006c6e00e0SHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 1016c6e00e0SHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 10242c57489SHong Zhang 10383445d95SHong Zhang /* first, compute symbolic AP = A_loc*P */ 104de0260b3SHong Zhang /*--------------------------------------*/ 10582412ba7SHong Zhang ierr = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr); 10682412ba7SHong Zhang ap->abi = api; 10782412ba7SHong Zhang api[0] = 0; 10883445d95SHong Zhang 10983445d95SHong Zhang /* create and initialize a linked list */ 11083445d95SHong Zhang nlnk = pN+1; 11183445d95SHong Zhang ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 11283445d95SHong Zhang 11382412ba7SHong Zhang /* Initial FreeSpace size is fill*nnz(A) */ 114de0260b3SHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(adi[am]+aoi[am])),&free_space);CHKERRQ(ierr); 115f4a743e1SHong Zhang current_space = free_space; 116f4a743e1SHong Zhang 117f4a743e1SHong Zhang for (i=0;i<am;i++) { 11882412ba7SHong Zhang apnz = 0; 119f4a743e1SHong Zhang /* diagonal portion of A */ 120ded4f561SHong Zhang nzi = adi[i+1] - adi[i]; 121ded4f561SHong Zhang for (j=0; j<nzi; j++){ 12282412ba7SHong Zhang row = *adj++; 12382412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 124ded4f561SHong Zhang Jptr = pj_loc + pi_loc[row]; 12583445d95SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 126ded4f561SHong Zhang ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 12782412ba7SHong Zhang apnz += nlnk; 128f4a743e1SHong Zhang } 129f4a743e1SHong Zhang /* off-diagonal portion of A */ 130ded4f561SHong Zhang nzi = aoi[i+1] - aoi[i]; 131ded4f561SHong Zhang for (j=0; j<nzi; j++){ 13282412ba7SHong Zhang row = *aoj++; 13382412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 134ded4f561SHong Zhang Jptr = pj_oth + pi_oth[row]; 135ded4f561SHong Zhang ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 13682412ba7SHong Zhang apnz += nlnk; 137f4a743e1SHong Zhang } 138f4a743e1SHong Zhang 13982412ba7SHong Zhang api[i+1] = api[i] + apnz; 14082412ba7SHong Zhang if (ap->abnz_max < apnz) ap->abnz_max = apnz; 141f4a743e1SHong Zhang 14283445d95SHong Zhang /* if free space is not available, double the total space in the list */ 14382412ba7SHong Zhang if (current_space->local_remaining<apnz) { 144f4a743e1SHong Zhang ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 145f4a743e1SHong Zhang } 146f4a743e1SHong Zhang 147f4a743e1SHong Zhang /* Copy data into free space, then initialize lnk */ 14882412ba7SHong Zhang ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 14982412ba7SHong Zhang current_space->array += apnz; 15082412ba7SHong Zhang current_space->local_used += apnz; 15182412ba7SHong Zhang current_space->local_remaining -= apnz; 152f4a743e1SHong Zhang } 15382412ba7SHong Zhang /* Allocate space for apj, initialize apj, and */ 154f4a743e1SHong Zhang /* destroy list of free space and other temporary array(s) */ 15582412ba7SHong Zhang ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ap->abj);CHKERRQ(ierr); 15682412ba7SHong Zhang apj = ap->abj; 157de0260b3SHong Zhang ierr = MakeSpaceContiguous(&free_space,ap->abj);CHKERRQ(ierr); 158f4a743e1SHong Zhang 159ded4f561SHong Zhang /* determine symbolic Co=(p->B)^T*AP - send to others */ 160ded4f561SHong Zhang /*----------------------------------------------------*/ 161de0260b3SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 16242c57489SHong Zhang 163ded4f561SHong Zhang /* then, compute symbolic Co = (p->B)^T*AP */ 164de0260b3SHong Zhang pon = (p->B)->n; /* total num of rows to be sent to other processors 16583445d95SHong Zhang >= (num of nonzero rows of C_seq) - pn */ 166de0260b3SHong Zhang ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr); 167de0260b3SHong Zhang coi[0] = 0; 16842c57489SHong Zhang 16982412ba7SHong Zhang /* set initial free space to be 3*pon*max( nnz(AP) per row) */ 17082412ba7SHong Zhang nnz = 3*pon*ap->abnz_max + 1; 171de0260b3SHong Zhang ierr = GetMoreSpace(nnz,&free_space); 17242c57489SHong Zhang current_space = free_space; 17342c57489SHong Zhang 174de0260b3SHong Zhang for (i=0; i<pon; i++) { 175ded4f561SHong Zhang nnz = 0; 17682412ba7SHong Zhang pnz = poti[i+1] - poti[i]; 17782412ba7SHong Zhang j = pnz; 178de0260b3SHong Zhang ptJ = potj + poti[i+1]; 17983445d95SHong Zhang while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 18083445d95SHong Zhang j--; ptJ--; 18182412ba7SHong Zhang row = *ptJ; /* row of AP == col of Pot */ 18282412ba7SHong Zhang apnz = api[row+1] - api[row]; 183ded4f561SHong Zhang Jptr = apj + api[row]; 18483445d95SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 185ded4f561SHong Zhang ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 186ded4f561SHong Zhang nnz += nlnk; 18742c57489SHong Zhang } 18842c57489SHong Zhang 18983445d95SHong Zhang /* If free space is not available, double the total space in the list */ 190ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 19142c57489SHong Zhang ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 19242c57489SHong Zhang } 19342c57489SHong Zhang 19442c57489SHong Zhang /* Copy data into free space, and zero out denserows */ 195ded4f561SHong Zhang ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 196ded4f561SHong Zhang current_space->array += nnz; 197ded4f561SHong Zhang current_space->local_used += nnz; 198ded4f561SHong Zhang current_space->local_remaining -= nnz; 199ded4f561SHong Zhang coi[i+1] = coi[i] + nnz; 20042c57489SHong Zhang } 201de0260b3SHong Zhang ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr); 202de0260b3SHong Zhang ierr = MakeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 203de0260b3SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 20442c57489SHong Zhang 205ded4f561SHong Zhang /* send j-array (coj) of Co to other processors */ 206ded4f561SHong Zhang /*----------------------------------------------*/ 20742c57489SHong Zhang /* determine row ownership */ 20883445d95SHong Zhang ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr); 20942c57489SHong Zhang ierr = PetscMapCreate(comm,&merge->rowmap);CHKERRQ(ierr); 21042c57489SHong Zhang ierr = PetscMapSetLocalSize(merge->rowmap,pn);CHKERRQ(ierr); 21142c57489SHong Zhang ierr = PetscMapSetType(merge->rowmap,MAP_MPI);CHKERRQ(ierr); 21242c57489SHong Zhang ierr = PetscMapGetGlobalRange(merge->rowmap,&owners);CHKERRQ(ierr); 21342c57489SHong Zhang 21442c57489SHong Zhang /* determine the number of messages to send, their lengths */ 21542c57489SHong Zhang ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr); 21683445d95SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 21742c57489SHong Zhang ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr); 21842c57489SHong Zhang len_s = merge->len_s; 21942c57489SHong Zhang merge->nsend = 0; 22083445d95SHong Zhang 221de0260b3SHong Zhang ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr); 222de0260b3SHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 223de0260b3SHong Zhang 22483445d95SHong Zhang proc = 0; 225de0260b3SHong Zhang for (i=0; i<pon; i++){ 226de0260b3SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 227de0260b3SHong Zhang len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 228de0260b3SHong Zhang len_s[proc] += coi[i+1] - coi[i]; 22983445d95SHong Zhang } 230de0260b3SHong Zhang 231de0260b3SHong Zhang len = 0; /* max length of buf_si[] */ 232de0260b3SHong Zhang owners_co[0] = 0; 23342c57489SHong Zhang for (proc=0; proc<size; proc++){ 234de0260b3SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 23583445d95SHong Zhang if (len_si[proc]){ 23642c57489SHong Zhang merge->nsend++; 23783445d95SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); 23842c57489SHong Zhang len += len_si[proc]; 23942c57489SHong Zhang } 24042c57489SHong Zhang } 24142c57489SHong Zhang 242ded4f561SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 24342c57489SHong Zhang ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 24442c57489SHong Zhang ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 24542c57489SHong Zhang 246ded4f561SHong Zhang /* post the Irecv and Isend of coj */ 247ded4f561SHong Zhang ierr = PetscObjectGetNewTag((PetscObject)merge->rowmap,&tag);CHKERRQ(ierr); 248ded4f561SHong Zhang ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 24942c57489SHong Zhang 250ded4f561SHong Zhang ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr); 25142c57489SHong Zhang 25242c57489SHong Zhang for (proc=0, k=0; proc<size; proc++){ 25342c57489SHong Zhang if (!len_s[proc]) continue; 254de0260b3SHong Zhang i = owners_co[proc]; 255ded4f561SHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr); 25642c57489SHong Zhang k++; 25742c57489SHong Zhang } 25842c57489SHong Zhang 259ded4f561SHong Zhang /* receives and sends of coj are complete */ 260ded4f561SHong Zhang ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr); 261ded4f561SHong Zhang i = merge->nrecv; 262ded4f561SHong Zhang while (i--) { 263ded4f561SHong Zhang ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 264ded4f561SHong Zhang } 265ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 266ded4f561SHong Zhang if (merge->nsend){ 267ded4f561SHong Zhang ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr); 26882412ba7SHong Zhang } 26982412ba7SHong Zhang 270ded4f561SHong Zhang /* send and recv coi */ 271ded4f561SHong Zhang /*-------------------*/ 272ded4f561SHong Zhang ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 27342c57489SHong Zhang 27442c57489SHong Zhang ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr); 27542c57489SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 27642c57489SHong Zhang for (proc=0,k=0; proc<size; proc++){ 27742c57489SHong Zhang if (!len_s[proc]) continue; 27842c57489SHong Zhang /* form outgoing message for i-structure: 27942c57489SHong Zhang buf_si[0]: nrows to be sent 28042c57489SHong Zhang [1:nrows]: row index (global) 28142c57489SHong Zhang [nrows+1:2*nrows+1]: i-structure index 28242c57489SHong Zhang */ 28342c57489SHong Zhang /*-------------------------------------------*/ 28442c57489SHong Zhang nrows = len_si[proc]/2 - 1; 28542c57489SHong Zhang buf_si_i = buf_si + nrows+1; 28642c57489SHong Zhang buf_si[0] = nrows; 28742c57489SHong Zhang buf_si_i[0] = 0; 28842c57489SHong Zhang nrows = 0; 289de0260b3SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++){ 290ded4f561SHong Zhang nzi = coi[i+1] - coi[i]; 291ded4f561SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 292de0260b3SHong Zhang buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */ 29342c57489SHong Zhang nrows++; 29442c57489SHong Zhang } 295ded4f561SHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr); 29642c57489SHong Zhang k++; 29742c57489SHong Zhang buf_si += len_si[proc]; 29842c57489SHong Zhang } 299ded4f561SHong Zhang i = merge->nrecv; 300ded4f561SHong Zhang while (i--) { 301ded4f561SHong Zhang ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 302ded4f561SHong Zhang } 303ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 304ded4f561SHong Zhang if (merge->nsend){ 305ded4f561SHong Zhang ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr); 306ded4f561SHong Zhang } 30742c57489SHong Zhang 30842c57489SHong Zhang ierr = PetscLogInfo((PetscObject)A,"MatMerge_SeqsToMPI: nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr); 30942c57489SHong Zhang for (i=0; i<merge->nrecv; i++){ 31042c57489SHong Zhang ierr = PetscLogInfo((PetscObject)A,"MatMerge_SeqsToMPI: recv len_ri=%d, len_rj=%d from [%d]\n",len_ri[i],merge->len_r[i],merge->id_r[i]);CHKERRQ(ierr); 31142c57489SHong Zhang } 31242c57489SHong Zhang 31342c57489SHong Zhang ierr = PetscFree(len_si);CHKERRQ(ierr); 31442c57489SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 315ded4f561SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 316ded4f561SHong Zhang ierr = PetscFree(sstatus);CHKERRQ(ierr); 31742c57489SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 31842c57489SHong Zhang 319de0260b3SHong Zhang /* compute the local portion of C (mpi mat) */ 320de0260b3SHong Zhang /*------------------------------------------*/ 321ded4f561SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 322ded4f561SHong Zhang 32342c57489SHong Zhang /* allocate bi array and free space for accumulating nonzero column info */ 32442c57489SHong Zhang ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 32542c57489SHong Zhang bi[0] = 0; 32642c57489SHong Zhang 327ded4f561SHong Zhang /* set initial free space to be 3*pn*max( nnz(AP) per row) */ 328ded4f561SHong Zhang nnz = 3*pn*ap->abnz_max + 1; 329ded4f561SHong Zhang ierr = GetMoreSpace(nnz,&free_space); 33042c57489SHong Zhang current_space = free_space; 33142c57489SHong Zhang 33242c57489SHong Zhang ierr = PetscMalloc((3*merge->nrecv+1)*sizeof(PetscInt**),&buf_ri_k);CHKERRQ(ierr); 33342c57489SHong Zhang nextrow = buf_ri_k + merge->nrecv; 33442c57489SHong Zhang nextci = nextrow + merge->nrecv; 33542c57489SHong Zhang for (k=0; k<merge->nrecv; k++){ 33642c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 33742c57489SHong Zhang nrows = *buf_ri_k[k]; 33842c57489SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 33942c57489SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 34042c57489SHong Zhang } 34142c57489SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 34242c57489SHong Zhang for (i=0; i<pn; i++) { 343ded4f561SHong Zhang /* add pdt[i,:]*AP into lnk */ 344ded4f561SHong Zhang nnz = 0; 345ded4f561SHong Zhang pnz = pdti[i+1] - pdti[i]; 346ded4f561SHong Zhang j = pnz; 347ded4f561SHong Zhang ptJ = pdtj + pdti[i+1]; 348ded4f561SHong Zhang while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 349ded4f561SHong Zhang j--; ptJ--; 350ded4f561SHong Zhang row = *ptJ; /* row of AP == col of Pt */ 351ded4f561SHong Zhang apnz = api[row+1] - api[row]; 352ded4f561SHong Zhang Jptr = apj + api[row]; 353ded4f561SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 354ded4f561SHong Zhang ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 355ded4f561SHong Zhang nnz += nlnk; 356ded4f561SHong Zhang } 35742c57489SHong Zhang /* add received col data into lnk */ 35842c57489SHong Zhang for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 35942c57489SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 360ded4f561SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 361ded4f561SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 362ded4f561SHong Zhang ierr = PetscLLAdd(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 363ded4f561SHong Zhang nnz += nlnk; 36442c57489SHong Zhang nextrow[k]++; nextci[k]++; 36542c57489SHong Zhang } 36642c57489SHong Zhang } 36742c57489SHong Zhang 36842c57489SHong Zhang /* if free space is not available, make more free space */ 369ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 37042c57489SHong Zhang ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 37142c57489SHong Zhang } 37242c57489SHong Zhang /* copy data into free space, then initialize lnk */ 373ded4f561SHong Zhang ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 374ded4f561SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 375ded4f561SHong Zhang current_space->array += nnz; 376ded4f561SHong Zhang current_space->local_used += nnz; 377ded4f561SHong Zhang current_space->local_remaining -= nnz; 378ded4f561SHong Zhang bi[i+1] = bi[i] + nnz; 37942c57489SHong Zhang } 380ded4f561SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 38142c57489SHong Zhang ierr = PetscFree(buf_ri_k);CHKERRQ(ierr); 38242c57489SHong Zhang 38342c57489SHong Zhang ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 38442c57489SHong Zhang ierr = MakeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 38542c57489SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 38642c57489SHong Zhang 38742c57489SHong Zhang /* create symbolic parallel matrix B_mpi */ 38842c57489SHong Zhang /*---------------------------------------*/ 38942c57489SHong Zhang ierr = MatCreate(comm,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE,&B_mpi);CHKERRQ(ierr); 39042c57489SHong Zhang ierr = MatSetType(B_mpi,MATMPIAIJ);CHKERRQ(ierr); 39142c57489SHong Zhang ierr = MatMPIAIJSetPreallocation(B_mpi,0,dnz,0,onz);CHKERRQ(ierr); 39242c57489SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 39342c57489SHong Zhang 39482412ba7SHong Zhang /* B_mpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 39542c57489SHong Zhang B_mpi->assembled = PETSC_FALSE; 39642c57489SHong Zhang B_mpi->ops->destroy = MatDestroy_MPIAIJ_SeqsToMPI; 39742c57489SHong Zhang merge->bi = bi; 39842c57489SHong Zhang merge->bj = bj; 399de0260b3SHong Zhang merge->coi = coi; 400de0260b3SHong Zhang merge->coj = coj; 40142c57489SHong Zhang merge->buf_ri = buf_ri; 40242c57489SHong Zhang merge->buf_rj = buf_rj; 403de0260b3SHong Zhang merge->owners_co = owners_co; 40442c57489SHong Zhang 40542c57489SHong Zhang /* attach the supporting struct to B_mpi for reuse */ 40642c57489SHong Zhang ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 40742c57489SHong Zhang ierr = PetscObjectContainerSetPointer(container,merge);CHKERRQ(ierr); 40842c57489SHong Zhang ierr = PetscObjectCompose((PetscObject)B_mpi,"MatMergeSeqsToMPI",(PetscObject)container);CHKERRQ(ierr); 40942c57489SHong Zhang *C = B_mpi; 41042c57489SHong Zhang PetscFunctionReturn(0); 41142c57489SHong Zhang } 41242c57489SHong Zhang 41342c57489SHong Zhang #undef __FUNCT__ 41442c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ" 41542c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 41642c57489SHong Zhang { 41742c57489SHong Zhang PetscErrorCode ierr; 41842c57489SHong Zhang Mat_Merge_SeqsToMPI *merge; 419de0260b3SHong Zhang Mat_MatMatMultMPI *ap; 42042c57489SHong Zhang PetscObjectContainer cont_merge,cont_ptap; 421de0260b3SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 42242c57489SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 423de0260b3SHong Zhang Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 42442c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 42582412ba7SHong Zhang PetscInt *adi=ad->i,*aoi=ao->i,*adj=ad->j,*aoj=ao->j,*apJ,nextp,flops=0; 42682412ba7SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj; 42782412ba7SHong Zhang PetscInt i,j,k,anz,pnz,apnz,nextap,row,*cj; 42882412ba7SHong Zhang MatScalar *ada=ad->a,*aoa=ao->a,*apa,*pa,*ca,*pa_loc,*pa_oth; 42982412ba7SHong Zhang PetscInt am=A->m,cm=C->m,pon=(p->B)->n; 43042c57489SHong Zhang MPI_Comm comm=C->comm; 43142c57489SHong Zhang PetscMPIInt size,rank,taga,*len_s; 43282412ba7SHong Zhang PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 43342c57489SHong Zhang PetscInt **buf_ri,**buf_rj; 43450f5bf86SHong Zhang PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 43542c57489SHong Zhang MPI_Request *s_waits,*r_waits; 43642c57489SHong Zhang MPI_Status *status; 43782412ba7SHong Zhang MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 43882412ba7SHong Zhang PetscInt *api,*apj,*coi,*coj; 43982412ba7SHong Zhang PetscInt *poJ=po->j,*pdJ=pd->j,pcstart=p->cstart,pcend=p->cend; 44042c57489SHong Zhang 44142c57489SHong Zhang PetscFunctionBegin; 44242c57489SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 44342c57489SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 44442c57489SHong Zhang 44542c57489SHong Zhang ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&cont_merge);CHKERRQ(ierr); 44642c57489SHong Zhang if (cont_merge) { 44742c57489SHong Zhang ierr = PetscObjectContainerGetPointer(cont_merge,(void **)&merge);CHKERRQ(ierr); 44842c57489SHong Zhang } else { 44942c57489SHong Zhang SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses an object container"); 45042c57489SHong Zhang } 4516c6e00e0SHong Zhang 452f33d1a9aSHong Zhang ierr = PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&cont_ptap);CHKERRQ(ierr); 45342c57489SHong Zhang if (cont_ptap) { 454de0260b3SHong Zhang ierr = PetscObjectContainerGetPointer(cont_ptap,(void **)&ap);CHKERRQ(ierr); 4556c6e00e0SHong Zhang if (ap->reuse == MAT_INITIAL_MATRIX){ 4566c6e00e0SHong Zhang ap->reuse = MAT_REUSE_MATRIX; 4576c6e00e0SHong Zhang } else { /* update numerical values of P_oth and P_loc */ 4586c6e00e0SHong Zhang ierr = MatGetBrowsOfAoCols(A,P,MAT_REUSE_MATRIX,&ap->startsj,&ap->bufa,&ap->B_oth);CHKERRQ(ierr); 4596c6e00e0SHong Zhang ierr = MatGetLocalMat(P,MAT_REUSE_MATRIX,&ap->B_loc);CHKERRQ(ierr); 4606c6e00e0SHong Zhang } 46142c57489SHong Zhang } else { 46242c57489SHong Zhang SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container"); 46342c57489SHong Zhang } 46442c57489SHong Zhang 46542c57489SHong Zhang /* get data from symbolic products */ 466de0260b3SHong Zhang p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data; 467de0260b3SHong Zhang p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data; 46882412ba7SHong Zhang pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a,pA=pa_loc; 46942c57489SHong Zhang pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 47042c57489SHong Zhang 471de0260b3SHong Zhang coi = merge->coi; coj = merge->coj; 472de0260b3SHong Zhang ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr); 473de0260b3SHong Zhang ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr); 474de0260b3SHong Zhang 475de0260b3SHong Zhang bi = merge->bi; bj = merge->bj; 476de0260b3SHong Zhang ierr = PetscMapGetGlobalRange(merge->rowmap,&owners);CHKERRQ(ierr); 477de0260b3SHong Zhang ierr = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr); 478de0260b3SHong Zhang ierr = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr); 47942c57489SHong Zhang 480f4a743e1SHong Zhang /* get data from symbolic A*P */ 481de0260b3SHong Zhang ierr = PetscMalloc((ap->abnz_max+1)*sizeof(MatScalar),&apa);CHKERRQ(ierr); 482de0260b3SHong Zhang ierr = PetscMemzero(apa,ap->abnz_max*sizeof(MatScalar));CHKERRQ(ierr); 48342c57489SHong Zhang 484f4a743e1SHong Zhang /* compute numeric C_seq=P_loc^T*A_loc*P */ 48582412ba7SHong Zhang api = ap->abi; apj = ap->abj; 48642c57489SHong Zhang for (i=0;i<am;i++) { 487f4a743e1SHong Zhang /* form i-th sparse row of A*P */ 48882412ba7SHong Zhang apnz = api[i+1] - api[i]; 48982412ba7SHong Zhang apJ = apj + api[i]; 49042c57489SHong Zhang /* diagonal portion of A */ 49182412ba7SHong Zhang anz = adi[i+1] - adi[i]; 49282412ba7SHong Zhang for (j=0;j<anz;j++) { 49382412ba7SHong Zhang row = *adj++; 49482412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 49582412ba7SHong Zhang pj = pj_loc + pi_loc[row]; 49682412ba7SHong Zhang pa = pa_loc + pi_loc[row]; 497f4a743e1SHong Zhang nextp = 0; 49882412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 49982412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 50082412ba7SHong Zhang apa[k] += (*ada)*pa[nextp++]; 50142c57489SHong Zhang } 50242c57489SHong Zhang } 50382412ba7SHong Zhang flops += 2*pnz; 50442c57489SHong Zhang ada++; 50542c57489SHong Zhang } 50642c57489SHong Zhang /* off-diagonal portion of A */ 50782412ba7SHong Zhang anz = aoi[i+1] - aoi[i]; 50882412ba7SHong Zhang for (j=0; j<anz; j++) { 50982412ba7SHong Zhang row = *aoj++; 51082412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 51182412ba7SHong Zhang pj = pj_oth + pi_oth[row]; 51282412ba7SHong Zhang pa = pa_oth + pi_oth[row]; 513f4a743e1SHong Zhang nextp = 0; 51482412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 51582412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 51682412ba7SHong Zhang apa[k] += (*aoa)*pa[nextp++]; 51742c57489SHong Zhang } 51842c57489SHong Zhang } 51982412ba7SHong Zhang flops += 2*pnz; 52042c57489SHong Zhang aoa++; 52142c57489SHong Zhang } 52242c57489SHong Zhang 52382412ba7SHong Zhang /* Compute P_loc[i,:]^T*AP[i,:] using outer product */ 52482412ba7SHong Zhang pnz = pi_loc[i+1] - pi_loc[i]; 52582412ba7SHong Zhang for (j=0; j<pnz; j++) { 52642c57489SHong Zhang nextap = 0; 52782412ba7SHong Zhang row = *pJ++; /* global index */ 52882412ba7SHong Zhang if (row < pcstart || row >=pcend) { /* put the value into Co */ 52982412ba7SHong Zhang cj = coj + coi[*poJ]; 53082412ba7SHong Zhang ca = coa + coi[*poJ++]; 53182412ba7SHong Zhang } else { /* put the value into Cd */ 53282412ba7SHong Zhang cj = bj + bi[*pdJ]; 53382412ba7SHong Zhang ca = ba + bi[*pdJ++]; 53442c57489SHong Zhang } 53582412ba7SHong Zhang for (k=0; nextap<apnz; k++) { 53682412ba7SHong Zhang if (cj[k]==apJ[nextap]) ca[k] += (*pA)*apa[nextap++]; 53742c57489SHong Zhang } 53882412ba7SHong Zhang flops += 2*apnz; 53982412ba7SHong Zhang pA++; 540de0260b3SHong Zhang } 541de0260b3SHong Zhang 54242c57489SHong Zhang /* zero the current row info for A*P */ 54382412ba7SHong Zhang ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr); 54442c57489SHong Zhang } 54542c57489SHong Zhang ierr = PetscFree(apa);CHKERRQ(ierr); 54642c57489SHong Zhang 547de0260b3SHong Zhang /* send and recv matrix values */ 548de0260b3SHong Zhang /*-----------------------------*/ 54942c57489SHong Zhang buf_ri = merge->buf_ri; 55042c57489SHong Zhang buf_rj = merge->buf_rj; 55142c57489SHong Zhang len_s = merge->len_s; 55242c57489SHong Zhang ierr = PetscObjectGetNewTag((PetscObject)merge->rowmap,&taga);CHKERRQ(ierr); 55342c57489SHong Zhang ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 55442c57489SHong Zhang 55542c57489SHong Zhang ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&s_waits);CHKERRQ(ierr); 55642c57489SHong Zhang for (proc=0,k=0; proc<size; proc++){ 55742c57489SHong Zhang if (!len_s[proc]) continue; 558de0260b3SHong Zhang i = merge->owners_co[proc]; 559de0260b3SHong Zhang ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 56042c57489SHong Zhang k++; 56142c57489SHong Zhang } 56282412ba7SHong Zhang ierr = PetscMalloc(size*sizeof(MPI_Status),&status);CHKERRQ(ierr); 56342c57489SHong Zhang ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr); 56442c57489SHong Zhang ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr); 56542c57489SHong Zhang ierr = PetscFree(status);CHKERRQ(ierr); 56642c57489SHong Zhang 56742c57489SHong Zhang ierr = PetscFree(s_waits);CHKERRQ(ierr); 56842c57489SHong Zhang ierr = PetscFree(r_waits);CHKERRQ(ierr); 56982412ba7SHong Zhang ierr = PetscFree(coa);CHKERRQ(ierr); 57042c57489SHong Zhang 57182412ba7SHong Zhang /* insert local and received values into C */ 57282412ba7SHong Zhang /*-----------------------------------------*/ 57342c57489SHong Zhang ierr = PetscMalloc((3*merge->nrecv+1)*sizeof(PetscInt**),&buf_ri_k);CHKERRQ(ierr); 57442c57489SHong Zhang nextrow = buf_ri_k + merge->nrecv; 57582412ba7SHong Zhang nextci = nextrow + merge->nrecv; 57642c57489SHong Zhang 57742c57489SHong Zhang for (k=0; k<merge->nrecv; k++){ 57842c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 57942c57489SHong Zhang nrows = *(buf_ri_k[k]); 58042c57489SHong Zhang nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 58182412ba7SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 58242c57489SHong Zhang } 58342c57489SHong Zhang 584de0260b3SHong Zhang for (i=0; i<cm; i++) { 58582412ba7SHong Zhang row = owners[rank] + i; /* global row index of C_seq */ 58642c57489SHong Zhang bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 587de0260b3SHong Zhang ba_i = ba + bi[i]; 58882412ba7SHong Zhang bnz = bi[i+1] - bi[i]; 58942c57489SHong Zhang /* add received vals into ba */ 59042c57489SHong Zhang for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 59142c57489SHong Zhang /* i-th row */ 59242c57489SHong Zhang if (i == *nextrow[k]) { 59382412ba7SHong Zhang cnz = *(nextci[k]+1) - *nextci[k]; 59482412ba7SHong Zhang cj = buf_rj[k] + *(nextci[k]); 59582412ba7SHong Zhang ca = abuf_r[k] + *(nextci[k]); 59682412ba7SHong Zhang nextcj = 0; 59782412ba7SHong Zhang for (j=0; nextcj<cnz; j++){ 59882412ba7SHong Zhang if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */ 59982412ba7SHong Zhang ba_i[j] += ca[nextcj++]; 60042c57489SHong Zhang } 60142c57489SHong Zhang } 60282412ba7SHong Zhang nextrow[k]++; nextci[k]++; 60342c57489SHong Zhang } 60442c57489SHong Zhang } 60582412ba7SHong Zhang ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 60682412ba7SHong Zhang flops += 2*cnz; 60742c57489SHong Zhang } 60882412ba7SHong Zhang ierr = MatSetOption(C,MAT_COLUMNS_SORTED);CHKERRQ(ierr); /* -cause delay? */ 60942c57489SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 61042c57489SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 61142c57489SHong Zhang 612de0260b3SHong Zhang ierr = PetscFree(ba);CHKERRQ(ierr); 61342c57489SHong Zhang ierr = PetscFree(abuf_r);CHKERRQ(ierr); 61442c57489SHong Zhang ierr = PetscFree(buf_ri_k);CHKERRQ(ierr); 61542c57489SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 61642c57489SHong Zhang PetscFunctionReturn(0); 61742c57489SHong Zhang } 618