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__ 1342c57489SHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ" 1442c57489SHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 1542c57489SHong Zhang { 1642c57489SHong Zhang PetscErrorCode ierr; 1742c57489SHong Zhang 1842c57489SHong Zhang PetscFunctionBegin; 1942c57489SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 2042c57489SHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 21*6c6e00e0SHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr); 2242c57489SHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 2342c57489SHong Zhang } 2442c57489SHong Zhang ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 2542c57489SHong Zhang ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr); 2642c57489SHong Zhang ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 2742c57489SHong Zhang PetscFunctionReturn(0); 2842c57489SHong Zhang } 2942c57489SHong Zhang 3042c57489SHong Zhang #undef __FUNCT__ 3142c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ" 3242c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 3342c57489SHong Zhang { 3442c57489SHong Zhang PetscErrorCode ierr; 35*6c6e00e0SHong Zhang Mat B_mpi; 36de0260b3SHong Zhang Mat_MatMatMultMPI *ap; 3742c57489SHong Zhang PetscObjectContainer container; 3842c57489SHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 3983445d95SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 4042c57489SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 4142c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 42ded4f561SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ; 4383445d95SHong Zhang PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz; 4413f74950SBarry Smith PetscInt nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row; 45a6b2eed2SHong Zhang PetscInt am=A->m,pN=P->N,pn=P->n; 4642c57489SHong Zhang PetscBT lnkbt; 4742c57489SHong Zhang MPI_Comm comm=A->comm; 48ded4f561SHong Zhang PetscMPIInt size,rank,tag,*len_si,*len_s,*len_ri; 4942c57489SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 5042c57489SHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 51ded4f561SHong Zhang PetscInt nzi,*bi,*bj; 5242c57489SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 53ded4f561SHong Zhang MPI_Request *swaits,*rwaits; 54ded4f561SHong Zhang MPI_Status *sstatus,rstatus; 5542c57489SHong Zhang Mat_Merge_SeqsToMPI *merge; 56ded4f561SHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,pon; 5713f74950SBarry Smith PetscMPIInt j; 5842c57489SHong Zhang 5942c57489SHong Zhang PetscFunctionBegin; 6083445d95SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 6183445d95SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 6283445d95SHong Zhang 63*6c6e00e0SHong Zhang /* destroy the container 'Mat_MatMatMultMPI' in case that P is attached to it */ 64*6c6e00e0SHong Zhang ierr = PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 65*6c6e00e0SHong Zhang if (container) { 66*6c6e00e0SHong Zhang ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr); 67*6c6e00e0SHong Zhang ierr = PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",0);CHKERRQ(ierr); 68*6c6e00e0SHong Zhang } 69*6c6e00e0SHong Zhang 70*6c6e00e0SHong Zhang /* create the container 'Mat_MatMatMultMPI' and attach it to P */ 71*6c6e00e0SHong Zhang ierr = PetscNew(Mat_MatMatMultMPI,&ap);CHKERRQ(ierr); 72*6c6e00e0SHong Zhang ap->abi=PETSC_NULL; ap->abj=PETSC_NULL; 73*6c6e00e0SHong Zhang ap->abnz_max = 0; 74*6c6e00e0SHong Zhang 75*6c6e00e0SHong Zhang ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 76*6c6e00e0SHong Zhang ierr = PetscObjectContainerSetPointer(container,ap);CHKERRQ(ierr); 77*6c6e00e0SHong Zhang ierr = PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr); 78*6c6e00e0SHong Zhang ap->MatDestroy = P->ops->destroy; 79*6c6e00e0SHong Zhang P->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 80*6c6e00e0SHong Zhang ap->reuse = MAT_INITIAL_MATRIX; 81*6c6e00e0SHong Zhang ierr = PetscObjectContainerSetUserDestroy(container,PetscObjectContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr); 82*6c6e00e0SHong Zhang 83*6c6e00e0SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 84*6c6e00e0SHong Zhang ierr = MatGetBrowsOfAoCols(A,P,MAT_INITIAL_MATRIX,&ap->startsj,&ap->bufa,&ap->B_oth);CHKERRQ(ierr); 85*6c6e00e0SHong Zhang /* get P_loc by taking all local rows of P */ 86*6c6e00e0SHong Zhang ierr = MatGetLocalMat(P,MAT_INITIAL_MATRIX,&ap->B_loc);CHKERRQ(ierr); 87*6c6e00e0SHong Zhang 88*6c6e00e0SHong Zhang p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data; 89*6c6e00e0SHong Zhang p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data; 90*6c6e00e0SHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 91*6c6e00e0SHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 9242c57489SHong Zhang 9383445d95SHong Zhang /* first, compute symbolic AP = A_loc*P */ 94de0260b3SHong Zhang /*--------------------------------------*/ 9582412ba7SHong Zhang ierr = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr); 9682412ba7SHong Zhang ap->abi = api; 9782412ba7SHong Zhang api[0] = 0; 9883445d95SHong Zhang 9983445d95SHong Zhang /* create and initialize a linked list */ 10083445d95SHong Zhang nlnk = pN+1; 10183445d95SHong Zhang ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 10283445d95SHong Zhang 10382412ba7SHong Zhang /* Initial FreeSpace size is fill*nnz(A) */ 104de0260b3SHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(adi[am]+aoi[am])),&free_space);CHKERRQ(ierr); 105f4a743e1SHong Zhang current_space = free_space; 106f4a743e1SHong Zhang 107f4a743e1SHong Zhang for (i=0;i<am;i++) { 10882412ba7SHong Zhang apnz = 0; 109f4a743e1SHong Zhang /* diagonal portion of A */ 110ded4f561SHong Zhang nzi = adi[i+1] - adi[i]; 111ded4f561SHong Zhang for (j=0; j<nzi; j++){ 11282412ba7SHong Zhang row = *adj++; 11382412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 114ded4f561SHong Zhang Jptr = pj_loc + pi_loc[row]; 11583445d95SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 116ded4f561SHong Zhang ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 11782412ba7SHong Zhang apnz += nlnk; 118f4a743e1SHong Zhang } 119f4a743e1SHong Zhang /* off-diagonal portion of A */ 120ded4f561SHong Zhang nzi = aoi[i+1] - aoi[i]; 121ded4f561SHong Zhang for (j=0; j<nzi; j++){ 12282412ba7SHong Zhang row = *aoj++; 12382412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 124ded4f561SHong Zhang Jptr = pj_oth + pi_oth[row]; 125ded4f561SHong Zhang ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 12682412ba7SHong Zhang apnz += nlnk; 127f4a743e1SHong Zhang } 128f4a743e1SHong Zhang 12982412ba7SHong Zhang api[i+1] = api[i] + apnz; 13082412ba7SHong Zhang if (ap->abnz_max < apnz) ap->abnz_max = apnz; 131f4a743e1SHong Zhang 13283445d95SHong Zhang /* if free space is not available, double the total space in the list */ 13382412ba7SHong Zhang if (current_space->local_remaining<apnz) { 134f4a743e1SHong Zhang ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 135f4a743e1SHong Zhang } 136f4a743e1SHong Zhang 137f4a743e1SHong Zhang /* Copy data into free space, then initialize lnk */ 13882412ba7SHong Zhang ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 13982412ba7SHong Zhang current_space->array += apnz; 14082412ba7SHong Zhang current_space->local_used += apnz; 14182412ba7SHong Zhang current_space->local_remaining -= apnz; 142f4a743e1SHong Zhang } 14382412ba7SHong Zhang /* Allocate space for apj, initialize apj, and */ 144f4a743e1SHong Zhang /* destroy list of free space and other temporary array(s) */ 14582412ba7SHong Zhang ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ap->abj);CHKERRQ(ierr); 14682412ba7SHong Zhang apj = ap->abj; 147de0260b3SHong Zhang ierr = MakeSpaceContiguous(&free_space,ap->abj);CHKERRQ(ierr); 148f4a743e1SHong Zhang 149ded4f561SHong Zhang /* determine symbolic Co=(p->B)^T*AP - send to others */ 150ded4f561SHong Zhang /*----------------------------------------------------*/ 151de0260b3SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 15242c57489SHong Zhang 153ded4f561SHong Zhang /* then, compute symbolic Co = (p->B)^T*AP */ 154de0260b3SHong Zhang pon = (p->B)->n; /* total num of rows to be sent to other processors 15583445d95SHong Zhang >= (num of nonzero rows of C_seq) - pn */ 156de0260b3SHong Zhang ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr); 157de0260b3SHong Zhang coi[0] = 0; 15842c57489SHong Zhang 15982412ba7SHong Zhang /* set initial free space to be 3*pon*max( nnz(AP) per row) */ 16082412ba7SHong Zhang nnz = 3*pon*ap->abnz_max + 1; 161de0260b3SHong Zhang ierr = GetMoreSpace(nnz,&free_space); 16242c57489SHong Zhang current_space = free_space; 16342c57489SHong Zhang 164de0260b3SHong Zhang for (i=0; i<pon; i++) { 165ded4f561SHong Zhang nnz = 0; 16682412ba7SHong Zhang pnz = poti[i+1] - poti[i]; 16782412ba7SHong Zhang j = pnz; 168de0260b3SHong Zhang ptJ = potj + poti[i+1]; 16983445d95SHong Zhang while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 17083445d95SHong Zhang j--; ptJ--; 17182412ba7SHong Zhang row = *ptJ; /* row of AP == col of Pot */ 17282412ba7SHong Zhang apnz = api[row+1] - api[row]; 173ded4f561SHong Zhang Jptr = apj + api[row]; 17483445d95SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 175ded4f561SHong Zhang ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 176ded4f561SHong Zhang nnz += nlnk; 17742c57489SHong Zhang } 17842c57489SHong Zhang 17983445d95SHong Zhang /* If free space is not available, double the total space in the list */ 180ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 18142c57489SHong Zhang ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 18242c57489SHong Zhang } 18342c57489SHong Zhang 18442c57489SHong Zhang /* Copy data into free space, and zero out denserows */ 185ded4f561SHong Zhang ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 186ded4f561SHong Zhang current_space->array += nnz; 187ded4f561SHong Zhang current_space->local_used += nnz; 188ded4f561SHong Zhang current_space->local_remaining -= nnz; 189ded4f561SHong Zhang coi[i+1] = coi[i] + nnz; 19042c57489SHong Zhang } 191de0260b3SHong Zhang ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr); 192de0260b3SHong Zhang ierr = MakeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 193de0260b3SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 19442c57489SHong Zhang 195ded4f561SHong Zhang /* send j-array (coj) of Co to other processors */ 196ded4f561SHong Zhang /*----------------------------------------------*/ 19742c57489SHong Zhang /* determine row ownership */ 19883445d95SHong Zhang ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr); 19942c57489SHong Zhang ierr = PetscMapCreate(comm,&merge->rowmap);CHKERRQ(ierr); 20042c57489SHong Zhang ierr = PetscMapSetLocalSize(merge->rowmap,pn);CHKERRQ(ierr); 20142c57489SHong Zhang ierr = PetscMapSetType(merge->rowmap,MAP_MPI);CHKERRQ(ierr); 20242c57489SHong Zhang ierr = PetscMapGetGlobalRange(merge->rowmap,&owners);CHKERRQ(ierr); 20342c57489SHong Zhang 20442c57489SHong Zhang /* determine the number of messages to send, their lengths */ 20542c57489SHong Zhang ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr); 20683445d95SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 20742c57489SHong Zhang ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr); 20842c57489SHong Zhang len_s = merge->len_s; 20942c57489SHong Zhang merge->nsend = 0; 21083445d95SHong Zhang 211de0260b3SHong Zhang ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr); 212de0260b3SHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 213de0260b3SHong Zhang 21483445d95SHong Zhang proc = 0; 215de0260b3SHong Zhang for (i=0; i<pon; i++){ 216de0260b3SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 217de0260b3SHong Zhang len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 218de0260b3SHong Zhang len_s[proc] += coi[i+1] - coi[i]; 21983445d95SHong Zhang } 220de0260b3SHong Zhang 221de0260b3SHong Zhang len = 0; /* max length of buf_si[] */ 222de0260b3SHong Zhang owners_co[0] = 0; 22342c57489SHong Zhang for (proc=0; proc<size; proc++){ 224de0260b3SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 22583445d95SHong Zhang if (len_si[proc]){ 22642c57489SHong Zhang merge->nsend++; 22783445d95SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); 22842c57489SHong Zhang len += len_si[proc]; 22942c57489SHong Zhang } 23042c57489SHong Zhang } 23142c57489SHong Zhang 232ded4f561SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 23342c57489SHong Zhang ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 23442c57489SHong Zhang ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 23542c57489SHong Zhang 236ded4f561SHong Zhang /* post the Irecv and Isend of coj */ 237ded4f561SHong Zhang ierr = PetscObjectGetNewTag((PetscObject)merge->rowmap,&tag);CHKERRQ(ierr); 238ded4f561SHong Zhang ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 23942c57489SHong Zhang 240ded4f561SHong Zhang ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr); 24142c57489SHong Zhang 24242c57489SHong Zhang for (proc=0, k=0; proc<size; proc++){ 24342c57489SHong Zhang if (!len_s[proc]) continue; 244de0260b3SHong Zhang i = owners_co[proc]; 245ded4f561SHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr); 24642c57489SHong Zhang k++; 24742c57489SHong Zhang } 24842c57489SHong Zhang 249ded4f561SHong Zhang /* receives and sends of coj are complete */ 250ded4f561SHong Zhang ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr); 251ded4f561SHong Zhang i = merge->nrecv; 252ded4f561SHong Zhang while (i--) { 253ded4f561SHong Zhang ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 254ded4f561SHong Zhang } 255ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 256ded4f561SHong Zhang if (merge->nsend){ 257ded4f561SHong Zhang ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr); 25882412ba7SHong Zhang } 25982412ba7SHong Zhang 260ded4f561SHong Zhang /* send and recv coi */ 261ded4f561SHong Zhang /*-------------------*/ 262ded4f561SHong Zhang ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 26342c57489SHong Zhang 26442c57489SHong Zhang ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr); 26542c57489SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 26642c57489SHong Zhang for (proc=0,k=0; proc<size; proc++){ 26742c57489SHong Zhang if (!len_s[proc]) continue; 26842c57489SHong Zhang /* form outgoing message for i-structure: 26942c57489SHong Zhang buf_si[0]: nrows to be sent 27042c57489SHong Zhang [1:nrows]: row index (global) 27142c57489SHong Zhang [nrows+1:2*nrows+1]: i-structure index 27242c57489SHong Zhang */ 27342c57489SHong Zhang /*-------------------------------------------*/ 27442c57489SHong Zhang nrows = len_si[proc]/2 - 1; 27542c57489SHong Zhang buf_si_i = buf_si + nrows+1; 27642c57489SHong Zhang buf_si[0] = nrows; 27742c57489SHong Zhang buf_si_i[0] = 0; 27842c57489SHong Zhang nrows = 0; 279de0260b3SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++){ 280ded4f561SHong Zhang nzi = coi[i+1] - coi[i]; 281ded4f561SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 282de0260b3SHong Zhang buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */ 28342c57489SHong Zhang nrows++; 28442c57489SHong Zhang } 285ded4f561SHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr); 28642c57489SHong Zhang k++; 28742c57489SHong Zhang buf_si += len_si[proc]; 28842c57489SHong Zhang } 289ded4f561SHong Zhang i = merge->nrecv; 290ded4f561SHong Zhang while (i--) { 291ded4f561SHong Zhang ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 292ded4f561SHong Zhang } 293ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 294ded4f561SHong Zhang if (merge->nsend){ 295ded4f561SHong Zhang ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr); 296ded4f561SHong Zhang } 29742c57489SHong Zhang 29842c57489SHong Zhang ierr = PetscLogInfo((PetscObject)A,"MatMerge_SeqsToMPI: nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr); 29942c57489SHong Zhang for (i=0; i<merge->nrecv; i++){ 30042c57489SHong 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); 30142c57489SHong Zhang } 30242c57489SHong Zhang 30342c57489SHong Zhang ierr = PetscFree(len_si);CHKERRQ(ierr); 30442c57489SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 305ded4f561SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 306ded4f561SHong Zhang ierr = PetscFree(sstatus);CHKERRQ(ierr); 30742c57489SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 30842c57489SHong Zhang 309de0260b3SHong Zhang /* compute the local portion of C (mpi mat) */ 310de0260b3SHong Zhang /*------------------------------------------*/ 311ded4f561SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 312ded4f561SHong Zhang 31342c57489SHong Zhang /* allocate bi array and free space for accumulating nonzero column info */ 31442c57489SHong Zhang ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 31542c57489SHong Zhang bi[0] = 0; 31642c57489SHong Zhang 317ded4f561SHong Zhang /* set initial free space to be 3*pn*max( nnz(AP) per row) */ 318ded4f561SHong Zhang nnz = 3*pn*ap->abnz_max + 1; 319ded4f561SHong Zhang ierr = GetMoreSpace(nnz,&free_space); 32042c57489SHong Zhang current_space = free_space; 32142c57489SHong Zhang 32242c57489SHong Zhang ierr = PetscMalloc((3*merge->nrecv+1)*sizeof(PetscInt**),&buf_ri_k);CHKERRQ(ierr); 32342c57489SHong Zhang nextrow = buf_ri_k + merge->nrecv; 32442c57489SHong Zhang nextci = nextrow + merge->nrecv; 32542c57489SHong Zhang for (k=0; k<merge->nrecv; k++){ 32642c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 32742c57489SHong Zhang nrows = *buf_ri_k[k]; 32842c57489SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 32942c57489SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 33042c57489SHong Zhang } 33142c57489SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 33242c57489SHong Zhang for (i=0; i<pn; i++) { 333ded4f561SHong Zhang /* add pdt[i,:]*AP into lnk */ 334ded4f561SHong Zhang nnz = 0; 335ded4f561SHong Zhang pnz = pdti[i+1] - pdti[i]; 336ded4f561SHong Zhang j = pnz; 337ded4f561SHong Zhang ptJ = pdtj + pdti[i+1]; 338ded4f561SHong Zhang while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 339ded4f561SHong Zhang j--; ptJ--; 340ded4f561SHong Zhang row = *ptJ; /* row of AP == col of Pt */ 341ded4f561SHong Zhang apnz = api[row+1] - api[row]; 342ded4f561SHong Zhang Jptr = apj + api[row]; 343ded4f561SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 344ded4f561SHong Zhang ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 345ded4f561SHong Zhang nnz += nlnk; 346ded4f561SHong Zhang } 34742c57489SHong Zhang /* add received col data into lnk */ 34842c57489SHong Zhang for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 34942c57489SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 350ded4f561SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 351ded4f561SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 352ded4f561SHong Zhang ierr = PetscLLAdd(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 353ded4f561SHong Zhang nnz += nlnk; 35442c57489SHong Zhang nextrow[k]++; nextci[k]++; 35542c57489SHong Zhang } 35642c57489SHong Zhang } 35742c57489SHong Zhang 35842c57489SHong Zhang /* if free space is not available, make more free space */ 359ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 36042c57489SHong Zhang ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 36142c57489SHong Zhang } 36242c57489SHong Zhang /* copy data into free space, then initialize lnk */ 363ded4f561SHong Zhang ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 364ded4f561SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 365ded4f561SHong Zhang current_space->array += nnz; 366ded4f561SHong Zhang current_space->local_used += nnz; 367ded4f561SHong Zhang current_space->local_remaining -= nnz; 368ded4f561SHong Zhang bi[i+1] = bi[i] + nnz; 36942c57489SHong Zhang } 370ded4f561SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 37142c57489SHong Zhang ierr = PetscFree(buf_ri_k);CHKERRQ(ierr); 37242c57489SHong Zhang 37342c57489SHong Zhang ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 37442c57489SHong Zhang ierr = MakeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 37542c57489SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 37642c57489SHong Zhang 37742c57489SHong Zhang /* create symbolic parallel matrix B_mpi */ 37842c57489SHong Zhang /*---------------------------------------*/ 37942c57489SHong Zhang ierr = MatCreate(comm,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE,&B_mpi);CHKERRQ(ierr); 38042c57489SHong Zhang ierr = MatSetType(B_mpi,MATMPIAIJ);CHKERRQ(ierr); 38142c57489SHong Zhang ierr = MatMPIAIJSetPreallocation(B_mpi,0,dnz,0,onz);CHKERRQ(ierr); 38242c57489SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 38342c57489SHong Zhang 38482412ba7SHong Zhang /* B_mpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 38542c57489SHong Zhang B_mpi->assembled = PETSC_FALSE; 38642c57489SHong Zhang B_mpi->ops->destroy = MatDestroy_MPIAIJ_SeqsToMPI; 38742c57489SHong Zhang merge->bi = bi; 38842c57489SHong Zhang merge->bj = bj; 389de0260b3SHong Zhang merge->coi = coi; 390de0260b3SHong Zhang merge->coj = coj; 39142c57489SHong Zhang merge->buf_ri = buf_ri; 39242c57489SHong Zhang merge->buf_rj = buf_rj; 393de0260b3SHong Zhang merge->owners_co = owners_co; 39442c57489SHong Zhang 39542c57489SHong Zhang /* attach the supporting struct to B_mpi for reuse */ 39642c57489SHong Zhang ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 39742c57489SHong Zhang ierr = PetscObjectContainerSetPointer(container,merge);CHKERRQ(ierr); 39842c57489SHong Zhang ierr = PetscObjectCompose((PetscObject)B_mpi,"MatMergeSeqsToMPI",(PetscObject)container);CHKERRQ(ierr); 39942c57489SHong Zhang *C = B_mpi; 40042c57489SHong Zhang PetscFunctionReturn(0); 40142c57489SHong Zhang } 40242c57489SHong Zhang 40342c57489SHong Zhang #undef __FUNCT__ 40442c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ" 40542c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 40642c57489SHong Zhang { 40742c57489SHong Zhang PetscErrorCode ierr; 40842c57489SHong Zhang Mat_Merge_SeqsToMPI *merge; 409de0260b3SHong Zhang Mat_MatMatMultMPI *ap; 41042c57489SHong Zhang PetscObjectContainer cont_merge,cont_ptap; 411de0260b3SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 41242c57489SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 413de0260b3SHong Zhang Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 41442c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 41582412ba7SHong Zhang PetscInt *adi=ad->i,*aoi=ao->i,*adj=ad->j,*aoj=ao->j,*apJ,nextp,flops=0; 41682412ba7SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj; 41782412ba7SHong Zhang PetscInt i,j,k,anz,pnz,apnz,nextap,row,*cj; 41882412ba7SHong Zhang MatScalar *ada=ad->a,*aoa=ao->a,*apa,*pa,*ca,*pa_loc,*pa_oth; 41982412ba7SHong Zhang PetscInt am=A->m,cm=C->m,pon=(p->B)->n; 42042c57489SHong Zhang MPI_Comm comm=C->comm; 42142c57489SHong Zhang PetscMPIInt size,rank,taga,*len_s; 42282412ba7SHong Zhang PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 42342c57489SHong Zhang PetscInt **buf_ri,**buf_rj; 42450f5bf86SHong Zhang PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 42542c57489SHong Zhang MPI_Request *s_waits,*r_waits; 42642c57489SHong Zhang MPI_Status *status; 42782412ba7SHong Zhang MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 42882412ba7SHong Zhang PetscInt *api,*apj,*coi,*coj; 42982412ba7SHong Zhang PetscInt *poJ=po->j,*pdJ=pd->j,pcstart=p->cstart,pcend=p->cend; 43042c57489SHong Zhang 43142c57489SHong Zhang PetscFunctionBegin; 43242c57489SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 43342c57489SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 43442c57489SHong Zhang 43542c57489SHong Zhang ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&cont_merge);CHKERRQ(ierr); 43642c57489SHong Zhang if (cont_merge) { 43742c57489SHong Zhang ierr = PetscObjectContainerGetPointer(cont_merge,(void **)&merge);CHKERRQ(ierr); 43842c57489SHong Zhang } else { 43942c57489SHong Zhang SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses an object container"); 44042c57489SHong Zhang } 441*6c6e00e0SHong Zhang 442f33d1a9aSHong Zhang ierr = PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&cont_ptap);CHKERRQ(ierr); 44342c57489SHong Zhang if (cont_ptap) { 444de0260b3SHong Zhang ierr = PetscObjectContainerGetPointer(cont_ptap,(void **)&ap);CHKERRQ(ierr); 445*6c6e00e0SHong Zhang if (ap->reuse == MAT_INITIAL_MATRIX){ 446*6c6e00e0SHong Zhang ap->reuse = MAT_REUSE_MATRIX; 447*6c6e00e0SHong Zhang } else { /* update numerical values of P_oth and P_loc */ 448*6c6e00e0SHong Zhang ierr = MatGetBrowsOfAoCols(A,P,MAT_REUSE_MATRIX,&ap->startsj,&ap->bufa,&ap->B_oth);CHKERRQ(ierr); 449*6c6e00e0SHong Zhang ierr = MatGetLocalMat(P,MAT_REUSE_MATRIX,&ap->B_loc);CHKERRQ(ierr); 450*6c6e00e0SHong Zhang } 45142c57489SHong Zhang } else { 45242c57489SHong Zhang SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container"); 45342c57489SHong Zhang } 45442c57489SHong Zhang 45542c57489SHong Zhang /* get data from symbolic products */ 456de0260b3SHong Zhang p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data; 457de0260b3SHong Zhang p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data; 45882412ba7SHong Zhang pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a,pA=pa_loc; 45942c57489SHong Zhang pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 46042c57489SHong Zhang 461de0260b3SHong Zhang coi = merge->coi; coj = merge->coj; 462de0260b3SHong Zhang ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr); 463de0260b3SHong Zhang ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr); 464de0260b3SHong Zhang 465de0260b3SHong Zhang bi = merge->bi; bj = merge->bj; 466de0260b3SHong Zhang ierr = PetscMapGetGlobalRange(merge->rowmap,&owners);CHKERRQ(ierr); 467de0260b3SHong Zhang ierr = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr); 468de0260b3SHong Zhang ierr = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr); 46942c57489SHong Zhang 470f4a743e1SHong Zhang /* get data from symbolic A*P */ 471de0260b3SHong Zhang ierr = PetscMalloc((ap->abnz_max+1)*sizeof(MatScalar),&apa);CHKERRQ(ierr); 472de0260b3SHong Zhang ierr = PetscMemzero(apa,ap->abnz_max*sizeof(MatScalar));CHKERRQ(ierr); 47342c57489SHong Zhang 474f4a743e1SHong Zhang /* compute numeric C_seq=P_loc^T*A_loc*P */ 47582412ba7SHong Zhang api = ap->abi; apj = ap->abj; 47642c57489SHong Zhang for (i=0;i<am;i++) { 477f4a743e1SHong Zhang /* form i-th sparse row of A*P */ 47882412ba7SHong Zhang apnz = api[i+1] - api[i]; 47982412ba7SHong Zhang apJ = apj + api[i]; 48042c57489SHong Zhang /* diagonal portion of A */ 48182412ba7SHong Zhang anz = adi[i+1] - adi[i]; 48282412ba7SHong Zhang for (j=0;j<anz;j++) { 48382412ba7SHong Zhang row = *adj++; 48482412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 48582412ba7SHong Zhang pj = pj_loc + pi_loc[row]; 48682412ba7SHong Zhang pa = pa_loc + pi_loc[row]; 487f4a743e1SHong Zhang nextp = 0; 48882412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 48982412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 49082412ba7SHong Zhang apa[k] += (*ada)*pa[nextp++]; 49142c57489SHong Zhang } 49242c57489SHong Zhang } 49382412ba7SHong Zhang flops += 2*pnz; 49442c57489SHong Zhang ada++; 49542c57489SHong Zhang } 49642c57489SHong Zhang /* off-diagonal portion of A */ 49782412ba7SHong Zhang anz = aoi[i+1] - aoi[i]; 49882412ba7SHong Zhang for (j=0; j<anz; j++) { 49982412ba7SHong Zhang row = *aoj++; 50082412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 50182412ba7SHong Zhang pj = pj_oth + pi_oth[row]; 50282412ba7SHong Zhang pa = pa_oth + pi_oth[row]; 503f4a743e1SHong Zhang nextp = 0; 50482412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 50582412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 50682412ba7SHong Zhang apa[k] += (*aoa)*pa[nextp++]; 50742c57489SHong Zhang } 50842c57489SHong Zhang } 50982412ba7SHong Zhang flops += 2*pnz; 51042c57489SHong Zhang aoa++; 51142c57489SHong Zhang } 51242c57489SHong Zhang 51382412ba7SHong Zhang /* Compute P_loc[i,:]^T*AP[i,:] using outer product */ 51482412ba7SHong Zhang pnz = pi_loc[i+1] - pi_loc[i]; 51582412ba7SHong Zhang for (j=0; j<pnz; j++) { 51642c57489SHong Zhang nextap = 0; 51782412ba7SHong Zhang row = *pJ++; /* global index */ 51882412ba7SHong Zhang if (row < pcstart || row >=pcend) { /* put the value into Co */ 51982412ba7SHong Zhang cj = coj + coi[*poJ]; 52082412ba7SHong Zhang ca = coa + coi[*poJ++]; 52182412ba7SHong Zhang } else { /* put the value into Cd */ 52282412ba7SHong Zhang cj = bj + bi[*pdJ]; 52382412ba7SHong Zhang ca = ba + bi[*pdJ++]; 52442c57489SHong Zhang } 52582412ba7SHong Zhang for (k=0; nextap<apnz; k++) { 52682412ba7SHong Zhang if (cj[k]==apJ[nextap]) ca[k] += (*pA)*apa[nextap++]; 52742c57489SHong Zhang } 52882412ba7SHong Zhang flops += 2*apnz; 52982412ba7SHong Zhang pA++; 530de0260b3SHong Zhang } 531de0260b3SHong Zhang 53242c57489SHong Zhang /* zero the current row info for A*P */ 53382412ba7SHong Zhang ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr); 53442c57489SHong Zhang } 53542c57489SHong Zhang ierr = PetscFree(apa);CHKERRQ(ierr); 53642c57489SHong Zhang 537de0260b3SHong Zhang /* send and recv matrix values */ 538de0260b3SHong Zhang /*-----------------------------*/ 53942c57489SHong Zhang buf_ri = merge->buf_ri; 54042c57489SHong Zhang buf_rj = merge->buf_rj; 54142c57489SHong Zhang len_s = merge->len_s; 54242c57489SHong Zhang ierr = PetscObjectGetNewTag((PetscObject)merge->rowmap,&taga);CHKERRQ(ierr); 54342c57489SHong Zhang ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 54442c57489SHong Zhang 54542c57489SHong Zhang ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&s_waits);CHKERRQ(ierr); 54642c57489SHong Zhang for (proc=0,k=0; proc<size; proc++){ 54742c57489SHong Zhang if (!len_s[proc]) continue; 548de0260b3SHong Zhang i = merge->owners_co[proc]; 549de0260b3SHong Zhang ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 55042c57489SHong Zhang k++; 55142c57489SHong Zhang } 55282412ba7SHong Zhang ierr = PetscMalloc(size*sizeof(MPI_Status),&status);CHKERRQ(ierr); 55342c57489SHong Zhang ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr); 55442c57489SHong Zhang ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr); 55542c57489SHong Zhang ierr = PetscFree(status);CHKERRQ(ierr); 55642c57489SHong Zhang 55742c57489SHong Zhang ierr = PetscFree(s_waits);CHKERRQ(ierr); 55842c57489SHong Zhang ierr = PetscFree(r_waits);CHKERRQ(ierr); 55982412ba7SHong Zhang ierr = PetscFree(coa);CHKERRQ(ierr); 56042c57489SHong Zhang 56182412ba7SHong Zhang /* insert local and received values into C */ 56282412ba7SHong Zhang /*-----------------------------------------*/ 56342c57489SHong Zhang ierr = PetscMalloc((3*merge->nrecv+1)*sizeof(PetscInt**),&buf_ri_k);CHKERRQ(ierr); 56442c57489SHong Zhang nextrow = buf_ri_k + merge->nrecv; 56582412ba7SHong Zhang nextci = nextrow + merge->nrecv; 56642c57489SHong Zhang 56742c57489SHong Zhang for (k=0; k<merge->nrecv; k++){ 56842c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 56942c57489SHong Zhang nrows = *(buf_ri_k[k]); 57042c57489SHong Zhang nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 57182412ba7SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 57242c57489SHong Zhang } 57342c57489SHong Zhang 574de0260b3SHong Zhang for (i=0; i<cm; i++) { 57582412ba7SHong Zhang row = owners[rank] + i; /* global row index of C_seq */ 57642c57489SHong Zhang bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 577de0260b3SHong Zhang ba_i = ba + bi[i]; 57882412ba7SHong Zhang bnz = bi[i+1] - bi[i]; 57942c57489SHong Zhang /* add received vals into ba */ 58042c57489SHong Zhang for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 58142c57489SHong Zhang /* i-th row */ 58242c57489SHong Zhang if (i == *nextrow[k]) { 58382412ba7SHong Zhang cnz = *(nextci[k]+1) - *nextci[k]; 58482412ba7SHong Zhang cj = buf_rj[k] + *(nextci[k]); 58582412ba7SHong Zhang ca = abuf_r[k] + *(nextci[k]); 58682412ba7SHong Zhang nextcj = 0; 58782412ba7SHong Zhang for (j=0; nextcj<cnz; j++){ 58882412ba7SHong Zhang if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */ 58982412ba7SHong Zhang ba_i[j] += ca[nextcj++]; 59042c57489SHong Zhang } 59142c57489SHong Zhang } 59282412ba7SHong Zhang nextrow[k]++; nextci[k]++; 59342c57489SHong Zhang } 59442c57489SHong Zhang } 59582412ba7SHong Zhang ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 59682412ba7SHong Zhang flops += 2*cnz; 59742c57489SHong Zhang } 59882412ba7SHong Zhang ierr = MatSetOption(C,MAT_COLUMNS_SORTED);CHKERRQ(ierr); /* -cause delay? */ 59942c57489SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 60042c57489SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 60142c57489SHong Zhang 602de0260b3SHong Zhang ierr = PetscFree(ba);CHKERRQ(ierr); 60342c57489SHong Zhang ierr = PetscFree(abuf_r);CHKERRQ(ierr); 60442c57489SHong Zhang ierr = PetscFree(buf_ri_k);CHKERRQ(ierr); 60542c57489SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 60642c57489SHong Zhang PetscFunctionReturn(0); 60742c57489SHong Zhang } 608