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 13*f671be37SHong 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() */ 320d3441aeSHong Zhang 3305d62848SHong Zhang ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr); 3405d62848SHong Zhang ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr); 35445158ffSHong Zhang 36445158ffSHong Zhang ierr = PetscFree(ptap->api);CHKERRQ(ierr); 37445158ffSHong Zhang ierr = PetscFree2(ptap->apj,ptap->apv);CHKERRQ(ierr); 380d3441aeSHong Zhang ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr); 392259aa2eSHong Zhang ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr); 402259aa2eSHong Zhang ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr); 410d3441aeSHong Zhang 42c5af039cSHong Zhang if (ptap->api) {ierr = PetscFree(ptap->api);CHKERRQ(ierr);} 43c5af039cSHong Zhang if (ptap->apj) {ierr = PetscFree(ptap->apj);CHKERRQ(ierr);} 44414894bdSHong Zhang if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);} 45c8b0795cSMark F. Adams if (merge) { 46cc6cb787SHong Zhang ierr = PetscFree(merge->id_r);CHKERRQ(ierr); 47cc6cb787SHong Zhang ierr = PetscFree(merge->len_s);CHKERRQ(ierr); 48cc6cb787SHong Zhang ierr = PetscFree(merge->len_r);CHKERRQ(ierr); 49cc6cb787SHong Zhang ierr = PetscFree(merge->bi);CHKERRQ(ierr); 50cc6cb787SHong Zhang ierr = PetscFree(merge->bj);CHKERRQ(ierr); 51c05d87d6SBarry Smith ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr); 52cc6cb787SHong Zhang ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr); 53c05d87d6SBarry Smith ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr); 54cc6cb787SHong Zhang ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr); 55445158ffSHong Zhang ierr = PetscFree(merge->coi);CHKERRQ(ierr); 56445158ffSHong Zhang ierr = PetscFree(merge->coj);CHKERRQ(ierr); 5705b42c5fSBarry Smith ierr = PetscFree(merge->owners_co);CHKERRQ(ierr); 586bf464f9SBarry Smith ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr); 59dce485f0SBarry Smith ierr = merge->destroy(A);CHKERRQ(ierr); 60f8487c73SHong Zhang ierr = PetscFree(ptap->merge);CHKERRQ(ierr); 61445158ffSHong Zhang } else { 62445158ffSHong Zhang ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr); 63bf0cc555SLisandro Dalcin } 64f8487c73SHong Zhang ierr = PetscFree(ptap);CHKERRQ(ierr); 65c8b0795cSMark F. Adams } 66cc6cb787SHong Zhang PetscFunctionReturn(0); 67cc6cb787SHong Zhang } 68cc6cb787SHong Zhang 69cc6cb787SHong Zhang #undef __FUNCT__ 70aa690a28SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP_old" 71aa690a28SHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP_old(Mat A, MatDuplicateOption op, Mat *M) 72f0c0a3a6SBarry Smith { 73cc6cb787SHong Zhang PetscErrorCode ierr; 7411a6856fSHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 7511a6856fSHong Zhang Mat_PtAPMPI *ptap = a->ptap; 7611a6856fSHong Zhang Mat_Merge_SeqsToMPI *merge = ptap->merge; 77cc6cb787SHong Zhang 78cc6cb787SHong Zhang PetscFunctionBegin; 79dce485f0SBarry Smith ierr = (*merge->duplicate)(A,op,M);CHKERRQ(ierr); 802205254eSKarl Rupp 8111a6856fSHong Zhang (*M)->ops->destroy = merge->destroy; 8211a6856fSHong Zhang (*M)->ops->duplicate = merge->duplicate; 83cc6cb787SHong Zhang PetscFunctionReturn(0); 84cc6cb787SHong Zhang } 85cc6cb787SHong Zhang 861a47ec54SHong Zhang 871a47ec54SHong Zhang #undef __FUNCT__ 88aa690a28SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP" 89aa690a28SHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M) 901a47ec54SHong Zhang { 911a47ec54SHong Zhang PetscErrorCode ierr; 921a47ec54SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 931a47ec54SHong Zhang Mat_PtAPMPI *ptap = a->ptap; 941a47ec54SHong Zhang 951a47ec54SHong Zhang PetscFunctionBegin; 961a47ec54SHong Zhang ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr); 971a47ec54SHong Zhang (*M)->ops->destroy = ptap->destroy; 981a47ec54SHong Zhang (*M)->ops->duplicate = ptap->duplicate; 991a47ec54SHong Zhang PetscFunctionReturn(0); 1001a47ec54SHong Zhang } 1011a47ec54SHong Zhang 10242c57489SHong Zhang #undef __FUNCT__ 103cf3ca8ceSHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ" 104cf3ca8ceSHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 10542c57489SHong Zhang { 10642c57489SHong Zhang PetscErrorCode ierr; 107aa690a28SHong Zhang PetscBool newalg=PETSC_TRUE; 10867a12041SHong Zhang MPI_Comm comm; 10942c57489SHong Zhang 11042c57489SHong Zhang PetscFunctionBegin; 11167a12041SHong Zhang /* check if matrix local sizes are compatible */ 11267a12041SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 11367a12041SHong Zhang if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) { 11467a12041SHong 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); 11567a12041SHong Zhang } 11667a12041SHong Zhang if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) { 11767a12041SHong 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); 11867a12041SHong Zhang } 11967a12041SHong Zhang 1200d3441aeSHong Zhang ierr = PetscOptionsGetBool(NULL,"-matptap_new",&newalg,NULL);CHKERRQ(ierr); 121cf3ca8ceSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 1223ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 1230d3441aeSHong Zhang if (newalg) { 124cf3ca8ceSHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr); 125aa690a28SHong Zhang } else { 126aa690a28SHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_old(A,P,fill,C);CHKERRQ(ierr); 1270d3441aeSHong Zhang } 1283ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 1297a7894deSKris Buschelman } 1303ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 1310d3441aeSHong Zhang if (newalg) { 132cf3ca8ceSHong Zhang ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr); 133aa690a28SHong Zhang } else { 134aa690a28SHong Zhang ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ_old(A,P,*C);CHKERRQ(ierr); 1350d3441aeSHong Zhang } 1363ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 13742c57489SHong Zhang PetscFunctionReturn(0); 13842c57489SHong Zhang } 13942c57489SHong Zhang 1408cb82516SHong Zhang /* requires array of size pN, thus nonscalable version */ 14142c57489SHong Zhang #undef __FUNCT__ 142aa690a28SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ" 143aa690a28SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 1440d3441aeSHong Zhang { 1450d3441aeSHong Zhang PetscErrorCode ierr; 1460d3441aeSHong Zhang Mat_PtAPMPI *ptap; 14767a12041SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 1482259aa2eSHong Zhang Mat_MPIAIJ *c; 1492259aa2eSHong Zhang MPI_Comm comm; 1502259aa2eSHong Zhang PetscMPIInt size,rank; 1518cb82516SHong Zhang #if defined(PTAP_PROFILE) 152445158ffSHong Zhang PetscLogDouble t0,t1,t11,t12,t2,t3,t4; 1538cb82516SHong Zhang #endif 15415a3b8e2SHong Zhang Mat Cmpi; 15515a3b8e2SHong Zhang PetscFreeSpaceList free_space=NULL,current_space=NULL; 15667a12041SHong Zhang PetscInt am=A->rmap->n,pm=P->rmap->n; 15767a12041SHong Zhang PetscInt *lnk,i,k,pnz,row,nnz; 158445158ffSHong Zhang PetscInt pN=P->cmap->N,pn=P->cmap->n; 15915a3b8e2SHong Zhang PetscBT lnkbt; 160*f671be37SHong Zhang PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 16115a3b8e2SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 16215a3b8e2SHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 16367a12041SHong Zhang PetscInt nzi,nspacedouble; 16415a3b8e2SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 16515a3b8e2SHong Zhang MPI_Request *swaits,*rwaits; 16615a3b8e2SHong Zhang MPI_Status *sstatus,rstatus; 167445158ffSHong Zhang PetscLayout rowmap; 168445158ffSHong Zhang PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 169*f671be37SHong Zhang PetscInt nsend; 170445158ffSHong Zhang PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 17167a12041SHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0; 17267a12041SHong Zhang PetscInt rmax,*aj,*ai,*pi; 17367a12041SHong Zhang Mat_SeqAIJ *p_loc,*p_oth,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*c_loc,*c_oth; 17467a12041SHong Zhang PetscScalar *apv; 175aa690a28SHong Zhang #if defined(PETSC_USE_INFO) 1768cb82516SHong Zhang PetscReal apfill; 177aa690a28SHong Zhang #endif 17867a12041SHong Zhang 17967a12041SHong Zhang PetscFunctionBegin; 18067a12041SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 18167a12041SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 18267a12041SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1838cb82516SHong Zhang #if defined(PTAP_PROFILE) 18480bb4639SHong Zhang ierr = PetscTime(&t0);CHKERRQ(ierr); 1858cb82516SHong Zhang #endif 1868cb82516SHong Zhang 18715a3b8e2SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 18815a3b8e2SHong Zhang ierr = PetscNew(&ptap);CHKERRQ(ierr); 18915a3b8e2SHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 19015a3b8e2SHong Zhang 19115a3b8e2SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 19215a3b8e2SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 19315a3b8e2SHong Zhang /* get P_loc by taking all local rows of P */ 19415a3b8e2SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 19515a3b8e2SHong Zhang 19667a12041SHong Zhang /* (0) compute Rd = Pd^T, Ro = Po^T */ 19767a12041SHong Zhang /* --------------------------------- */ 198de817e96SHong Zhang ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 199de817e96SHong Zhang ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 2008cb82516SHong Zhang #if defined(PTAP_PROFILE) 201445158ffSHong Zhang ierr = PetscTime(&t11);CHKERRQ(ierr); 2028cb82516SHong Zhang #endif 20315a3b8e2SHong Zhang 20467a12041SHong Zhang /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 20567a12041SHong Zhang /* ----------------------------------------------------------------- */ 20667a12041SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 20767a12041SHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 208445158ffSHong Zhang 2098cb82516SHong Zhang /* create and initialize a linked list - nonscalable version */ 210445158ffSHong Zhang ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 211445158ffSHong Zhang 2128cb82516SHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 21367a12041SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ad->i[am]+ao->i[am]+p_loc->i[pm])),&free_space);CHKERRQ(ierr); 214445158ffSHong Zhang current_space = free_space; 21567a12041SHong Zhang nspacedouble = 0; 21667a12041SHong Zhang 21767a12041SHong Zhang ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 21867a12041SHong Zhang api[0] = 0; 219445158ffSHong Zhang for (i=0; i<am; i++) { 22067a12041SHong Zhang /* diagonal portion: Ad[i,:]*P */ 22167a12041SHong Zhang ai = ad->i; pi = p_loc->i; 22267a12041SHong Zhang nzi = ai[i+1] - ai[i]; 22367a12041SHong Zhang aj = ad->j + ai[i]; 224445158ffSHong Zhang for (j=0; j<nzi; j++) { 225445158ffSHong Zhang row = aj[j]; 22667a12041SHong Zhang pnz = pi[row+1] - pi[row]; 22767a12041SHong Zhang Jptr = p_loc->j + pi[row]; 228445158ffSHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 229445158ffSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 230445158ffSHong Zhang } 23167a12041SHong Zhang /* off-diagonal portion: Ao[i,:]*P */ 23267a12041SHong Zhang ai = ao->i; pi = p_oth->i; 23367a12041SHong Zhang nzi = ai[i+1] - ai[i]; 23467a12041SHong Zhang aj = ao->j + ai[i]; 235445158ffSHong Zhang for (j=0; j<nzi; j++) { 236445158ffSHong Zhang row = aj[j]; 23767a12041SHong Zhang pnz = pi[row+1] - pi[row]; 23867a12041SHong Zhang Jptr = p_oth->j + pi[row]; 239445158ffSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 240445158ffSHong Zhang } 241445158ffSHong Zhang apnz = lnk[0]; 242445158ffSHong Zhang api[i+1] = api[i] + apnz; 243445158ffSHong Zhang if (ap_rmax < apnz) ap_rmax = apnz; 244445158ffSHong Zhang 245445158ffSHong Zhang /* if free space is not available, double the total space in the list */ 246445158ffSHong Zhang if (current_space->local_remaining<apnz) { 247445158ffSHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 248445158ffSHong Zhang nspacedouble++; 249445158ffSHong Zhang } 250445158ffSHong Zhang 251445158ffSHong Zhang /* Copy data into free space, then initialize lnk */ 252445158ffSHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 253445158ffSHong Zhang 254445158ffSHong Zhang current_space->array += apnz; 255445158ffSHong Zhang current_space->local_used += apnz; 256445158ffSHong Zhang current_space->local_remaining -= apnz; 257445158ffSHong Zhang } 258445158ffSHong Zhang /* Allocate space for apj, initialize apj, and */ 259445158ffSHong Zhang /* destroy list of free space and other temporary array(s) */ 260445158ffSHong Zhang ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 261445158ffSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 262445158ffSHong Zhang 263aa690a28SHong Zhang /* Create AP_loc for reuse */ 264445158ffSHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 265aa690a28SHong Zhang 266aa690a28SHong Zhang #if defined(PETSC_USE_INFO) 267aa690a28SHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 268aa690a28SHong Zhang ptap->AP_loc->info.mallocs = nspacedouble; 269aa690a28SHong Zhang ptap->AP_loc->info.fill_ratio_given = fill; 270aa690a28SHong Zhang ptap->AP_loc->info.fill_ratio_needed = apfill; 271aa690a28SHong Zhang 272aa690a28SHong Zhang if (api[am]) { 273aa690a28SHong Zhang ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr); 274aa690a28SHong Zhang ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 275aa690a28SHong Zhang } else { 276aa690a28SHong Zhang ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr); 277aa690a28SHong Zhang } 278aa690a28SHong Zhang #endif 279aa690a28SHong Zhang 2808cb82516SHong Zhang #if defined(PTAP_PROFILE) 281445158ffSHong Zhang ierr = PetscTime(&t12);CHKERRQ(ierr); 2828cb82516SHong Zhang #endif 28315a3b8e2SHong Zhang 2848cb82516SHong Zhang /* (2) compute symbolic C_loc = Rd*AP_loc, Co = Ro*AP_loc */ 2858cb82516SHong Zhang /* ------------------------------------------------------- */ 28667a12041SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 28767a12041SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 28867a12041SHong Zhang c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 28967a12041SHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 2908cb82516SHong Zhang #if defined(PTAP_PROFILE) 29180bb4639SHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 2928cb82516SHong Zhang #endif 29315a3b8e2SHong Zhang 29467a12041SHong Zhang /* (3) send coj of C_oth to other processors */ 29567a12041SHong Zhang /* ------------------------------------------ */ 29615a3b8e2SHong Zhang /* determine row ownership */ 297445158ffSHong Zhang ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 298445158ffSHong Zhang rowmap->n = pn; 299445158ffSHong Zhang rowmap->bs = 1; 300445158ffSHong Zhang ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 301445158ffSHong Zhang owners = rowmap->range; 30215a3b8e2SHong Zhang 30315a3b8e2SHong Zhang /* determine the number of messages to send, their lengths */ 3048cb82516SHong Zhang ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 3058cb82516SHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 30615a3b8e2SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 30715a3b8e2SHong Zhang 30867a12041SHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 30967a12041SHong Zhang coi = c_oth->i; coj = c_oth->j; 31067a12041SHong Zhang con = ptap->C_oth->rmap->n; 31115a3b8e2SHong Zhang proc = 0; 31267a12041SHong Zhang for (i=0; i<con; i++) { 31315a3b8e2SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 31415a3b8e2SHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 31515a3b8e2SHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 31615a3b8e2SHong Zhang } 31715a3b8e2SHong Zhang 31815a3b8e2SHong Zhang len = 0; /* max length of buf_si[], see (4) */ 31915a3b8e2SHong Zhang owners_co[0] = 0; 32067a12041SHong Zhang nsend = 0; 32115a3b8e2SHong Zhang for (proc=0; proc<size; proc++) { 32215a3b8e2SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 32315a3b8e2SHong Zhang if (len_s[proc]) { 324445158ffSHong Zhang nsend++; 32515a3b8e2SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 32615a3b8e2SHong Zhang len += len_si[proc]; 32715a3b8e2SHong Zhang } 32815a3b8e2SHong Zhang } 32915a3b8e2SHong Zhang 33015a3b8e2SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 331445158ffSHong Zhang ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 332445158ffSHong Zhang ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 33315a3b8e2SHong Zhang 33415a3b8e2SHong Zhang /* post the Irecv and Isend of coj */ 33515a3b8e2SHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 336445158ffSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 337445158ffSHong Zhang ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 33815a3b8e2SHong Zhang for (proc=0, k=0; proc<size; proc++) { 33915a3b8e2SHong Zhang if (!len_s[proc]) continue; 34015a3b8e2SHong Zhang i = owners_co[proc]; 34115a3b8e2SHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 34215a3b8e2SHong Zhang k++; 34315a3b8e2SHong Zhang } 34415a3b8e2SHong Zhang 34515a3b8e2SHong Zhang /* receives and sends of coj are complete */ 346445158ffSHong Zhang for (i=0; i<nrecv; i++) { 347445158ffSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 34815a3b8e2SHong Zhang } 34915a3b8e2SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 350445158ffSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 35115a3b8e2SHong 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 } 381445158ffSHong Zhang i = nrecv; 38215a3b8e2SHong Zhang while (i--) { 383445158ffSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 38415a3b8e2SHong Zhang } 38515a3b8e2SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 386445158ffSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 38715a3b8e2SHong Zhang 3888cb82516SHong Zhang ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 38915a3b8e2SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 39015a3b8e2SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 39115a3b8e2SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 3928cb82516SHong Zhang #if defined(PTAP_PROFILE) 39380bb4639SHong Zhang ierr = PetscTime(&t2);CHKERRQ(ierr); 3948cb82516SHong Zhang #endif 39567a12041SHong Zhang /* (5) compute the local portion of Cmpi */ 39667a12041SHong Zhang /* ------------------------------------------ */ 39715a3b8e2SHong Zhang /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 3988cb82516SHong Zhang ierr = PetscFreeSpaceGet(pN,&free_space);CHKERRQ(ierr); /* non-scalable version */ 39915a3b8e2SHong Zhang current_space = free_space; 40015a3b8e2SHong Zhang 401445158ffSHong Zhang ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 402445158ffSHong Zhang for (k=0; k<nrecv; k++) { 40315a3b8e2SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 40415a3b8e2SHong Zhang nrows = *buf_ri_k[k]; 40515a3b8e2SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 40615a3b8e2SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 40715a3b8e2SHong Zhang } 40815a3b8e2SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 409445158ffSHong Zhang 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; 461445158ffSHong Zhang ptap->api = api; 462445158ffSHong Zhang ptap->apj = apj; 463445158ffSHong Zhang ptap->apv = apv; 46415a3b8e2SHong Zhang ptap->rmax = ap_rmax; 4651a47ec54SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 4661a47ec54SHong Zhang ptap->destroy = Cmpi->ops->destroy; 46715a3b8e2SHong Zhang 4688cb82516SHong Zhang /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ_new() */ 4698cb82516SHong Zhang ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 4702259aa2eSHong Zhang 4711a47ec54SHong Zhang 4721a47ec54SHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 4731a47ec54SHong Zhang Cmpi->assembled = PETSC_FALSE; 4741a47ec54SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 475aa690a28SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 476aa690a28SHong Zhang Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ; 4771a47ec54SHong Zhang *C = Cmpi; 4781a47ec54SHong Zhang 4798cb82516SHong Zhang #if defined(PTAP_PROFILE) 48080bb4639SHong Zhang ierr = PetscTime(&t4);CHKERRQ(ierr); 481e9c1f85fSHong Zhang if (rank == 1) { 482445158ffSHong 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); 483e9c1f85fSHong Zhang } 4848cb82516SHong Zhang #endif 4850d3441aeSHong Zhang PetscFunctionReturn(0); 4860d3441aeSHong Zhang } 4870d3441aeSHong Zhang 4880d3441aeSHong Zhang #undef __FUNCT__ 489aa690a28SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ" 490aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 4910d3441aeSHong Zhang { 4920d3441aeSHong Zhang PetscErrorCode ierr; 4930d3441aeSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 4940d3441aeSHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 4958cb82516SHong Zhang Mat_SeqAIJ *ap,*p_loc,*p_oth,*c_seq; 4960d3441aeSHong Zhang Mat_PtAPMPI *ptap = c->ptap; 4979ce11a7cSHong Zhang Mat AP_loc,C_loc,C_oth; 4980d3441aeSHong Zhang PetscInt i,rstart,rend,cm,ncols,row; 4998cb82516SHong Zhang PetscInt *api,*apj,am = A->rmap->n,j,col,apnz; 5008cb82516SHong Zhang PetscScalar *apa; 5010d3441aeSHong Zhang const PetscInt *cols; 5020d3441aeSHong Zhang const PetscScalar *vals; 5038cb82516SHong Zhang #if defined(PTAP_PROFILE) 5048cb82516SHong Zhang PetscMPIInt rank; 5058cb82516SHong Zhang MPI_Comm comm; 5060d3441aeSHong Zhang PetscLogDouble t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi; 5078cb82516SHong Zhang #endif 5080d3441aeSHong Zhang 5090d3441aeSHong Zhang PetscFunctionBegin; 5100d3441aeSHong Zhang ierr = MatZeroEntries(C);CHKERRQ(ierr); 5110d3441aeSHong Zhang 512e497d3c8SHong Zhang /* 1) get R = Pd^T,Ro = Po^T */ 5138cb82516SHong Zhang #if defined(PTAP_PROFILE) 5140d3441aeSHong Zhang ierr = PetscTime(&t0);CHKERRQ(ierr); 5158cb82516SHong Zhang #endif 516748c7196SHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 5170d3441aeSHong Zhang ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 5180d3441aeSHong Zhang ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 519748c7196SHong Zhang } 5208cb82516SHong Zhang #if defined(PTAP_PROFILE) 5210d3441aeSHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 5220d3441aeSHong Zhang eR = t1 - t0; 5238cb82516SHong Zhang #endif 5240d3441aeSHong Zhang 525e497d3c8SHong Zhang /* 2) get AP_loc */ 5260d3441aeSHong Zhang AP_loc = ptap->AP_loc; 5278cb82516SHong Zhang ap = (Mat_SeqAIJ*)AP_loc->data; 5280d3441aeSHong Zhang 529e497d3c8SHong Zhang /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 5300d3441aeSHong Zhang /*-----------------------------------------------------*/ 531748c7196SHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 532748c7196SHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 5330d3441aeSHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 5340d3441aeSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 535e497d3c8SHong Zhang } 5360d3441aeSHong Zhang 5378cb82516SHong Zhang /* 2-2) compute numeric A_loc*P - dominating part */ 5388cb82516SHong Zhang /* ---------------------------------------------- */ 5390d3441aeSHong Zhang /* get data from symbolic products */ 5408cb82516SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 5418cb82516SHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 5428cb82516SHong Zhang apa = ptap->apa; 5438cb82516SHong Zhang api = ptap->api; 5448cb82516SHong Zhang apj = ptap->apj; 545e497d3c8SHong Zhang for (i=0; i<am; i++) { 546445158ffSHong Zhang /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 547e497d3c8SHong Zhang AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa); 548e497d3c8SHong Zhang apnz = api[i+1] - api[i]; 549e497d3c8SHong Zhang for (j=0; j<apnz; j++) { 550e497d3c8SHong Zhang col = apj[j+api[i]]; 551e497d3c8SHong Zhang ap->a[j+ap->i[i]] = apa[col]; 552e497d3c8SHong Zhang apa[col] = 0.0; 5530d3441aeSHong Zhang } 554aa690a28SHong Zhang ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 5550d3441aeSHong Zhang } 5568cb82516SHong Zhang #if defined(PTAP_PROFILE) 5570d3441aeSHong Zhang ierr = PetscTime(&t2);CHKERRQ(ierr); 5580d3441aeSHong Zhang eAP = t2 - t1; 5598cb82516SHong Zhang #endif 5600d3441aeSHong Zhang 5618cb82516SHong Zhang /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 562445158ffSHong Zhang ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 563445158ffSHong Zhang ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 5649ce11a7cSHong Zhang C_loc = ptap->C_loc; 5659ce11a7cSHong Zhang C_oth = ptap->C_oth; 5668cb82516SHong Zhang #if defined(PTAP_PROFILE) 5670d3441aeSHong Zhang ierr = PetscTime(&t3);CHKERRQ(ierr); 5680d3441aeSHong Zhang eCseq = t3 - t2; 5698cb82516SHong Zhang #endif 5700d3441aeSHong Zhang 5710d3441aeSHong Zhang /* add C_loc and Co to to C */ 5720d3441aeSHong Zhang ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 5730d3441aeSHong Zhang 5740d3441aeSHong Zhang /* C_loc -> C */ 5750d3441aeSHong Zhang cm = C_loc->rmap->N; 5769ce11a7cSHong Zhang c_seq = (Mat_SeqAIJ*)C_loc->data; 5778cb82516SHong Zhang cols = c_seq->j; 5788cb82516SHong Zhang vals = c_seq->a; 5790d3441aeSHong Zhang for (i=0; i<cm; i++) { 5809ce11a7cSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 5810d3441aeSHong Zhang row = rstart + i; 5820d3441aeSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 5838cb82516SHong Zhang cols += ncols; vals += ncols; 5840d3441aeSHong Zhang } 5850d3441aeSHong Zhang 5860d3441aeSHong Zhang /* Co -> C, off-processor part */ 5879ce11a7cSHong Zhang cm = C_oth->rmap->N; 5889ce11a7cSHong Zhang c_seq = (Mat_SeqAIJ*)C_oth->data; 5898cb82516SHong Zhang cols = c_seq->j; 5908cb82516SHong Zhang vals = c_seq->a; 5919ce11a7cSHong Zhang for (i=0; i<cm; i++) { 5929ce11a7cSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 5930d3441aeSHong Zhang row = p->garray[i]; 5940d3441aeSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 5958cb82516SHong Zhang cols += ncols; vals += ncols; 5960d3441aeSHong Zhang } 5970d3441aeSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5980d3441aeSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5998cb82516SHong Zhang #if defined(PTAP_PROFILE) 6000d3441aeSHong Zhang ierr = PetscTime(&t4);CHKERRQ(ierr); 6010d3441aeSHong Zhang eCmpi = t4 - t3; 6020d3441aeSHong Zhang 6038cb82516SHong Zhang ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 6048cb82516SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 6050d3441aeSHong Zhang if (rank==1) { 6060d3441aeSHong 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); 6070d3441aeSHong Zhang } 6088cb82516SHong Zhang #endif 609748c7196SHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 6100d3441aeSHong Zhang PetscFunctionReturn(0); 6110d3441aeSHong Zhang } 6120d3441aeSHong Zhang 6130d3441aeSHong Zhang #undef __FUNCT__ 614aa690a28SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ_old" 615aa690a28SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_old(Mat A,Mat P,PetscReal fill,Mat *C) 61642c57489SHong Zhang { 61742c57489SHong Zhang PetscErrorCode ierr; 61877584fe6SHong Zhang Mat Cmpi; 619f8487c73SHong Zhang Mat_PtAPMPI *ptap; 6200298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 621f8487c73SHong Zhang Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 62242c57489SHong Zhang Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 62342c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 624ded4f561SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ; 62577584fe6SHong Zhang PetscInt *adi=ad->i,*aj,*aoi=ao->i,nnz; 626a914927fSHong Zhang PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 627d16ca5f0SHong Zhang PetscInt am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n; 62842c57489SHong Zhang PetscBT lnkbt; 629ce94432eSBarry Smith MPI_Comm comm; 630a914927fSHong Zhang PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0; 63142c57489SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 63242c57489SHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 63324ecddacSHong Zhang PetscInt nzi,*pti,*ptj; 63442c57489SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 635ded4f561SHong Zhang MPI_Request *swaits,*rwaits; 636ded4f561SHong Zhang MPI_Status *sstatus,rstatus; 63742c57489SHong Zhang Mat_Merge_SeqsToMPI *merge; 63877584fe6SHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0; 639d16ca5f0SHong Zhang PetscReal afill=1.0,afill_tmp; 640a914927fSHong Zhang PetscInt rmax; 64142c57489SHong Zhang 64242c57489SHong Zhang PetscFunctionBegin; 643ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 64483445d95SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 64583445d95SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 64683445d95SHong Zhang 647f8487c73SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 648b00a9115SJed Brown ierr = PetscNew(&ptap);CHKERRQ(ierr); 649b00a9115SJed Brown ierr = PetscNew(&merge);CHKERRQ(ierr); 6509f4155fbSHong Zhang ptap->merge = merge; 651f8487c73SHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 6526c6e00e0SHong Zhang 6536c6e00e0SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 654b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 655fe615159SHong Zhang 6566c6e00e0SHong Zhang /* get P_loc by taking all local rows of P */ 657a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 6586c6e00e0SHong Zhang 659a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 660a1a4e84aSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 6616c6e00e0SHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 6626c6e00e0SHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 66342c57489SHong Zhang 6642addb229SHong Zhang /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */ 6652addb229SHong Zhang /*--------------------------------------------------------------------------*/ 666854ce69bSBarry Smith ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 66782412ba7SHong Zhang api[0] = 0; 66883445d95SHong Zhang 66983445d95SHong Zhang /* create and initialize a linked list */ 670a914927fSHong Zhang ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 67183445d95SHong Zhang 672a914927fSHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */ 673d16ca5f0SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 674f4a743e1SHong Zhang current_space = free_space; 675f4a743e1SHong Zhang 676f4a743e1SHong Zhang for (i=0; i<am; i++) { 677f4a743e1SHong Zhang /* diagonal portion of A */ 678ded4f561SHong Zhang nzi = adi[i+1] - adi[i]; 67977584fe6SHong Zhang aj = ad->j + adi[i]; 680ded4f561SHong Zhang for (j=0; j<nzi; j++) { 68177584fe6SHong Zhang row = aj[j]; 68282412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 683ded4f561SHong Zhang Jptr = pj_loc + pi_loc[row]; 68483445d95SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 685a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 686f4a743e1SHong Zhang } 687f4a743e1SHong Zhang /* off-diagonal portion of A */ 688ded4f561SHong Zhang nzi = aoi[i+1] - aoi[i]; 68977584fe6SHong Zhang aj = ao->j + aoi[i]; 690ded4f561SHong Zhang for (j=0; j<nzi; j++) { 69177584fe6SHong Zhang row = aj[j]; 69282412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 693ded4f561SHong Zhang Jptr = pj_oth + pi_oth[row]; 694a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 695f4a743e1SHong Zhang } 696a914927fSHong Zhang apnz = lnk[0]; 69782412ba7SHong Zhang api[i+1] = api[i] + apnz; 69877584fe6SHong Zhang if (ap_rmax < apnz) ap_rmax = apnz; 699f4a743e1SHong Zhang 70083445d95SHong Zhang /* if free space is not available, double the total space in the list */ 70182412ba7SHong Zhang if (current_space->local_remaining<apnz) { 7022ba1acd5SHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 703f2b054eeSHong Zhang nspacedouble++; 704f4a743e1SHong Zhang } 705f4a743e1SHong Zhang 706f4a743e1SHong Zhang /* Copy data into free space, then initialize lnk */ 707a914927fSHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 7082205254eSKarl Rupp 70982412ba7SHong Zhang current_space->array += apnz; 71082412ba7SHong Zhang current_space->local_used += apnz; 71182412ba7SHong Zhang current_space->local_remaining -= apnz; 712f4a743e1SHong Zhang } 713a914927fSHong Zhang 71482412ba7SHong Zhang /* Allocate space for apj, initialize apj, and */ 715f4a743e1SHong Zhang /* destroy list of free space and other temporary array(s) */ 716854ce69bSBarry Smith ierr = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr); 71777584fe6SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 718118727c9SMark F. Adams afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1); 719d16ca5f0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 720f4a743e1SHong Zhang 7212addb229SHong Zhang /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/ 7222addb229SHong Zhang /*-----------------------------------------------------------------*/ 723de0260b3SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 72442c57489SHong Zhang 725ded4f561SHong Zhang /* then, compute symbolic Co = (p->B)^T*AP */ 726d0f46423SBarry Smith pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 72783445d95SHong Zhang >= (num of nonzero rows of C_seq) - pn */ 728854ce69bSBarry Smith ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 729de0260b3SHong Zhang coi[0] = 0; 73042c57489SHong Zhang 731d16ca5f0SHong Zhang /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */ 732d16ca5f0SHong Zhang nnz = fill*(poti[pon] + api[am]); 73322d28d08SBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 73442c57489SHong Zhang current_space = free_space; 73542c57489SHong Zhang 736de0260b3SHong Zhang for (i=0; i<pon; i++) { 73782412ba7SHong Zhang pnz = poti[i+1] - poti[i]; 73877584fe6SHong Zhang ptJ = potj + poti[i]; 73977584fe6SHong Zhang for (j=0; j<pnz; j++) { 74077584fe6SHong Zhang row = ptJ[j]; /* row of AP == col of Pot */ 74182412ba7SHong Zhang apnz = api[row+1] - api[row]; 742ded4f561SHong Zhang Jptr = apj + api[row]; 74383445d95SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 744a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 74542c57489SHong Zhang } 746a914927fSHong Zhang nnz = lnk[0]; 74742c57489SHong Zhang 74883445d95SHong Zhang /* If free space is not available, double the total space in the list */ 749ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 7504238b7adSHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 751d16ca5f0SHong Zhang nspacedouble++; 75242c57489SHong Zhang } 75342c57489SHong Zhang 75442c57489SHong Zhang /* Copy data into free space, and zero out denserows */ 755a914927fSHong Zhang ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 7562205254eSKarl Rupp 757ded4f561SHong Zhang current_space->array += nnz; 758ded4f561SHong Zhang current_space->local_used += nnz; 759ded4f561SHong Zhang current_space->local_remaining -= nnz; 7602205254eSKarl Rupp 761ded4f561SHong Zhang coi[i+1] = coi[i] + nnz; 76242c57489SHong Zhang } 7636b911d16SHong Zhang 7646b911d16SHong Zhang ierr = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr); 765a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 766118727c9SMark F. Adams afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1); 767d16ca5f0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 768de0260b3SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 76942c57489SHong Zhang 7706b911d16SHong Zhang /* (3) send j-array (coj) of Co to other processors */ 7716b911d16SHong Zhang /*--------------------------------------------------*/ 7726b911d16SHong Zhang ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr); 7736b911d16SHong Zhang len_s = merge->len_s; 7746b911d16SHong Zhang merge->nsend = 0; 7756b911d16SHong Zhang 7766b911d16SHong Zhang 77742c57489SHong Zhang /* determine row ownership */ 77826283091SBarry Smith ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 7797a2fc3feSBarry Smith merge->rowmap->n = pn; 7807a2fc3feSBarry Smith merge->rowmap->bs = 1; 7812205254eSKarl Rupp 78226283091SBarry Smith ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 7837a2fc3feSBarry Smith owners = merge->rowmap->range; 78442c57489SHong Zhang 78542c57489SHong Zhang /* determine the number of messages to send, their lengths */ 786dcca6d9dSJed Brown ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr); 78783445d95SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 788854ce69bSBarry Smith ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 789de0260b3SHong Zhang 79083445d95SHong Zhang proc = 0; 791de0260b3SHong Zhang for (i=0; i<pon; i++) { 792de0260b3SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 7932addb229SHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 794ee6e4b5bSHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 79583445d95SHong Zhang } 796de0260b3SHong Zhang 7972addb229SHong Zhang len = 0; /* max length of buf_si[], see (4) */ 798de0260b3SHong Zhang owners_co[0] = 0; 79942c57489SHong Zhang for (proc=0; proc<size; proc++) { 800de0260b3SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 801ee6e4b5bSHong Zhang if (len_s[proc]) { 80242c57489SHong Zhang merge->nsend++; 8032addb229SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 80442c57489SHong Zhang len += len_si[proc]; 80542c57489SHong Zhang } 80642c57489SHong Zhang } 80742c57489SHong Zhang 808ded4f561SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 8090298fd71SBarry Smith ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 81042c57489SHong Zhang ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 81142c57489SHong Zhang 812ded4f561SHong Zhang /* post the Irecv and Isend of coj */ 813529bc5dcSHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 814529bc5dcSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 815854ce69bSBarry Smith ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 81642c57489SHong Zhang for (proc=0, k=0; proc<size; proc++) { 81742c57489SHong Zhang if (!len_s[proc]) continue; 818de0260b3SHong Zhang i = owners_co[proc]; 819529bc5dcSHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 82042c57489SHong Zhang k++; 82142c57489SHong Zhang } 82242c57489SHong Zhang 823ded4f561SHong Zhang /* receives and sends of coj are complete */ 82477584fe6SHong Zhang for (i=0; i<merge->nrecv; i++) { 825c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 826ded4f561SHong Zhang } 827ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 8280c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 82982412ba7SHong Zhang 8306b911d16SHong Zhang /* (4) send and recv coi */ 8316b911d16SHong Zhang /*-----------------------*/ 832529bc5dcSHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 833529bc5dcSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 834854ce69bSBarry Smith ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 83542c57489SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 83642c57489SHong Zhang for (proc=0,k=0; proc<size; proc++) { 83742c57489SHong Zhang if (!len_s[proc]) continue; 83842c57489SHong Zhang /* form outgoing message for i-structure: 83942c57489SHong Zhang buf_si[0]: nrows to be sent 84042c57489SHong Zhang [1:nrows]: row index (global) 84142c57489SHong Zhang [nrows+1:2*nrows+1]: i-structure index 84242c57489SHong Zhang */ 84342c57489SHong Zhang /*-------------------------------------------*/ 8442addb229SHong Zhang nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 84542c57489SHong Zhang buf_si_i = buf_si + nrows+1; 84642c57489SHong Zhang buf_si[0] = nrows; 84742c57489SHong Zhang buf_si_i[0] = 0; 84842c57489SHong Zhang nrows = 0; 849de0260b3SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 850ded4f561SHong Zhang nzi = coi[i+1] - coi[i]; 851ded4f561SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 852de0260b3SHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 85342c57489SHong Zhang nrows++; 85442c57489SHong Zhang } 855529bc5dcSHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 85642c57489SHong Zhang k++; 85742c57489SHong Zhang buf_si += len_si[proc]; 85842c57489SHong Zhang } 859ded4f561SHong Zhang i = merge->nrecv; 860ded4f561SHong Zhang while (i--) { 861c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 862ded4f561SHong Zhang } 863ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 8640c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 865a914927fSHong Zhang 86624ecddacSHong Zhang ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr); 86742c57489SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 868ded4f561SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 86942c57489SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 87042c57489SHong Zhang 8716b911d16SHong Zhang /* (5) compute the local portion of C (mpi mat) */ 8726b911d16SHong Zhang /*----------------------------------------------*/ 873ded4f561SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 874ded4f561SHong Zhang 87524ecddacSHong Zhang /* allocate pti array and free space for accumulating nonzero column info */ 876854ce69bSBarry Smith ierr = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr); 87724ecddacSHong Zhang pti[0] = 0; 87842c57489SHong Zhang 879d16ca5f0SHong Zhang /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 880d16ca5f0SHong Zhang nnz = fill*(pi_loc[pm] + api[am]); 88122d28d08SBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 88242c57489SHong Zhang current_space = free_space; 88342c57489SHong Zhang 884dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 88542c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { 88642c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 88742c57489SHong Zhang nrows = *buf_ri_k[k]; 88842c57489SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 88942c57489SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 89042c57489SHong Zhang } 89142c57489SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 89208cb4532SHong Zhang rmax = 0; 89342c57489SHong Zhang for (i=0; i<pn; i++) { 894ded4f561SHong Zhang /* add pdt[i,:]*AP into lnk */ 895ded4f561SHong Zhang pnz = pdti[i+1] - pdti[i]; 89677584fe6SHong Zhang ptJ = pdtj + pdti[i]; 89777584fe6SHong Zhang for (j=0; j<pnz; j++) { 89877584fe6SHong Zhang row = ptJ[j]; /* row of AP == col of Pt */ 899ded4f561SHong Zhang apnz = api[row+1] - api[row]; 900ded4f561SHong Zhang Jptr = apj + api[row]; 901ded4f561SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 902a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 903ded4f561SHong Zhang } 904a914927fSHong Zhang 90542c57489SHong Zhang /* add received col data into lnk */ 90642c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 90742c57489SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 908ded4f561SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 909ded4f561SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 910a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 91142c57489SHong Zhang nextrow[k]++; nextci[k]++; 91242c57489SHong Zhang } 91342c57489SHong Zhang } 914a914927fSHong Zhang nnz = lnk[0]; 91542c57489SHong Zhang 91642c57489SHong Zhang /* if free space is not available, make more free space */ 917ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 9184238b7adSHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 919d16ca5f0SHong Zhang nspacedouble++; 92042c57489SHong Zhang } 92142c57489SHong Zhang /* copy data into free space, then initialize lnk */ 922a914927fSHong Zhang ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 923ded4f561SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 9242205254eSKarl Rupp 925ded4f561SHong Zhang current_space->array += nnz; 926ded4f561SHong Zhang current_space->local_used += nnz; 927ded4f561SHong Zhang current_space->local_remaining -= nnz; 9282205254eSKarl Rupp 92924ecddacSHong Zhang pti[i+1] = pti[i] + nnz; 93008cb4532SHong Zhang if (nnz > rmax) rmax = nnz; 93142c57489SHong Zhang } 932ded4f561SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 9330572522cSBarry Smith ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 93442c57489SHong Zhang 935854ce69bSBarry Smith ierr = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr); 93624ecddacSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr); 93724ecddacSHong Zhang afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1); 938d16ca5f0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 93942c57489SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 94042c57489SHong Zhang 9416b911d16SHong Zhang /* (6) create symbolic parallel matrix Cmpi */ 9426b911d16SHong Zhang /*------------------------------------------*/ 94377584fe6SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 94477584fe6SHong Zhang ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 94533d57670SJed Brown ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 94677584fe6SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 94777584fe6SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 94842c57489SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 949a2f3521dSMark F. Adams 950b091e509SHong Zhang merge->bi = pti; /* Cseq->i */ 951b091e509SHong Zhang merge->bj = ptj; /* Cseq->j */ 952b091e509SHong Zhang merge->coi = coi; /* Co->i */ 953b091e509SHong Zhang merge->coj = coj; /* Co->j */ 95442c57489SHong Zhang merge->buf_ri = buf_ri; 95542c57489SHong Zhang merge->buf_rj = buf_rj; 956de0260b3SHong Zhang merge->owners_co = owners_co; 95777584fe6SHong Zhang merge->destroy = Cmpi->ops->destroy; 95877584fe6SHong Zhang merge->duplicate = Cmpi->ops->duplicate; 959cc6cb787SHong Zhang 96077584fe6SHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 961a914927fSHong Zhang Cmpi->assembled = PETSC_FALSE; 96277584fe6SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 963aa690a28SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP_old; 964aa690a28SHong Zhang Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_old; 96542c57489SHong Zhang 96677584fe6SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 96777584fe6SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 968f8487c73SHong Zhang c->ptap = ptap; 96977584fe6SHong Zhang ptap->api = api; 97077584fe6SHong Zhang ptap->apj = apj; 971d6ab1dc0SHong Zhang ptap->rmax = ap_rmax; 97277584fe6SHong Zhang *C = Cmpi; 973414894bdSHong Zhang 974414894bdSHong Zhang /* flag 'scalable' determines which implementations to be used: 975414894bdSHong Zhang 0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa; 976414894bdSHong Zhang 1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */ 977414894bdSHong Zhang /* set default scalable */ 9789516a85cSHong Zhang ptap->scalable = PETSC_FALSE; //PETSC_TRUE; 9792205254eSKarl Rupp 9800298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr); 981414894bdSHong Zhang if (!ptap->scalable) { /* Do dense axpy */ 9821795a4d1SJed Brown ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 983414894bdSHong Zhang } else { 9841795a4d1SJed Brown ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr); 985414894bdSHong Zhang } 986414894bdSHong Zhang 987f2b054eeSHong Zhang #if defined(PETSC_USE_INFO) 98824ecddacSHong Zhang if (pti[pn] != 0) { 98957622a8eSBarry Smith ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 99057622a8eSBarry Smith ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 991f2b054eeSHong Zhang } else { 99277584fe6SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 993f2b054eeSHong Zhang } 994f2b054eeSHong Zhang #endif 99542c57489SHong Zhang PetscFunctionReturn(0); 99642c57489SHong Zhang } 99742c57489SHong Zhang 99842c57489SHong Zhang #undef __FUNCT__ 999aa690a28SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ_old" 1000aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_old(Mat A,Mat P,Mat C) 100142c57489SHong Zhang { 100242c57489SHong Zhang PetscErrorCode ierr; 1003f8487c73SHong Zhang Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 100442c57489SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 1005de0260b3SHong Zhang Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 100642c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 1007f8487c73SHong Zhang Mat_PtAPMPI *ptap; 10089f4155fbSHong Zhang Mat_Merge_SeqsToMPI *merge; 10091d633602SHong Zhang PetscInt *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp; 101082412ba7SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj; 101182412ba7SHong Zhang PetscInt i,j,k,anz,pnz,apnz,nextap,row,*cj; 1012e900a757SHong Zhang MatScalar *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp; 1013d0f46423SBarry Smith PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1014ce94432eSBarry Smith MPI_Comm comm; 101542c57489SHong Zhang PetscMPIInt size,rank,taga,*len_s; 101682412ba7SHong Zhang PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 101742c57489SHong Zhang PetscInt **buf_ri,**buf_rj; 101850f5bf86SHong Zhang PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 101942c57489SHong Zhang MPI_Request *s_waits,*r_waits; 102042c57489SHong Zhang MPI_Status *status; 102182412ba7SHong Zhang MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 102282412ba7SHong Zhang PetscInt *api,*apj,*coi,*coj; 1023d0f46423SBarry Smith PetscInt *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend; 10246ce94e8fSHong Zhang PetscBool scalable; 102538571addSHong Zhang #if defined(PTAP_PROFILE) 1026e497d3c8SHong Zhang PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2; 102738571addSHong Zhang #endif 102842c57489SHong Zhang 102942c57489SHong Zhang PetscFunctionBegin; 1030ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 103138571addSHong Zhang #if defined(PTAP_PROFILE) 10328563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 103338571addSHong Zhang #endif 103442c57489SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 103542c57489SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 103642c57489SHong Zhang 1037f8487c73SHong Zhang ptap = c->ptap; 1038ce94432eSBarry 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"); 1039f8487c73SHong Zhang merge = ptap->merge; 1040414894bdSHong Zhang apa = ptap->apa; 10416ce94e8fSHong Zhang scalable = ptap->scalable; 10426c6e00e0SHong Zhang 1043a1a4e84aSHong Zhang /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 10446b911d16SHong Zhang /*-----------------------------------------------------*/ 1045f8487c73SHong Zhang if (ptap->reuse == MAT_INITIAL_MATRIX) { 10469f4155fbSHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */ 1047f8487c73SHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 10486c6e00e0SHong Zhang } else { /* update numerical values of P_oth and P_loc */ 1049b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1050a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 10516c6e00e0SHong Zhang } 105238571addSHong Zhang #if defined(PTAP_PROFILE) 10538563dfccSBarry Smith ierr = PetscTime(&t1);CHKERRQ(ierr); 105405d62848SHong Zhang eP = t1-t0; 105538571addSHong Zhang #endif 105605d62848SHong Zhang /* 105705d62848SHong Zhang printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank, 105805d62848SHong Zhang a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N, 105905d62848SHong Zhang ptap->P_loc->rmap->N,ptap->P_loc->cmap->N, 106005d62848SHong Zhang ptap->P_oth->rmap->N,ptap->P_oth->cmap->N); 106105d62848SHong Zhang */ 1062f8487c73SHong Zhang 1063f9473cf7SHong Zhang /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */ 1064f9473cf7SHong Zhang /*--------------------------------------------------------------*/ 106542c57489SHong Zhang /* get data from symbolic products */ 1066a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1067a1a4e84aSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 1068a9d06231SHong Zhang pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a; 106942c57489SHong Zhang pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 107042c57489SHong Zhang 1071de0260b3SHong Zhang coi = merge->coi; coj = merge->coj; 10721795a4d1SJed Brown ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 1073de0260b3SHong Zhang 1074de0260b3SHong Zhang bi = merge->bi; bj = merge->bj; 10757a2fc3feSBarry Smith owners = merge->rowmap->range; 10761795a4d1SJed Brown ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); /* ba: Cseq->a */ 107742c57489SHong Zhang 1078a1a4e84aSHong Zhang api = ptap->api; apj = ptap->apj; 1079d9824c15SHong Zhang 10809516a85cSHong Zhang if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */ 1081b38d1a2dSJed Brown ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr); 108205d62848SHong Zhang #if 0 10830d3441aeSHong Zhang /* ------ 10x slower -------------- */ 10840d3441aeSHong Zhang /*==================================*/ 10859516a85cSHong Zhang Mat R = ptap->R; 10869516a85cSHong Zhang Mat_SeqAIJ *r = (Mat_SeqAIJ*)R->data; 10879516a85cSHong Zhang PetscInt *ri=r->i,*rj=r->j,rnz,arow,l,prow,pcol,pN=P->cmap->N; 10889516a85cSHong Zhang PetscScalar *ra=r->a,tmp,cdense[pN]; 10899516a85cSHong Zhang 10909516a85cSHong Zhang ierr = PetscMemzero(cdense,pN*sizeof(PetscScalar));CHKERRQ(ierr); 10919516a85cSHong Zhang for (i=0; i<cm; i++) { /* each row of C or R */ 10929516a85cSHong Zhang rnz = ri[i+1] - ri[i]; 10939516a85cSHong Zhang 10949516a85cSHong Zhang for (j=0; j<rnz; j++) { /* each nz of R */ 10959516a85cSHong Zhang arow = rj[ri[i] + j]; 10969516a85cSHong Zhang 10979516a85cSHong Zhang /* diagonal portion of A */ 10989516a85cSHong Zhang anz = ad->i[arow+1] - ad->i[arow]; 10999516a85cSHong Zhang for (k=0; k<anz; k++) { /* each nz of Ad */ 11009516a85cSHong Zhang tmp = ra[ri[i] + j]*ad->a[ad->i[arow] + k]; 11019516a85cSHong Zhang prow = ad->j[ad->i[arow] + k]; 11029516a85cSHong Zhang pnz = pi_loc[prow+1] - pi_loc[prow]; 11039516a85cSHong Zhang 11049516a85cSHong Zhang for (l=0; l<pnz; l++) { /* each nz of P_loc */ 11059516a85cSHong Zhang pcol = pj_loc[pi_loc[prow] + l]; 11069516a85cSHong Zhang cdense[pcol] += tmp*pa_loc[pi_loc[prow] + l]; 11079516a85cSHong Zhang } 11089516a85cSHong Zhang } 11099516a85cSHong Zhang 11109516a85cSHong Zhang /* off-diagonal portion of A */ 11119516a85cSHong Zhang anz = ao->i[arow+1] - ao->i[arow]; 11129516a85cSHong Zhang for (k=0; k<anz; k++) { /* each nz of Ao */ 11139516a85cSHong Zhang tmp = ra[ri[i] + j]*ao->a[ao->i[arow] + k]; 11149516a85cSHong Zhang prow = ao->j[ao->i[arow] + k]; 11159516a85cSHong Zhang pnz = pi_oth[prow+1] - pi_oth[prow]; 11169516a85cSHong Zhang 11179516a85cSHong Zhang for (l=0; l<pnz; l++) { /* each nz of P_oth */ 11189516a85cSHong Zhang pcol = pj_oth[pi_oth[prow] + l]; 11199516a85cSHong Zhang cdense[pcol] += tmp*pa_oth[pi_oth[prow] + l]; 11209516a85cSHong Zhang } 11219516a85cSHong Zhang } 11229516a85cSHong Zhang 11239516a85cSHong Zhang } //for (j=0; j<rnz; j++) 11249516a85cSHong Zhang 11259516a85cSHong Zhang /* copy cdense[] into ca; zero cdense[] */ 11269516a85cSHong Zhang cnz = bi[i+1] - bi[i]; 11279516a85cSHong Zhang cj = bj + bi[i]; 11289516a85cSHong Zhang ca = ba + bi[i]; 11299516a85cSHong Zhang for (j=0; j<cnz; j++) { 11309516a85cSHong Zhang ca[j] += cdense[cj[j]]; 11319516a85cSHong Zhang cdense[cj[j]] = 0.0; 11329516a85cSHong Zhang } 11339516a85cSHong Zhang #if 0 11349516a85cSHong Zhang if (rank == 0) { 11359516a85cSHong Zhang printf("[%d] row %d: ",rank,i); 11369516a85cSHong Zhang for (j=0; j<pN; j++) printf(" %g,",cdense[j]); 11379516a85cSHong Zhang printf("\n"); 11389516a85cSHong Zhang } 11399516a85cSHong Zhang for (j=0; j<pN; j++) cdense[j]=0.0; // zero cdnese[] 11409516a85cSHong Zhang #endif 11419516a85cSHong Zhang } //for (i=0; i<cm; i++) { 114205d62848SHong Zhang #endif 11439516a85cSHong Zhang 1144*f671be37SHong Zhang /* ========================================== */ 11458cb82516SHong Zhang #if defined(PTAP_PROFILE) 114605d62848SHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 11478cb82516SHong Zhang #endif 1148a9d06231SHong Zhang for (i=0; i<am; i++) { 1149b091e509SHong Zhang #if defined(PTAP_PROFILE) 11508563dfccSBarry Smith ierr = PetscTime(&t2_0);CHKERRQ(ierr); 1151b091e509SHong Zhang #endif 1152a9d06231SHong Zhang /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */ 1153a9d06231SHong Zhang /*------------------------------------------------------------*/ 1154a9d06231SHong Zhang apJ = apj + api[i]; 1155a9d06231SHong Zhang 1156a9d06231SHong Zhang /* diagonal portion of A */ 1157a9d06231SHong Zhang anz = adi[i+1] - adi[i]; 1158a9d06231SHong Zhang adj = ad->j + adi[i]; 1159a9d06231SHong Zhang ada = ad->a + adi[i]; 1160a9d06231SHong Zhang for (j=0; j<anz; j++) { 1161a9d06231SHong Zhang row = adj[j]; 1162a9d06231SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 1163a9d06231SHong Zhang pj = pj_loc + pi_loc[row]; 1164a9d06231SHong Zhang pa = pa_loc + pi_loc[row]; 1165a9d06231SHong Zhang 1166a9d06231SHong Zhang /* perform dense axpy */ 1167e900a757SHong Zhang valtmp = ada[j]; 1168a9d06231SHong Zhang for (k=0; k<pnz; k++) { 1169e900a757SHong Zhang apa[pj[k]] += valtmp*pa[k]; 1170a9d06231SHong Zhang } 1171a9d06231SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1172a9d06231SHong Zhang } 1173a9d06231SHong Zhang 1174a9d06231SHong Zhang /* off-diagonal portion of A */ 1175a9d06231SHong Zhang anz = aoi[i+1] - aoi[i]; 1176a9d06231SHong Zhang aoj = ao->j + aoi[i]; 1177a9d06231SHong Zhang aoa = ao->a + aoi[i]; 1178a9d06231SHong Zhang for (j=0; j<anz; j++) { 1179a9d06231SHong Zhang row = aoj[j]; 1180a9d06231SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 1181a9d06231SHong Zhang pj = pj_oth + pi_oth[row]; 1182a9d06231SHong Zhang pa = pa_oth + pi_oth[row]; 1183a9d06231SHong Zhang 1184a9d06231SHong Zhang /* perform dense axpy */ 1185e900a757SHong Zhang valtmp = aoa[j]; 1186a9d06231SHong Zhang for (k=0; k<pnz; k++) { 1187e900a757SHong Zhang apa[pj[k]] += valtmp*pa[k]; 1188a9d06231SHong Zhang } 1189a9d06231SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1190a9d06231SHong Zhang } 1191b091e509SHong Zhang #if defined(PTAP_PROFILE) 11928563dfccSBarry Smith ierr = PetscTime(&t2_1);CHKERRQ(ierr); 1193b091e509SHong Zhang et2_AP += t2_1 - t2_0; 1194b091e509SHong Zhang #endif 1195a9d06231SHong Zhang 1196a9d06231SHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 1197a9d06231SHong Zhang /*--------------------------------------------------------------*/ 1198a9d06231SHong Zhang apnz = api[i+1] - api[i]; 1199a9d06231SHong Zhang /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */ 1200a9d06231SHong Zhang pnz = po->i[i+1] - po->i[i]; 1201e900a757SHong Zhang poJ = po->j + po->i[i]; 1202a9d06231SHong Zhang pA = po->a + po->i[i]; 1203a9d06231SHong Zhang for (j=0; j<pnz; j++) { 1204e900a757SHong Zhang row = poJ[j]; 1205e900a757SHong Zhang cnz = coi[row+1] - coi[row]; 1206e900a757SHong Zhang cj = coj + coi[row]; 1207e900a757SHong Zhang ca = coa + coi[row]; 1208a9d06231SHong Zhang /* perform dense axpy */ 1209e900a757SHong Zhang valtmp = pA[j]; 1210a9d06231SHong Zhang for (k=0; k<cnz; k++) { 1211e900a757SHong Zhang ca[k] += valtmp*apa[cj[k]]; 1212a9d06231SHong Zhang } 1213a9d06231SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1214a9d06231SHong Zhang } 121505d62848SHong Zhang #if 1 1216a9d06231SHong Zhang /* put the value into Cd (diagonal part) */ 1217a9d06231SHong Zhang pnz = pd->i[i+1] - pd->i[i]; 1218e900a757SHong Zhang pdJ = pd->j + pd->i[i]; 1219a9d06231SHong Zhang pA = pd->a + pd->i[i]; 1220a9d06231SHong Zhang for (j=0; j<pnz; j++) { 1221e900a757SHong Zhang row = pdJ[j]; 1222e900a757SHong Zhang cnz = bi[row+1] - bi[row]; 1223e900a757SHong Zhang cj = bj + bi[row]; 1224e900a757SHong Zhang ca = ba + bi[row]; 1225a9d06231SHong Zhang /* perform dense axpy */ 1226e900a757SHong Zhang valtmp = pA[j]; 1227a9d06231SHong Zhang for (k=0; k<cnz; k++) { 1228e900a757SHong Zhang ca[k] += valtmp*apa[cj[k]]; 1229a9d06231SHong Zhang } 1230a9d06231SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1231a9d06231SHong Zhang } 12329516a85cSHong Zhang #endif 1233a9d06231SHong Zhang /* zero the current row of A*P */ 1234a9d06231SHong Zhang for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0; 1235b091e509SHong Zhang #if defined(PTAP_PROFILE) 12368563dfccSBarry Smith ierr = PetscTime(&t2_2);CHKERRQ(ierr); 123705d62848SHong Zhang ePtAP += t2_2 - t2_1; 1238b091e509SHong Zhang #endif 1239a9d06231SHong Zhang } 124038571addSHong Zhang } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */ 1241b38d1a2dSJed Brown ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr); 124238571addSHong Zhang /*-----------------------------------------------------------------------------------------*/ 1243a9d06231SHong Zhang pA=pa_loc; 124442c57489SHong Zhang for (i=0; i<am; i++) { 1245b091e509SHong Zhang #if defined(PTAP_PROFILE) 12468563dfccSBarry Smith ierr = PetscTime(&t2_0);CHKERRQ(ierr); 1247b091e509SHong Zhang #endif 1248f4a743e1SHong Zhang /* form i-th sparse row of A*P */ 124982412ba7SHong Zhang apnz = api[i+1] - api[i]; 125082412ba7SHong Zhang apJ = apj + api[i]; 125142c57489SHong Zhang /* diagonal portion of A */ 125282412ba7SHong Zhang anz = adi[i+1] - adi[i]; 12531d633602SHong Zhang adj = ad->j + adi[i]; 12541d633602SHong Zhang ada = ad->a + adi[i]; 125582412ba7SHong Zhang for (j=0; j<anz; j++) { 12561d633602SHong Zhang row = adj[j]; 125782412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 125882412ba7SHong Zhang pj = pj_loc + pi_loc[row]; 125982412ba7SHong Zhang pa = pa_loc + pi_loc[row]; 1260e900a757SHong Zhang valtmp = ada[j]; 1261f4a743e1SHong Zhang nextp = 0; 126282412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 126382412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 1264e900a757SHong Zhang apa[k] += valtmp*pa[nextp++]; 126542c57489SHong Zhang } 126642c57489SHong Zhang } 1267dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 126842c57489SHong Zhang } 126942c57489SHong Zhang /* off-diagonal portion of A */ 127082412ba7SHong Zhang anz = aoi[i+1] - aoi[i]; 12711d633602SHong Zhang aoj = ao->j + aoi[i]; 12721d633602SHong Zhang aoa = ao->a + aoi[i]; 127382412ba7SHong Zhang for (j=0; j<anz; j++) { 12741d633602SHong Zhang row = aoj[j]; 127582412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 127682412ba7SHong Zhang pj = pj_oth + pi_oth[row]; 127782412ba7SHong Zhang pa = pa_oth + pi_oth[row]; 1278e900a757SHong Zhang valtmp = aoa[j]; 1279f4a743e1SHong Zhang nextp = 0; 128082412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 128182412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 1282e900a757SHong Zhang apa[k] += valtmp*pa[nextp++]; 128342c57489SHong Zhang } 128442c57489SHong Zhang } 1285dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 128642c57489SHong Zhang } 1287b091e509SHong Zhang #if defined(PTAP_PROFILE) 12888563dfccSBarry Smith ierr = PetscTime(&t2_1);CHKERRQ(ierr); 1289b091e509SHong Zhang et2_AP += t2_1 - t2_0; 1290b091e509SHong Zhang #endif 129142c57489SHong Zhang 1292a9d06231SHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 1293a9d06231SHong Zhang /*--------------------------------------------------------------*/ 129482412ba7SHong Zhang pnz = pi_loc[i+1] - pi_loc[i]; 1295f9473cf7SHong Zhang pJ = pj_loc + pi_loc[i]; 129682412ba7SHong Zhang for (j=0; j<pnz; j++) { 129742c57489SHong Zhang nextap = 0; 1298f9473cf7SHong Zhang row = pJ[j]; /* global index */ 129982412ba7SHong Zhang if (row < pcstart || row >=pcend) { /* put the value into Co */ 1300e900a757SHong Zhang row = *poJ; 1301e900a757SHong Zhang cj = coj + coi[row]; 1302e900a757SHong Zhang ca = coa + coi[row]; poJ++; 130382412ba7SHong Zhang } else { /* put the value into Cd */ 1304e900a757SHong Zhang row = *pdJ; 1305e900a757SHong Zhang cj = bj + bi[row]; 1306e900a757SHong Zhang ca = ba + bi[row]; pdJ++; 130742c57489SHong Zhang } 1308e900a757SHong Zhang valtmp = pA[j]; 130982412ba7SHong Zhang for (k=0; nextap<apnz; k++) { 1310e900a757SHong Zhang if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++]; 131142c57489SHong Zhang } 1312dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 1313de0260b3SHong Zhang } 13140496f32cSHong Zhang pA += pnz; 131542c57489SHong Zhang /* zero the current row info for A*P */ 131682412ba7SHong Zhang ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr); 1317b091e509SHong Zhang #if defined(PTAP_PROFILE) 13188563dfccSBarry Smith ierr = PetscTime(&t2_2);CHKERRQ(ierr); 131905d62848SHong Zhang ePtAP += t2_2 - t2_1; 1320b091e509SHong Zhang #endif 132142c57489SHong Zhang } 1322d9824c15SHong Zhang } 132338571addSHong Zhang #if defined(PTAP_PROFILE) 13248563dfccSBarry Smith ierr = PetscTime(&t2);CHKERRQ(ierr); 132538571addSHong Zhang #endif 132642c57489SHong Zhang 1327a9d06231SHong Zhang /* 3) send and recv matrix values coa */ 1328a9d06231SHong Zhang /*------------------------------------*/ 132942c57489SHong Zhang buf_ri = merge->buf_ri; 133042c57489SHong Zhang buf_rj = merge->buf_rj; 133142c57489SHong Zhang len_s = merge->len_s; 1332fc42d0c8SSatish Balay ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 133342c57489SHong Zhang ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 133442c57489SHong Zhang 1335dcca6d9dSJed Brown ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 133642c57489SHong Zhang for (proc=0,k=0; proc<size; proc++) { 133742c57489SHong Zhang if (!len_s[proc]) continue; 1338de0260b3SHong Zhang i = merge->owners_co[proc]; 1339de0260b3SHong Zhang ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 134042c57489SHong Zhang k++; 134142c57489SHong Zhang } 13420c468ba9SBarry Smith if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 13430c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 134442c57489SHong Zhang 1345a9d06231SHong Zhang ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 134642c57489SHong Zhang ierr = PetscFree(r_waits);CHKERRQ(ierr); 134782412ba7SHong Zhang ierr = PetscFree(coa);CHKERRQ(ierr); 134838571addSHong Zhang #if defined(PTAP_PROFILE) 13498563dfccSBarry Smith ierr = PetscTime(&t3);CHKERRQ(ierr); 135038571addSHong Zhang #endif 135142c57489SHong Zhang 1352a9d06231SHong Zhang /* 4) insert local Cseq and received values into Cmpi */ 1353a9d06231SHong Zhang /*------------------------------------------------------*/ 1354dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 135542c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { 135642c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 135742c57489SHong Zhang nrows = *(buf_ri_k[k]); 135842c57489SHong Zhang nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 135982412ba7SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 136042c57489SHong Zhang } 136142c57489SHong Zhang 1362de0260b3SHong Zhang for (i=0; i<cm; i++) { 136382412ba7SHong Zhang row = owners[rank] + i; /* global row index of C_seq */ 136442c57489SHong Zhang bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1365de0260b3SHong Zhang ba_i = ba + bi[i]; 136682412ba7SHong Zhang bnz = bi[i+1] - bi[i]; 136742c57489SHong Zhang /* add received vals into ba */ 136842c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 136942c57489SHong Zhang /* i-th row */ 137042c57489SHong Zhang if (i == *nextrow[k]) { 137182412ba7SHong Zhang cnz = *(nextci[k]+1) - *nextci[k]; 137282412ba7SHong Zhang cj = buf_rj[k] + *(nextci[k]); 137382412ba7SHong Zhang ca = abuf_r[k] + *(nextci[k]); 137482412ba7SHong Zhang nextcj = 0; 137582412ba7SHong Zhang for (j=0; nextcj<cnz; j++) { 137682412ba7SHong Zhang if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 137782412ba7SHong Zhang ba_i[j] += ca[nextcj++]; 137842c57489SHong Zhang } 137942c57489SHong Zhang } 138082412ba7SHong Zhang nextrow[k]++; nextci[k]++; 1381c5af039cSHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 138242c57489SHong Zhang } 138342c57489SHong Zhang } 138482412ba7SHong Zhang ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 138542c57489SHong Zhang } 138642c57489SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 138742c57489SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 138842c57489SHong Zhang 1389de0260b3SHong Zhang ierr = PetscFree(ba);CHKERRQ(ierr); 1390c05d87d6SBarry Smith ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 139142c57489SHong Zhang ierr = PetscFree(abuf_r);CHKERRQ(ierr); 13920572522cSBarry Smith ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 139338571addSHong Zhang #if defined(PTAP_PROFILE) 13948563dfccSBarry Smith ierr = PetscTime(&t4);CHKERRQ(ierr); 13959516a85cSHong Zhang if (rank==1) { 139605d62848SHong 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); 13979516a85cSHong Zhang } 139838571addSHong Zhang #endif 139942c57489SHong Zhang PetscFunctionReturn(0); 140042c57489SHong Zhang } 1401