182412ba7SHong Zhang 242c57489SHong Zhang /* 342c57489SHong Zhang Defines projective product routines where A is a MPIAIJ matrix 442c57489SHong Zhang C = P^T * A * P 542c57489SHong Zhang */ 642c57489SHong Zhang 7c6db04a5SJed Brown #include <../src/mat/impls/aij/seq/aij.h> /*I "petscmat.h" I*/ 8c6db04a5SJed Brown #include <../src/mat/utils/freespace.h> 9c6db04a5SJed Brown #include <../src/mat/impls/aij/mpi/mpiaij.h> 10c6db04a5SJed Brown #include <petscbt.h> 118563dfccSBarry Smith #include <petsctime.h> 1242c57489SHong Zhang 13f671be37SHong Zhang /* #define PTAP_PROFILE */ 1424ecddacSHong Zhang 1509573ac7SBarry Smith extern PetscErrorCode MatDestroy_MPIAIJ(Mat); 16cc6cb787SHong Zhang #undef __FUNCT__ 17f8487c73SHong Zhang #define __FUNCT__ "MatDestroy_MPIAIJ_PtAP" 18f8487c73SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A) 19cc6cb787SHong Zhang { 20cc6cb787SHong Zhang PetscErrorCode ierr; 21f8487c73SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 22f8487c73SHong Zhang Mat_PtAPMPI *ptap=a->ptap; 23cc6cb787SHong Zhang 24cc6cb787SHong Zhang PetscFunctionBegin; 25f8487c73SHong Zhang if (ptap) { 26c8b0795cSMark F. Adams Mat_Merge_SeqsToMPI *merge=ptap->merge; 27b7f45c76SHong Zhang ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr); 28f8487c73SHong Zhang ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 29a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 30a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 31c5af039cSHong Zhang ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */ 3205d62848SHong Zhang ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr); 3305d62848SHong Zhang ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr); 348403a639SHong Zhang if (ptap->AP_loc) { /* used by alg_rap */ 35681d504bSHong Zhang Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data; 36681d504bSHong Zhang ierr = PetscFree(ap->i);CHKERRQ(ierr); 37681d504bSHong Zhang ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr); 380d3441aeSHong Zhang ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr); 398403a639SHong Zhang } else { /* used by alg_ptap */ 408403a639SHong Zhang ierr = PetscFree(ptap->api);CHKERRQ(ierr); 418403a639SHong Zhang ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 42681d504bSHong Zhang } 432259aa2eSHong Zhang ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr); 442259aa2eSHong Zhang ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr); 45414894bdSHong Zhang if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);} 468403a639SHong Zhang if (merge) { /* used by alg_ptap */ 47cc6cb787SHong Zhang ierr = PetscFree(merge->id_r);CHKERRQ(ierr); 48cc6cb787SHong Zhang ierr = PetscFree(merge->len_s);CHKERRQ(ierr); 49cc6cb787SHong Zhang ierr = PetscFree(merge->len_r);CHKERRQ(ierr); 50cc6cb787SHong Zhang ierr = PetscFree(merge->bi);CHKERRQ(ierr); 51cc6cb787SHong Zhang ierr = PetscFree(merge->bj);CHKERRQ(ierr); 52c05d87d6SBarry Smith ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr); 53cc6cb787SHong Zhang ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr); 54c05d87d6SBarry Smith ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr); 55cc6cb787SHong Zhang ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr); 56445158ffSHong Zhang ierr = PetscFree(merge->coi);CHKERRQ(ierr); 57445158ffSHong Zhang ierr = PetscFree(merge->coj);CHKERRQ(ierr); 5805b42c5fSBarry Smith ierr = PetscFree(merge->owners_co);CHKERRQ(ierr); 596bf464f9SBarry Smith ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr); 60dce485f0SBarry Smith ierr = merge->destroy(A);CHKERRQ(ierr); 61f8487c73SHong Zhang ierr = PetscFree(ptap->merge);CHKERRQ(ierr); 62445158ffSHong Zhang } else { 63445158ffSHong Zhang ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr); 64bf0cc555SLisandro Dalcin } 65f8487c73SHong Zhang ierr = PetscFree(ptap);CHKERRQ(ierr); 66c8b0795cSMark F. Adams } 67cc6cb787SHong Zhang PetscFunctionReturn(0); 68cc6cb787SHong Zhang } 69cc6cb787SHong Zhang 70cc6cb787SHong Zhang #undef __FUNCT__ 71aa690a28SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP" 72aa690a28SHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M) 731a47ec54SHong Zhang { 741a47ec54SHong Zhang PetscErrorCode ierr; 751a47ec54SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 761a47ec54SHong Zhang Mat_PtAPMPI *ptap = a->ptap; 771a47ec54SHong Zhang 781a47ec54SHong Zhang PetscFunctionBegin; 791a47ec54SHong Zhang ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr); 801a47ec54SHong Zhang (*M)->ops->destroy = ptap->destroy; 811a47ec54SHong Zhang (*M)->ops->duplicate = ptap->duplicate; 821a47ec54SHong Zhang PetscFunctionReturn(0); 831a47ec54SHong Zhang } 841a47ec54SHong Zhang 8542c57489SHong Zhang #undef __FUNCT__ 86cf3ca8ceSHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ" 87cf3ca8ceSHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 8842c57489SHong Zhang { 8942c57489SHong Zhang PetscErrorCode ierr; 908403a639SHong Zhang PetscBool rap=PETSC_TRUE; /* do R=P^T locally, then C=R*A*P */ 9167a12041SHong Zhang MPI_Comm comm; 9242c57489SHong Zhang 9342c57489SHong Zhang PetscFunctionBegin; 9467a12041SHong Zhang /* check if matrix local sizes are compatible */ 9567a12041SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 9667a12041SHong Zhang if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) { 9767a12041SHong Zhang SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Arow (%D, %D) != Prow (%D,%D)",A->rmap->rstart,A->rmap->rend,P->rmap->rstart,P->rmap->rend); 9867a12041SHong Zhang } 9967a12041SHong Zhang if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) { 10067a12041SHong Zhang SETERRQ4(comm,PETSC_ERR_ARG_SIZ,"Matrix local dimensions are incompatible, Acol (%D, %D) != Prow (%D,%D)",A->cmap->rstart,A->cmap->rend,P->rmap->rstart,P->rmap->rend); 10167a12041SHong Zhang } 10267a12041SHong Zhang 1038403a639SHong Zhang ierr = PetscOptionsGetBool(NULL,"-matrap",&rap,NULL);CHKERRQ(ierr); 104cf3ca8ceSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 1053ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 1068403a639SHong Zhang if (rap) { /* do R=P^T locally, then C=R*A*P */ 107cf3ca8ceSHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr); 1088403a639SHong Zhang } else { /* do P^T*A*P */ 1098403a639SHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(A,P,fill,C);CHKERRQ(ierr); 1100d3441aeSHong Zhang } 1113ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 1127a7894deSKris Buschelman } 1133ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 1148403a639SHong Zhang ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr); 1153ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 11642c57489SHong Zhang PetscFunctionReturn(0); 11742c57489SHong Zhang } 11842c57489SHong Zhang 11942c57489SHong Zhang #undef __FUNCT__ 120aa690a28SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ" 121aa690a28SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 1220d3441aeSHong Zhang { 1230d3441aeSHong Zhang PetscErrorCode ierr; 1240d3441aeSHong Zhang Mat_PtAPMPI *ptap; 125681d504bSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 1262259aa2eSHong Zhang MPI_Comm comm; 1272259aa2eSHong Zhang PetscMPIInt size,rank; 12815a3b8e2SHong Zhang Mat Cmpi; 12915a3b8e2SHong Zhang PetscFreeSpaceList free_space=NULL,current_space=NULL; 130*d0e9a2f7SHong Zhang PetscInt am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n; 131*d0e9a2f7SHong Zhang PetscInt *lnk,i,k,pnz,row,nsend; 13215a3b8e2SHong Zhang PetscBT lnkbt; 133f671be37SHong Zhang PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 13415a3b8e2SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 135*d0e9a2f7SHong Zhang PetscInt len,proc,*dnz,*onz,*owners,nzi,nspacedouble; 13615a3b8e2SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 13715a3b8e2SHong Zhang MPI_Request *swaits,*rwaits; 13815a3b8e2SHong Zhang MPI_Status *sstatus,rstatus; 139445158ffSHong Zhang PetscLayout rowmap; 140445158ffSHong Zhang PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 141445158ffSHong Zhang PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 142*d0e9a2f7SHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,rmax,*aj,*ai,*pi; 14367a12041SHong Zhang Mat_SeqAIJ *p_loc,*p_oth,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*c_loc,*c_oth; 14467a12041SHong Zhang PetscScalar *apv; 14578d30b94SHong Zhang PetscTable ta; 14678d30b94SHong Zhang const char *algTypes[2] = {"scalable","nonscalable"}; 14778d30b94SHong Zhang PetscInt alg=1; /* set default algorithm */ 148aa690a28SHong Zhang #if defined(PETSC_USE_INFO) 1498cb82516SHong Zhang PetscReal apfill; 150aa690a28SHong Zhang #endif 151681d504bSHong Zhang #if defined(PTAP_PROFILE) 152681d504bSHong Zhang PetscLogDouble t0,t1,t11,t12,t2,t3,t4; 153681d504bSHong Zhang #endif 15467a12041SHong Zhang 15567a12041SHong Zhang PetscFunctionBegin; 15667a12041SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 15767a12041SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 15867a12041SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1598cb82516SHong Zhang #if defined(PTAP_PROFILE) 16080bb4639SHong Zhang ierr = PetscTime(&t0);CHKERRQ(ierr); 1618cb82516SHong Zhang #endif 1628cb82516SHong Zhang 16315a3b8e2SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 16415a3b8e2SHong Zhang ierr = PetscNew(&ptap);CHKERRQ(ierr); 16515a3b8e2SHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 16615a3b8e2SHong Zhang 16715a3b8e2SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 16815a3b8e2SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 16915a3b8e2SHong Zhang /* get P_loc by taking all local rows of P */ 17015a3b8e2SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 17115a3b8e2SHong Zhang 17267a12041SHong Zhang /* (0) compute Rd = Pd^T, Ro = Po^T */ 17367a12041SHong Zhang /* --------------------------------- */ 174de817e96SHong Zhang ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 175de817e96SHong Zhang ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 1768cb82516SHong Zhang #if defined(PTAP_PROFILE) 177445158ffSHong Zhang ierr = PetscTime(&t11);CHKERRQ(ierr); 1788cb82516SHong Zhang #endif 17915a3b8e2SHong Zhang 18067a12041SHong Zhang /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 18167a12041SHong Zhang /* ----------------------------------------------------------------- */ 18267a12041SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 18367a12041SHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 184445158ffSHong Zhang 1859a6dcab7SHong Zhang /* create and initialize a linked list */ 186*d0e9a2f7SHong Zhang Crmax = 5*(p_loc->rmax+p_oth->rmax + (PetscInt)(1.e-2*pN)); /* expected Crmax */ 187*d0e9a2f7SHong Zhang if (Crmax > pN) Crmax=pN; 188*d0e9a2f7SHong Zhang ierr = PetscTableCreate(Crmax,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */ 1894b38b95cSHong Zhang MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta); 1904b38b95cSHong Zhang MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta); 19178d30b94SHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */ 192*d0e9a2f7SHong Zhang /* printf("[%d] est %d, Crmax %d; pN %d\n",rank,5*(p_loc->rmax+p_oth->rmax + (PetscInt)(1.e-2*pN)),Crmax,pN); */ 19378d30b94SHong Zhang 19478d30b94SHong Zhang ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 195445158ffSHong Zhang 1968cb82516SHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 19767a12041SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ad->i[am]+ao->i[am]+p_loc->i[pm])),&free_space);CHKERRQ(ierr); 198445158ffSHong Zhang current_space = free_space; 19967a12041SHong Zhang nspacedouble = 0; 20067a12041SHong Zhang 20167a12041SHong Zhang ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 20267a12041SHong Zhang api[0] = 0; 203445158ffSHong Zhang for (i=0; i<am; i++) { 20467a12041SHong Zhang /* diagonal portion: Ad[i,:]*P */ 20567a12041SHong Zhang ai = ad->i; pi = p_loc->i; 20667a12041SHong Zhang nzi = ai[i+1] - ai[i]; 20767a12041SHong Zhang aj = ad->j + ai[i]; 208445158ffSHong Zhang for (j=0; j<nzi; j++) { 209445158ffSHong Zhang row = aj[j]; 21067a12041SHong Zhang pnz = pi[row+1] - pi[row]; 21167a12041SHong Zhang Jptr = p_loc->j + pi[row]; 212445158ffSHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 213445158ffSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 214445158ffSHong Zhang } 21567a12041SHong Zhang /* off-diagonal portion: Ao[i,:]*P */ 21667a12041SHong Zhang ai = ao->i; pi = p_oth->i; 21767a12041SHong Zhang nzi = ai[i+1] - ai[i]; 21867a12041SHong Zhang aj = ao->j + ai[i]; 219445158ffSHong Zhang for (j=0; j<nzi; j++) { 220445158ffSHong Zhang row = aj[j]; 22167a12041SHong Zhang pnz = pi[row+1] - pi[row]; 22267a12041SHong Zhang Jptr = p_oth->j + pi[row]; 223445158ffSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 224445158ffSHong Zhang } 225445158ffSHong Zhang apnz = lnk[0]; 226445158ffSHong Zhang api[i+1] = api[i] + apnz; 227445158ffSHong Zhang if (ap_rmax < apnz) ap_rmax = apnz; 228445158ffSHong Zhang 229445158ffSHong Zhang /* if free space is not available, double the total space in the list */ 230445158ffSHong Zhang if (current_space->local_remaining<apnz) { 231445158ffSHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 232445158ffSHong Zhang nspacedouble++; 233445158ffSHong Zhang } 234445158ffSHong Zhang 235445158ffSHong Zhang /* Copy data into free space, then initialize lnk */ 236445158ffSHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 237445158ffSHong Zhang 238445158ffSHong Zhang current_space->array += apnz; 239445158ffSHong Zhang current_space->local_used += apnz; 240445158ffSHong Zhang current_space->local_remaining -= apnz; 241445158ffSHong Zhang } 242681d504bSHong Zhang /* Allocate space for apj and apv, initialize apj, and */ 243445158ffSHong Zhang /* destroy list of free space and other temporary array(s) */ 244445158ffSHong Zhang ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 245445158ffSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 2469a6dcab7SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 2479a6dcab7SHong Zhang 248aa690a28SHong Zhang /* Create AP_loc for reuse */ 249445158ffSHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 250aa690a28SHong Zhang 251aa690a28SHong Zhang #if defined(PETSC_USE_INFO) 252aa690a28SHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 253aa690a28SHong Zhang ptap->AP_loc->info.mallocs = nspacedouble; 254aa690a28SHong Zhang ptap->AP_loc->info.fill_ratio_given = fill; 255aa690a28SHong Zhang ptap->AP_loc->info.fill_ratio_needed = apfill; 256aa690a28SHong Zhang 257aa690a28SHong Zhang if (api[am]) { 258aa690a28SHong Zhang ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr); 259aa690a28SHong Zhang ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 260aa690a28SHong Zhang } else { 261aa690a28SHong Zhang ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr); 262aa690a28SHong Zhang } 263aa690a28SHong Zhang #endif 264aa690a28SHong Zhang 2658cb82516SHong Zhang #if defined(PTAP_PROFILE) 266445158ffSHong Zhang ierr = PetscTime(&t12);CHKERRQ(ierr); 2678cb82516SHong Zhang #endif 26815a3b8e2SHong Zhang 269681d504bSHong Zhang /* (2-1) compute symbolic Co = Ro*AP_loc */ 270681d504bSHong Zhang /* ------------------------------------ */ 27167a12041SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 27267a12041SHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 2738cb82516SHong Zhang #if defined(PTAP_PROFILE) 27480bb4639SHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 2758cb82516SHong Zhang #endif 27615a3b8e2SHong Zhang 27767a12041SHong Zhang /* (3) send coj of C_oth to other processors */ 27867a12041SHong Zhang /* ------------------------------------------ */ 27915a3b8e2SHong Zhang /* determine row ownership */ 280445158ffSHong Zhang ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 281445158ffSHong Zhang rowmap->n = pn; 282445158ffSHong Zhang rowmap->bs = 1; 283445158ffSHong Zhang ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 284445158ffSHong Zhang owners = rowmap->range; 28515a3b8e2SHong Zhang 28615a3b8e2SHong Zhang /* determine the number of messages to send, their lengths */ 2878cb82516SHong Zhang ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 2888cb82516SHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 28915a3b8e2SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 29015a3b8e2SHong Zhang 29167a12041SHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 29267a12041SHong Zhang coi = c_oth->i; coj = c_oth->j; 29367a12041SHong Zhang con = ptap->C_oth->rmap->n; 29415a3b8e2SHong Zhang proc = 0; 29567a12041SHong Zhang for (i=0; i<con; i++) { 29615a3b8e2SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 29715a3b8e2SHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 29815a3b8e2SHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 29915a3b8e2SHong Zhang } 30015a3b8e2SHong Zhang 30115a3b8e2SHong Zhang len = 0; /* max length of buf_si[], see (4) */ 30215a3b8e2SHong Zhang owners_co[0] = 0; 30367a12041SHong Zhang nsend = 0; 30415a3b8e2SHong Zhang for (proc=0; proc<size; proc++) { 30515a3b8e2SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 30615a3b8e2SHong Zhang if (len_s[proc]) { 307445158ffSHong Zhang nsend++; 30815a3b8e2SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 30915a3b8e2SHong Zhang len += len_si[proc]; 31015a3b8e2SHong Zhang } 31115a3b8e2SHong Zhang } 31215a3b8e2SHong Zhang 31315a3b8e2SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 314445158ffSHong Zhang ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 315445158ffSHong Zhang ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 31615a3b8e2SHong Zhang 31715a3b8e2SHong Zhang /* post the Irecv and Isend of coj */ 31815a3b8e2SHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 319445158ffSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 320445158ffSHong Zhang ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 32115a3b8e2SHong Zhang for (proc=0, k=0; proc<size; proc++) { 32215a3b8e2SHong Zhang if (!len_s[proc]) continue; 32315a3b8e2SHong Zhang i = owners_co[proc]; 32415a3b8e2SHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 32515a3b8e2SHong Zhang k++; 32615a3b8e2SHong Zhang } 32715a3b8e2SHong Zhang 328681d504bSHong Zhang /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 329681d504bSHong Zhang /* ---------------------------------------- */ 330681d504bSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 331681d504bSHong Zhang c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 332681d504bSHong Zhang 333681d504bSHong Zhang /* receives coj are complete */ 334445158ffSHong Zhang for (i=0; i<nrecv; i++) { 335445158ffSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 33615a3b8e2SHong Zhang } 33715a3b8e2SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 338445158ffSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 33915a3b8e2SHong Zhang 34078d30b94SHong Zhang /* add received column indices into ta to update Crmax */ 34178d30b94SHong Zhang for (k=0; k<nrecv; k++) {/* k-th received message */ 34278d30b94SHong Zhang Jptr = buf_rj[k]; 34378d30b94SHong Zhang for (j=0; j<len_r[k]; j++) { 34478d30b94SHong Zhang ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 34578d30b94SHong Zhang } 34678d30b94SHong Zhang } 34778d30b94SHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 34878d30b94SHong Zhang ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 3499a6dcab7SHong Zhang 35015a3b8e2SHong Zhang /* (4) send and recv coi */ 35115a3b8e2SHong Zhang /*-----------------------*/ 35215a3b8e2SHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 353445158ffSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 35415a3b8e2SHong Zhang ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 35515a3b8e2SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 35615a3b8e2SHong Zhang for (proc=0,k=0; proc<size; proc++) { 35715a3b8e2SHong Zhang if (!len_s[proc]) continue; 35815a3b8e2SHong Zhang /* form outgoing message for i-structure: 35915a3b8e2SHong Zhang buf_si[0]: nrows to be sent 36015a3b8e2SHong Zhang [1:nrows]: row index (global) 36115a3b8e2SHong Zhang [nrows+1:2*nrows+1]: i-structure index 36215a3b8e2SHong Zhang */ 36315a3b8e2SHong Zhang /*-------------------------------------------*/ 36415a3b8e2SHong Zhang nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 36515a3b8e2SHong Zhang buf_si_i = buf_si + nrows+1; 36615a3b8e2SHong Zhang buf_si[0] = nrows; 36715a3b8e2SHong Zhang buf_si_i[0] = 0; 36815a3b8e2SHong Zhang nrows = 0; 36915a3b8e2SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 37015a3b8e2SHong Zhang nzi = coi[i+1] - coi[i]; 37115a3b8e2SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 37215a3b8e2SHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 37315a3b8e2SHong Zhang nrows++; 37415a3b8e2SHong Zhang } 37515a3b8e2SHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 37615a3b8e2SHong Zhang k++; 37715a3b8e2SHong Zhang buf_si += len_si[proc]; 37815a3b8e2SHong Zhang } 379681d504bSHong Zhang for (i=0; i<nrecv; i++) { 380445158ffSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 38115a3b8e2SHong Zhang } 38215a3b8e2SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 383445158ffSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 38415a3b8e2SHong Zhang 3858cb82516SHong Zhang ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 38615a3b8e2SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 38715a3b8e2SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 38815a3b8e2SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 3898cb82516SHong Zhang #if defined(PTAP_PROFILE) 39080bb4639SHong Zhang ierr = PetscTime(&t2);CHKERRQ(ierr); 3918cb82516SHong Zhang #endif 39267a12041SHong Zhang /* (5) compute the local portion of Cmpi */ 39367a12041SHong Zhang /* ------------------------------------------ */ 39478d30b94SHong Zhang /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 39578d30b94SHong Zhang ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 39615a3b8e2SHong Zhang current_space = free_space; 39715a3b8e2SHong Zhang 398445158ffSHong Zhang ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 399445158ffSHong Zhang for (k=0; k<nrecv; k++) { 40015a3b8e2SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 40115a3b8e2SHong Zhang nrows = *buf_ri_k[k]; 40215a3b8e2SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 40315a3b8e2SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 40415a3b8e2SHong Zhang } 405445158ffSHong Zhang 40678d30b94SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 40778d30b94SHong Zhang ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 40815a3b8e2SHong Zhang rmax = 0; 40915a3b8e2SHong Zhang for (i=0; i<pn; i++) { 41067a12041SHong Zhang /* add C_loc into Cmpi */ 41167a12041SHong Zhang nzi = c_loc->i[i+1] - c_loc->i[i]; 41267a12041SHong Zhang Jptr = c_loc->j + c_loc->i[i]; 41367a12041SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 41415a3b8e2SHong Zhang 41515a3b8e2SHong Zhang /* add received col data into lnk */ 416445158ffSHong Zhang for (k=0; k<nrecv; k++) { /* k-th received message */ 41715a3b8e2SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 41815a3b8e2SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 41915a3b8e2SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 42015a3b8e2SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 42115a3b8e2SHong Zhang nextrow[k]++; nextci[k]++; 42215a3b8e2SHong Zhang } 42315a3b8e2SHong Zhang } 424*d0e9a2f7SHong Zhang nzi = lnk[0]; 4258cb82516SHong Zhang 42615a3b8e2SHong Zhang /* copy data into free space, then initialize lnk */ 427*d0e9a2f7SHong Zhang ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 428*d0e9a2f7SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 429*d0e9a2f7SHong Zhang if (nzi > rmax) rmax = nzi; 43015a3b8e2SHong Zhang } 43115a3b8e2SHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 43215a3b8e2SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 433445158ffSHong Zhang ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 4348cb82516SHong Zhang #if defined(PTAP_PROFILE) 43580bb4639SHong Zhang ierr = PetscTime(&t3);CHKERRQ(ierr); 4368cb82516SHong Zhang #endif 43780bb4639SHong Zhang 43815a3b8e2SHong Zhang /* (6) create symbolic parallel matrix Cmpi */ 43915a3b8e2SHong Zhang /*------------------------------------------*/ 44015a3b8e2SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 44115a3b8e2SHong Zhang ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 44215a3b8e2SHong Zhang ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 44315a3b8e2SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 44415a3b8e2SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 44515a3b8e2SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 44615a3b8e2SHong Zhang 447445158ffSHong Zhang /* members in merge */ 448445158ffSHong Zhang ierr = PetscFree(id_r);CHKERRQ(ierr); 449445158ffSHong Zhang ierr = PetscFree(len_r);CHKERRQ(ierr); 450445158ffSHong Zhang ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 451445158ffSHong Zhang ierr = PetscFree(buf_ri);CHKERRQ(ierr); 452445158ffSHong Zhang ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 453445158ffSHong Zhang ierr = PetscFree(buf_rj);CHKERRQ(ierr); 454445158ffSHong Zhang ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 45515a3b8e2SHong Zhang 45615a3b8e2SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 45715a3b8e2SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 45815a3b8e2SHong Zhang c->ptap = ptap; 45915a3b8e2SHong Zhang ptap->rmax = ap_rmax; 4601a47ec54SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 4611a47ec54SHong Zhang ptap->destroy = Cmpi->ops->destroy; 46215a3b8e2SHong Zhang 46378d30b94SHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 46478d30b94SHong Zhang ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,2,algTypes[1],&alg,NULL);CHKERRQ(ierr); 46578d30b94SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 46678d30b94SHong Zhang 46778d30b94SHong Zhang if (alg == 1) { 46878d30b94SHong Zhang /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */ 4698cb82516SHong Zhang ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 47078d30b94SHong Zhang Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ; 47178d30b94SHong Zhang } else { 472e55be12dSHong Zhang Cmpi->ops->ptapnumeric = 0; /* not done yet */ 47378d30b94SHong Zhang SETERRQ(comm,PETSC_ERR_ARG_SIZ,"Not done yet"); 47478d30b94SHong Zhang } 4752259aa2eSHong Zhang 4761a47ec54SHong Zhang 4771a47ec54SHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 4781a47ec54SHong Zhang Cmpi->assembled = PETSC_FALSE; 4791a47ec54SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 480aa690a28SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 4811a47ec54SHong Zhang *C = Cmpi; 4821a47ec54SHong Zhang 4838cb82516SHong Zhang #if defined(PTAP_PROFILE) 48480bb4639SHong Zhang ierr = PetscTime(&t4);CHKERRQ(ierr); 485e9c1f85fSHong Zhang if (rank == 1) { 486445158ffSHong Zhang printf("PtAPSym: %g + %g + %g + %g + %g + %g = %g\n",t11-t0,t1-t11,t12-t11,t2-t2,t3-t2,t4-t3,t4-t0); 487e9c1f85fSHong Zhang } 4888cb82516SHong Zhang #endif 4890d3441aeSHong Zhang PetscFunctionReturn(0); 4900d3441aeSHong Zhang } 4910d3441aeSHong Zhang 4920d3441aeSHong Zhang #undef __FUNCT__ 493aa690a28SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ" 494aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 4950d3441aeSHong Zhang { 4960d3441aeSHong Zhang PetscErrorCode ierr; 4970d3441aeSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 4980d3441aeSHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 4998cb82516SHong Zhang Mat_SeqAIJ *ap,*p_loc,*p_oth,*c_seq; 5000d3441aeSHong Zhang Mat_PtAPMPI *ptap = c->ptap; 5019ce11a7cSHong Zhang Mat AP_loc,C_loc,C_oth; 5020d3441aeSHong Zhang PetscInt i,rstart,rend,cm,ncols,row; 5038cb82516SHong Zhang PetscInt *api,*apj,am = A->rmap->n,j,col,apnz; 5048cb82516SHong Zhang PetscScalar *apa; 5050d3441aeSHong Zhang const PetscInt *cols; 5060d3441aeSHong Zhang const PetscScalar *vals; 5078cb82516SHong Zhang #if defined(PTAP_PROFILE) 5088cb82516SHong Zhang PetscMPIInt rank; 5098cb82516SHong Zhang MPI_Comm comm; 5100d3441aeSHong Zhang PetscLogDouble t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi; 5118cb82516SHong Zhang #endif 5120d3441aeSHong Zhang 5130d3441aeSHong Zhang PetscFunctionBegin; 5140d3441aeSHong Zhang ierr = MatZeroEntries(C);CHKERRQ(ierr); 5150d3441aeSHong Zhang 516e497d3c8SHong Zhang /* 1) get R = Pd^T,Ro = Po^T */ 5178cb82516SHong Zhang #if defined(PTAP_PROFILE) 5180d3441aeSHong Zhang ierr = PetscTime(&t0);CHKERRQ(ierr); 5198cb82516SHong Zhang #endif 520748c7196SHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 5210d3441aeSHong Zhang ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 5220d3441aeSHong Zhang ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 523748c7196SHong Zhang } 5248cb82516SHong Zhang #if defined(PTAP_PROFILE) 5250d3441aeSHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 5260d3441aeSHong Zhang eR = t1 - t0; 5278cb82516SHong Zhang #endif 5280d3441aeSHong Zhang 529e497d3c8SHong Zhang /* 2) get AP_loc */ 5300d3441aeSHong Zhang AP_loc = ptap->AP_loc; 5318cb82516SHong Zhang ap = (Mat_SeqAIJ*)AP_loc->data; 5320d3441aeSHong Zhang 533e497d3c8SHong Zhang /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 5340d3441aeSHong Zhang /*-----------------------------------------------------*/ 535748c7196SHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 536748c7196SHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 5370d3441aeSHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 5380d3441aeSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 539e497d3c8SHong Zhang } 5400d3441aeSHong Zhang 5418cb82516SHong Zhang /* 2-2) compute numeric A_loc*P - dominating part */ 5428cb82516SHong Zhang /* ---------------------------------------------- */ 5430d3441aeSHong Zhang /* get data from symbolic products */ 5448cb82516SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 5458cb82516SHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 5468cb82516SHong Zhang apa = ptap->apa; 547681d504bSHong Zhang api = ap->i; 548681d504bSHong Zhang apj = ap->j; 549e497d3c8SHong Zhang for (i=0; i<am; i++) { 550445158ffSHong Zhang /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 551e497d3c8SHong Zhang AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa); 552e497d3c8SHong Zhang apnz = api[i+1] - api[i]; 553e497d3c8SHong Zhang for (j=0; j<apnz; j++) { 554e497d3c8SHong Zhang col = apj[j+api[i]]; 555e497d3c8SHong Zhang ap->a[j+ap->i[i]] = apa[col]; 556e497d3c8SHong Zhang apa[col] = 0.0; 5570d3441aeSHong Zhang } 558aa690a28SHong Zhang ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 5590d3441aeSHong Zhang } 5608cb82516SHong Zhang #if defined(PTAP_PROFILE) 5610d3441aeSHong Zhang ierr = PetscTime(&t2);CHKERRQ(ierr); 5620d3441aeSHong Zhang eAP = t2 - t1; 5638cb82516SHong Zhang #endif 5640d3441aeSHong Zhang 5658cb82516SHong Zhang /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 566445158ffSHong Zhang ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 567445158ffSHong Zhang ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 5689ce11a7cSHong Zhang C_loc = ptap->C_loc; 5699ce11a7cSHong Zhang C_oth = ptap->C_oth; 5708cb82516SHong Zhang #if defined(PTAP_PROFILE) 5710d3441aeSHong Zhang ierr = PetscTime(&t3);CHKERRQ(ierr); 5720d3441aeSHong Zhang eCseq = t3 - t2; 5738cb82516SHong Zhang #endif 5740d3441aeSHong Zhang 5750d3441aeSHong Zhang /* add C_loc and Co to to C */ 5760d3441aeSHong Zhang ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 5770d3441aeSHong Zhang 5780d3441aeSHong Zhang /* C_loc -> C */ 5790d3441aeSHong Zhang cm = C_loc->rmap->N; 5809ce11a7cSHong Zhang c_seq = (Mat_SeqAIJ*)C_loc->data; 5818cb82516SHong Zhang cols = c_seq->j; 5828cb82516SHong Zhang vals = c_seq->a; 5830d3441aeSHong Zhang for (i=0; i<cm; i++) { 5849ce11a7cSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 5850d3441aeSHong Zhang row = rstart + i; 5860d3441aeSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 5878cb82516SHong Zhang cols += ncols; vals += ncols; 5880d3441aeSHong Zhang } 5890d3441aeSHong Zhang 5900d3441aeSHong Zhang /* Co -> C, off-processor part */ 5919ce11a7cSHong Zhang cm = C_oth->rmap->N; 5929ce11a7cSHong Zhang c_seq = (Mat_SeqAIJ*)C_oth->data; 5938cb82516SHong Zhang cols = c_seq->j; 5948cb82516SHong Zhang vals = c_seq->a; 5959ce11a7cSHong Zhang for (i=0; i<cm; i++) { 5969ce11a7cSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 5970d3441aeSHong Zhang row = p->garray[i]; 5980d3441aeSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 5998cb82516SHong Zhang cols += ncols; vals += ncols; 6000d3441aeSHong Zhang } 6010d3441aeSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6020d3441aeSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6038cb82516SHong Zhang #if defined(PTAP_PROFILE) 6040d3441aeSHong Zhang ierr = PetscTime(&t4);CHKERRQ(ierr); 6050d3441aeSHong Zhang eCmpi = t4 - t3; 6060d3441aeSHong Zhang 6078cb82516SHong Zhang ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 6088cb82516SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 6090d3441aeSHong Zhang if (rank==1) { 6100d3441aeSHong Zhang ierr = PetscPrintf(MPI_COMM_SELF," R %g, AP %g, Cseq %g, Cmpi %g = %g\n", eR,eAP,eCseq,eCmpi,eR+eAP+eCseq+eCmpi);CHKERRQ(ierr); 6110d3441aeSHong Zhang } 6128cb82516SHong Zhang #endif 613748c7196SHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 6140d3441aeSHong Zhang PetscFunctionReturn(0); 6150d3441aeSHong Zhang } 6160d3441aeSHong Zhang 6170d3441aeSHong Zhang #undef __FUNCT__ 6188403a639SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap" 6198403a639SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,PetscReal fill,Mat *C) 62042c57489SHong Zhang { 62142c57489SHong Zhang PetscErrorCode ierr; 62277584fe6SHong Zhang Mat Cmpi; 623f8487c73SHong Zhang Mat_PtAPMPI *ptap; 6240298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 625f8487c73SHong Zhang Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 62642c57489SHong Zhang Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 62742c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 628ded4f561SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ; 62977584fe6SHong Zhang PetscInt *adi=ad->i,*aj,*aoi=ao->i,nnz; 630a914927fSHong Zhang PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 631d16ca5f0SHong Zhang PetscInt am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n; 63242c57489SHong Zhang PetscBT lnkbt; 633ce94432eSBarry Smith MPI_Comm comm; 634a914927fSHong Zhang PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0; 63542c57489SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 63642c57489SHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 63724ecddacSHong Zhang PetscInt nzi,*pti,*ptj; 63842c57489SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 639ded4f561SHong Zhang MPI_Request *swaits,*rwaits; 640ded4f561SHong Zhang MPI_Status *sstatus,rstatus; 64142c57489SHong Zhang Mat_Merge_SeqsToMPI *merge; 64277584fe6SHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0; 643d16ca5f0SHong Zhang PetscReal afill=1.0,afill_tmp; 644a914927fSHong Zhang PetscInt rmax; 64542c57489SHong Zhang 64642c57489SHong Zhang PetscFunctionBegin; 647ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 64883445d95SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 64983445d95SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 65083445d95SHong Zhang 651f8487c73SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 652b00a9115SJed Brown ierr = PetscNew(&ptap);CHKERRQ(ierr); 653b00a9115SJed Brown ierr = PetscNew(&merge);CHKERRQ(ierr); 6549f4155fbSHong Zhang ptap->merge = merge; 655f8487c73SHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 6566c6e00e0SHong Zhang 6576c6e00e0SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 658b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 659fe615159SHong Zhang 6606c6e00e0SHong Zhang /* get P_loc by taking all local rows of P */ 661a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 6626c6e00e0SHong Zhang 663a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 664a1a4e84aSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 6656c6e00e0SHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 6666c6e00e0SHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 66742c57489SHong Zhang 6682addb229SHong Zhang /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */ 6692addb229SHong Zhang /*--------------------------------------------------------------------------*/ 670854ce69bSBarry Smith ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 67182412ba7SHong Zhang api[0] = 0; 67283445d95SHong Zhang 67383445d95SHong Zhang /* create and initialize a linked list */ 674a914927fSHong Zhang ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 67583445d95SHong Zhang 676a914927fSHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */ 677d16ca5f0SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 678f4a743e1SHong Zhang current_space = free_space; 679f4a743e1SHong Zhang 680f4a743e1SHong Zhang for (i=0; i<am; i++) { 681f4a743e1SHong Zhang /* diagonal portion of A */ 682ded4f561SHong Zhang nzi = adi[i+1] - adi[i]; 68377584fe6SHong Zhang aj = ad->j + adi[i]; 684ded4f561SHong Zhang for (j=0; j<nzi; j++) { 68577584fe6SHong Zhang row = aj[j]; 68682412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 687ded4f561SHong Zhang Jptr = pj_loc + pi_loc[row]; 68883445d95SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 689a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 690f4a743e1SHong Zhang } 691f4a743e1SHong Zhang /* off-diagonal portion of A */ 692ded4f561SHong Zhang nzi = aoi[i+1] - aoi[i]; 69377584fe6SHong Zhang aj = ao->j + aoi[i]; 694ded4f561SHong Zhang for (j=0; j<nzi; j++) { 69577584fe6SHong Zhang row = aj[j]; 69682412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 697ded4f561SHong Zhang Jptr = pj_oth + pi_oth[row]; 698a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 699f4a743e1SHong Zhang } 700a914927fSHong Zhang apnz = lnk[0]; 70182412ba7SHong Zhang api[i+1] = api[i] + apnz; 70277584fe6SHong Zhang if (ap_rmax < apnz) ap_rmax = apnz; 703f4a743e1SHong Zhang 70483445d95SHong Zhang /* if free space is not available, double the total space in the list */ 70582412ba7SHong Zhang if (current_space->local_remaining<apnz) { 7062ba1acd5SHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 707f2b054eeSHong Zhang nspacedouble++; 708f4a743e1SHong Zhang } 709f4a743e1SHong Zhang 710f4a743e1SHong Zhang /* Copy data into free space, then initialize lnk */ 711a914927fSHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 7122205254eSKarl Rupp 71382412ba7SHong Zhang current_space->array += apnz; 71482412ba7SHong Zhang current_space->local_used += apnz; 71582412ba7SHong Zhang current_space->local_remaining -= apnz; 716f4a743e1SHong Zhang } 717a914927fSHong Zhang 71882412ba7SHong Zhang /* Allocate space for apj, initialize apj, and */ 719f4a743e1SHong Zhang /* destroy list of free space and other temporary array(s) */ 720854ce69bSBarry Smith ierr = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr); 72177584fe6SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 722118727c9SMark F. Adams afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1); 723d16ca5f0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 724f4a743e1SHong Zhang 7252addb229SHong Zhang /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/ 7262addb229SHong Zhang /*-----------------------------------------------------------------*/ 727de0260b3SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 72842c57489SHong Zhang 729ded4f561SHong Zhang /* then, compute symbolic Co = (p->B)^T*AP */ 730d0f46423SBarry Smith pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 73183445d95SHong Zhang >= (num of nonzero rows of C_seq) - pn */ 732854ce69bSBarry Smith ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 733de0260b3SHong Zhang coi[0] = 0; 73442c57489SHong Zhang 735d16ca5f0SHong Zhang /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */ 736d16ca5f0SHong Zhang nnz = fill*(poti[pon] + api[am]); 73722d28d08SBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 73842c57489SHong Zhang current_space = free_space; 73942c57489SHong Zhang 740de0260b3SHong Zhang for (i=0; i<pon; i++) { 74182412ba7SHong Zhang pnz = poti[i+1] - poti[i]; 74277584fe6SHong Zhang ptJ = potj + poti[i]; 74377584fe6SHong Zhang for (j=0; j<pnz; j++) { 74477584fe6SHong Zhang row = ptJ[j]; /* row of AP == col of Pot */ 74582412ba7SHong Zhang apnz = api[row+1] - api[row]; 746ded4f561SHong Zhang Jptr = apj + api[row]; 74783445d95SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 748a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 74942c57489SHong Zhang } 750a914927fSHong Zhang nnz = lnk[0]; 75142c57489SHong Zhang 75283445d95SHong Zhang /* If free space is not available, double the total space in the list */ 753ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 7544238b7adSHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 755d16ca5f0SHong Zhang nspacedouble++; 75642c57489SHong Zhang } 75742c57489SHong Zhang 75842c57489SHong Zhang /* Copy data into free space, and zero out denserows */ 759a914927fSHong Zhang ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 7602205254eSKarl Rupp 761ded4f561SHong Zhang current_space->array += nnz; 762ded4f561SHong Zhang current_space->local_used += nnz; 763ded4f561SHong Zhang current_space->local_remaining -= nnz; 7642205254eSKarl Rupp 765ded4f561SHong Zhang coi[i+1] = coi[i] + nnz; 76642c57489SHong Zhang } 7676b911d16SHong Zhang 7686b911d16SHong Zhang ierr = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr); 769a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 770118727c9SMark F. Adams afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1); 771d16ca5f0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 772de0260b3SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 77342c57489SHong Zhang 7746b911d16SHong Zhang /* (3) send j-array (coj) of Co to other processors */ 7756b911d16SHong Zhang /*--------------------------------------------------*/ 7766b911d16SHong Zhang ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr); 7776b911d16SHong Zhang len_s = merge->len_s; 7786b911d16SHong Zhang merge->nsend = 0; 7796b911d16SHong Zhang 7806b911d16SHong Zhang 78142c57489SHong Zhang /* determine row ownership */ 78226283091SBarry Smith ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 7837a2fc3feSBarry Smith merge->rowmap->n = pn; 7847a2fc3feSBarry Smith merge->rowmap->bs = 1; 7852205254eSKarl Rupp 78626283091SBarry Smith ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 7877a2fc3feSBarry Smith owners = merge->rowmap->range; 78842c57489SHong Zhang 78942c57489SHong Zhang /* determine the number of messages to send, their lengths */ 790dcca6d9dSJed Brown ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr); 79183445d95SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 792854ce69bSBarry Smith ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 793de0260b3SHong Zhang 79483445d95SHong Zhang proc = 0; 795de0260b3SHong Zhang for (i=0; i<pon; i++) { 796de0260b3SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 7972addb229SHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 798ee6e4b5bSHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 79983445d95SHong Zhang } 800de0260b3SHong Zhang 8012addb229SHong Zhang len = 0; /* max length of buf_si[], see (4) */ 802de0260b3SHong Zhang owners_co[0] = 0; 80342c57489SHong Zhang for (proc=0; proc<size; proc++) { 804de0260b3SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 805ee6e4b5bSHong Zhang if (len_s[proc]) { 80642c57489SHong Zhang merge->nsend++; 8072addb229SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 80842c57489SHong Zhang len += len_si[proc]; 80942c57489SHong Zhang } 81042c57489SHong Zhang } 81142c57489SHong Zhang 812ded4f561SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 8130298fd71SBarry Smith ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 81442c57489SHong Zhang ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 81542c57489SHong Zhang 816ded4f561SHong Zhang /* post the Irecv and Isend of coj */ 817529bc5dcSHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 818529bc5dcSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 819854ce69bSBarry Smith ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 82042c57489SHong Zhang for (proc=0, k=0; proc<size; proc++) { 82142c57489SHong Zhang if (!len_s[proc]) continue; 822de0260b3SHong Zhang i = owners_co[proc]; 823529bc5dcSHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 82442c57489SHong Zhang k++; 82542c57489SHong Zhang } 82642c57489SHong Zhang 827ded4f561SHong Zhang /* receives and sends of coj are complete */ 82877584fe6SHong Zhang for (i=0; i<merge->nrecv; i++) { 829c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 830ded4f561SHong Zhang } 831ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 8320c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 83382412ba7SHong Zhang 8346b911d16SHong Zhang /* (4) send and recv coi */ 8356b911d16SHong Zhang /*-----------------------*/ 836529bc5dcSHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 837529bc5dcSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 838854ce69bSBarry Smith ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 83942c57489SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 84042c57489SHong Zhang for (proc=0,k=0; proc<size; proc++) { 84142c57489SHong Zhang if (!len_s[proc]) continue; 84242c57489SHong Zhang /* form outgoing message for i-structure: 84342c57489SHong Zhang buf_si[0]: nrows to be sent 84442c57489SHong Zhang [1:nrows]: row index (global) 84542c57489SHong Zhang [nrows+1:2*nrows+1]: i-structure index 84642c57489SHong Zhang */ 84742c57489SHong Zhang /*-------------------------------------------*/ 8482addb229SHong Zhang nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 84942c57489SHong Zhang buf_si_i = buf_si + nrows+1; 85042c57489SHong Zhang buf_si[0] = nrows; 85142c57489SHong Zhang buf_si_i[0] = 0; 85242c57489SHong Zhang nrows = 0; 853de0260b3SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 854ded4f561SHong Zhang nzi = coi[i+1] - coi[i]; 855ded4f561SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 856de0260b3SHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 85742c57489SHong Zhang nrows++; 85842c57489SHong Zhang } 859529bc5dcSHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 86042c57489SHong Zhang k++; 86142c57489SHong Zhang buf_si += len_si[proc]; 86242c57489SHong Zhang } 863ded4f561SHong Zhang i = merge->nrecv; 864ded4f561SHong Zhang while (i--) { 865c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 866ded4f561SHong Zhang } 867ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 8680c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 869a914927fSHong Zhang 87024ecddacSHong Zhang ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr); 87142c57489SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 872ded4f561SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 87342c57489SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 87442c57489SHong Zhang 8756b911d16SHong Zhang /* (5) compute the local portion of C (mpi mat) */ 8766b911d16SHong Zhang /*----------------------------------------------*/ 877ded4f561SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 878ded4f561SHong Zhang 87924ecddacSHong Zhang /* allocate pti array and free space for accumulating nonzero column info */ 880854ce69bSBarry Smith ierr = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr); 88124ecddacSHong Zhang pti[0] = 0; 88242c57489SHong Zhang 883d16ca5f0SHong Zhang /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 884d16ca5f0SHong Zhang nnz = fill*(pi_loc[pm] + api[am]); 88522d28d08SBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 88642c57489SHong Zhang current_space = free_space; 88742c57489SHong Zhang 888dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 88942c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { 89042c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 89142c57489SHong Zhang nrows = *buf_ri_k[k]; 89242c57489SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 89342c57489SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 89442c57489SHong Zhang } 89542c57489SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 89608cb4532SHong Zhang rmax = 0; 89742c57489SHong Zhang for (i=0; i<pn; i++) { 898ded4f561SHong Zhang /* add pdt[i,:]*AP into lnk */ 899ded4f561SHong Zhang pnz = pdti[i+1] - pdti[i]; 90077584fe6SHong Zhang ptJ = pdtj + pdti[i]; 90177584fe6SHong Zhang for (j=0; j<pnz; j++) { 90277584fe6SHong Zhang row = ptJ[j]; /* row of AP == col of Pt */ 903ded4f561SHong Zhang apnz = api[row+1] - api[row]; 904ded4f561SHong Zhang Jptr = apj + api[row]; 905ded4f561SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 906a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 907ded4f561SHong Zhang } 908a914927fSHong Zhang 90942c57489SHong Zhang /* add received col data into lnk */ 91042c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 91142c57489SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 912ded4f561SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 913ded4f561SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 914a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 91542c57489SHong Zhang nextrow[k]++; nextci[k]++; 91642c57489SHong Zhang } 91742c57489SHong Zhang } 918a914927fSHong Zhang nnz = lnk[0]; 91942c57489SHong Zhang 92042c57489SHong Zhang /* if free space is not available, make more free space */ 921ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 9224238b7adSHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 923d16ca5f0SHong Zhang nspacedouble++; 92442c57489SHong Zhang } 92542c57489SHong Zhang /* copy data into free space, then initialize lnk */ 926a914927fSHong Zhang ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 927ded4f561SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 9282205254eSKarl Rupp 929ded4f561SHong Zhang current_space->array += nnz; 930ded4f561SHong Zhang current_space->local_used += nnz; 931ded4f561SHong Zhang current_space->local_remaining -= nnz; 9322205254eSKarl Rupp 93324ecddacSHong Zhang pti[i+1] = pti[i] + nnz; 93408cb4532SHong Zhang if (nnz > rmax) rmax = nnz; 93542c57489SHong Zhang } 936ded4f561SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 9370572522cSBarry Smith ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 93842c57489SHong Zhang 939854ce69bSBarry Smith ierr = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr); 94024ecddacSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr); 94124ecddacSHong Zhang afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1); 942d16ca5f0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 94342c57489SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 94442c57489SHong Zhang 9456b911d16SHong Zhang /* (6) create symbolic parallel matrix Cmpi */ 9466b911d16SHong Zhang /*------------------------------------------*/ 94777584fe6SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 94877584fe6SHong Zhang ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 94933d57670SJed Brown ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 95077584fe6SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 95177584fe6SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 95242c57489SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 953a2f3521dSMark F. Adams 954b091e509SHong Zhang merge->bi = pti; /* Cseq->i */ 955b091e509SHong Zhang merge->bj = ptj; /* Cseq->j */ 956b091e509SHong Zhang merge->coi = coi; /* Co->i */ 957b091e509SHong Zhang merge->coj = coj; /* Co->j */ 95842c57489SHong Zhang merge->buf_ri = buf_ri; 95942c57489SHong Zhang merge->buf_rj = buf_rj; 960de0260b3SHong Zhang merge->owners_co = owners_co; 96177584fe6SHong Zhang merge->destroy = Cmpi->ops->destroy; 96242c57489SHong Zhang 96377584fe6SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 96477584fe6SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 965f8487c73SHong Zhang c->ptap = ptap; 96677584fe6SHong Zhang ptap->api = api; 96777584fe6SHong Zhang ptap->apj = apj; 968d6ab1dc0SHong Zhang ptap->rmax = ap_rmax; 9698403a639SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 9708403a639SHong Zhang ptap->destroy = Cmpi->ops->destroy; 9718403a639SHong Zhang 9728403a639SHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 9738403a639SHong Zhang Cmpi->assembled = PETSC_FALSE; 9748403a639SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 9758403a639SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 9768403a639SHong Zhang Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap; 97777584fe6SHong Zhang *C = Cmpi; 978414894bdSHong Zhang 979414894bdSHong Zhang /* flag 'scalable' determines which implementations to be used: 980414894bdSHong Zhang 0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa; 981414894bdSHong Zhang 1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */ 982414894bdSHong Zhang /* set default scalable */ 983da0a95b2SSatish Balay ptap->scalable = PETSC_FALSE; /* PETSC_TRUE; */ 9842205254eSKarl Rupp 9850298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr); 986414894bdSHong Zhang if (!ptap->scalable) { /* Do dense axpy */ 9871795a4d1SJed Brown ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 988414894bdSHong Zhang } else { 9891795a4d1SJed Brown ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr); 990414894bdSHong Zhang } 991414894bdSHong Zhang 992f2b054eeSHong Zhang #if defined(PETSC_USE_INFO) 99324ecddacSHong Zhang if (pti[pn] != 0) { 99457622a8eSBarry Smith ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 99557622a8eSBarry Smith ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 996f2b054eeSHong Zhang } else { 99777584fe6SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 998f2b054eeSHong Zhang } 999f2b054eeSHong Zhang #endif 100042c57489SHong Zhang PetscFunctionReturn(0); 100142c57489SHong Zhang } 100242c57489SHong Zhang 100342c57489SHong Zhang #undef __FUNCT__ 10048403a639SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap" 10058403a639SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,Mat C) 100642c57489SHong Zhang { 100742c57489SHong Zhang PetscErrorCode ierr; 1008f8487c73SHong Zhang Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 100942c57489SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 1010de0260b3SHong Zhang Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 101142c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 1012f8487c73SHong Zhang Mat_PtAPMPI *ptap; 10139f4155fbSHong Zhang Mat_Merge_SeqsToMPI *merge; 10141d633602SHong Zhang PetscInt *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp; 101582412ba7SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj; 101682412ba7SHong Zhang PetscInt i,j,k,anz,pnz,apnz,nextap,row,*cj; 1017e900a757SHong Zhang MatScalar *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp; 1018d0f46423SBarry Smith PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1019ce94432eSBarry Smith MPI_Comm comm; 102042c57489SHong Zhang PetscMPIInt size,rank,taga,*len_s; 102182412ba7SHong Zhang PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 102242c57489SHong Zhang PetscInt **buf_ri,**buf_rj; 102350f5bf86SHong Zhang PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 102442c57489SHong Zhang MPI_Request *s_waits,*r_waits; 102542c57489SHong Zhang MPI_Status *status; 102682412ba7SHong Zhang MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 102782412ba7SHong Zhang PetscInt *api,*apj,*coi,*coj; 1028d0f46423SBarry Smith PetscInt *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend; 10296ce94e8fSHong Zhang PetscBool scalable; 103038571addSHong Zhang #if defined(PTAP_PROFILE) 1031e497d3c8SHong Zhang PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2; 103238571addSHong Zhang #endif 103342c57489SHong Zhang 103442c57489SHong Zhang PetscFunctionBegin; 1035ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 103638571addSHong Zhang #if defined(PTAP_PROFILE) 10378563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 103838571addSHong Zhang #endif 103942c57489SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 104042c57489SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 104142c57489SHong Zhang 1042f8487c73SHong Zhang ptap = c->ptap; 1043ce94432eSBarry Smith if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_INCOMP,"MatPtAP() has not been called to create matrix C yet, cannot use MAT_REUSE_MATRIX"); 1044f8487c73SHong Zhang merge = ptap->merge; 1045414894bdSHong Zhang apa = ptap->apa; 10466ce94e8fSHong Zhang scalable = ptap->scalable; 10476c6e00e0SHong Zhang 1048a1a4e84aSHong Zhang /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 10496b911d16SHong Zhang /*-----------------------------------------------------*/ 1050f8487c73SHong Zhang if (ptap->reuse == MAT_INITIAL_MATRIX) { 10519f4155fbSHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */ 1052f8487c73SHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 10536c6e00e0SHong Zhang } else { /* update numerical values of P_oth and P_loc */ 1054b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1055a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 10566c6e00e0SHong Zhang } 105738571addSHong Zhang #if defined(PTAP_PROFILE) 10588563dfccSBarry Smith ierr = PetscTime(&t1);CHKERRQ(ierr); 105905d62848SHong Zhang eP = t1-t0; 106038571addSHong Zhang #endif 106105d62848SHong Zhang /* 106205d62848SHong Zhang printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank, 106305d62848SHong Zhang a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N, 106405d62848SHong Zhang ptap->P_loc->rmap->N,ptap->P_loc->cmap->N, 106505d62848SHong Zhang ptap->P_oth->rmap->N,ptap->P_oth->cmap->N); 106605d62848SHong Zhang */ 1067f8487c73SHong Zhang 1068f9473cf7SHong Zhang /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */ 1069f9473cf7SHong Zhang /*--------------------------------------------------------------*/ 107042c57489SHong Zhang /* get data from symbolic products */ 1071a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1072a1a4e84aSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 1073a9d06231SHong Zhang pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a; 107442c57489SHong Zhang pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 107542c57489SHong Zhang 1076de0260b3SHong Zhang coi = merge->coi; coj = merge->coj; 10771795a4d1SJed Brown ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 1078de0260b3SHong Zhang 1079de0260b3SHong Zhang bi = merge->bi; bj = merge->bj; 10807a2fc3feSBarry Smith owners = merge->rowmap->range; 10811795a4d1SJed Brown ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); /* ba: Cseq->a */ 108242c57489SHong Zhang 1083a1a4e84aSHong Zhang api = ptap->api; apj = ptap->apj; 1084d9824c15SHong Zhang 10859516a85cSHong Zhang if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */ 1086b38d1a2dSJed Brown ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr); 10878cb82516SHong Zhang #if defined(PTAP_PROFILE) 108805d62848SHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 10898cb82516SHong Zhang #endif 1090a9d06231SHong Zhang for (i=0; i<am; i++) { 1091b091e509SHong Zhang #if defined(PTAP_PROFILE) 10928563dfccSBarry Smith ierr = PetscTime(&t2_0);CHKERRQ(ierr); 1093b091e509SHong Zhang #endif 1094a9d06231SHong Zhang /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */ 1095a9d06231SHong Zhang /*------------------------------------------------------------*/ 1096a9d06231SHong Zhang apJ = apj + api[i]; 1097a9d06231SHong Zhang 1098a9d06231SHong Zhang /* diagonal portion of A */ 1099a9d06231SHong Zhang anz = adi[i+1] - adi[i]; 1100a9d06231SHong Zhang adj = ad->j + adi[i]; 1101a9d06231SHong Zhang ada = ad->a + adi[i]; 1102a9d06231SHong Zhang for (j=0; j<anz; j++) { 1103a9d06231SHong Zhang row = adj[j]; 1104a9d06231SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 1105a9d06231SHong Zhang pj = pj_loc + pi_loc[row]; 1106a9d06231SHong Zhang pa = pa_loc + pi_loc[row]; 1107a9d06231SHong Zhang 1108a9d06231SHong Zhang /* perform dense axpy */ 1109e900a757SHong Zhang valtmp = ada[j]; 1110a9d06231SHong Zhang for (k=0; k<pnz; k++) { 1111e900a757SHong Zhang apa[pj[k]] += valtmp*pa[k]; 1112a9d06231SHong Zhang } 1113a9d06231SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1114a9d06231SHong Zhang } 1115a9d06231SHong Zhang 1116a9d06231SHong Zhang /* off-diagonal portion of A */ 1117a9d06231SHong Zhang anz = aoi[i+1] - aoi[i]; 1118a9d06231SHong Zhang aoj = ao->j + aoi[i]; 1119a9d06231SHong Zhang aoa = ao->a + aoi[i]; 1120a9d06231SHong Zhang for (j=0; j<anz; j++) { 1121a9d06231SHong Zhang row = aoj[j]; 1122a9d06231SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 1123a9d06231SHong Zhang pj = pj_oth + pi_oth[row]; 1124a9d06231SHong Zhang pa = pa_oth + pi_oth[row]; 1125a9d06231SHong Zhang 1126a9d06231SHong Zhang /* perform dense axpy */ 1127e900a757SHong Zhang valtmp = aoa[j]; 1128a9d06231SHong Zhang for (k=0; k<pnz; k++) { 1129e900a757SHong Zhang apa[pj[k]] += valtmp*pa[k]; 1130a9d06231SHong Zhang } 1131a9d06231SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1132a9d06231SHong Zhang } 1133b091e509SHong Zhang #if defined(PTAP_PROFILE) 11348563dfccSBarry Smith ierr = PetscTime(&t2_1);CHKERRQ(ierr); 1135b091e509SHong Zhang et2_AP += t2_1 - t2_0; 1136b091e509SHong Zhang #endif 1137a9d06231SHong Zhang 1138a9d06231SHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 1139a9d06231SHong Zhang /*--------------------------------------------------------------*/ 1140a9d06231SHong Zhang apnz = api[i+1] - api[i]; 1141a9d06231SHong Zhang /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */ 1142a9d06231SHong Zhang pnz = po->i[i+1] - po->i[i]; 1143e900a757SHong Zhang poJ = po->j + po->i[i]; 1144a9d06231SHong Zhang pA = po->a + po->i[i]; 1145a9d06231SHong Zhang for (j=0; j<pnz; j++) { 1146e900a757SHong Zhang row = poJ[j]; 1147e900a757SHong Zhang cnz = coi[row+1] - coi[row]; 1148e900a757SHong Zhang cj = coj + coi[row]; 1149e900a757SHong Zhang ca = coa + coi[row]; 1150a9d06231SHong Zhang /* perform dense axpy */ 1151e900a757SHong Zhang valtmp = pA[j]; 1152a9d06231SHong Zhang for (k=0; k<cnz; k++) { 1153e900a757SHong Zhang ca[k] += valtmp*apa[cj[k]]; 1154a9d06231SHong Zhang } 1155a9d06231SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1156a9d06231SHong Zhang } 1157a9d06231SHong Zhang /* put the value into Cd (diagonal part) */ 1158a9d06231SHong Zhang pnz = pd->i[i+1] - pd->i[i]; 1159e900a757SHong Zhang pdJ = pd->j + pd->i[i]; 1160a9d06231SHong Zhang pA = pd->a + pd->i[i]; 1161a9d06231SHong Zhang for (j=0; j<pnz; j++) { 1162e900a757SHong Zhang row = pdJ[j]; 1163e900a757SHong Zhang cnz = bi[row+1] - bi[row]; 1164e900a757SHong Zhang cj = bj + bi[row]; 1165e900a757SHong Zhang ca = ba + bi[row]; 1166a9d06231SHong Zhang /* perform dense axpy */ 1167e900a757SHong Zhang valtmp = pA[j]; 1168a9d06231SHong Zhang for (k=0; k<cnz; k++) { 1169e900a757SHong Zhang ca[k] += valtmp*apa[cj[k]]; 1170a9d06231SHong Zhang } 1171a9d06231SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1172a9d06231SHong Zhang } 11738403a639SHong Zhang 1174a9d06231SHong Zhang /* zero the current row of A*P */ 1175a9d06231SHong Zhang for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0; 1176b091e509SHong Zhang #if defined(PTAP_PROFILE) 11778563dfccSBarry Smith ierr = PetscTime(&t2_2);CHKERRQ(ierr); 117805d62848SHong Zhang ePtAP += t2_2 - t2_1; 1179b091e509SHong Zhang #endif 1180a9d06231SHong Zhang } 118138571addSHong Zhang } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */ 1182b38d1a2dSJed Brown ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr); 118338571addSHong Zhang /*-----------------------------------------------------------------------------------------*/ 1184a9d06231SHong Zhang pA=pa_loc; 118542c57489SHong Zhang for (i=0; i<am; i++) { 1186b091e509SHong Zhang #if defined(PTAP_PROFILE) 11878563dfccSBarry Smith ierr = PetscTime(&t2_0);CHKERRQ(ierr); 1188b091e509SHong Zhang #endif 1189f4a743e1SHong Zhang /* form i-th sparse row of A*P */ 119082412ba7SHong Zhang apnz = api[i+1] - api[i]; 119182412ba7SHong Zhang apJ = apj + api[i]; 119242c57489SHong Zhang /* diagonal portion of A */ 119382412ba7SHong Zhang anz = adi[i+1] - adi[i]; 11941d633602SHong Zhang adj = ad->j + adi[i]; 11951d633602SHong Zhang ada = ad->a + adi[i]; 119682412ba7SHong Zhang for (j=0; j<anz; j++) { 11971d633602SHong Zhang row = adj[j]; 119882412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 119982412ba7SHong Zhang pj = pj_loc + pi_loc[row]; 120082412ba7SHong Zhang pa = pa_loc + pi_loc[row]; 1201e900a757SHong Zhang valtmp = ada[j]; 1202f4a743e1SHong Zhang nextp = 0; 120382412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 120482412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 1205e900a757SHong Zhang apa[k] += valtmp*pa[nextp++]; 120642c57489SHong Zhang } 120742c57489SHong Zhang } 1208dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 120942c57489SHong Zhang } 121042c57489SHong Zhang /* off-diagonal portion of A */ 121182412ba7SHong Zhang anz = aoi[i+1] - aoi[i]; 12121d633602SHong Zhang aoj = ao->j + aoi[i]; 12131d633602SHong Zhang aoa = ao->a + aoi[i]; 121482412ba7SHong Zhang for (j=0; j<anz; j++) { 12151d633602SHong Zhang row = aoj[j]; 121682412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 121782412ba7SHong Zhang pj = pj_oth + pi_oth[row]; 121882412ba7SHong Zhang pa = pa_oth + pi_oth[row]; 1219e900a757SHong Zhang valtmp = aoa[j]; 1220f4a743e1SHong Zhang nextp = 0; 122182412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 122282412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 1223e900a757SHong Zhang apa[k] += valtmp*pa[nextp++]; 122442c57489SHong Zhang } 122542c57489SHong Zhang } 1226dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 122742c57489SHong Zhang } 1228b091e509SHong Zhang #if defined(PTAP_PROFILE) 12298563dfccSBarry Smith ierr = PetscTime(&t2_1);CHKERRQ(ierr); 1230b091e509SHong Zhang et2_AP += t2_1 - t2_0; 1231b091e509SHong Zhang #endif 123242c57489SHong Zhang 1233a9d06231SHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 1234a9d06231SHong Zhang /*--------------------------------------------------------------*/ 123582412ba7SHong Zhang pnz = pi_loc[i+1] - pi_loc[i]; 1236f9473cf7SHong Zhang pJ = pj_loc + pi_loc[i]; 123782412ba7SHong Zhang for (j=0; j<pnz; j++) { 123842c57489SHong Zhang nextap = 0; 1239f9473cf7SHong Zhang row = pJ[j]; /* global index */ 124082412ba7SHong Zhang if (row < pcstart || row >=pcend) { /* put the value into Co */ 1241e900a757SHong Zhang row = *poJ; 1242e900a757SHong Zhang cj = coj + coi[row]; 1243e900a757SHong Zhang ca = coa + coi[row]; poJ++; 124482412ba7SHong Zhang } else { /* put the value into Cd */ 1245e900a757SHong Zhang row = *pdJ; 1246e900a757SHong Zhang cj = bj + bi[row]; 1247e900a757SHong Zhang ca = ba + bi[row]; pdJ++; 124842c57489SHong Zhang } 1249e900a757SHong Zhang valtmp = pA[j]; 125082412ba7SHong Zhang for (k=0; nextap<apnz; k++) { 1251e900a757SHong Zhang if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++]; 125242c57489SHong Zhang } 1253dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 1254de0260b3SHong Zhang } 12550496f32cSHong Zhang pA += pnz; 125642c57489SHong Zhang /* zero the current row info for A*P */ 125782412ba7SHong Zhang ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr); 1258b091e509SHong Zhang #if defined(PTAP_PROFILE) 12598563dfccSBarry Smith ierr = PetscTime(&t2_2);CHKERRQ(ierr); 126005d62848SHong Zhang ePtAP += t2_2 - t2_1; 1261b091e509SHong Zhang #endif 126242c57489SHong Zhang } 1263d9824c15SHong Zhang } 126438571addSHong Zhang #if defined(PTAP_PROFILE) 12658563dfccSBarry Smith ierr = PetscTime(&t2);CHKERRQ(ierr); 126638571addSHong Zhang #endif 126742c57489SHong Zhang 1268a9d06231SHong Zhang /* 3) send and recv matrix values coa */ 1269a9d06231SHong Zhang /*------------------------------------*/ 127042c57489SHong Zhang buf_ri = merge->buf_ri; 127142c57489SHong Zhang buf_rj = merge->buf_rj; 127242c57489SHong Zhang len_s = merge->len_s; 1273fc42d0c8SSatish Balay ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 127442c57489SHong Zhang ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 127542c57489SHong Zhang 1276dcca6d9dSJed Brown ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 127742c57489SHong Zhang for (proc=0,k=0; proc<size; proc++) { 127842c57489SHong Zhang if (!len_s[proc]) continue; 1279de0260b3SHong Zhang i = merge->owners_co[proc]; 1280de0260b3SHong Zhang ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 128142c57489SHong Zhang k++; 128242c57489SHong Zhang } 12830c468ba9SBarry Smith if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 12840c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 128542c57489SHong Zhang 1286a9d06231SHong Zhang ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 128742c57489SHong Zhang ierr = PetscFree(r_waits);CHKERRQ(ierr); 128882412ba7SHong Zhang ierr = PetscFree(coa);CHKERRQ(ierr); 128938571addSHong Zhang #if defined(PTAP_PROFILE) 12908563dfccSBarry Smith ierr = PetscTime(&t3);CHKERRQ(ierr); 129138571addSHong Zhang #endif 129242c57489SHong Zhang 1293a9d06231SHong Zhang /* 4) insert local Cseq and received values into Cmpi */ 1294a9d06231SHong Zhang /*------------------------------------------------------*/ 1295dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 129642c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { 129742c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 129842c57489SHong Zhang nrows = *(buf_ri_k[k]); 129942c57489SHong Zhang nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 130082412ba7SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 130142c57489SHong Zhang } 130242c57489SHong Zhang 1303de0260b3SHong Zhang for (i=0; i<cm; i++) { 130482412ba7SHong Zhang row = owners[rank] + i; /* global row index of C_seq */ 130542c57489SHong Zhang bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1306de0260b3SHong Zhang ba_i = ba + bi[i]; 130782412ba7SHong Zhang bnz = bi[i+1] - bi[i]; 130842c57489SHong Zhang /* add received vals into ba */ 130942c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 131042c57489SHong Zhang /* i-th row */ 131142c57489SHong Zhang if (i == *nextrow[k]) { 131282412ba7SHong Zhang cnz = *(nextci[k]+1) - *nextci[k]; 131382412ba7SHong Zhang cj = buf_rj[k] + *(nextci[k]); 131482412ba7SHong Zhang ca = abuf_r[k] + *(nextci[k]); 131582412ba7SHong Zhang nextcj = 0; 131682412ba7SHong Zhang for (j=0; nextcj<cnz; j++) { 131782412ba7SHong Zhang if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 131882412ba7SHong Zhang ba_i[j] += ca[nextcj++]; 131942c57489SHong Zhang } 132042c57489SHong Zhang } 132182412ba7SHong Zhang nextrow[k]++; nextci[k]++; 1322c5af039cSHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 132342c57489SHong Zhang } 132442c57489SHong Zhang } 132582412ba7SHong Zhang ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 132642c57489SHong Zhang } 132742c57489SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 132842c57489SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 132942c57489SHong Zhang 1330de0260b3SHong Zhang ierr = PetscFree(ba);CHKERRQ(ierr); 1331c05d87d6SBarry Smith ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 133242c57489SHong Zhang ierr = PetscFree(abuf_r);CHKERRQ(ierr); 13330572522cSBarry Smith ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 133438571addSHong Zhang #if defined(PTAP_PROFILE) 13358563dfccSBarry Smith ierr = PetscTime(&t4);CHKERRQ(ierr); 13369516a85cSHong Zhang if (rank==1) { 133705d62848SHong Zhang ierr = PetscPrintf(MPI_COMM_SELF," [%d] PtAPNum %g/P + %g/PtAP( %g/A*P + %g/Pt*AP ) + %g/comm + %g/Cloc = %g\n\n",rank,eP,t2-t1,et2_AP,ePtAP,t3-t2,t4-t3,t4-t0);CHKERRQ(ierr); 13389516a85cSHong Zhang } 133938571addSHong Zhang #endif 134042c57489SHong Zhang PetscFunctionReturn(0); 134142c57489SHong Zhang } 1342