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 EXTERN PetscErrorCode MatDestroy_MPIAIJ(Mat); 1342c57489SHong Zhang #undef __FUNCT__ 1442c57489SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP" 1542c57489SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A) 1642c57489SHong Zhang { 1742c57489SHong Zhang PetscErrorCode ierr; 18a6b2eed2SHong Zhang Mat_MatMatMultMPI *mult; 1942c57489SHong Zhang PetscObjectContainer container; 2042c57489SHong Zhang 2142c57489SHong Zhang PetscFunctionBegin; 2242c57489SHong Zhang ierr = PetscObjectQuery((PetscObject)A,"MatPtAPMPI",(PetscObject *)&container);CHKERRQ(ierr); 2342c57489SHong Zhang if (container) { 24a6b2eed2SHong Zhang ierr = PetscObjectContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 25a6b2eed2SHong Zhang ierr = MatDestroy(mult->B_loc);CHKERRQ(ierr); 26a6b2eed2SHong Zhang ierr = MatDestroy(mult->B_oth);CHKERRQ(ierr); 27a6b2eed2SHong Zhang if (mult->abi) ierr = PetscFree(mult->abi);CHKERRQ(ierr); 28a6b2eed2SHong Zhang if (mult->abj) ierr = PetscFree(mult->abj);CHKERRQ(ierr); 29a6b2eed2SHong Zhang if (mult->startsj) ierr = PetscFree(mult->startsj);CHKERRQ(ierr); 30a6b2eed2SHong Zhang if (mult->bufa) ierr = PetscFree(mult->bufa);CHKERRQ(ierr); 3142c57489SHong Zhang ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr); 3242c57489SHong Zhang ierr = PetscObjectCompose((PetscObject)A,"MatPtAPMPI",0);CHKERRQ(ierr); 3342c57489SHong Zhang } 34a6b2eed2SHong Zhang ierr = PetscFree(mult);CHKERRQ(ierr); 3542c57489SHong Zhang 3642c57489SHong Zhang ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr); 3742c57489SHong Zhang PetscFunctionReturn(0); 3842c57489SHong Zhang } 3942c57489SHong Zhang 4042c57489SHong Zhang #undef __FUNCT__ 4142c57489SHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ" 4242c57489SHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 4342c57489SHong Zhang { 4442c57489SHong Zhang PetscErrorCode ierr; 45a6b2eed2SHong Zhang Mat_MatMatMultMPI *mult; 4642c57489SHong Zhang PetscObjectContainer container; 4742c57489SHong Zhang 4842c57489SHong Zhang PetscFunctionBegin; 4942c57489SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 5042c57489SHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 51a6b2eed2SHong Zhang ierr = PetscNew(Mat_MatMatMultMPI,&mult);CHKERRQ(ierr); 52a6b2eed2SHong Zhang mult->B_loc=PETSC_NULL; mult->B_oth=PETSC_NULL; 53a6b2eed2SHong Zhang mult->abi=PETSC_NULL; mult->abj=PETSC_NULL; 54a6b2eed2SHong Zhang mult->abnz_max = 0; /* symbolic A*P is not done yet */ 5542c57489SHong Zhang 5642c57489SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 57a6b2eed2SHong Zhang ierr = MatGetBrowsOfAoCols(A,P,scall,&mult->startsj,&mult->bufa,&mult->B_oth);CHKERRQ(ierr); 5842c57489SHong Zhang 5942c57489SHong Zhang /* get P_loc by taking all local rows of P */ 60a6b2eed2SHong Zhang ierr = MatGetLocalMat(P,scall,&mult->B_loc);CHKERRQ(ierr); 6142c57489SHong Zhang 6242c57489SHong Zhang /* attach the supporting struct to P for reuse */ 6342c57489SHong Zhang P->ops->destroy = MatDestroy_MPIAIJ_PtAP; 6442c57489SHong Zhang ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 65a6b2eed2SHong Zhang ierr = PetscObjectContainerSetPointer(container,mult);CHKERRQ(ierr); 6642c57489SHong Zhang ierr = PetscObjectCompose((PetscObject)P,"MatPtAPMPI",(PetscObject)container);CHKERRQ(ierr); 6742c57489SHong Zhang 6842c57489SHong Zhang /* now, compute symbolic local P^T*A*P */ 6942c57489SHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr);/* numeric product is computed as well */ 7042c57489SHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 7142c57489SHong Zhang } else if (scall == MAT_REUSE_MATRIX){ 7242c57489SHong Zhang ierr = PetscObjectQuery((PetscObject)P,"MatPtAPMPI",(PetscObject *)&container);CHKERRQ(ierr); 7342c57489SHong Zhang if (container) { 74a6b2eed2SHong Zhang ierr = PetscObjectContainerGetPointer(container,(void **)&mult);CHKERRQ(ierr); 7542c57489SHong Zhang } else { 7642c57489SHong Zhang SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container"); 7742c57489SHong Zhang } 7842c57489SHong Zhang 7942c57489SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 80a6b2eed2SHong Zhang ierr = MatGetBrowsOfAoCols(A,P,scall,&mult->startsj,&mult->bufa,&mult->B_oth);CHKERRQ(ierr); 8142c57489SHong Zhang 8242c57489SHong Zhang /* get P_loc by taking all local rows of P */ 83a6b2eed2SHong Zhang ierr = MatGetLocalMat(P,scall,&mult->B_loc);CHKERRQ(ierr); 8442c57489SHong Zhang 8542c57489SHong Zhang } else { 8642c57489SHong Zhang SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",scall); 8742c57489SHong Zhang } 8842c57489SHong Zhang /* now, compute numeric local P^T*A*P */ 8942c57489SHong Zhang ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 9042c57489SHong Zhang ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr); 9142c57489SHong Zhang ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 9242c57489SHong Zhang 9342c57489SHong Zhang PetscFunctionReturn(0); 9442c57489SHong Zhang } 9542c57489SHong Zhang 9642c57489SHong Zhang #undef __FUNCT__ 9742c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ" 9842c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 9942c57489SHong Zhang { 10042c57489SHong Zhang PetscErrorCode ierr; 101de0260b3SHong Zhang Mat P_loc,P_oth,B_mpi; 102de0260b3SHong Zhang Mat_MatMatMultMPI *ap; 10342c57489SHong Zhang PetscObjectContainer container; 10442c57489SHong Zhang FreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 10583445d95SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 10642c57489SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 10742c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 108*ded4f561SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ; 10983445d95SHong Zhang PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz; 110*ded4f561SHong Zhang PetscInt nlnk,*lnk,*owners_co,*coi,*coj,i,j,k,pnz,row; 111a6b2eed2SHong Zhang PetscInt am=A->m,pN=P->N,pn=P->n; 11242c57489SHong Zhang PetscBT lnkbt; 11342c57489SHong Zhang MPI_Comm comm=A->comm; 114*ded4f561SHong Zhang PetscMPIInt size,rank,tag,*len_si,*len_s,*len_ri; 11542c57489SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 11642c57489SHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 117*ded4f561SHong Zhang PetscInt nzi,*bi,*bj; 11842c57489SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 119*ded4f561SHong Zhang MPI_Request *swaits,*rwaits; 120*ded4f561SHong Zhang MPI_Status *sstatus,rstatus; 12142c57489SHong Zhang Mat_Merge_SeqsToMPI *merge; 122*ded4f561SHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,pon; 12342c57489SHong Zhang 12442c57489SHong Zhang PetscFunctionBegin; 12542c57489SHong Zhang ierr = PetscObjectQuery((PetscObject)P,"MatPtAPMPI",(PetscObject *)&container);CHKERRQ(ierr); 12642c57489SHong Zhang if (container) { 127de0260b3SHong Zhang ierr = PetscObjectContainerGetPointer(container,(void **)&ap);CHKERRQ(ierr); 12842c57489SHong Zhang } else { 12942c57489SHong Zhang SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container"); 13042c57489SHong Zhang } 13183445d95SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 13283445d95SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 13383445d95SHong Zhang 13442c57489SHong Zhang /* get data from P_loc and P_oth */ 135a6b2eed2SHong Zhang P_loc=ap->B_loc; P_oth=ap->B_oth; 13642c57489SHong Zhang p_loc=(Mat_SeqAIJ*)P_loc->data; p_oth=(Mat_SeqAIJ*)P_oth->data; 13742c57489SHong Zhang pi_loc=p_loc->i; pj_loc=p_loc->j; pi_oth=p_oth->i; pj_oth=p_oth->j; 13842c57489SHong Zhang 13983445d95SHong Zhang /* first, compute symbolic AP = A_loc*P */ 140de0260b3SHong Zhang /*--------------------------------------*/ 14182412ba7SHong Zhang ierr = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr); 14282412ba7SHong Zhang ap->abi = api; 14382412ba7SHong Zhang api[0] = 0; 14483445d95SHong Zhang 14583445d95SHong Zhang /* create and initialize a linked list */ 14683445d95SHong Zhang nlnk = pN+1; 14783445d95SHong Zhang ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 14883445d95SHong Zhang 14982412ba7SHong Zhang /* Initial FreeSpace size is fill*nnz(A) */ 150de0260b3SHong Zhang ierr = GetMoreSpace((PetscInt)(fill*(adi[am]+aoi[am])),&free_space);CHKERRQ(ierr); 151f4a743e1SHong Zhang current_space = free_space; 152f4a743e1SHong Zhang 153f4a743e1SHong Zhang for (i=0;i<am;i++) { 15482412ba7SHong Zhang apnz = 0; 155f4a743e1SHong Zhang /* diagonal portion of A */ 156*ded4f561SHong Zhang nzi = adi[i+1] - adi[i]; 157*ded4f561SHong Zhang for (j=0; j<nzi; j++){ 15882412ba7SHong Zhang row = *adj++; 15982412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 160*ded4f561SHong Zhang Jptr = pj_loc + pi_loc[row]; 16183445d95SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 162*ded4f561SHong Zhang ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 16382412ba7SHong Zhang apnz += nlnk; 164f4a743e1SHong Zhang } 165f4a743e1SHong Zhang /* off-diagonal portion of A */ 166*ded4f561SHong Zhang nzi = aoi[i+1] - aoi[i]; 167*ded4f561SHong Zhang for (j=0; j<nzi; j++){ 16882412ba7SHong Zhang row = *aoj++; 16982412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 170*ded4f561SHong Zhang Jptr = pj_oth + pi_oth[row]; 171*ded4f561SHong Zhang ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 17282412ba7SHong Zhang apnz += nlnk; 173f4a743e1SHong Zhang } 174f4a743e1SHong Zhang 17582412ba7SHong Zhang api[i+1] = api[i] + apnz; 17682412ba7SHong Zhang if (ap->abnz_max < apnz) ap->abnz_max = apnz; 177f4a743e1SHong Zhang 17883445d95SHong Zhang /* if free space is not available, double the total space in the list */ 17982412ba7SHong Zhang if (current_space->local_remaining<apnz) { 180f4a743e1SHong Zhang ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 181f4a743e1SHong Zhang } 182f4a743e1SHong Zhang 183f4a743e1SHong Zhang /* Copy data into free space, then initialize lnk */ 18482412ba7SHong Zhang ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 18582412ba7SHong Zhang current_space->array += apnz; 18682412ba7SHong Zhang current_space->local_used += apnz; 18782412ba7SHong Zhang current_space->local_remaining -= apnz; 188f4a743e1SHong Zhang } 18982412ba7SHong Zhang /* Allocate space for apj, initialize apj, and */ 190f4a743e1SHong Zhang /* destroy list of free space and other temporary array(s) */ 19182412ba7SHong Zhang ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ap->abj);CHKERRQ(ierr); 19282412ba7SHong Zhang apj = ap->abj; 193de0260b3SHong Zhang ierr = MakeSpaceContiguous(&free_space,ap->abj);CHKERRQ(ierr); 194f4a743e1SHong Zhang 195*ded4f561SHong Zhang /* determine symbolic Co=(p->B)^T*AP - send to others */ 196*ded4f561SHong Zhang /*----------------------------------------------------*/ 197de0260b3SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 19842c57489SHong Zhang 199*ded4f561SHong Zhang /* then, compute symbolic Co = (p->B)^T*AP */ 200de0260b3SHong Zhang pon = (p->B)->n; /* total num of rows to be sent to other processors 20183445d95SHong Zhang >= (num of nonzero rows of C_seq) - pn */ 202de0260b3SHong Zhang ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr); 203de0260b3SHong Zhang coi[0] = 0; 20442c57489SHong Zhang 20582412ba7SHong Zhang /* set initial free space to be 3*pon*max( nnz(AP) per row) */ 20682412ba7SHong Zhang nnz = 3*pon*ap->abnz_max + 1; 207de0260b3SHong Zhang ierr = GetMoreSpace(nnz,&free_space); 20842c57489SHong Zhang current_space = free_space; 20942c57489SHong Zhang 210de0260b3SHong Zhang for (i=0; i<pon; i++) { 211*ded4f561SHong Zhang nnz = 0; 21282412ba7SHong Zhang pnz = poti[i+1] - poti[i]; 21382412ba7SHong Zhang j = pnz; 214de0260b3SHong Zhang ptJ = potj + poti[i+1]; 21583445d95SHong Zhang while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 21683445d95SHong Zhang j--; ptJ--; 21782412ba7SHong Zhang row = *ptJ; /* row of AP == col of Pot */ 21882412ba7SHong Zhang apnz = api[row+1] - api[row]; 219*ded4f561SHong Zhang Jptr = apj + api[row]; 22083445d95SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 221*ded4f561SHong Zhang ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 222*ded4f561SHong Zhang nnz += nlnk; 22342c57489SHong Zhang } 22442c57489SHong Zhang 22583445d95SHong Zhang /* If free space is not available, double the total space in the list */ 226*ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 22742c57489SHong Zhang ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 22842c57489SHong Zhang } 22942c57489SHong Zhang 23042c57489SHong Zhang /* Copy data into free space, and zero out denserows */ 231*ded4f561SHong Zhang ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 232*ded4f561SHong Zhang current_space->array += nnz; 233*ded4f561SHong Zhang current_space->local_used += nnz; 234*ded4f561SHong Zhang current_space->local_remaining -= nnz; 235*ded4f561SHong Zhang coi[i+1] = coi[i] + nnz; 23642c57489SHong Zhang } 237de0260b3SHong Zhang ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr); 238de0260b3SHong Zhang ierr = MakeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 239de0260b3SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 24042c57489SHong Zhang 241*ded4f561SHong Zhang /* send j-array (coj) of Co to other processors */ 242*ded4f561SHong Zhang /*----------------------------------------------*/ 24342c57489SHong Zhang /* determine row ownership */ 24483445d95SHong Zhang ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr); 24542c57489SHong Zhang ierr = PetscMapCreate(comm,&merge->rowmap);CHKERRQ(ierr); 24642c57489SHong Zhang ierr = PetscMapSetLocalSize(merge->rowmap,pn);CHKERRQ(ierr); 24742c57489SHong Zhang ierr = PetscMapSetType(merge->rowmap,MAP_MPI);CHKERRQ(ierr); 24842c57489SHong Zhang ierr = PetscMapGetGlobalRange(merge->rowmap,&owners);CHKERRQ(ierr); 24942c57489SHong Zhang 25042c57489SHong Zhang /* determine the number of messages to send, their lengths */ 25142c57489SHong Zhang ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr); 25283445d95SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 25342c57489SHong Zhang ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr); 25442c57489SHong Zhang len_s = merge->len_s; 25542c57489SHong Zhang merge->nsend = 0; 25683445d95SHong Zhang 257de0260b3SHong Zhang ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr); 258de0260b3SHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 259de0260b3SHong Zhang 26083445d95SHong Zhang proc = 0; 261de0260b3SHong Zhang for (i=0; i<pon; i++){ 262de0260b3SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 263de0260b3SHong Zhang len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 264de0260b3SHong Zhang len_s[proc] += coi[i+1] - coi[i]; 26583445d95SHong Zhang } 266de0260b3SHong Zhang 267de0260b3SHong Zhang len = 0; /* max length of buf_si[] */ 268de0260b3SHong Zhang owners_co[0] = 0; 26942c57489SHong Zhang for (proc=0; proc<size; proc++){ 270de0260b3SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 27183445d95SHong Zhang if (len_si[proc]){ 27242c57489SHong Zhang merge->nsend++; 27383445d95SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); 27442c57489SHong Zhang len += len_si[proc]; 27542c57489SHong Zhang } 27642c57489SHong Zhang } 27742c57489SHong Zhang 278*ded4f561SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 27942c57489SHong Zhang ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 28042c57489SHong Zhang ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 28142c57489SHong Zhang 282*ded4f561SHong Zhang /* post the Irecv and Isend of coj */ 283*ded4f561SHong Zhang ierr = PetscObjectGetNewTag((PetscObject)merge->rowmap,&tag);CHKERRQ(ierr); 284*ded4f561SHong Zhang ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 28542c57489SHong Zhang 286*ded4f561SHong Zhang ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr); 28742c57489SHong Zhang 28842c57489SHong Zhang for (proc=0, k=0; proc<size; proc++){ 28942c57489SHong Zhang if (!len_s[proc]) continue; 290de0260b3SHong Zhang i = owners_co[proc]; 291*ded4f561SHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr); 29242c57489SHong Zhang k++; 29342c57489SHong Zhang } 29442c57489SHong Zhang 295*ded4f561SHong Zhang /* receives and sends of coj are complete */ 296*ded4f561SHong Zhang ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr); 297*ded4f561SHong Zhang i = merge->nrecv; 298*ded4f561SHong Zhang while (i--) { 299*ded4f561SHong Zhang ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 300*ded4f561SHong Zhang } 301*ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 302*ded4f561SHong Zhang if (merge->nsend){ 303*ded4f561SHong Zhang ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr); 30482412ba7SHong Zhang } 30582412ba7SHong Zhang 306*ded4f561SHong Zhang /* send and recv coi */ 307*ded4f561SHong Zhang /*-------------------*/ 308*ded4f561SHong Zhang ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 30942c57489SHong Zhang 31042c57489SHong Zhang ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr); 31142c57489SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 31242c57489SHong Zhang for (proc=0,k=0; proc<size; proc++){ 31342c57489SHong Zhang if (!len_s[proc]) continue; 31442c57489SHong Zhang /* form outgoing message for i-structure: 31542c57489SHong Zhang buf_si[0]: nrows to be sent 31642c57489SHong Zhang [1:nrows]: row index (global) 31742c57489SHong Zhang [nrows+1:2*nrows+1]: i-structure index 31842c57489SHong Zhang */ 31942c57489SHong Zhang /*-------------------------------------------*/ 32042c57489SHong Zhang nrows = len_si[proc]/2 - 1; 32142c57489SHong Zhang buf_si_i = buf_si + nrows+1; 32242c57489SHong Zhang buf_si[0] = nrows; 32342c57489SHong Zhang buf_si_i[0] = 0; 32442c57489SHong Zhang nrows = 0; 325de0260b3SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++){ 326*ded4f561SHong Zhang nzi = coi[i+1] - coi[i]; 327*ded4f561SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 328de0260b3SHong Zhang buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */ 32942c57489SHong Zhang nrows++; 33042c57489SHong Zhang } 331*ded4f561SHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr); 33242c57489SHong Zhang k++; 33342c57489SHong Zhang buf_si += len_si[proc]; 33442c57489SHong Zhang } 335*ded4f561SHong Zhang i = merge->nrecv; 336*ded4f561SHong Zhang while (i--) { 337*ded4f561SHong Zhang ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 338*ded4f561SHong Zhang } 339*ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 340*ded4f561SHong Zhang if (merge->nsend){ 341*ded4f561SHong Zhang ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr); 342*ded4f561SHong Zhang } 34342c57489SHong Zhang 34442c57489SHong Zhang ierr = PetscLogInfo((PetscObject)A,"MatMerge_SeqsToMPI: nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr); 34542c57489SHong Zhang for (i=0; i<merge->nrecv; i++){ 34642c57489SHong 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); 34742c57489SHong Zhang } 34842c57489SHong Zhang 34942c57489SHong Zhang ierr = PetscFree(len_si);CHKERRQ(ierr); 35042c57489SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 351*ded4f561SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 352*ded4f561SHong Zhang ierr = PetscFree(sstatus);CHKERRQ(ierr); 35342c57489SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 35442c57489SHong Zhang 355de0260b3SHong Zhang /* compute the local portion of C (mpi mat) */ 356de0260b3SHong Zhang /*------------------------------------------*/ 357*ded4f561SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 358*ded4f561SHong Zhang 35942c57489SHong Zhang /* allocate bi array and free space for accumulating nonzero column info */ 36042c57489SHong Zhang ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 36142c57489SHong Zhang bi[0] = 0; 36242c57489SHong Zhang 363*ded4f561SHong Zhang /* set initial free space to be 3*pn*max( nnz(AP) per row) */ 364*ded4f561SHong Zhang nnz = 3*pn*ap->abnz_max + 1; 365*ded4f561SHong Zhang ierr = GetMoreSpace(nnz,&free_space); 36642c57489SHong Zhang current_space = free_space; 36742c57489SHong Zhang 36842c57489SHong Zhang ierr = PetscMalloc((3*merge->nrecv+1)*sizeof(PetscInt**),&buf_ri_k);CHKERRQ(ierr); 36942c57489SHong Zhang nextrow = buf_ri_k + merge->nrecv; 37042c57489SHong Zhang nextci = nextrow + merge->nrecv; 37142c57489SHong Zhang for (k=0; k<merge->nrecv; k++){ 37242c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 37342c57489SHong Zhang nrows = *buf_ri_k[k]; 37442c57489SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 37542c57489SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 37642c57489SHong Zhang } 37742c57489SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 37842c57489SHong Zhang for (i=0; i<pn; i++) { 379*ded4f561SHong Zhang /* add pdt[i,:]*AP into lnk */ 380*ded4f561SHong Zhang nnz = 0; 381*ded4f561SHong Zhang pnz = pdti[i+1] - pdti[i]; 382*ded4f561SHong Zhang j = pnz; 383*ded4f561SHong Zhang ptJ = pdtj + pdti[i+1]; 384*ded4f561SHong Zhang while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 385*ded4f561SHong Zhang j--; ptJ--; 386*ded4f561SHong Zhang row = *ptJ; /* row of AP == col of Pt */ 387*ded4f561SHong Zhang apnz = api[row+1] - api[row]; 388*ded4f561SHong Zhang Jptr = apj + api[row]; 389*ded4f561SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 390*ded4f561SHong Zhang ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 391*ded4f561SHong Zhang nnz += nlnk; 392*ded4f561SHong Zhang } 39342c57489SHong Zhang /* add received col data into lnk */ 39442c57489SHong Zhang for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 39542c57489SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 396*ded4f561SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 397*ded4f561SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 398*ded4f561SHong Zhang ierr = PetscLLAdd(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 399*ded4f561SHong Zhang nnz += nlnk; 40042c57489SHong Zhang nextrow[k]++; nextci[k]++; 40142c57489SHong Zhang } 40242c57489SHong Zhang } 40342c57489SHong Zhang 40442c57489SHong Zhang /* if free space is not available, make more free space */ 405*ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 40642c57489SHong Zhang ierr = GetMoreSpace(current_space->total_array_size,¤t_space);CHKERRQ(ierr); 40742c57489SHong Zhang } 40842c57489SHong Zhang /* copy data into free space, then initialize lnk */ 409*ded4f561SHong Zhang ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 410*ded4f561SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 411*ded4f561SHong Zhang current_space->array += nnz; 412*ded4f561SHong Zhang current_space->local_used += nnz; 413*ded4f561SHong Zhang current_space->local_remaining -= nnz; 414*ded4f561SHong Zhang bi[i+1] = bi[i] + nnz; 41542c57489SHong Zhang } 416*ded4f561SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 41742c57489SHong Zhang ierr = PetscFree(buf_ri_k);CHKERRQ(ierr); 41842c57489SHong Zhang 41942c57489SHong Zhang ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 42042c57489SHong Zhang ierr = MakeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 42142c57489SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 42242c57489SHong Zhang 42342c57489SHong Zhang /* create symbolic parallel matrix B_mpi */ 42442c57489SHong Zhang /*---------------------------------------*/ 42542c57489SHong Zhang ierr = MatCreate(comm,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE,&B_mpi);CHKERRQ(ierr); 42642c57489SHong Zhang ierr = MatSetType(B_mpi,MATMPIAIJ);CHKERRQ(ierr); 42742c57489SHong Zhang ierr = MatMPIAIJSetPreallocation(B_mpi,0,dnz,0,onz);CHKERRQ(ierr); 42842c57489SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 42942c57489SHong Zhang 43082412ba7SHong Zhang /* B_mpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 43142c57489SHong Zhang B_mpi->assembled = PETSC_FALSE; 43242c57489SHong Zhang B_mpi->ops->destroy = MatDestroy_MPIAIJ_SeqsToMPI; 43342c57489SHong Zhang merge->bi = bi; 43442c57489SHong Zhang merge->bj = bj; 435de0260b3SHong Zhang merge->coi = coi; 436de0260b3SHong Zhang merge->coj = coj; 43742c57489SHong Zhang merge->buf_ri = buf_ri; 43842c57489SHong Zhang merge->buf_rj = buf_rj; 439de0260b3SHong Zhang merge->owners_co = owners_co; 44042c57489SHong Zhang 44142c57489SHong Zhang /* attach the supporting struct to B_mpi for reuse */ 44242c57489SHong Zhang ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 44342c57489SHong Zhang ierr = PetscObjectContainerSetPointer(container,merge);CHKERRQ(ierr); 44442c57489SHong Zhang ierr = PetscObjectCompose((PetscObject)B_mpi,"MatMergeSeqsToMPI",(PetscObject)container);CHKERRQ(ierr); 44542c57489SHong Zhang *C = B_mpi; 44642c57489SHong Zhang 44742c57489SHong Zhang PetscFunctionReturn(0); 44842c57489SHong Zhang } 44942c57489SHong Zhang 45042c57489SHong Zhang #undef __FUNCT__ 45142c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ" 45242c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 45342c57489SHong Zhang { 45442c57489SHong Zhang PetscErrorCode ierr; 45542c57489SHong Zhang Mat_Merge_SeqsToMPI *merge; 456de0260b3SHong Zhang Mat_MatMatMultMPI *ap; 45742c57489SHong Zhang PetscObjectContainer cont_merge,cont_ptap; 458de0260b3SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 45942c57489SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 460de0260b3SHong Zhang Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 46142c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 46282412ba7SHong Zhang PetscInt *adi=ad->i,*aoi=ao->i,*adj=ad->j,*aoj=ao->j,*apJ,nextp,flops=0; 46382412ba7SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj; 46482412ba7SHong Zhang PetscInt i,j,k,anz,pnz,apnz,nextap,row,*cj; 46582412ba7SHong Zhang MatScalar *ada=ad->a,*aoa=ao->a,*apa,*pa,*ca,*pa_loc,*pa_oth; 46682412ba7SHong Zhang PetscInt am=A->m,cm=C->m,pon=(p->B)->n; 46742c57489SHong Zhang MPI_Comm comm=C->comm; 46842c57489SHong Zhang PetscMPIInt size,rank,taga,*len_s; 46982412ba7SHong Zhang PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 47042c57489SHong Zhang PetscInt **buf_ri,**buf_rj; 47150f5bf86SHong Zhang PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 47242c57489SHong Zhang MPI_Request *s_waits,*r_waits; 47342c57489SHong Zhang MPI_Status *status; 47482412ba7SHong Zhang MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 47582412ba7SHong Zhang PetscInt *api,*apj,*coi,*coj; 47682412ba7SHong Zhang PetscInt *poJ=po->j,*pdJ=pd->j,pcstart=p->cstart,pcend=p->cend; 47742c57489SHong Zhang 47842c57489SHong Zhang PetscFunctionBegin; 47942c57489SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 48042c57489SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 48142c57489SHong Zhang 48242c57489SHong Zhang ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&cont_merge);CHKERRQ(ierr); 48342c57489SHong Zhang if (cont_merge) { 48442c57489SHong Zhang ierr = PetscObjectContainerGetPointer(cont_merge,(void **)&merge);CHKERRQ(ierr); 48542c57489SHong Zhang } else { 48642c57489SHong Zhang SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses an object container"); 48742c57489SHong Zhang } 48842c57489SHong Zhang ierr = PetscObjectQuery((PetscObject)P,"MatPtAPMPI",(PetscObject *)&cont_ptap);CHKERRQ(ierr); 48942c57489SHong Zhang if (cont_ptap) { 490de0260b3SHong Zhang ierr = PetscObjectContainerGetPointer(cont_ptap,(void **)&ap);CHKERRQ(ierr); 49142c57489SHong Zhang } else { 49242c57489SHong Zhang SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container"); 49342c57489SHong Zhang } 49442c57489SHong Zhang 49542c57489SHong Zhang /* get data from symbolic products */ 496de0260b3SHong Zhang p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data; 497de0260b3SHong Zhang p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data; 49882412ba7SHong Zhang pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a,pA=pa_loc; 49942c57489SHong Zhang pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 50042c57489SHong Zhang 501de0260b3SHong Zhang coi = merge->coi; coj = merge->coj; 502de0260b3SHong Zhang ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr); 503de0260b3SHong Zhang ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr); 504de0260b3SHong Zhang 505de0260b3SHong Zhang bi = merge->bi; bj = merge->bj; 506de0260b3SHong Zhang ierr = PetscMapGetGlobalRange(merge->rowmap,&owners);CHKERRQ(ierr); 507de0260b3SHong Zhang ierr = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr); 508de0260b3SHong Zhang ierr = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr); 50942c57489SHong Zhang 510f4a743e1SHong Zhang /* get data from symbolic A*P */ 511de0260b3SHong Zhang ierr = PetscMalloc((ap->abnz_max+1)*sizeof(MatScalar),&apa);CHKERRQ(ierr); 512de0260b3SHong Zhang ierr = PetscMemzero(apa,ap->abnz_max*sizeof(MatScalar));CHKERRQ(ierr); 51342c57489SHong Zhang 514f4a743e1SHong Zhang /* compute numeric C_seq=P_loc^T*A_loc*P */ 51582412ba7SHong Zhang api = ap->abi; apj = ap->abj; 51642c57489SHong Zhang for (i=0;i<am;i++) { 517f4a743e1SHong Zhang /* form i-th sparse row of A*P */ 51882412ba7SHong Zhang apnz = api[i+1] - api[i]; 51982412ba7SHong Zhang apJ = apj + api[i]; 52042c57489SHong Zhang /* diagonal portion of A */ 52182412ba7SHong Zhang anz = adi[i+1] - adi[i]; 52282412ba7SHong Zhang for (j=0;j<anz;j++) { 52382412ba7SHong Zhang row = *adj++; 52482412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 52582412ba7SHong Zhang pj = pj_loc + pi_loc[row]; 52682412ba7SHong Zhang pa = pa_loc + pi_loc[row]; 527f4a743e1SHong Zhang nextp = 0; 52882412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 52982412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 53082412ba7SHong Zhang apa[k] += (*ada)*pa[nextp++]; 53142c57489SHong Zhang } 53242c57489SHong Zhang } 53382412ba7SHong Zhang flops += 2*pnz; 53442c57489SHong Zhang ada++; 53542c57489SHong Zhang } 53642c57489SHong Zhang /* off-diagonal portion of A */ 53782412ba7SHong Zhang anz = aoi[i+1] - aoi[i]; 53882412ba7SHong Zhang for (j=0; j<anz; j++) { 53982412ba7SHong Zhang row = *aoj++; 54082412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 54182412ba7SHong Zhang pj = pj_oth + pi_oth[row]; 54282412ba7SHong Zhang pa = pa_oth + pi_oth[row]; 543f4a743e1SHong Zhang nextp = 0; 54482412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 54582412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 54682412ba7SHong Zhang apa[k] += (*aoa)*pa[nextp++]; 54742c57489SHong Zhang } 54842c57489SHong Zhang } 54982412ba7SHong Zhang flops += 2*pnz; 55042c57489SHong Zhang aoa++; 55142c57489SHong Zhang } 55242c57489SHong Zhang 55382412ba7SHong Zhang /* Compute P_loc[i,:]^T*AP[i,:] using outer product */ 55482412ba7SHong Zhang pnz = pi_loc[i+1] - pi_loc[i]; 55582412ba7SHong Zhang for (j=0; j<pnz; j++) { 55642c57489SHong Zhang nextap = 0; 55782412ba7SHong Zhang row = *pJ++; /* global index */ 55882412ba7SHong Zhang if (row < pcstart || row >=pcend) { /* put the value into Co */ 55982412ba7SHong Zhang cj = coj + coi[*poJ]; 56082412ba7SHong Zhang ca = coa + coi[*poJ++]; 56182412ba7SHong Zhang } else { /* put the value into Cd */ 56282412ba7SHong Zhang cj = bj + bi[*pdJ]; 56382412ba7SHong Zhang ca = ba + bi[*pdJ++]; 56442c57489SHong Zhang } 56582412ba7SHong Zhang for (k=0; nextap<apnz; k++) { 56682412ba7SHong Zhang if (cj[k]==apJ[nextap]) ca[k] += (*pA)*apa[nextap++]; 56742c57489SHong Zhang } 56882412ba7SHong Zhang flops += 2*apnz; 56982412ba7SHong Zhang pA++; 570de0260b3SHong Zhang } 571de0260b3SHong Zhang 57242c57489SHong Zhang /* zero the current row info for A*P */ 57382412ba7SHong Zhang ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr); 57442c57489SHong Zhang } 57542c57489SHong Zhang ierr = PetscFree(apa);CHKERRQ(ierr); 57642c57489SHong Zhang 577de0260b3SHong Zhang /* send and recv matrix values */ 578de0260b3SHong Zhang /*-----------------------------*/ 57942c57489SHong Zhang buf_ri = merge->buf_ri; 58042c57489SHong Zhang buf_rj = merge->buf_rj; 58142c57489SHong Zhang len_s = merge->len_s; 58242c57489SHong Zhang ierr = PetscObjectGetNewTag((PetscObject)merge->rowmap,&taga);CHKERRQ(ierr); 58342c57489SHong Zhang ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 58442c57489SHong Zhang 58542c57489SHong Zhang ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&s_waits);CHKERRQ(ierr); 58642c57489SHong Zhang for (proc=0,k=0; proc<size; proc++){ 58742c57489SHong Zhang if (!len_s[proc]) continue; 588de0260b3SHong Zhang i = merge->owners_co[proc]; 589de0260b3SHong Zhang ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 59042c57489SHong Zhang k++; 59142c57489SHong Zhang } 59282412ba7SHong Zhang ierr = PetscMalloc(size*sizeof(MPI_Status),&status);CHKERRQ(ierr); 59342c57489SHong Zhang ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr); 59442c57489SHong Zhang ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr); 59542c57489SHong Zhang ierr = PetscFree(status);CHKERRQ(ierr); 59642c57489SHong Zhang 59742c57489SHong Zhang ierr = PetscFree(s_waits);CHKERRQ(ierr); 59842c57489SHong Zhang ierr = PetscFree(r_waits);CHKERRQ(ierr); 59982412ba7SHong Zhang ierr = PetscFree(coa);CHKERRQ(ierr); 60042c57489SHong Zhang 60182412ba7SHong Zhang /* insert local and received values into C */ 60282412ba7SHong Zhang /*-----------------------------------------*/ 60342c57489SHong Zhang ierr = PetscMalloc((3*merge->nrecv+1)*sizeof(PetscInt**),&buf_ri_k);CHKERRQ(ierr); 60442c57489SHong Zhang nextrow = buf_ri_k + merge->nrecv; 60582412ba7SHong Zhang nextci = nextrow + merge->nrecv; 60642c57489SHong Zhang 60742c57489SHong Zhang for (k=0; k<merge->nrecv; k++){ 60842c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 60942c57489SHong Zhang nrows = *(buf_ri_k[k]); 61042c57489SHong Zhang nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 61182412ba7SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 61242c57489SHong Zhang } 61342c57489SHong Zhang 614de0260b3SHong Zhang for (i=0; i<cm; i++) { 61582412ba7SHong Zhang row = owners[rank] + i; /* global row index of C_seq */ 61642c57489SHong Zhang bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 617de0260b3SHong Zhang ba_i = ba + bi[i]; 61882412ba7SHong Zhang bnz = bi[i+1] - bi[i]; 61942c57489SHong Zhang /* add received vals into ba */ 62042c57489SHong Zhang for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 62142c57489SHong Zhang /* i-th row */ 62242c57489SHong Zhang if (i == *nextrow[k]) { 62382412ba7SHong Zhang cnz = *(nextci[k]+1) - *nextci[k]; 62482412ba7SHong Zhang cj = buf_rj[k] + *(nextci[k]); 62582412ba7SHong Zhang ca = abuf_r[k] + *(nextci[k]); 62682412ba7SHong Zhang nextcj = 0; 62782412ba7SHong Zhang for (j=0; nextcj<cnz; j++){ 62882412ba7SHong Zhang if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */ 62982412ba7SHong Zhang ba_i[j] += ca[nextcj++]; 63042c57489SHong Zhang } 63142c57489SHong Zhang } 63282412ba7SHong Zhang nextrow[k]++; nextci[k]++; 63342c57489SHong Zhang } 63442c57489SHong Zhang } 63582412ba7SHong Zhang ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 63682412ba7SHong Zhang flops += 2*cnz; 63742c57489SHong Zhang } 63882412ba7SHong Zhang ierr = MatSetOption(C,MAT_COLUMNS_SORTED);CHKERRQ(ierr); /* -cause delay? */ 63942c57489SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 64042c57489SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 64142c57489SHong Zhang 642de0260b3SHong Zhang ierr = PetscFree(ba);CHKERRQ(ierr); 64342c57489SHong Zhang ierr = PetscFree(abuf_r);CHKERRQ(ierr); 64442c57489SHong Zhang ierr = PetscFree(buf_ri_k);CHKERRQ(ierr); 64542c57489SHong Zhang ierr = PetscLogFlops(flops);CHKERRQ(ierr); 64642c57489SHong Zhang 64742c57489SHong Zhang PetscFunctionReturn(0); 64842c57489SHong Zhang } 649