1be1d678aSKris Buschelman #define PETSCMAT_DLL 282412ba7SHong Zhang 342c57489SHong Zhang /* 442c57489SHong Zhang Defines projective product routines where A is a MPIAIJ matrix 542c57489SHong Zhang C = P^T * A * P 642c57489SHong Zhang */ 742c57489SHong Zhang 87c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h" /*I "petscmat.h" I*/ 97c4f633dSBarry Smith #include "../src/mat/utils/freespace.h" 107c4f633dSBarry Smith #include "../src/mat/impls/aij/mpi/mpiaij.h" 1142c57489SHong Zhang #include "petscbt.h" 1242c57489SHong Zhang 13cc6cb787SHong Zhang EXTERN PetscErrorCode MatDestroy_MPIAIJ(Mat); 14cc6cb787SHong Zhang #undef __FUNCT__ 15cc6cb787SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_MatPtAP" 16cc6cb787SHong Zhang PetscErrorCode PETSCMAT_DLLEXPORT MatDestroy_MPIAIJ_MatPtAP(Mat A) 17cc6cb787SHong Zhang { 18cc6cb787SHong Zhang PetscErrorCode ierr; 19cc6cb787SHong Zhang Mat_Merge_SeqsToMPI *merge; 20776b82aeSLisandro Dalcin PetscContainer container; 21cc6cb787SHong Zhang 22cc6cb787SHong Zhang PetscFunctionBegin; 23cc6cb787SHong Zhang ierr = PetscObjectQuery((PetscObject)A,"MatMergeSeqsToMPI",(PetscObject *)&container);CHKERRQ(ierr); 24cc6cb787SHong Zhang if (container) { 25776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&merge);CHKERRQ(ierr); 26cc6cb787SHong Zhang ierr = PetscFree(merge->id_r);CHKERRQ(ierr); 27cc6cb787SHong Zhang ierr = PetscFree(merge->len_s);CHKERRQ(ierr); 28cc6cb787SHong Zhang ierr = PetscFree(merge->len_r);CHKERRQ(ierr); 29cc6cb787SHong Zhang ierr = PetscFree(merge->bi);CHKERRQ(ierr); 30cc6cb787SHong Zhang ierr = PetscFree(merge->bj);CHKERRQ(ierr); 31c05d87d6SBarry Smith ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr); 32cc6cb787SHong Zhang ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr); 33c05d87d6SBarry Smith ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr); 34cc6cb787SHong Zhang ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr); 3505b42c5fSBarry Smith ierr = PetscFree(merge->coi);CHKERRQ(ierr); 3605b42c5fSBarry Smith ierr = PetscFree(merge->coj);CHKERRQ(ierr); 3705b42c5fSBarry Smith ierr = PetscFree(merge->owners_co);CHKERRQ(ierr); 3826283091SBarry Smith ierr = PetscLayoutDestroy(merge->rowmap);CHKERRQ(ierr); 39cc6cb787SHong Zhang 40776b82aeSLisandro Dalcin ierr = PetscContainerDestroy(container);CHKERRQ(ierr); 41cc6cb787SHong Zhang ierr = PetscObjectCompose((PetscObject)A,"MatMergeSeqsToMPI",0);CHKERRQ(ierr); 42cc6cb787SHong Zhang } 43cc6cb787SHong Zhang ierr = merge->MatDestroy(A);CHKERRQ(ierr); 44cc6cb787SHong Zhang ierr = PetscFree(merge);CHKERRQ(ierr); 45cc6cb787SHong Zhang PetscFunctionReturn(0); 46cc6cb787SHong Zhang } 47cc6cb787SHong Zhang 48cc6cb787SHong Zhang #undef __FUNCT__ 49cc6cb787SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP" 50f0c0a3a6SBarry Smith PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M) 51f0c0a3a6SBarry Smith { 52cc6cb787SHong Zhang PetscErrorCode ierr; 53cc6cb787SHong Zhang Mat_Merge_SeqsToMPI *merge; 54776b82aeSLisandro Dalcin PetscContainer container; 55cc6cb787SHong Zhang 56cc6cb787SHong Zhang PetscFunctionBegin; 57cc6cb787SHong Zhang ierr = PetscObjectQuery((PetscObject)A,"MatMergeSeqsToMPI",(PetscObject *)&container);CHKERRQ(ierr); 58cc6cb787SHong Zhang if (container) { 59776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&merge);CHKERRQ(ierr); 60*e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Container does not exit"); 61cc6cb787SHong Zhang ierr = (*merge->MatDuplicate)(A,op,M);CHKERRQ(ierr); 62cc6cb787SHong Zhang (*M)->ops->destroy = merge->MatDestroy; /* =MatDestroy_MPIAIJ, *M doesn't duplicate A's container! */ 63cc6cb787SHong Zhang (*M)->ops->duplicate = merge->MatDuplicate; /* =MatDuplicate_ MPIAIJ */ 64cc6cb787SHong Zhang PetscFunctionReturn(0); 65cc6cb787SHong Zhang } 66cc6cb787SHong Zhang 6742c57489SHong Zhang #undef __FUNCT__ 687a7894deSKris Buschelman #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ" 697a7894deSKris Buschelman PetscErrorCode MatPtAPSymbolic_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 7042c57489SHong Zhang { 7142c57489SHong Zhang PetscErrorCode ierr; 7242c57489SHong Zhang 7342c57489SHong Zhang PetscFunctionBegin; 747a7894deSKris Buschelman if (!P->ops->ptapsymbolic_mpiaij) { 75e32f2f54SBarry Smith SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",((PetscObject)A)->type_name,((PetscObject)P)->type_name); 7642c57489SHong Zhang } 777a7894deSKris Buschelman ierr = (*P->ops->ptapsymbolic_mpiaij)(A,P,fill,C);CHKERRQ(ierr); 787a7894deSKris Buschelman PetscFunctionReturn(0); 797a7894deSKris Buschelman } 807a7894deSKris Buschelman 817a7894deSKris Buschelman #undef __FUNCT__ 827a7894deSKris Buschelman #define __FUNCT__ "MatPtAPNumeric_MPIAIJ" 837a7894deSKris Buschelman PetscErrorCode MatPtAPNumeric_MPIAIJ(Mat A,Mat P,Mat C) 847a7894deSKris Buschelman { 857a7894deSKris Buschelman PetscErrorCode ierr; 867a7894deSKris Buschelman 877a7894deSKris Buschelman PetscFunctionBegin; 887a7894deSKris Buschelman if (!P->ops->ptapnumeric_mpiaij) { 89e32f2f54SBarry Smith SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_SUP,"Not implemented for A=%s and P=%s",((PetscObject)A)->type_name,((PetscObject)P)->type_name); 907a7894deSKris Buschelman } 917a7894deSKris Buschelman ierr = (*P->ops->ptapnumeric_mpiaij)(A,P,C);CHKERRQ(ierr); 9242c57489SHong Zhang PetscFunctionReturn(0); 9342c57489SHong Zhang } 9442c57489SHong Zhang 9542c57489SHong Zhang #undef __FUNCT__ 9642c57489SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ" 9742c57489SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 9842c57489SHong Zhang { 9942c57489SHong Zhang PetscErrorCode ierr; 1006c6e00e0SHong Zhang Mat B_mpi; 101de0260b3SHong Zhang Mat_MatMatMultMPI *ap; 102776b82aeSLisandro Dalcin PetscContainer container; 103a1a86e44SBarry Smith PetscFreeSpaceList free_space=PETSC_NULL,current_space=PETSC_NULL; 10483445d95SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 10542c57489SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 10642c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 107ded4f561SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ; 10883445d95SHong Zhang PetscInt *adi=ad->i,*adj=ad->j,*aoi=ao->i,*aoj=ao->j,nnz; 10913f74950SBarry Smith PetscInt nlnk,*lnk,*owners_co,*coi,*coj,i,k,pnz,row; 110d0f46423SBarry Smith PetscInt am=A->rmap->n,pN=P->cmap->N,pn=P->cmap->n; 11142c57489SHong Zhang PetscBT lnkbt; 1127adad957SLisandro Dalcin MPI_Comm comm=((PetscObject)A)->comm; 113ded4f561SHong Zhang PetscMPIInt size,rank,tag,*len_si,*len_s,*len_ri; 11442c57489SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 11542c57489SHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 116ded4f561SHong Zhang PetscInt nzi,*bi,*bj; 11742c57489SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 118ded4f561SHong Zhang MPI_Request *swaits,*rwaits; 119ded4f561SHong Zhang MPI_Status *sstatus,rstatus; 12042c57489SHong Zhang Mat_Merge_SeqsToMPI *merge; 121f2b054eeSHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0; 12213f74950SBarry Smith PetscMPIInt j; 12342c57489SHong Zhang 12442c57489SHong Zhang PetscFunctionBegin; 12583445d95SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 12683445d95SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 12783445d95SHong Zhang 1286c6e00e0SHong Zhang /* destroy the container 'Mat_MatMatMultMPI' in case that P is attached to it */ 1296c6e00e0SHong Zhang ierr = PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&container);CHKERRQ(ierr); 1306c6e00e0SHong Zhang if (container) { 1314067f9b5SMatthew Knepley /* reset functions */ 1324067f9b5SMatthew Knepley ierr = PetscContainerGetPointer(container,(void **)&ap);CHKERRQ(ierr); 1334067f9b5SMatthew Knepley P->ops->destroy = ap->MatDestroy; 1344067f9b5SMatthew Knepley P->ops->duplicate = ap->MatDuplicate; 1354067f9b5SMatthew Knepley /* destroy container and contents */ 136776b82aeSLisandro Dalcin ierr = PetscContainerDestroy(container);CHKERRQ(ierr); 1376c6e00e0SHong Zhang ierr = PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",0);CHKERRQ(ierr); 1386c6e00e0SHong Zhang } 1396c6e00e0SHong Zhang 1406c6e00e0SHong Zhang /* create the container 'Mat_MatMatMultMPI' and attach it to P */ 1416c6e00e0SHong Zhang ierr = PetscNew(Mat_MatMatMultMPI,&ap);CHKERRQ(ierr); 1426c6e00e0SHong Zhang ap->abi=PETSC_NULL; ap->abj=PETSC_NULL; 1436c6e00e0SHong Zhang ap->abnz_max = 0; 1446c6e00e0SHong Zhang 145776b82aeSLisandro Dalcin ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 146776b82aeSLisandro Dalcin ierr = PetscContainerSetPointer(container,ap);CHKERRQ(ierr); 1476c6e00e0SHong Zhang ierr = PetscObjectCompose((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject)container);CHKERRQ(ierr); 1486c6e00e0SHong Zhang ap->MatDestroy = P->ops->destroy; 1496c6e00e0SHong Zhang P->ops->destroy = MatDestroy_MPIAIJ_MatMatMult; 1506c6e00e0SHong Zhang ap->reuse = MAT_INITIAL_MATRIX; 151776b82aeSLisandro Dalcin ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_Mat_MatMatMultMPI);CHKERRQ(ierr); 1526c6e00e0SHong Zhang 1536c6e00e0SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 1541d79065fSBarry Smith ierr = MatGetBrowsOfAoCols(A,P,MAT_INITIAL_MATRIX,&ap->startsj,&ap->startsj_r,&ap->bufa,&ap->B_oth);CHKERRQ(ierr); 1556c6e00e0SHong Zhang /* get P_loc by taking all local rows of P */ 1566c6e00e0SHong Zhang ierr = MatGetLocalMat(P,MAT_INITIAL_MATRIX,&ap->B_loc);CHKERRQ(ierr); 1576c6e00e0SHong Zhang 1586c6e00e0SHong Zhang p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data; 1596c6e00e0SHong Zhang p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data; 1606c6e00e0SHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 1616c6e00e0SHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 16242c57489SHong Zhang 163fd0ff01cSHong Zhang /* first, compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth */ 164fd0ff01cSHong Zhang /*-------------------------------------------------------------------*/ 16582412ba7SHong Zhang ierr = PetscMalloc((am+2)*sizeof(PetscInt),&api);CHKERRQ(ierr); 16682412ba7SHong Zhang ap->abi = api; 16782412ba7SHong Zhang api[0] = 0; 16883445d95SHong Zhang 16983445d95SHong Zhang /* create and initialize a linked list */ 17083445d95SHong Zhang nlnk = pN+1; 17183445d95SHong Zhang ierr = PetscLLCreate(pN,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 17283445d95SHong Zhang 17382412ba7SHong Zhang /* Initial FreeSpace size is fill*nnz(A) */ 174a1a86e44SBarry Smith ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am])),&free_space);CHKERRQ(ierr); 175f4a743e1SHong Zhang current_space = free_space; 176f4a743e1SHong Zhang 177f4a743e1SHong Zhang for (i=0;i<am;i++) { 17882412ba7SHong Zhang apnz = 0; 179f4a743e1SHong Zhang /* diagonal portion of A */ 180ded4f561SHong Zhang nzi = adi[i+1] - adi[i]; 181ded4f561SHong Zhang for (j=0; j<nzi; j++){ 18282412ba7SHong Zhang row = *adj++; 18382412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 184ded4f561SHong Zhang Jptr = pj_loc + pi_loc[row]; 18583445d95SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 186ded4f561SHong Zhang ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 18782412ba7SHong Zhang apnz += nlnk; 188f4a743e1SHong Zhang } 189f4a743e1SHong Zhang /* off-diagonal portion of A */ 190ded4f561SHong Zhang nzi = aoi[i+1] - aoi[i]; 191ded4f561SHong Zhang for (j=0; j<nzi; j++){ 19282412ba7SHong Zhang row = *aoj++; 19382412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 194ded4f561SHong Zhang Jptr = pj_oth + pi_oth[row]; 195ded4f561SHong Zhang ierr = PetscLLAdd(pnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 19682412ba7SHong Zhang apnz += nlnk; 197f4a743e1SHong Zhang } 198f4a743e1SHong Zhang 19982412ba7SHong Zhang api[i+1] = api[i] + apnz; 20082412ba7SHong Zhang if (ap->abnz_max < apnz) ap->abnz_max = apnz; 201f4a743e1SHong Zhang 20283445d95SHong Zhang /* if free space is not available, double the total space in the list */ 20382412ba7SHong Zhang if (current_space->local_remaining<apnz) { 2042ba1acd5SHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 205f2b054eeSHong Zhang nspacedouble++; 206f4a743e1SHong Zhang } 207f4a743e1SHong Zhang 208f4a743e1SHong Zhang /* Copy data into free space, then initialize lnk */ 20982412ba7SHong Zhang ierr = PetscLLClean(pN,pN,apnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 21082412ba7SHong Zhang current_space->array += apnz; 21182412ba7SHong Zhang current_space->local_used += apnz; 21282412ba7SHong Zhang current_space->local_remaining -= apnz; 213f4a743e1SHong Zhang } 21482412ba7SHong Zhang /* Allocate space for apj, initialize apj, and */ 215f4a743e1SHong Zhang /* destroy list of free space and other temporary array(s) */ 21682412ba7SHong Zhang ierr = PetscMalloc((api[am]+1)*sizeof(PetscInt),&ap->abj);CHKERRQ(ierr); 21782412ba7SHong Zhang apj = ap->abj; 218a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,ap->abj);CHKERRQ(ierr); 219f4a743e1SHong Zhang 220ded4f561SHong Zhang /* determine symbolic Co=(p->B)^T*AP - send to others */ 221ded4f561SHong Zhang /*----------------------------------------------------*/ 222de0260b3SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 22342c57489SHong Zhang 224ded4f561SHong Zhang /* then, compute symbolic Co = (p->B)^T*AP */ 225d0f46423SBarry Smith pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 22683445d95SHong Zhang >= (num of nonzero rows of C_seq) - pn */ 227de0260b3SHong Zhang ierr = PetscMalloc((pon+1)*sizeof(PetscInt),&coi);CHKERRQ(ierr); 228de0260b3SHong Zhang coi[0] = 0; 22942c57489SHong Zhang 23082412ba7SHong Zhang /* set initial free space to be 3*pon*max( nnz(AP) per row) */ 23182412ba7SHong Zhang nnz = 3*pon*ap->abnz_max + 1; 232a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space); 23342c57489SHong Zhang current_space = free_space; 23442c57489SHong Zhang 235de0260b3SHong Zhang for (i=0; i<pon; i++) { 236ded4f561SHong Zhang nnz = 0; 23782412ba7SHong Zhang pnz = poti[i+1] - poti[i]; 23882412ba7SHong Zhang j = pnz; 239de0260b3SHong Zhang ptJ = potj + poti[i+1]; 24083445d95SHong Zhang while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 24183445d95SHong Zhang j--; ptJ--; 24282412ba7SHong Zhang row = *ptJ; /* row of AP == col of Pot */ 24382412ba7SHong Zhang apnz = api[row+1] - api[row]; 244ded4f561SHong Zhang Jptr = apj + api[row]; 24583445d95SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 246ded4f561SHong Zhang ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 247ded4f561SHong Zhang nnz += nlnk; 24842c57489SHong Zhang } 24942c57489SHong Zhang 25083445d95SHong Zhang /* If free space is not available, double the total space in the list */ 251ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 2524238b7adSHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 25342c57489SHong Zhang } 25442c57489SHong Zhang 25542c57489SHong Zhang /* Copy data into free space, and zero out denserows */ 256ded4f561SHong Zhang ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 257ded4f561SHong Zhang current_space->array += nnz; 258ded4f561SHong Zhang current_space->local_used += nnz; 259ded4f561SHong Zhang current_space->local_remaining -= nnz; 260ded4f561SHong Zhang coi[i+1] = coi[i] + nnz; 26142c57489SHong Zhang } 262de0260b3SHong Zhang ierr = PetscMalloc((coi[pon]+1)*sizeof(PetscInt),&coj);CHKERRQ(ierr); 263a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 264de0260b3SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 26542c57489SHong Zhang 266ded4f561SHong Zhang /* send j-array (coj) of Co to other processors */ 267ded4f561SHong Zhang /*----------------------------------------------*/ 26842c57489SHong Zhang /* determine row ownership */ 26983445d95SHong Zhang ierr = PetscNew(Mat_Merge_SeqsToMPI,&merge);CHKERRQ(ierr); 27026283091SBarry Smith ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 2717a2fc3feSBarry Smith merge->rowmap->n = pn; 2727a2fc3feSBarry Smith merge->rowmap->bs = 1; 27326283091SBarry Smith ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 2747a2fc3feSBarry Smith owners = merge->rowmap->range; 27542c57489SHong Zhang 27642c57489SHong Zhang /* determine the number of messages to send, their lengths */ 27742c57489SHong Zhang ierr = PetscMalloc(size*sizeof(PetscMPIInt),&len_si);CHKERRQ(ierr); 27883445d95SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 27942c57489SHong Zhang ierr = PetscMalloc(size*sizeof(PetscMPIInt),&merge->len_s);CHKERRQ(ierr); 28042c57489SHong Zhang len_s = merge->len_s; 28142c57489SHong Zhang merge->nsend = 0; 28283445d95SHong Zhang 283de0260b3SHong Zhang ierr = PetscMalloc((size+2)*sizeof(PetscInt),&owners_co);CHKERRQ(ierr); 284de0260b3SHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 285de0260b3SHong Zhang 28683445d95SHong Zhang proc = 0; 287de0260b3SHong Zhang for (i=0; i<pon; i++){ 288de0260b3SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 289de0260b3SHong Zhang len_si[proc]++; /* num of rows in Co to be sent to [proc] */ 290de0260b3SHong Zhang len_s[proc] += coi[i+1] - coi[i]; 29183445d95SHong Zhang } 292de0260b3SHong Zhang 293de0260b3SHong Zhang len = 0; /* max length of buf_si[] */ 294de0260b3SHong Zhang owners_co[0] = 0; 29542c57489SHong Zhang for (proc=0; proc<size; proc++){ 296de0260b3SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 29783445d95SHong Zhang if (len_si[proc]){ 29842c57489SHong Zhang merge->nsend++; 29983445d95SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); 30042c57489SHong Zhang len += len_si[proc]; 30142c57489SHong Zhang } 30242c57489SHong Zhang } 30342c57489SHong Zhang 304ded4f561SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 30542c57489SHong Zhang ierr = PetscGatherNumberOfMessages(comm,PETSC_NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 30642c57489SHong Zhang ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 30742c57489SHong Zhang 308ded4f561SHong Zhang /* post the Irecv and Isend of coj */ 309fc42d0c8SSatish Balay ierr = PetscCommGetNewTag(comm,&tag);CHKERRQ(ierr); 310ded4f561SHong Zhang ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 31142c57489SHong Zhang 312ded4f561SHong Zhang ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&swaits);CHKERRQ(ierr); 31342c57489SHong Zhang 31442c57489SHong Zhang for (proc=0, k=0; proc<size; proc++){ 31542c57489SHong Zhang if (!len_s[proc]) continue; 316de0260b3SHong Zhang i = owners_co[proc]; 317ded4f561SHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr); 31842c57489SHong Zhang k++; 31942c57489SHong Zhang } 32042c57489SHong Zhang 321ded4f561SHong Zhang /* receives and sends of coj are complete */ 322ded4f561SHong Zhang ierr = PetscMalloc(size*sizeof(MPI_Status),&sstatus);CHKERRQ(ierr); 323ded4f561SHong Zhang i = merge->nrecv; 324ded4f561SHong Zhang while (i--) { 325ded4f561SHong Zhang ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 326ded4f561SHong Zhang } 327ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 3280c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 32982412ba7SHong Zhang 330ded4f561SHong Zhang /* send and recv coi */ 331ded4f561SHong Zhang /*-------------------*/ 332ded4f561SHong Zhang ierr = PetscPostIrecvInt(comm,tag,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 33342c57489SHong Zhang 33442c57489SHong Zhang ierr = PetscMalloc((len+1)*sizeof(PetscInt),&buf_s);CHKERRQ(ierr); 33542c57489SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 33642c57489SHong Zhang for (proc=0,k=0; proc<size; proc++){ 33742c57489SHong Zhang if (!len_s[proc]) continue; 33842c57489SHong Zhang /* form outgoing message for i-structure: 33942c57489SHong Zhang buf_si[0]: nrows to be sent 34042c57489SHong Zhang [1:nrows]: row index (global) 34142c57489SHong Zhang [nrows+1:2*nrows+1]: i-structure index 34242c57489SHong Zhang */ 34342c57489SHong Zhang /*-------------------------------------------*/ 34442c57489SHong Zhang nrows = len_si[proc]/2 - 1; 34542c57489SHong Zhang buf_si_i = buf_si + nrows+1; 34642c57489SHong Zhang buf_si[0] = nrows; 34742c57489SHong Zhang buf_si_i[0] = 0; 34842c57489SHong Zhang nrows = 0; 349de0260b3SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++){ 350ded4f561SHong Zhang nzi = coi[i+1] - coi[i]; 351ded4f561SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 352de0260b3SHong Zhang buf_si[nrows+1] =prmap[i] -owners[proc]; /* local row index */ 35342c57489SHong Zhang nrows++; 35442c57489SHong Zhang } 355ded4f561SHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tag,comm,swaits+k);CHKERRQ(ierr); 35642c57489SHong Zhang k++; 35742c57489SHong Zhang buf_si += len_si[proc]; 35842c57489SHong Zhang } 359ded4f561SHong Zhang i = merge->nrecv; 360ded4f561SHong Zhang while (i--) { 361ded4f561SHong Zhang ierr = MPI_Waitany(merge->nrecv,rwaits,&j,&rstatus);CHKERRQ(ierr); 362ded4f561SHong Zhang } 363ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 3640c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 365f2b054eeSHong Zhang /* 366ae15b995SBarry Smith ierr = PetscInfo2(A,"nsend: %d, nrecv: %d\n",merge->nsend,merge->nrecv);CHKERRQ(ierr); 36742c57489SHong Zhang for (i=0; i<merge->nrecv; i++){ 368ae15b995SBarry Smith ierr = PetscInfo3(A,"recv len_ri=%d, len_rj=%d from [%d]\n",len_ri[i],merge->len_r[i],merge->id_r[i]);CHKERRQ(ierr); 36942c57489SHong Zhang } 370f2b054eeSHong Zhang */ 37142c57489SHong Zhang ierr = PetscFree(len_si);CHKERRQ(ierr); 37242c57489SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 373ded4f561SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 374ded4f561SHong Zhang ierr = PetscFree(sstatus);CHKERRQ(ierr); 37542c57489SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 37642c57489SHong Zhang 377de0260b3SHong Zhang /* compute the local portion of C (mpi mat) */ 378de0260b3SHong Zhang /*------------------------------------------*/ 379ded4f561SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 380ded4f561SHong Zhang 38142c57489SHong Zhang /* allocate bi array and free space for accumulating nonzero column info */ 38242c57489SHong Zhang ierr = PetscMalloc((pn+1)*sizeof(PetscInt),&bi);CHKERRQ(ierr); 38342c57489SHong Zhang bi[0] = 0; 38442c57489SHong Zhang 385ded4f561SHong Zhang /* set initial free space to be 3*pn*max( nnz(AP) per row) */ 386ded4f561SHong Zhang nnz = 3*pn*ap->abnz_max + 1; 387a1a86e44SBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space); 38842c57489SHong Zhang current_space = free_space; 38942c57489SHong Zhang 390687f1162SBarry Smith ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr); 39142c57489SHong Zhang for (k=0; k<merge->nrecv; k++){ 39242c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 39342c57489SHong Zhang nrows = *buf_ri_k[k]; 39442c57489SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 39542c57489SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 39642c57489SHong Zhang } 39742c57489SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 39842c57489SHong Zhang for (i=0; i<pn; i++) { 399ded4f561SHong Zhang /* add pdt[i,:]*AP into lnk */ 400ded4f561SHong Zhang nnz = 0; 401ded4f561SHong Zhang pnz = pdti[i+1] - pdti[i]; 402ded4f561SHong Zhang j = pnz; 403ded4f561SHong Zhang ptJ = pdtj + pdti[i+1]; 404ded4f561SHong Zhang while (j){/* assume cols are almost in increasing order, starting from its end saves computation */ 405ded4f561SHong Zhang j--; ptJ--; 406ded4f561SHong Zhang row = *ptJ; /* row of AP == col of Pt */ 407ded4f561SHong Zhang apnz = api[row+1] - api[row]; 408ded4f561SHong Zhang Jptr = apj + api[row]; 409ded4f561SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 410ded4f561SHong Zhang ierr = PetscLLAdd(apnz,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 411ded4f561SHong Zhang nnz += nlnk; 412ded4f561SHong Zhang } 41342c57489SHong Zhang /* add received col data into lnk */ 41442c57489SHong Zhang for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 41542c57489SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 416ded4f561SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 417ded4f561SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 418ded4f561SHong Zhang ierr = PetscLLAdd(nzi,Jptr,pN,nlnk,lnk,lnkbt);CHKERRQ(ierr); 419ded4f561SHong Zhang nnz += nlnk; 42042c57489SHong Zhang nextrow[k]++; nextci[k]++; 42142c57489SHong Zhang } 42242c57489SHong Zhang } 42342c57489SHong Zhang 42442c57489SHong Zhang /* if free space is not available, make more free space */ 425ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 4264238b7adSHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 42742c57489SHong Zhang } 42842c57489SHong Zhang /* copy data into free space, then initialize lnk */ 429ded4f561SHong Zhang ierr = PetscLLClean(pN,pN,nnz,lnk,current_space->array,lnkbt);CHKERRQ(ierr); 430ded4f561SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 431ded4f561SHong Zhang current_space->array += nnz; 432ded4f561SHong Zhang current_space->local_used += nnz; 433ded4f561SHong Zhang current_space->local_remaining -= nnz; 434ded4f561SHong Zhang bi[i+1] = bi[i] + nnz; 43542c57489SHong Zhang } 436ded4f561SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 4370572522cSBarry Smith ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 43842c57489SHong Zhang 43942c57489SHong Zhang ierr = PetscMalloc((bi[pn]+1)*sizeof(PetscInt),&bj);CHKERRQ(ierr); 440a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,bj);CHKERRQ(ierr); 44142c57489SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 44242c57489SHong Zhang 44342c57489SHong Zhang /* create symbolic parallel matrix B_mpi */ 44442c57489SHong Zhang /*---------------------------------------*/ 445f69a0ea3SMatthew Knepley ierr = MatCreate(comm,&B_mpi);CHKERRQ(ierr); 446f69a0ea3SMatthew Knepley ierr = MatSetSizes(B_mpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 44742c57489SHong Zhang ierr = MatSetType(B_mpi,MATMPIAIJ);CHKERRQ(ierr); 44842c57489SHong Zhang ierr = MatMPIAIJSetPreallocation(B_mpi,0,dnz,0,onz);CHKERRQ(ierr); 44942c57489SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 45042c57489SHong Zhang 45142c57489SHong Zhang merge->bi = bi; 45242c57489SHong Zhang merge->bj = bj; 453de0260b3SHong Zhang merge->coi = coi; 454de0260b3SHong Zhang merge->coj = coj; 45542c57489SHong Zhang merge->buf_ri = buf_ri; 45642c57489SHong Zhang merge->buf_rj = buf_rj; 457de0260b3SHong Zhang merge->owners_co = owners_co; 458cc6cb787SHong Zhang merge->MatDestroy = B_mpi->ops->destroy; 459cc6cb787SHong Zhang merge->MatDuplicate = B_mpi->ops->duplicate; 460cc6cb787SHong Zhang 461cc6cb787SHong Zhang /* B_mpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 462cc6cb787SHong Zhang B_mpi->assembled = PETSC_FALSE; 463cc6cb787SHong Zhang B_mpi->ops->destroy = MatDestroy_MPIAIJ_MatPtAP; 464cc6cb787SHong Zhang B_mpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 46542c57489SHong Zhang 46642c57489SHong Zhang /* attach the supporting struct to B_mpi for reuse */ 467776b82aeSLisandro Dalcin ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 468776b82aeSLisandro Dalcin ierr = PetscContainerSetPointer(container,merge);CHKERRQ(ierr); 46942c57489SHong Zhang ierr = PetscObjectCompose((PetscObject)B_mpi,"MatMergeSeqsToMPI",(PetscObject)container);CHKERRQ(ierr); 47042c57489SHong Zhang *C = B_mpi; 471f2b054eeSHong Zhang #if defined(PETSC_USE_INFO) 472f2b054eeSHong Zhang if (bi[pn] != 0) { 473f2b054eeSHong Zhang PetscReal afill = ((PetscReal)bi[pn])/(adi[am]+aoi[am]); 474f2b054eeSHong Zhang if (afill < 1.0) afill = 1.0; 475f2b054eeSHong Zhang ierr = PetscInfo3(B_mpi,"Reallocs %D; Fill ratio: given %G needed %G when computing A*P.\n",nspacedouble,fill,afill);CHKERRQ(ierr); 476f2b054eeSHong Zhang ierr = PetscInfo1(B_mpi,"Use MatPtAP(A,P,MatReuse,%G,&C) for best performance.\n",afill);CHKERRQ(ierr); 477f2b054eeSHong Zhang } else { 478f2b054eeSHong Zhang ierr = PetscInfo(B_mpi,"Empty matrix product\n");CHKERRQ(ierr); 479f2b054eeSHong Zhang } 480f2b054eeSHong Zhang #endif 48142c57489SHong Zhang PetscFunctionReturn(0); 48242c57489SHong Zhang } 48342c57489SHong Zhang 48442c57489SHong Zhang #undef __FUNCT__ 48542c57489SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ" 48642c57489SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 48742c57489SHong Zhang { 48842c57489SHong Zhang PetscErrorCode ierr; 48942c57489SHong Zhang Mat_Merge_SeqsToMPI *merge; 490de0260b3SHong Zhang Mat_MatMatMultMPI *ap; 491776b82aeSLisandro Dalcin PetscContainer cont_merge,cont_ptap; 492de0260b3SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 49342c57489SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 494de0260b3SHong Zhang Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 49542c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 496bf3909cdSBarry Smith PetscInt *adi=ad->i,*aoi=ao->i,*adj=ad->j,*aoj=ao->j,*apJ,nextp; 49782412ba7SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj; 49882412ba7SHong Zhang PetscInt i,j,k,anz,pnz,apnz,nextap,row,*cj; 49982412ba7SHong Zhang MatScalar *ada=ad->a,*aoa=ao->a,*apa,*pa,*ca,*pa_loc,*pa_oth; 500d0f46423SBarry Smith PetscInt am=A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 5017adad957SLisandro Dalcin MPI_Comm comm=((PetscObject)C)->comm; 50242c57489SHong Zhang PetscMPIInt size,rank,taga,*len_s; 50382412ba7SHong Zhang PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 50442c57489SHong Zhang PetscInt **buf_ri,**buf_rj; 50550f5bf86SHong Zhang PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 50642c57489SHong Zhang MPI_Request *s_waits,*r_waits; 50742c57489SHong Zhang MPI_Status *status; 50882412ba7SHong Zhang MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 50982412ba7SHong Zhang PetscInt *api,*apj,*coi,*coj; 510d0f46423SBarry Smith PetscInt *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend; 51142c57489SHong Zhang 51242c57489SHong Zhang PetscFunctionBegin; 51342c57489SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 51442c57489SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 51542c57489SHong Zhang 51642c57489SHong Zhang ierr = PetscObjectQuery((PetscObject)C,"MatMergeSeqsToMPI",(PetscObject *)&cont_merge);CHKERRQ(ierr); 51742c57489SHong Zhang if (cont_merge) { 518776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(cont_merge,(void **)&merge);CHKERRQ(ierr); 519*e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Matrix C does not posses an object container"); 5206c6e00e0SHong Zhang 521f33d1a9aSHong Zhang ierr = PetscObjectQuery((PetscObject)P,"Mat_MatMatMultMPI",(PetscObject *)&cont_ptap);CHKERRQ(ierr); 52242c57489SHong Zhang if (cont_ptap) { 523776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(cont_ptap,(void **)&ap);CHKERRQ(ierr); 5246c6e00e0SHong Zhang if (ap->reuse == MAT_INITIAL_MATRIX){ 5256c6e00e0SHong Zhang ap->reuse = MAT_REUSE_MATRIX; 5266c6e00e0SHong Zhang } else { /* update numerical values of P_oth and P_loc */ 5271d79065fSBarry Smith ierr = MatGetBrowsOfAoCols(A,P,MAT_REUSE_MATRIX,&ap->startsj,&ap->startsj_r,&ap->bufa,&ap->B_oth);CHKERRQ(ierr); 5286c6e00e0SHong Zhang ierr = MatGetLocalMat(P,MAT_REUSE_MATRIX,&ap->B_loc);CHKERRQ(ierr); 5296c6e00e0SHong Zhang } 530*e7e72b3dSBarry Smith } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE, "Matrix P does not posses an object container"); 53142c57489SHong Zhang 53242c57489SHong Zhang /* get data from symbolic products */ 533de0260b3SHong Zhang p_loc = (Mat_SeqAIJ*)(ap->B_loc)->data; 534de0260b3SHong Zhang p_oth = (Mat_SeqAIJ*)(ap->B_oth)->data; 53582412ba7SHong Zhang pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a,pA=pa_loc; 53642c57489SHong Zhang pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 53742c57489SHong Zhang 538de0260b3SHong Zhang coi = merge->coi; coj = merge->coj; 539de0260b3SHong Zhang ierr = PetscMalloc((coi[pon]+1)*sizeof(MatScalar),&coa);CHKERRQ(ierr); 540de0260b3SHong Zhang ierr = PetscMemzero(coa,coi[pon]*sizeof(MatScalar));CHKERRQ(ierr); 541de0260b3SHong Zhang 542de0260b3SHong Zhang bi = merge->bi; bj = merge->bj; 5437a2fc3feSBarry Smith owners = merge->rowmap->range; 544de0260b3SHong Zhang ierr = PetscMalloc((bi[cm]+1)*sizeof(MatScalar),&ba);CHKERRQ(ierr); 545de0260b3SHong Zhang ierr = PetscMemzero(ba,bi[cm]*sizeof(MatScalar));CHKERRQ(ierr); 54642c57489SHong Zhang 547f4a743e1SHong Zhang /* get data from symbolic A*P */ 548de0260b3SHong Zhang ierr = PetscMalloc((ap->abnz_max+1)*sizeof(MatScalar),&apa);CHKERRQ(ierr); 549de0260b3SHong Zhang ierr = PetscMemzero(apa,ap->abnz_max*sizeof(MatScalar));CHKERRQ(ierr); 55042c57489SHong Zhang 551f4a743e1SHong Zhang /* compute numeric C_seq=P_loc^T*A_loc*P */ 55282412ba7SHong Zhang api = ap->abi; apj = ap->abj; 55342c57489SHong Zhang for (i=0;i<am;i++) { 554f4a743e1SHong Zhang /* form i-th sparse row of A*P */ 55582412ba7SHong Zhang apnz = api[i+1] - api[i]; 55682412ba7SHong Zhang apJ = apj + api[i]; 55742c57489SHong Zhang /* diagonal portion of A */ 55882412ba7SHong Zhang anz = adi[i+1] - adi[i]; 55982412ba7SHong Zhang for (j=0;j<anz;j++) { 56082412ba7SHong Zhang row = *adj++; 56182412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 56282412ba7SHong Zhang pj = pj_loc + pi_loc[row]; 56382412ba7SHong Zhang pa = pa_loc + pi_loc[row]; 564f4a743e1SHong Zhang nextp = 0; 56582412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 56682412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 56782412ba7SHong Zhang apa[k] += (*ada)*pa[nextp++]; 56842c57489SHong Zhang } 56942c57489SHong Zhang } 570dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 57142c57489SHong Zhang ada++; 57242c57489SHong Zhang } 57342c57489SHong Zhang /* off-diagonal portion of A */ 57482412ba7SHong Zhang anz = aoi[i+1] - aoi[i]; 57582412ba7SHong Zhang for (j=0; j<anz; j++) { 57682412ba7SHong Zhang row = *aoj++; 57782412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 57882412ba7SHong Zhang pj = pj_oth + pi_oth[row]; 57982412ba7SHong Zhang pa = pa_oth + pi_oth[row]; 580f4a743e1SHong Zhang nextp = 0; 58182412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 58282412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 58382412ba7SHong Zhang apa[k] += (*aoa)*pa[nextp++]; 58442c57489SHong Zhang } 58542c57489SHong Zhang } 586dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 58742c57489SHong Zhang aoa++; 58842c57489SHong Zhang } 58942c57489SHong Zhang 59082412ba7SHong Zhang /* Compute P_loc[i,:]^T*AP[i,:] using outer product */ 59182412ba7SHong Zhang pnz = pi_loc[i+1] - pi_loc[i]; 59282412ba7SHong Zhang for (j=0; j<pnz; j++) { 59342c57489SHong Zhang nextap = 0; 59482412ba7SHong Zhang row = *pJ++; /* global index */ 59582412ba7SHong Zhang if (row < pcstart || row >=pcend) { /* put the value into Co */ 59682412ba7SHong Zhang cj = coj + coi[*poJ]; 59782412ba7SHong Zhang ca = coa + coi[*poJ++]; 59882412ba7SHong Zhang } else { /* put the value into Cd */ 59982412ba7SHong Zhang cj = bj + bi[*pdJ]; 60082412ba7SHong Zhang ca = ba + bi[*pdJ++]; 60142c57489SHong Zhang } 60282412ba7SHong Zhang for (k=0; nextap<apnz; k++) { 60382412ba7SHong Zhang if (cj[k]==apJ[nextap]) ca[k] += (*pA)*apa[nextap++]; 60442c57489SHong Zhang } 605dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 60682412ba7SHong Zhang pA++; 607de0260b3SHong Zhang } 608de0260b3SHong Zhang 60942c57489SHong Zhang /* zero the current row info for A*P */ 61082412ba7SHong Zhang ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr); 61142c57489SHong Zhang } 61242c57489SHong Zhang ierr = PetscFree(apa);CHKERRQ(ierr); 61342c57489SHong Zhang 614de0260b3SHong Zhang /* send and recv matrix values */ 615de0260b3SHong Zhang /*-----------------------------*/ 61642c57489SHong Zhang buf_ri = merge->buf_ri; 61742c57489SHong Zhang buf_rj = merge->buf_rj; 61842c57489SHong Zhang len_s = merge->len_s; 619fc42d0c8SSatish Balay ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 62042c57489SHong Zhang ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 62142c57489SHong Zhang 62242c57489SHong Zhang ierr = PetscMalloc((merge->nsend+1)*sizeof(MPI_Request),&s_waits);CHKERRQ(ierr); 62342c57489SHong Zhang for (proc=0,k=0; proc<size; proc++){ 62442c57489SHong Zhang if (!len_s[proc]) continue; 625de0260b3SHong Zhang i = merge->owners_co[proc]; 626de0260b3SHong Zhang ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 62742c57489SHong Zhang k++; 62842c57489SHong Zhang } 62982412ba7SHong Zhang ierr = PetscMalloc(size*sizeof(MPI_Status),&status);CHKERRQ(ierr); 6300c468ba9SBarry Smith if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 6310c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 63242c57489SHong Zhang ierr = PetscFree(status);CHKERRQ(ierr); 63342c57489SHong Zhang 63442c57489SHong Zhang ierr = PetscFree(s_waits);CHKERRQ(ierr); 63542c57489SHong Zhang ierr = PetscFree(r_waits);CHKERRQ(ierr); 63682412ba7SHong Zhang ierr = PetscFree(coa);CHKERRQ(ierr); 63742c57489SHong Zhang 63882412ba7SHong Zhang /* insert local and received values into C */ 63982412ba7SHong Zhang /*-----------------------------------------*/ 640687f1162SBarry Smith ierr = PetscMalloc3(merge->nrecv,PetscInt**,&buf_ri_k,merge->nrecv,PetscInt*,&nextrow,merge->nrecv,PetscInt*,&nextci);CHKERRQ(ierr); 64142c57489SHong Zhang 64242c57489SHong Zhang for (k=0; k<merge->nrecv; k++){ 64342c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 64442c57489SHong Zhang nrows = *(buf_ri_k[k]); 64542c57489SHong Zhang nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 64682412ba7SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1);/* poins to the next i-structure of k-th recved i-structure */ 64742c57489SHong Zhang } 64842c57489SHong Zhang 649de0260b3SHong Zhang for (i=0; i<cm; i++) { 65082412ba7SHong Zhang row = owners[rank] + i; /* global row index of C_seq */ 65142c57489SHong Zhang bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 652de0260b3SHong Zhang ba_i = ba + bi[i]; 65382412ba7SHong Zhang bnz = bi[i+1] - bi[i]; 65442c57489SHong Zhang /* add received vals into ba */ 65542c57489SHong Zhang for (k=0; k<merge->nrecv; k++){ /* k-th received message */ 65642c57489SHong Zhang /* i-th row */ 65742c57489SHong Zhang if (i == *nextrow[k]) { 65882412ba7SHong Zhang cnz = *(nextci[k]+1) - *nextci[k]; 65982412ba7SHong Zhang cj = buf_rj[k] + *(nextci[k]); 66082412ba7SHong Zhang ca = abuf_r[k] + *(nextci[k]); 66182412ba7SHong Zhang nextcj = 0; 66282412ba7SHong Zhang for (j=0; nextcj<cnz; j++){ 66382412ba7SHong Zhang if (bj_i[j] == cj[nextcj]){ /* bcol == ccol */ 66482412ba7SHong Zhang ba_i[j] += ca[nextcj++]; 66542c57489SHong Zhang } 66642c57489SHong Zhang } 66782412ba7SHong Zhang nextrow[k]++; nextci[k]++; 66842c57489SHong Zhang } 66942c57489SHong Zhang } 67082412ba7SHong Zhang ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 671dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 67242c57489SHong Zhang } 673fd0ff01cSHong Zhang ierr = MatSetBlockSize(C,1);CHKERRQ(ierr); 67442c57489SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 67542c57489SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 67642c57489SHong Zhang 677de0260b3SHong Zhang ierr = PetscFree(ba);CHKERRQ(ierr); 678c05d87d6SBarry Smith ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 67942c57489SHong Zhang ierr = PetscFree(abuf_r);CHKERRQ(ierr); 6800572522cSBarry Smith ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 68142c57489SHong Zhang PetscFunctionReturn(0); 68242c57489SHong Zhang } 683