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; 13067a12041SHong Zhang PetscInt am=A->rmap->n,pm=P->rmap->n; 13167a12041SHong Zhang PetscInt *lnk,i,k,pnz,row,nnz; 132445158ffSHong Zhang PetscInt pN=P->cmap->N,pn=P->cmap->n; 13315a3b8e2SHong Zhang PetscBT lnkbt; 134f671be37SHong Zhang PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 13515a3b8e2SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 13615a3b8e2SHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 13767a12041SHong Zhang PetscInt nzi,nspacedouble; 13815a3b8e2SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 13915a3b8e2SHong Zhang MPI_Request *swaits,*rwaits; 14015a3b8e2SHong Zhang MPI_Status *sstatus,rstatus; 141445158ffSHong Zhang PetscLayout rowmap; 142445158ffSHong Zhang PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 143f671be37SHong Zhang PetscInt nsend; 144445158ffSHong Zhang PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 14578d30b94SHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax; 14667a12041SHong Zhang PetscInt rmax,*aj,*ai,*pi; 14767a12041SHong Zhang Mat_SeqAIJ *p_loc,*p_oth,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*c_loc,*c_oth; 14867a12041SHong Zhang PetscScalar *apv; 14978d30b94SHong Zhang PetscTable ta; 15078d30b94SHong Zhang const char *algTypes[2] = {"scalable","nonscalable"}; 15178d30b94SHong Zhang PetscInt alg=1; /* set default algorithm */ 152aa690a28SHong Zhang #if defined(PETSC_USE_INFO) 1538cb82516SHong Zhang PetscReal apfill; 154aa690a28SHong Zhang #endif 155681d504bSHong Zhang #if defined(PTAP_PROFILE) 156681d504bSHong Zhang PetscLogDouble t0,t1,t11,t12,t2,t3,t4; 157681d504bSHong Zhang #endif 15867a12041SHong Zhang 15967a12041SHong Zhang PetscFunctionBegin; 16067a12041SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 16167a12041SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 16267a12041SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1638cb82516SHong Zhang #if defined(PTAP_PROFILE) 16480bb4639SHong Zhang ierr = PetscTime(&t0);CHKERRQ(ierr); 1658cb82516SHong Zhang #endif 1668cb82516SHong Zhang 16715a3b8e2SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 16815a3b8e2SHong Zhang ierr = PetscNew(&ptap);CHKERRQ(ierr); 16915a3b8e2SHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 17015a3b8e2SHong Zhang 17115a3b8e2SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 17215a3b8e2SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 17315a3b8e2SHong Zhang /* get P_loc by taking all local rows of P */ 17415a3b8e2SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 17515a3b8e2SHong Zhang 17667a12041SHong Zhang /* (0) compute Rd = Pd^T, Ro = Po^T */ 17767a12041SHong Zhang /* --------------------------------- */ 178de817e96SHong Zhang ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 179de817e96SHong Zhang ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 1808cb82516SHong Zhang #if defined(PTAP_PROFILE) 181445158ffSHong Zhang ierr = PetscTime(&t11);CHKERRQ(ierr); 1828cb82516SHong Zhang #endif 18315a3b8e2SHong Zhang 18467a12041SHong Zhang /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 18567a12041SHong Zhang /* ----------------------------------------------------------------- */ 18667a12041SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 18767a12041SHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 188445158ffSHong Zhang 1899a6dcab7SHong Zhang /* create and initialize a linked list */ 19078d30b94SHong Zhang Crmax = p_loc->rmax+p_oth->rmax; 19178d30b94SHong Zhang ierr = PetscTableCreate(2*Crmax,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */ 192*4b38b95cSHong Zhang MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta); 193*4b38b95cSHong Zhang MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta); 19478d30b94SHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */ 19578d30b94SHong Zhang 19678d30b94SHong Zhang ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 197445158ffSHong Zhang 1988cb82516SHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 19967a12041SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ad->i[am]+ao->i[am]+p_loc->i[pm])),&free_space);CHKERRQ(ierr); 200445158ffSHong Zhang current_space = free_space; 20167a12041SHong Zhang nspacedouble = 0; 20267a12041SHong Zhang 20367a12041SHong Zhang ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 20467a12041SHong Zhang api[0] = 0; 205445158ffSHong Zhang for (i=0; i<am; i++) { 20667a12041SHong Zhang /* diagonal portion: Ad[i,:]*P */ 20767a12041SHong Zhang ai = ad->i; pi = p_loc->i; 20867a12041SHong Zhang nzi = ai[i+1] - ai[i]; 20967a12041SHong Zhang aj = ad->j + ai[i]; 210445158ffSHong Zhang for (j=0; j<nzi; j++) { 211445158ffSHong Zhang row = aj[j]; 21267a12041SHong Zhang pnz = pi[row+1] - pi[row]; 21367a12041SHong Zhang Jptr = p_loc->j + pi[row]; 214445158ffSHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 215445158ffSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 216445158ffSHong Zhang } 21767a12041SHong Zhang /* off-diagonal portion: Ao[i,:]*P */ 21867a12041SHong Zhang ai = ao->i; pi = p_oth->i; 21967a12041SHong Zhang nzi = ai[i+1] - ai[i]; 22067a12041SHong Zhang aj = ao->j + ai[i]; 221445158ffSHong Zhang for (j=0; j<nzi; j++) { 222445158ffSHong Zhang row = aj[j]; 22367a12041SHong Zhang pnz = pi[row+1] - pi[row]; 22467a12041SHong Zhang Jptr = p_oth->j + pi[row]; 225445158ffSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 226445158ffSHong Zhang } 227445158ffSHong Zhang apnz = lnk[0]; 228445158ffSHong Zhang api[i+1] = api[i] + apnz; 229445158ffSHong Zhang if (ap_rmax < apnz) ap_rmax = apnz; 230445158ffSHong Zhang 231445158ffSHong Zhang /* if free space is not available, double the total space in the list */ 232445158ffSHong Zhang if (current_space->local_remaining<apnz) { 233445158ffSHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 234445158ffSHong Zhang nspacedouble++; 235445158ffSHong Zhang } 236445158ffSHong Zhang 237445158ffSHong Zhang /* Copy data into free space, then initialize lnk */ 238445158ffSHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 239445158ffSHong Zhang 240445158ffSHong Zhang current_space->array += apnz; 241445158ffSHong Zhang current_space->local_used += apnz; 242445158ffSHong Zhang current_space->local_remaining -= apnz; 243445158ffSHong Zhang } 244681d504bSHong Zhang /* Allocate space for apj and apv, initialize apj, and */ 245445158ffSHong Zhang /* destroy list of free space and other temporary array(s) */ 246445158ffSHong Zhang ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 247445158ffSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 2489a6dcab7SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 2499a6dcab7SHong Zhang 250aa690a28SHong Zhang /* Create AP_loc for reuse */ 251445158ffSHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 252aa690a28SHong Zhang 253aa690a28SHong Zhang #if defined(PETSC_USE_INFO) 254aa690a28SHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 255aa690a28SHong Zhang ptap->AP_loc->info.mallocs = nspacedouble; 256aa690a28SHong Zhang ptap->AP_loc->info.fill_ratio_given = fill; 257aa690a28SHong Zhang ptap->AP_loc->info.fill_ratio_needed = apfill; 258aa690a28SHong Zhang 259aa690a28SHong Zhang if (api[am]) { 260aa690a28SHong Zhang ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr); 261aa690a28SHong Zhang ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 262aa690a28SHong Zhang } else { 263aa690a28SHong Zhang ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr); 264aa690a28SHong Zhang } 265aa690a28SHong Zhang #endif 266aa690a28SHong Zhang 2678cb82516SHong Zhang #if defined(PTAP_PROFILE) 268445158ffSHong Zhang ierr = PetscTime(&t12);CHKERRQ(ierr); 2698cb82516SHong Zhang #endif 27015a3b8e2SHong Zhang 271681d504bSHong Zhang /* (2-1) compute symbolic Co = Ro*AP_loc */ 272681d504bSHong Zhang /* ------------------------------------ */ 27367a12041SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 27467a12041SHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 2758cb82516SHong Zhang #if defined(PTAP_PROFILE) 27680bb4639SHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 2778cb82516SHong Zhang #endif 27815a3b8e2SHong Zhang 27967a12041SHong Zhang /* (3) send coj of C_oth to other processors */ 28067a12041SHong Zhang /* ------------------------------------------ */ 28115a3b8e2SHong Zhang /* determine row ownership */ 282445158ffSHong Zhang ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 283445158ffSHong Zhang rowmap->n = pn; 284445158ffSHong Zhang rowmap->bs = 1; 285445158ffSHong Zhang ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 286445158ffSHong Zhang owners = rowmap->range; 28715a3b8e2SHong Zhang 28815a3b8e2SHong Zhang /* determine the number of messages to send, their lengths */ 2898cb82516SHong Zhang ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 2908cb82516SHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 29115a3b8e2SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 29215a3b8e2SHong Zhang 29367a12041SHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 29467a12041SHong Zhang coi = c_oth->i; coj = c_oth->j; 29567a12041SHong Zhang con = ptap->C_oth->rmap->n; 29615a3b8e2SHong Zhang proc = 0; 29767a12041SHong Zhang for (i=0; i<con; i++) { 29815a3b8e2SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 29915a3b8e2SHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 30015a3b8e2SHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 30115a3b8e2SHong Zhang } 30215a3b8e2SHong Zhang 30315a3b8e2SHong Zhang len = 0; /* max length of buf_si[], see (4) */ 30415a3b8e2SHong Zhang owners_co[0] = 0; 30567a12041SHong Zhang nsend = 0; 30615a3b8e2SHong Zhang for (proc=0; proc<size; proc++) { 30715a3b8e2SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 30815a3b8e2SHong Zhang if (len_s[proc]) { 309445158ffSHong Zhang nsend++; 31015a3b8e2SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 31115a3b8e2SHong Zhang len += len_si[proc]; 31215a3b8e2SHong Zhang } 31315a3b8e2SHong Zhang } 31415a3b8e2SHong Zhang 31515a3b8e2SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 316445158ffSHong Zhang ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 317445158ffSHong Zhang ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 31815a3b8e2SHong Zhang 31915a3b8e2SHong Zhang /* post the Irecv and Isend of coj */ 32015a3b8e2SHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 321445158ffSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 322445158ffSHong Zhang ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 32315a3b8e2SHong Zhang for (proc=0, k=0; proc<size; proc++) { 32415a3b8e2SHong Zhang if (!len_s[proc]) continue; 32515a3b8e2SHong Zhang i = owners_co[proc]; 32615a3b8e2SHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 32715a3b8e2SHong Zhang k++; 32815a3b8e2SHong Zhang } 32915a3b8e2SHong Zhang 330681d504bSHong Zhang /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 331681d504bSHong Zhang /* ---------------------------------------- */ 332681d504bSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 333681d504bSHong Zhang c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 334681d504bSHong Zhang 335681d504bSHong Zhang /* receives coj are complete */ 336445158ffSHong Zhang for (i=0; i<nrecv; i++) { 337445158ffSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 33815a3b8e2SHong Zhang } 33915a3b8e2SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 340445158ffSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 34115a3b8e2SHong Zhang 34278d30b94SHong Zhang /* add received column indices into ta to update Crmax */ 34378d30b94SHong Zhang for (k=0; k<nrecv; k++) {/* k-th received message */ 34478d30b94SHong Zhang Jptr = buf_rj[k]; 34578d30b94SHong Zhang for (j=0; j<len_r[k]; j++) { 34678d30b94SHong Zhang ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 34778d30b94SHong Zhang } 34878d30b94SHong Zhang } 34978d30b94SHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 35078d30b94SHong Zhang ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 3519a6dcab7SHong Zhang 35215a3b8e2SHong Zhang /* (4) send and recv coi */ 35315a3b8e2SHong Zhang /*-----------------------*/ 35415a3b8e2SHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 355445158ffSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 35615a3b8e2SHong Zhang ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 35715a3b8e2SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 35815a3b8e2SHong Zhang for (proc=0,k=0; proc<size; proc++) { 35915a3b8e2SHong Zhang if (!len_s[proc]) continue; 36015a3b8e2SHong Zhang /* form outgoing message for i-structure: 36115a3b8e2SHong Zhang buf_si[0]: nrows to be sent 36215a3b8e2SHong Zhang [1:nrows]: row index (global) 36315a3b8e2SHong Zhang [nrows+1:2*nrows+1]: i-structure index 36415a3b8e2SHong Zhang */ 36515a3b8e2SHong Zhang /*-------------------------------------------*/ 36615a3b8e2SHong Zhang nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 36715a3b8e2SHong Zhang buf_si_i = buf_si + nrows+1; 36815a3b8e2SHong Zhang buf_si[0] = nrows; 36915a3b8e2SHong Zhang buf_si_i[0] = 0; 37015a3b8e2SHong Zhang nrows = 0; 37115a3b8e2SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 37215a3b8e2SHong Zhang nzi = coi[i+1] - coi[i]; 37315a3b8e2SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 37415a3b8e2SHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 37515a3b8e2SHong Zhang nrows++; 37615a3b8e2SHong Zhang } 37715a3b8e2SHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 37815a3b8e2SHong Zhang k++; 37915a3b8e2SHong Zhang buf_si += len_si[proc]; 38015a3b8e2SHong Zhang } 381681d504bSHong Zhang for (i=0; i<nrecv; i++) { 382445158ffSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 38315a3b8e2SHong Zhang } 38415a3b8e2SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 385445158ffSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 38615a3b8e2SHong Zhang 3878cb82516SHong Zhang ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 38815a3b8e2SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 38915a3b8e2SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 39015a3b8e2SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 3918cb82516SHong Zhang #if defined(PTAP_PROFILE) 39280bb4639SHong Zhang ierr = PetscTime(&t2);CHKERRQ(ierr); 3938cb82516SHong Zhang #endif 39467a12041SHong Zhang /* (5) compute the local portion of Cmpi */ 39567a12041SHong Zhang /* ------------------------------------------ */ 39678d30b94SHong Zhang /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 39778d30b94SHong Zhang ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 39815a3b8e2SHong Zhang current_space = free_space; 39915a3b8e2SHong Zhang 400445158ffSHong Zhang ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 401445158ffSHong Zhang for (k=0; k<nrecv; k++) { 40215a3b8e2SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 40315a3b8e2SHong Zhang nrows = *buf_ri_k[k]; 40415a3b8e2SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 40515a3b8e2SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 40615a3b8e2SHong Zhang } 407445158ffSHong Zhang 40878d30b94SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 40978d30b94SHong Zhang ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 41015a3b8e2SHong Zhang rmax = 0; 41115a3b8e2SHong Zhang for (i=0; i<pn; i++) { 41267a12041SHong Zhang /* add C_loc into Cmpi */ 41367a12041SHong Zhang nzi = c_loc->i[i+1] - c_loc->i[i]; 41467a12041SHong Zhang Jptr = c_loc->j + c_loc->i[i]; 41567a12041SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 41615a3b8e2SHong Zhang 41715a3b8e2SHong Zhang /* add received col data into lnk */ 418445158ffSHong Zhang for (k=0; k<nrecv; k++) { /* k-th received message */ 41915a3b8e2SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 42015a3b8e2SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 42115a3b8e2SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 42215a3b8e2SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 42315a3b8e2SHong Zhang nextrow[k]++; nextci[k]++; 42415a3b8e2SHong Zhang } 42515a3b8e2SHong Zhang } 42615a3b8e2SHong Zhang nnz = lnk[0]; 4278cb82516SHong Zhang 42815a3b8e2SHong Zhang /* copy data into free space, then initialize lnk */ 42915a3b8e2SHong Zhang ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 43015a3b8e2SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 43115a3b8e2SHong Zhang if (nnz > rmax) rmax = nnz; 43215a3b8e2SHong Zhang } 43315a3b8e2SHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 43415a3b8e2SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 435445158ffSHong Zhang ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 4368cb82516SHong Zhang #if defined(PTAP_PROFILE) 43780bb4639SHong Zhang ierr = PetscTime(&t3);CHKERRQ(ierr); 4388cb82516SHong Zhang #endif 43980bb4639SHong Zhang 44015a3b8e2SHong Zhang /* (6) create symbolic parallel matrix Cmpi */ 44115a3b8e2SHong Zhang /*------------------------------------------*/ 44215a3b8e2SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 44315a3b8e2SHong Zhang ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 44415a3b8e2SHong Zhang ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 44515a3b8e2SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 44615a3b8e2SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 44715a3b8e2SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 44815a3b8e2SHong Zhang 449445158ffSHong Zhang /* members in merge */ 450445158ffSHong Zhang ierr = PetscFree(id_r);CHKERRQ(ierr); 451445158ffSHong Zhang ierr = PetscFree(len_r);CHKERRQ(ierr); 452445158ffSHong Zhang ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 453445158ffSHong Zhang ierr = PetscFree(buf_ri);CHKERRQ(ierr); 454445158ffSHong Zhang ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 455445158ffSHong Zhang ierr = PetscFree(buf_rj);CHKERRQ(ierr); 456445158ffSHong Zhang ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 45715a3b8e2SHong Zhang 45815a3b8e2SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 45915a3b8e2SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 46015a3b8e2SHong Zhang c->ptap = ptap; 46115a3b8e2SHong Zhang ptap->rmax = ap_rmax; 4621a47ec54SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 4631a47ec54SHong Zhang ptap->destroy = Cmpi->ops->destroy; 46415a3b8e2SHong Zhang 46578d30b94SHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 46678d30b94SHong Zhang ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,2,algTypes[1],&alg,NULL);CHKERRQ(ierr); 46778d30b94SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 46878d30b94SHong Zhang 46978d30b94SHong Zhang if (alg == 1) { 47078d30b94SHong Zhang /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */ 4718cb82516SHong Zhang ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 47278d30b94SHong Zhang Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ; 47378d30b94SHong Zhang } else { 474e55be12dSHong Zhang Cmpi->ops->ptapnumeric = 0; /* not done yet */ 47578d30b94SHong Zhang SETERRQ(comm,PETSC_ERR_ARG_SIZ,"Not done yet"); 47678d30b94SHong Zhang } 4772259aa2eSHong Zhang 4781a47ec54SHong Zhang 4791a47ec54SHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 4801a47ec54SHong Zhang Cmpi->assembled = PETSC_FALSE; 4811a47ec54SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 482aa690a28SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 4831a47ec54SHong Zhang *C = Cmpi; 4841a47ec54SHong Zhang 4858cb82516SHong Zhang #if defined(PTAP_PROFILE) 48680bb4639SHong Zhang ierr = PetscTime(&t4);CHKERRQ(ierr); 487e9c1f85fSHong Zhang if (rank == 1) { 488445158ffSHong 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); 489e9c1f85fSHong Zhang } 4908cb82516SHong Zhang #endif 4910d3441aeSHong Zhang PetscFunctionReturn(0); 4920d3441aeSHong Zhang } 4930d3441aeSHong Zhang 4940d3441aeSHong Zhang #undef __FUNCT__ 495aa690a28SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ" 496aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 4970d3441aeSHong Zhang { 4980d3441aeSHong Zhang PetscErrorCode ierr; 4990d3441aeSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 5000d3441aeSHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 5018cb82516SHong Zhang Mat_SeqAIJ *ap,*p_loc,*p_oth,*c_seq; 5020d3441aeSHong Zhang Mat_PtAPMPI *ptap = c->ptap; 5039ce11a7cSHong Zhang Mat AP_loc,C_loc,C_oth; 5040d3441aeSHong Zhang PetscInt i,rstart,rend,cm,ncols,row; 5058cb82516SHong Zhang PetscInt *api,*apj,am = A->rmap->n,j,col,apnz; 5068cb82516SHong Zhang PetscScalar *apa; 5070d3441aeSHong Zhang const PetscInt *cols; 5080d3441aeSHong Zhang const PetscScalar *vals; 5098cb82516SHong Zhang #if defined(PTAP_PROFILE) 5108cb82516SHong Zhang PetscMPIInt rank; 5118cb82516SHong Zhang MPI_Comm comm; 5120d3441aeSHong Zhang PetscLogDouble t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi; 5138cb82516SHong Zhang #endif 5140d3441aeSHong Zhang 5150d3441aeSHong Zhang PetscFunctionBegin; 5160d3441aeSHong Zhang ierr = MatZeroEntries(C);CHKERRQ(ierr); 5170d3441aeSHong Zhang 518e497d3c8SHong Zhang /* 1) get R = Pd^T,Ro = Po^T */ 5198cb82516SHong Zhang #if defined(PTAP_PROFILE) 5200d3441aeSHong Zhang ierr = PetscTime(&t0);CHKERRQ(ierr); 5218cb82516SHong Zhang #endif 522748c7196SHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 5230d3441aeSHong Zhang ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 5240d3441aeSHong Zhang ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 525748c7196SHong Zhang } 5268cb82516SHong Zhang #if defined(PTAP_PROFILE) 5270d3441aeSHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 5280d3441aeSHong Zhang eR = t1 - t0; 5298cb82516SHong Zhang #endif 5300d3441aeSHong Zhang 531e497d3c8SHong Zhang /* 2) get AP_loc */ 5320d3441aeSHong Zhang AP_loc = ptap->AP_loc; 5338cb82516SHong Zhang ap = (Mat_SeqAIJ*)AP_loc->data; 5340d3441aeSHong Zhang 535e497d3c8SHong Zhang /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 5360d3441aeSHong Zhang /*-----------------------------------------------------*/ 537748c7196SHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 538748c7196SHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 5390d3441aeSHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 5400d3441aeSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 541e497d3c8SHong Zhang } 5420d3441aeSHong Zhang 5438cb82516SHong Zhang /* 2-2) compute numeric A_loc*P - dominating part */ 5448cb82516SHong Zhang /* ---------------------------------------------- */ 5450d3441aeSHong Zhang /* get data from symbolic products */ 5468cb82516SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 5478cb82516SHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 5488cb82516SHong Zhang apa = ptap->apa; 549681d504bSHong Zhang api = ap->i; 550681d504bSHong Zhang apj = ap->j; 551e497d3c8SHong Zhang for (i=0; i<am; i++) { 552445158ffSHong Zhang /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 553e497d3c8SHong Zhang AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa); 554e497d3c8SHong Zhang apnz = api[i+1] - api[i]; 555e497d3c8SHong Zhang for (j=0; j<apnz; j++) { 556e497d3c8SHong Zhang col = apj[j+api[i]]; 557e497d3c8SHong Zhang ap->a[j+ap->i[i]] = apa[col]; 558e497d3c8SHong Zhang apa[col] = 0.0; 5590d3441aeSHong Zhang } 560aa690a28SHong Zhang ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 5610d3441aeSHong Zhang } 5628cb82516SHong Zhang #if defined(PTAP_PROFILE) 5630d3441aeSHong Zhang ierr = PetscTime(&t2);CHKERRQ(ierr); 5640d3441aeSHong Zhang eAP = t2 - t1; 5658cb82516SHong Zhang #endif 5660d3441aeSHong Zhang 5678cb82516SHong Zhang /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 568445158ffSHong Zhang ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 569445158ffSHong Zhang ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 5709ce11a7cSHong Zhang C_loc = ptap->C_loc; 5719ce11a7cSHong Zhang C_oth = ptap->C_oth; 5728cb82516SHong Zhang #if defined(PTAP_PROFILE) 5730d3441aeSHong Zhang ierr = PetscTime(&t3);CHKERRQ(ierr); 5740d3441aeSHong Zhang eCseq = t3 - t2; 5758cb82516SHong Zhang #endif 5760d3441aeSHong Zhang 5770d3441aeSHong Zhang /* add C_loc and Co to to C */ 5780d3441aeSHong Zhang ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 5790d3441aeSHong Zhang 5800d3441aeSHong Zhang /* C_loc -> C */ 5810d3441aeSHong Zhang cm = C_loc->rmap->N; 5829ce11a7cSHong Zhang c_seq = (Mat_SeqAIJ*)C_loc->data; 5838cb82516SHong Zhang cols = c_seq->j; 5848cb82516SHong Zhang vals = c_seq->a; 5850d3441aeSHong Zhang for (i=0; i<cm; i++) { 5869ce11a7cSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 5870d3441aeSHong Zhang row = rstart + i; 5880d3441aeSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 5898cb82516SHong Zhang cols += ncols; vals += ncols; 5900d3441aeSHong Zhang } 5910d3441aeSHong Zhang 5920d3441aeSHong Zhang /* Co -> C, off-processor part */ 5939ce11a7cSHong Zhang cm = C_oth->rmap->N; 5949ce11a7cSHong Zhang c_seq = (Mat_SeqAIJ*)C_oth->data; 5958cb82516SHong Zhang cols = c_seq->j; 5968cb82516SHong Zhang vals = c_seq->a; 5979ce11a7cSHong Zhang for (i=0; i<cm; i++) { 5989ce11a7cSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 5990d3441aeSHong Zhang row = p->garray[i]; 6000d3441aeSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 6018cb82516SHong Zhang cols += ncols; vals += ncols; 6020d3441aeSHong Zhang } 6030d3441aeSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6040d3441aeSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6058cb82516SHong Zhang #if defined(PTAP_PROFILE) 6060d3441aeSHong Zhang ierr = PetscTime(&t4);CHKERRQ(ierr); 6070d3441aeSHong Zhang eCmpi = t4 - t3; 6080d3441aeSHong Zhang 6098cb82516SHong Zhang ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 6108cb82516SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 6110d3441aeSHong Zhang if (rank==1) { 6120d3441aeSHong 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); 6130d3441aeSHong Zhang } 6148cb82516SHong Zhang #endif 615748c7196SHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 6160d3441aeSHong Zhang PetscFunctionReturn(0); 6170d3441aeSHong Zhang } 6180d3441aeSHong Zhang 6190d3441aeSHong Zhang #undef __FUNCT__ 6208403a639SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap" 6218403a639SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,PetscReal fill,Mat *C) 62242c57489SHong Zhang { 62342c57489SHong Zhang PetscErrorCode ierr; 62477584fe6SHong Zhang Mat Cmpi; 625f8487c73SHong Zhang Mat_PtAPMPI *ptap; 6260298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 627f8487c73SHong Zhang Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 62842c57489SHong Zhang Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 62942c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 630ded4f561SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ; 63177584fe6SHong Zhang PetscInt *adi=ad->i,*aj,*aoi=ao->i,nnz; 632a914927fSHong Zhang PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 633d16ca5f0SHong Zhang PetscInt am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n; 63442c57489SHong Zhang PetscBT lnkbt; 635ce94432eSBarry Smith MPI_Comm comm; 636a914927fSHong Zhang PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0; 63742c57489SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 63842c57489SHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 63924ecddacSHong Zhang PetscInt nzi,*pti,*ptj; 64042c57489SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 641ded4f561SHong Zhang MPI_Request *swaits,*rwaits; 642ded4f561SHong Zhang MPI_Status *sstatus,rstatus; 64342c57489SHong Zhang Mat_Merge_SeqsToMPI *merge; 64477584fe6SHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0; 645d16ca5f0SHong Zhang PetscReal afill=1.0,afill_tmp; 646a914927fSHong Zhang PetscInt rmax; 64742c57489SHong Zhang 64842c57489SHong Zhang PetscFunctionBegin; 649ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 65083445d95SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 65183445d95SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 65283445d95SHong Zhang 653f8487c73SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 654b00a9115SJed Brown ierr = PetscNew(&ptap);CHKERRQ(ierr); 655b00a9115SJed Brown ierr = PetscNew(&merge);CHKERRQ(ierr); 6569f4155fbSHong Zhang ptap->merge = merge; 657f8487c73SHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 6586c6e00e0SHong Zhang 6596c6e00e0SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 660b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 661fe615159SHong Zhang 6626c6e00e0SHong Zhang /* get P_loc by taking all local rows of P */ 663a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 6646c6e00e0SHong Zhang 665a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 666a1a4e84aSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 6676c6e00e0SHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 6686c6e00e0SHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 66942c57489SHong Zhang 6702addb229SHong Zhang /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */ 6712addb229SHong Zhang /*--------------------------------------------------------------------------*/ 672854ce69bSBarry Smith ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 67382412ba7SHong Zhang api[0] = 0; 67483445d95SHong Zhang 67583445d95SHong Zhang /* create and initialize a linked list */ 676a914927fSHong Zhang ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 67783445d95SHong Zhang 678a914927fSHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */ 679d16ca5f0SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 680f4a743e1SHong Zhang current_space = free_space; 681f4a743e1SHong Zhang 682f4a743e1SHong Zhang for (i=0; i<am; i++) { 683f4a743e1SHong Zhang /* diagonal portion of A */ 684ded4f561SHong Zhang nzi = adi[i+1] - adi[i]; 68577584fe6SHong Zhang aj = ad->j + adi[i]; 686ded4f561SHong Zhang for (j=0; j<nzi; j++) { 68777584fe6SHong Zhang row = aj[j]; 68882412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 689ded4f561SHong Zhang Jptr = pj_loc + pi_loc[row]; 69083445d95SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 691a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 692f4a743e1SHong Zhang } 693f4a743e1SHong Zhang /* off-diagonal portion of A */ 694ded4f561SHong Zhang nzi = aoi[i+1] - aoi[i]; 69577584fe6SHong Zhang aj = ao->j + aoi[i]; 696ded4f561SHong Zhang for (j=0; j<nzi; j++) { 69777584fe6SHong Zhang row = aj[j]; 69882412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 699ded4f561SHong Zhang Jptr = pj_oth + pi_oth[row]; 700a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 701f4a743e1SHong Zhang } 702a914927fSHong Zhang apnz = lnk[0]; 70382412ba7SHong Zhang api[i+1] = api[i] + apnz; 70477584fe6SHong Zhang if (ap_rmax < apnz) ap_rmax = apnz; 705f4a743e1SHong Zhang 70683445d95SHong Zhang /* if free space is not available, double the total space in the list */ 70782412ba7SHong Zhang if (current_space->local_remaining<apnz) { 7082ba1acd5SHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 709f2b054eeSHong Zhang nspacedouble++; 710f4a743e1SHong Zhang } 711f4a743e1SHong Zhang 712f4a743e1SHong Zhang /* Copy data into free space, then initialize lnk */ 713a914927fSHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 7142205254eSKarl Rupp 71582412ba7SHong Zhang current_space->array += apnz; 71682412ba7SHong Zhang current_space->local_used += apnz; 71782412ba7SHong Zhang current_space->local_remaining -= apnz; 718f4a743e1SHong Zhang } 719a914927fSHong Zhang 72082412ba7SHong Zhang /* Allocate space for apj, initialize apj, and */ 721f4a743e1SHong Zhang /* destroy list of free space and other temporary array(s) */ 722854ce69bSBarry Smith ierr = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr); 72377584fe6SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 724118727c9SMark F. Adams afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1); 725d16ca5f0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 726f4a743e1SHong Zhang 7272addb229SHong Zhang /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/ 7282addb229SHong Zhang /*-----------------------------------------------------------------*/ 729de0260b3SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 73042c57489SHong Zhang 731ded4f561SHong Zhang /* then, compute symbolic Co = (p->B)^T*AP */ 732d0f46423SBarry Smith pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 73383445d95SHong Zhang >= (num of nonzero rows of C_seq) - pn */ 734854ce69bSBarry Smith ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 735de0260b3SHong Zhang coi[0] = 0; 73642c57489SHong Zhang 737d16ca5f0SHong Zhang /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */ 738d16ca5f0SHong Zhang nnz = fill*(poti[pon] + api[am]); 73922d28d08SBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 74042c57489SHong Zhang current_space = free_space; 74142c57489SHong Zhang 742de0260b3SHong Zhang for (i=0; i<pon; i++) { 74382412ba7SHong Zhang pnz = poti[i+1] - poti[i]; 74477584fe6SHong Zhang ptJ = potj + poti[i]; 74577584fe6SHong Zhang for (j=0; j<pnz; j++) { 74677584fe6SHong Zhang row = ptJ[j]; /* row of AP == col of Pot */ 74782412ba7SHong Zhang apnz = api[row+1] - api[row]; 748ded4f561SHong Zhang Jptr = apj + api[row]; 74983445d95SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 750a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 75142c57489SHong Zhang } 752a914927fSHong Zhang nnz = lnk[0]; 75342c57489SHong Zhang 75483445d95SHong Zhang /* If free space is not available, double the total space in the list */ 755ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 7564238b7adSHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 757d16ca5f0SHong Zhang nspacedouble++; 75842c57489SHong Zhang } 75942c57489SHong Zhang 76042c57489SHong Zhang /* Copy data into free space, and zero out denserows */ 761a914927fSHong Zhang ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 7622205254eSKarl Rupp 763ded4f561SHong Zhang current_space->array += nnz; 764ded4f561SHong Zhang current_space->local_used += nnz; 765ded4f561SHong Zhang current_space->local_remaining -= nnz; 7662205254eSKarl Rupp 767ded4f561SHong Zhang coi[i+1] = coi[i] + nnz; 76842c57489SHong Zhang } 7696b911d16SHong Zhang 7706b911d16SHong Zhang ierr = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr); 771a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 772118727c9SMark F. Adams afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1); 773d16ca5f0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 774de0260b3SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 77542c57489SHong Zhang 7766b911d16SHong Zhang /* (3) send j-array (coj) of Co to other processors */ 7776b911d16SHong Zhang /*--------------------------------------------------*/ 7786b911d16SHong Zhang ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr); 7796b911d16SHong Zhang len_s = merge->len_s; 7806b911d16SHong Zhang merge->nsend = 0; 7816b911d16SHong Zhang 7826b911d16SHong Zhang 78342c57489SHong Zhang /* determine row ownership */ 78426283091SBarry Smith ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 7857a2fc3feSBarry Smith merge->rowmap->n = pn; 7867a2fc3feSBarry Smith merge->rowmap->bs = 1; 7872205254eSKarl Rupp 78826283091SBarry Smith ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 7897a2fc3feSBarry Smith owners = merge->rowmap->range; 79042c57489SHong Zhang 79142c57489SHong Zhang /* determine the number of messages to send, their lengths */ 792dcca6d9dSJed Brown ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr); 79383445d95SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 794854ce69bSBarry Smith ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 795de0260b3SHong Zhang 79683445d95SHong Zhang proc = 0; 797de0260b3SHong Zhang for (i=0; i<pon; i++) { 798de0260b3SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 7992addb229SHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 800ee6e4b5bSHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 80183445d95SHong Zhang } 802de0260b3SHong Zhang 8032addb229SHong Zhang len = 0; /* max length of buf_si[], see (4) */ 804de0260b3SHong Zhang owners_co[0] = 0; 80542c57489SHong Zhang for (proc=0; proc<size; proc++) { 806de0260b3SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 807ee6e4b5bSHong Zhang if (len_s[proc]) { 80842c57489SHong Zhang merge->nsend++; 8092addb229SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 81042c57489SHong Zhang len += len_si[proc]; 81142c57489SHong Zhang } 81242c57489SHong Zhang } 81342c57489SHong Zhang 814ded4f561SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 8150298fd71SBarry Smith ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 81642c57489SHong Zhang ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 81742c57489SHong Zhang 818ded4f561SHong Zhang /* post the Irecv and Isend of coj */ 819529bc5dcSHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 820529bc5dcSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 821854ce69bSBarry Smith ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 82242c57489SHong Zhang for (proc=0, k=0; proc<size; proc++) { 82342c57489SHong Zhang if (!len_s[proc]) continue; 824de0260b3SHong Zhang i = owners_co[proc]; 825529bc5dcSHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 82642c57489SHong Zhang k++; 82742c57489SHong Zhang } 82842c57489SHong Zhang 829ded4f561SHong Zhang /* receives and sends of coj are complete */ 83077584fe6SHong Zhang for (i=0; i<merge->nrecv; i++) { 831c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 832ded4f561SHong Zhang } 833ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 8340c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 83582412ba7SHong Zhang 8366b911d16SHong Zhang /* (4) send and recv coi */ 8376b911d16SHong Zhang /*-----------------------*/ 838529bc5dcSHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 839529bc5dcSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 840854ce69bSBarry Smith ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 84142c57489SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 84242c57489SHong Zhang for (proc=0,k=0; proc<size; proc++) { 84342c57489SHong Zhang if (!len_s[proc]) continue; 84442c57489SHong Zhang /* form outgoing message for i-structure: 84542c57489SHong Zhang buf_si[0]: nrows to be sent 84642c57489SHong Zhang [1:nrows]: row index (global) 84742c57489SHong Zhang [nrows+1:2*nrows+1]: i-structure index 84842c57489SHong Zhang */ 84942c57489SHong Zhang /*-------------------------------------------*/ 8502addb229SHong Zhang nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 85142c57489SHong Zhang buf_si_i = buf_si + nrows+1; 85242c57489SHong Zhang buf_si[0] = nrows; 85342c57489SHong Zhang buf_si_i[0] = 0; 85442c57489SHong Zhang nrows = 0; 855de0260b3SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 856ded4f561SHong Zhang nzi = coi[i+1] - coi[i]; 857ded4f561SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 858de0260b3SHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 85942c57489SHong Zhang nrows++; 86042c57489SHong Zhang } 861529bc5dcSHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 86242c57489SHong Zhang k++; 86342c57489SHong Zhang buf_si += len_si[proc]; 86442c57489SHong Zhang } 865ded4f561SHong Zhang i = merge->nrecv; 866ded4f561SHong Zhang while (i--) { 867c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 868ded4f561SHong Zhang } 869ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 8700c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 871a914927fSHong Zhang 87224ecddacSHong Zhang ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr); 87342c57489SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 874ded4f561SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 87542c57489SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 87642c57489SHong Zhang 8776b911d16SHong Zhang /* (5) compute the local portion of C (mpi mat) */ 8786b911d16SHong Zhang /*----------------------------------------------*/ 879ded4f561SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 880ded4f561SHong Zhang 88124ecddacSHong Zhang /* allocate pti array and free space for accumulating nonzero column info */ 882854ce69bSBarry Smith ierr = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr); 88324ecddacSHong Zhang pti[0] = 0; 88442c57489SHong Zhang 885d16ca5f0SHong Zhang /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 886d16ca5f0SHong Zhang nnz = fill*(pi_loc[pm] + api[am]); 88722d28d08SBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 88842c57489SHong Zhang current_space = free_space; 88942c57489SHong Zhang 890dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 89142c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { 89242c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 89342c57489SHong Zhang nrows = *buf_ri_k[k]; 89442c57489SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 89542c57489SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 89642c57489SHong Zhang } 89742c57489SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 89808cb4532SHong Zhang rmax = 0; 89942c57489SHong Zhang for (i=0; i<pn; i++) { 900ded4f561SHong Zhang /* add pdt[i,:]*AP into lnk */ 901ded4f561SHong Zhang pnz = pdti[i+1] - pdti[i]; 90277584fe6SHong Zhang ptJ = pdtj + pdti[i]; 90377584fe6SHong Zhang for (j=0; j<pnz; j++) { 90477584fe6SHong Zhang row = ptJ[j]; /* row of AP == col of Pt */ 905ded4f561SHong Zhang apnz = api[row+1] - api[row]; 906ded4f561SHong Zhang Jptr = apj + api[row]; 907ded4f561SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 908a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 909ded4f561SHong Zhang } 910a914927fSHong Zhang 91142c57489SHong Zhang /* add received col data into lnk */ 91242c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 91342c57489SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 914ded4f561SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 915ded4f561SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 916a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 91742c57489SHong Zhang nextrow[k]++; nextci[k]++; 91842c57489SHong Zhang } 91942c57489SHong Zhang } 920a914927fSHong Zhang nnz = lnk[0]; 92142c57489SHong Zhang 92242c57489SHong Zhang /* if free space is not available, make more free space */ 923ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 9244238b7adSHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 925d16ca5f0SHong Zhang nspacedouble++; 92642c57489SHong Zhang } 92742c57489SHong Zhang /* copy data into free space, then initialize lnk */ 928a914927fSHong Zhang ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 929ded4f561SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 9302205254eSKarl Rupp 931ded4f561SHong Zhang current_space->array += nnz; 932ded4f561SHong Zhang current_space->local_used += nnz; 933ded4f561SHong Zhang current_space->local_remaining -= nnz; 9342205254eSKarl Rupp 93524ecddacSHong Zhang pti[i+1] = pti[i] + nnz; 93608cb4532SHong Zhang if (nnz > rmax) rmax = nnz; 93742c57489SHong Zhang } 938ded4f561SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 9390572522cSBarry Smith ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 94042c57489SHong Zhang 941854ce69bSBarry Smith ierr = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr); 94224ecddacSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr); 94324ecddacSHong Zhang afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1); 944d16ca5f0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 94542c57489SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 94642c57489SHong Zhang 9476b911d16SHong Zhang /* (6) create symbolic parallel matrix Cmpi */ 9486b911d16SHong Zhang /*------------------------------------------*/ 94977584fe6SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 95077584fe6SHong Zhang ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 95133d57670SJed Brown ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 95277584fe6SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 95377584fe6SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 95442c57489SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 955a2f3521dSMark F. Adams 956b091e509SHong Zhang merge->bi = pti; /* Cseq->i */ 957b091e509SHong Zhang merge->bj = ptj; /* Cseq->j */ 958b091e509SHong Zhang merge->coi = coi; /* Co->i */ 959b091e509SHong Zhang merge->coj = coj; /* Co->j */ 96042c57489SHong Zhang merge->buf_ri = buf_ri; 96142c57489SHong Zhang merge->buf_rj = buf_rj; 962de0260b3SHong Zhang merge->owners_co = owners_co; 96377584fe6SHong Zhang merge->destroy = Cmpi->ops->destroy; 96442c57489SHong Zhang 96577584fe6SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 96677584fe6SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 967f8487c73SHong Zhang c->ptap = ptap; 96877584fe6SHong Zhang ptap->api = api; 96977584fe6SHong Zhang ptap->apj = apj; 970d6ab1dc0SHong Zhang ptap->rmax = ap_rmax; 9718403a639SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 9728403a639SHong Zhang ptap->destroy = Cmpi->ops->destroy; 9738403a639SHong Zhang 9748403a639SHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 9758403a639SHong Zhang Cmpi->assembled = PETSC_FALSE; 9768403a639SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 9778403a639SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 9788403a639SHong Zhang Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap; 97977584fe6SHong Zhang *C = Cmpi; 980414894bdSHong Zhang 981414894bdSHong Zhang /* flag 'scalable' determines which implementations to be used: 982414894bdSHong Zhang 0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa; 983414894bdSHong Zhang 1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */ 984414894bdSHong Zhang /* set default scalable */ 985da0a95b2SSatish Balay ptap->scalable = PETSC_FALSE; /* PETSC_TRUE; */ 9862205254eSKarl Rupp 9870298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr); 988414894bdSHong Zhang if (!ptap->scalable) { /* Do dense axpy */ 9891795a4d1SJed Brown ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 990414894bdSHong Zhang } else { 9911795a4d1SJed Brown ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr); 992414894bdSHong Zhang } 993414894bdSHong Zhang 994f2b054eeSHong Zhang #if defined(PETSC_USE_INFO) 99524ecddacSHong Zhang if (pti[pn] != 0) { 99657622a8eSBarry Smith ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 99757622a8eSBarry Smith ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 998f2b054eeSHong Zhang } else { 99977584fe6SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1000f2b054eeSHong Zhang } 1001f2b054eeSHong Zhang #endif 100242c57489SHong Zhang PetscFunctionReturn(0); 100342c57489SHong Zhang } 100442c57489SHong Zhang 100542c57489SHong Zhang #undef __FUNCT__ 10068403a639SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap" 10078403a639SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,Mat C) 100842c57489SHong Zhang { 100942c57489SHong Zhang PetscErrorCode ierr; 1010f8487c73SHong Zhang Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 101142c57489SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 1012de0260b3SHong Zhang Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 101342c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 1014f8487c73SHong Zhang Mat_PtAPMPI *ptap; 10159f4155fbSHong Zhang Mat_Merge_SeqsToMPI *merge; 10161d633602SHong Zhang PetscInt *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp; 101782412ba7SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj; 101882412ba7SHong Zhang PetscInt i,j,k,anz,pnz,apnz,nextap,row,*cj; 1019e900a757SHong Zhang MatScalar *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp; 1020d0f46423SBarry Smith PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1021ce94432eSBarry Smith MPI_Comm comm; 102242c57489SHong Zhang PetscMPIInt size,rank,taga,*len_s; 102382412ba7SHong Zhang PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 102442c57489SHong Zhang PetscInt **buf_ri,**buf_rj; 102550f5bf86SHong Zhang PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 102642c57489SHong Zhang MPI_Request *s_waits,*r_waits; 102742c57489SHong Zhang MPI_Status *status; 102882412ba7SHong Zhang MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 102982412ba7SHong Zhang PetscInt *api,*apj,*coi,*coj; 1030d0f46423SBarry Smith PetscInt *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend; 10316ce94e8fSHong Zhang PetscBool scalable; 103238571addSHong Zhang #if defined(PTAP_PROFILE) 1033e497d3c8SHong Zhang PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2; 103438571addSHong Zhang #endif 103542c57489SHong Zhang 103642c57489SHong Zhang PetscFunctionBegin; 1037ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 103838571addSHong Zhang #if defined(PTAP_PROFILE) 10398563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 104038571addSHong Zhang #endif 104142c57489SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 104242c57489SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 104342c57489SHong Zhang 1044f8487c73SHong Zhang ptap = c->ptap; 1045ce94432eSBarry 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"); 1046f8487c73SHong Zhang merge = ptap->merge; 1047414894bdSHong Zhang apa = ptap->apa; 10486ce94e8fSHong Zhang scalable = ptap->scalable; 10496c6e00e0SHong Zhang 1050a1a4e84aSHong Zhang /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 10516b911d16SHong Zhang /*-----------------------------------------------------*/ 1052f8487c73SHong Zhang if (ptap->reuse == MAT_INITIAL_MATRIX) { 10539f4155fbSHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */ 1054f8487c73SHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 10556c6e00e0SHong Zhang } else { /* update numerical values of P_oth and P_loc */ 1056b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1057a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 10586c6e00e0SHong Zhang } 105938571addSHong Zhang #if defined(PTAP_PROFILE) 10608563dfccSBarry Smith ierr = PetscTime(&t1);CHKERRQ(ierr); 106105d62848SHong Zhang eP = t1-t0; 106238571addSHong Zhang #endif 106305d62848SHong Zhang /* 106405d62848SHong Zhang printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank, 106505d62848SHong Zhang a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N, 106605d62848SHong Zhang ptap->P_loc->rmap->N,ptap->P_loc->cmap->N, 106705d62848SHong Zhang ptap->P_oth->rmap->N,ptap->P_oth->cmap->N); 106805d62848SHong Zhang */ 1069f8487c73SHong Zhang 1070f9473cf7SHong Zhang /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */ 1071f9473cf7SHong Zhang /*--------------------------------------------------------------*/ 107242c57489SHong Zhang /* get data from symbolic products */ 1073a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1074a1a4e84aSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 1075a9d06231SHong Zhang pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a; 107642c57489SHong Zhang pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 107742c57489SHong Zhang 1078de0260b3SHong Zhang coi = merge->coi; coj = merge->coj; 10791795a4d1SJed Brown ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 1080de0260b3SHong Zhang 1081de0260b3SHong Zhang bi = merge->bi; bj = merge->bj; 10827a2fc3feSBarry Smith owners = merge->rowmap->range; 10831795a4d1SJed Brown ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); /* ba: Cseq->a */ 108442c57489SHong Zhang 1085a1a4e84aSHong Zhang api = ptap->api; apj = ptap->apj; 1086d9824c15SHong Zhang 10879516a85cSHong Zhang if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */ 1088b38d1a2dSJed Brown ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr); 10898cb82516SHong Zhang #if defined(PTAP_PROFILE) 109005d62848SHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 10918cb82516SHong Zhang #endif 1092a9d06231SHong Zhang for (i=0; i<am; i++) { 1093b091e509SHong Zhang #if defined(PTAP_PROFILE) 10948563dfccSBarry Smith ierr = PetscTime(&t2_0);CHKERRQ(ierr); 1095b091e509SHong Zhang #endif 1096a9d06231SHong Zhang /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */ 1097a9d06231SHong Zhang /*------------------------------------------------------------*/ 1098a9d06231SHong Zhang apJ = apj + api[i]; 1099a9d06231SHong Zhang 1100a9d06231SHong Zhang /* diagonal portion of A */ 1101a9d06231SHong Zhang anz = adi[i+1] - adi[i]; 1102a9d06231SHong Zhang adj = ad->j + adi[i]; 1103a9d06231SHong Zhang ada = ad->a + adi[i]; 1104a9d06231SHong Zhang for (j=0; j<anz; j++) { 1105a9d06231SHong Zhang row = adj[j]; 1106a9d06231SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 1107a9d06231SHong Zhang pj = pj_loc + pi_loc[row]; 1108a9d06231SHong Zhang pa = pa_loc + pi_loc[row]; 1109a9d06231SHong Zhang 1110a9d06231SHong Zhang /* perform dense axpy */ 1111e900a757SHong Zhang valtmp = ada[j]; 1112a9d06231SHong Zhang for (k=0; k<pnz; k++) { 1113e900a757SHong Zhang apa[pj[k]] += valtmp*pa[k]; 1114a9d06231SHong Zhang } 1115a9d06231SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1116a9d06231SHong Zhang } 1117a9d06231SHong Zhang 1118a9d06231SHong Zhang /* off-diagonal portion of A */ 1119a9d06231SHong Zhang anz = aoi[i+1] - aoi[i]; 1120a9d06231SHong Zhang aoj = ao->j + aoi[i]; 1121a9d06231SHong Zhang aoa = ao->a + aoi[i]; 1122a9d06231SHong Zhang for (j=0; j<anz; j++) { 1123a9d06231SHong Zhang row = aoj[j]; 1124a9d06231SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 1125a9d06231SHong Zhang pj = pj_oth + pi_oth[row]; 1126a9d06231SHong Zhang pa = pa_oth + pi_oth[row]; 1127a9d06231SHong Zhang 1128a9d06231SHong Zhang /* perform dense axpy */ 1129e900a757SHong Zhang valtmp = aoa[j]; 1130a9d06231SHong Zhang for (k=0; k<pnz; k++) { 1131e900a757SHong Zhang apa[pj[k]] += valtmp*pa[k]; 1132a9d06231SHong Zhang } 1133a9d06231SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1134a9d06231SHong Zhang } 1135b091e509SHong Zhang #if defined(PTAP_PROFILE) 11368563dfccSBarry Smith ierr = PetscTime(&t2_1);CHKERRQ(ierr); 1137b091e509SHong Zhang et2_AP += t2_1 - t2_0; 1138b091e509SHong Zhang #endif 1139a9d06231SHong Zhang 1140a9d06231SHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 1141a9d06231SHong Zhang /*--------------------------------------------------------------*/ 1142a9d06231SHong Zhang apnz = api[i+1] - api[i]; 1143a9d06231SHong Zhang /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */ 1144a9d06231SHong Zhang pnz = po->i[i+1] - po->i[i]; 1145e900a757SHong Zhang poJ = po->j + po->i[i]; 1146a9d06231SHong Zhang pA = po->a + po->i[i]; 1147a9d06231SHong Zhang for (j=0; j<pnz; j++) { 1148e900a757SHong Zhang row = poJ[j]; 1149e900a757SHong Zhang cnz = coi[row+1] - coi[row]; 1150e900a757SHong Zhang cj = coj + coi[row]; 1151e900a757SHong Zhang ca = coa + coi[row]; 1152a9d06231SHong Zhang /* perform dense axpy */ 1153e900a757SHong Zhang valtmp = pA[j]; 1154a9d06231SHong Zhang for (k=0; k<cnz; k++) { 1155e900a757SHong Zhang ca[k] += valtmp*apa[cj[k]]; 1156a9d06231SHong Zhang } 1157a9d06231SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1158a9d06231SHong Zhang } 1159a9d06231SHong Zhang /* put the value into Cd (diagonal part) */ 1160a9d06231SHong Zhang pnz = pd->i[i+1] - pd->i[i]; 1161e900a757SHong Zhang pdJ = pd->j + pd->i[i]; 1162a9d06231SHong Zhang pA = pd->a + pd->i[i]; 1163a9d06231SHong Zhang for (j=0; j<pnz; j++) { 1164e900a757SHong Zhang row = pdJ[j]; 1165e900a757SHong Zhang cnz = bi[row+1] - bi[row]; 1166e900a757SHong Zhang cj = bj + bi[row]; 1167e900a757SHong Zhang ca = ba + bi[row]; 1168a9d06231SHong Zhang /* perform dense axpy */ 1169e900a757SHong Zhang valtmp = pA[j]; 1170a9d06231SHong Zhang for (k=0; k<cnz; k++) { 1171e900a757SHong Zhang ca[k] += valtmp*apa[cj[k]]; 1172a9d06231SHong Zhang } 1173a9d06231SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1174a9d06231SHong Zhang } 11758403a639SHong Zhang 1176a9d06231SHong Zhang /* zero the current row of A*P */ 1177a9d06231SHong Zhang for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0; 1178b091e509SHong Zhang #if defined(PTAP_PROFILE) 11798563dfccSBarry Smith ierr = PetscTime(&t2_2);CHKERRQ(ierr); 118005d62848SHong Zhang ePtAP += t2_2 - t2_1; 1181b091e509SHong Zhang #endif 1182a9d06231SHong Zhang } 118338571addSHong Zhang } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */ 1184b38d1a2dSJed Brown ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr); 118538571addSHong Zhang /*-----------------------------------------------------------------------------------------*/ 1186a9d06231SHong Zhang pA=pa_loc; 118742c57489SHong Zhang for (i=0; i<am; i++) { 1188b091e509SHong Zhang #if defined(PTAP_PROFILE) 11898563dfccSBarry Smith ierr = PetscTime(&t2_0);CHKERRQ(ierr); 1190b091e509SHong Zhang #endif 1191f4a743e1SHong Zhang /* form i-th sparse row of A*P */ 119282412ba7SHong Zhang apnz = api[i+1] - api[i]; 119382412ba7SHong Zhang apJ = apj + api[i]; 119442c57489SHong Zhang /* diagonal portion of A */ 119582412ba7SHong Zhang anz = adi[i+1] - adi[i]; 11961d633602SHong Zhang adj = ad->j + adi[i]; 11971d633602SHong Zhang ada = ad->a + adi[i]; 119882412ba7SHong Zhang for (j=0; j<anz; j++) { 11991d633602SHong Zhang row = adj[j]; 120082412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 120182412ba7SHong Zhang pj = pj_loc + pi_loc[row]; 120282412ba7SHong Zhang pa = pa_loc + pi_loc[row]; 1203e900a757SHong Zhang valtmp = ada[j]; 1204f4a743e1SHong Zhang nextp = 0; 120582412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 120682412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 1207e900a757SHong Zhang apa[k] += valtmp*pa[nextp++]; 120842c57489SHong Zhang } 120942c57489SHong Zhang } 1210dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 121142c57489SHong Zhang } 121242c57489SHong Zhang /* off-diagonal portion of A */ 121382412ba7SHong Zhang anz = aoi[i+1] - aoi[i]; 12141d633602SHong Zhang aoj = ao->j + aoi[i]; 12151d633602SHong Zhang aoa = ao->a + aoi[i]; 121682412ba7SHong Zhang for (j=0; j<anz; j++) { 12171d633602SHong Zhang row = aoj[j]; 121882412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 121982412ba7SHong Zhang pj = pj_oth + pi_oth[row]; 122082412ba7SHong Zhang pa = pa_oth + pi_oth[row]; 1221e900a757SHong Zhang valtmp = aoa[j]; 1222f4a743e1SHong Zhang nextp = 0; 122382412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 122482412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 1225e900a757SHong Zhang apa[k] += valtmp*pa[nextp++]; 122642c57489SHong Zhang } 122742c57489SHong Zhang } 1228dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 122942c57489SHong Zhang } 1230b091e509SHong Zhang #if defined(PTAP_PROFILE) 12318563dfccSBarry Smith ierr = PetscTime(&t2_1);CHKERRQ(ierr); 1232b091e509SHong Zhang et2_AP += t2_1 - t2_0; 1233b091e509SHong Zhang #endif 123442c57489SHong Zhang 1235a9d06231SHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 1236a9d06231SHong Zhang /*--------------------------------------------------------------*/ 123782412ba7SHong Zhang pnz = pi_loc[i+1] - pi_loc[i]; 1238f9473cf7SHong Zhang pJ = pj_loc + pi_loc[i]; 123982412ba7SHong Zhang for (j=0; j<pnz; j++) { 124042c57489SHong Zhang nextap = 0; 1241f9473cf7SHong Zhang row = pJ[j]; /* global index */ 124282412ba7SHong Zhang if (row < pcstart || row >=pcend) { /* put the value into Co */ 1243e900a757SHong Zhang row = *poJ; 1244e900a757SHong Zhang cj = coj + coi[row]; 1245e900a757SHong Zhang ca = coa + coi[row]; poJ++; 124682412ba7SHong Zhang } else { /* put the value into Cd */ 1247e900a757SHong Zhang row = *pdJ; 1248e900a757SHong Zhang cj = bj + bi[row]; 1249e900a757SHong Zhang ca = ba + bi[row]; pdJ++; 125042c57489SHong Zhang } 1251e900a757SHong Zhang valtmp = pA[j]; 125282412ba7SHong Zhang for (k=0; nextap<apnz; k++) { 1253e900a757SHong Zhang if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++]; 125442c57489SHong Zhang } 1255dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 1256de0260b3SHong Zhang } 12570496f32cSHong Zhang pA += pnz; 125842c57489SHong Zhang /* zero the current row info for A*P */ 125982412ba7SHong Zhang ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr); 1260b091e509SHong Zhang #if defined(PTAP_PROFILE) 12618563dfccSBarry Smith ierr = PetscTime(&t2_2);CHKERRQ(ierr); 126205d62848SHong Zhang ePtAP += t2_2 - t2_1; 1263b091e509SHong Zhang #endif 126442c57489SHong Zhang } 1265d9824c15SHong Zhang } 126638571addSHong Zhang #if defined(PTAP_PROFILE) 12678563dfccSBarry Smith ierr = PetscTime(&t2);CHKERRQ(ierr); 126838571addSHong Zhang #endif 126942c57489SHong Zhang 1270a9d06231SHong Zhang /* 3) send and recv matrix values coa */ 1271a9d06231SHong Zhang /*------------------------------------*/ 127242c57489SHong Zhang buf_ri = merge->buf_ri; 127342c57489SHong Zhang buf_rj = merge->buf_rj; 127442c57489SHong Zhang len_s = merge->len_s; 1275fc42d0c8SSatish Balay ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 127642c57489SHong Zhang ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 127742c57489SHong Zhang 1278dcca6d9dSJed Brown ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 127942c57489SHong Zhang for (proc=0,k=0; proc<size; proc++) { 128042c57489SHong Zhang if (!len_s[proc]) continue; 1281de0260b3SHong Zhang i = merge->owners_co[proc]; 1282de0260b3SHong Zhang ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 128342c57489SHong Zhang k++; 128442c57489SHong Zhang } 12850c468ba9SBarry Smith if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 12860c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 128742c57489SHong Zhang 1288a9d06231SHong Zhang ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 128942c57489SHong Zhang ierr = PetscFree(r_waits);CHKERRQ(ierr); 129082412ba7SHong Zhang ierr = PetscFree(coa);CHKERRQ(ierr); 129138571addSHong Zhang #if defined(PTAP_PROFILE) 12928563dfccSBarry Smith ierr = PetscTime(&t3);CHKERRQ(ierr); 129338571addSHong Zhang #endif 129442c57489SHong Zhang 1295a9d06231SHong Zhang /* 4) insert local Cseq and received values into Cmpi */ 1296a9d06231SHong Zhang /*------------------------------------------------------*/ 1297dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 129842c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { 129942c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 130042c57489SHong Zhang nrows = *(buf_ri_k[k]); 130142c57489SHong Zhang nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 130282412ba7SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 130342c57489SHong Zhang } 130442c57489SHong Zhang 1305de0260b3SHong Zhang for (i=0; i<cm; i++) { 130682412ba7SHong Zhang row = owners[rank] + i; /* global row index of C_seq */ 130742c57489SHong Zhang bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1308de0260b3SHong Zhang ba_i = ba + bi[i]; 130982412ba7SHong Zhang bnz = bi[i+1] - bi[i]; 131042c57489SHong Zhang /* add received vals into ba */ 131142c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 131242c57489SHong Zhang /* i-th row */ 131342c57489SHong Zhang if (i == *nextrow[k]) { 131482412ba7SHong Zhang cnz = *(nextci[k]+1) - *nextci[k]; 131582412ba7SHong Zhang cj = buf_rj[k] + *(nextci[k]); 131682412ba7SHong Zhang ca = abuf_r[k] + *(nextci[k]); 131782412ba7SHong Zhang nextcj = 0; 131882412ba7SHong Zhang for (j=0; nextcj<cnz; j++) { 131982412ba7SHong Zhang if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 132082412ba7SHong Zhang ba_i[j] += ca[nextcj++]; 132142c57489SHong Zhang } 132242c57489SHong Zhang } 132382412ba7SHong Zhang nextrow[k]++; nextci[k]++; 1324c5af039cSHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 132542c57489SHong Zhang } 132642c57489SHong Zhang } 132782412ba7SHong Zhang ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 132842c57489SHong Zhang } 132942c57489SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 133042c57489SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 133142c57489SHong Zhang 1332de0260b3SHong Zhang ierr = PetscFree(ba);CHKERRQ(ierr); 1333c05d87d6SBarry Smith ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 133442c57489SHong Zhang ierr = PetscFree(abuf_r);CHKERRQ(ierr); 13350572522cSBarry Smith ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 133638571addSHong Zhang #if defined(PTAP_PROFILE) 13378563dfccSBarry Smith ierr = PetscTime(&t4);CHKERRQ(ierr); 13389516a85cSHong Zhang if (rank==1) { 133905d62848SHong 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); 13409516a85cSHong Zhang } 134138571addSHong Zhang #endif 134242c57489SHong Zhang PetscFunctionReturn(0); 134342c57489SHong Zhang } 1344