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); 34*681d504bSHong Zhang if (ptap->AP_loc) { 35*681d504bSHong Zhang Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data; 36*681d504bSHong Zhang ierr = PetscFree(ap->i);CHKERRQ(ierr); 37*681d504bSHong Zhang ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr); 380d3441aeSHong Zhang ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr); 39*681d504bSHong Zhang } 402259aa2eSHong Zhang ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr); 412259aa2eSHong Zhang ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr); 42414894bdSHong Zhang if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);} 43c8b0795cSMark F. Adams if (merge) { 44cc6cb787SHong Zhang ierr = PetscFree(merge->id_r);CHKERRQ(ierr); 45cc6cb787SHong Zhang ierr = PetscFree(merge->len_s);CHKERRQ(ierr); 46cc6cb787SHong Zhang ierr = PetscFree(merge->len_r);CHKERRQ(ierr); 47cc6cb787SHong Zhang ierr = PetscFree(merge->bi);CHKERRQ(ierr); 48cc6cb787SHong Zhang ierr = PetscFree(merge->bj);CHKERRQ(ierr); 49c05d87d6SBarry Smith ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr); 50cc6cb787SHong Zhang ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr); 51c05d87d6SBarry Smith ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr); 52cc6cb787SHong Zhang ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr); 53445158ffSHong Zhang ierr = PetscFree(merge->coi);CHKERRQ(ierr); 54445158ffSHong Zhang ierr = PetscFree(merge->coj);CHKERRQ(ierr); 5505b42c5fSBarry Smith ierr = PetscFree(merge->owners_co);CHKERRQ(ierr); 566bf464f9SBarry Smith ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr); 57dce485f0SBarry Smith ierr = merge->destroy(A);CHKERRQ(ierr); 58f8487c73SHong Zhang ierr = PetscFree(ptap->merge);CHKERRQ(ierr); 59445158ffSHong Zhang } else { 60445158ffSHong Zhang ierr = MatDestroy_MPIAIJ(A);CHKERRQ(ierr); 61bf0cc555SLisandro Dalcin } 62f8487c73SHong Zhang ierr = PetscFree(ptap);CHKERRQ(ierr); 63c8b0795cSMark F. Adams } 64cc6cb787SHong Zhang PetscFunctionReturn(0); 65cc6cb787SHong Zhang } 66cc6cb787SHong Zhang 67cc6cb787SHong Zhang #undef __FUNCT__ 68aa690a28SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP_old" 69aa690a28SHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP_old(Mat A, MatDuplicateOption op, Mat *M) 70f0c0a3a6SBarry Smith { 71cc6cb787SHong Zhang PetscErrorCode ierr; 7211a6856fSHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 7311a6856fSHong Zhang Mat_PtAPMPI *ptap = a->ptap; 7411a6856fSHong Zhang Mat_Merge_SeqsToMPI *merge = ptap->merge; 75cc6cb787SHong Zhang 76cc6cb787SHong Zhang PetscFunctionBegin; 77dce485f0SBarry Smith ierr = (*merge->duplicate)(A,op,M);CHKERRQ(ierr); 782205254eSKarl Rupp 7911a6856fSHong Zhang (*M)->ops->destroy = merge->destroy; 8011a6856fSHong Zhang (*M)->ops->duplicate = merge->duplicate; 81cc6cb787SHong Zhang PetscFunctionReturn(0); 82cc6cb787SHong Zhang } 83cc6cb787SHong Zhang 841a47ec54SHong Zhang 851a47ec54SHong Zhang #undef __FUNCT__ 86aa690a28SHong Zhang #define __FUNCT__ "MatDuplicate_MPIAIJ_MatPtAP" 87aa690a28SHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M) 881a47ec54SHong Zhang { 891a47ec54SHong Zhang PetscErrorCode ierr; 901a47ec54SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 911a47ec54SHong Zhang Mat_PtAPMPI *ptap = a->ptap; 921a47ec54SHong Zhang 931a47ec54SHong Zhang PetscFunctionBegin; 941a47ec54SHong Zhang ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr); 951a47ec54SHong Zhang (*M)->ops->destroy = ptap->destroy; 961a47ec54SHong Zhang (*M)->ops->duplicate = ptap->duplicate; 971a47ec54SHong Zhang PetscFunctionReturn(0); 981a47ec54SHong Zhang } 991a47ec54SHong Zhang 10042c57489SHong Zhang #undef __FUNCT__ 101cf3ca8ceSHong Zhang #define __FUNCT__ "MatPtAP_MPIAIJ_MPIAIJ" 102cf3ca8ceSHong Zhang PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 10342c57489SHong Zhang { 10442c57489SHong Zhang PetscErrorCode ierr; 105aa690a28SHong Zhang PetscBool newalg=PETSC_TRUE; 10667a12041SHong Zhang MPI_Comm comm; 10742c57489SHong Zhang 10842c57489SHong Zhang PetscFunctionBegin; 10967a12041SHong Zhang /* check if matrix local sizes are compatible */ 11067a12041SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 11167a12041SHong Zhang if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) { 11267a12041SHong 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); 11367a12041SHong Zhang } 11467a12041SHong Zhang if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) { 11567a12041SHong 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); 11667a12041SHong Zhang } 11767a12041SHong Zhang 1180d3441aeSHong Zhang ierr = PetscOptionsGetBool(NULL,"-matptap_new",&newalg,NULL);CHKERRQ(ierr); 119cf3ca8ceSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 1203ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 1210d3441aeSHong Zhang if (newalg) { 122cf3ca8ceSHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr); 123aa690a28SHong Zhang } else { 124aa690a28SHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_old(A,P,fill,C);CHKERRQ(ierr); 1250d3441aeSHong Zhang } 1263ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 1277a7894deSKris Buschelman } 1283ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 1290d3441aeSHong Zhang if (newalg) { 130cf3ca8ceSHong Zhang ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ(A,P,*C);CHKERRQ(ierr); 131aa690a28SHong Zhang } else { 132aa690a28SHong Zhang ierr = MatPtAPNumeric_MPIAIJ_MPIAIJ_old(A,P,*C);CHKERRQ(ierr); 1330d3441aeSHong Zhang } 1343ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 13542c57489SHong Zhang PetscFunctionReturn(0); 13642c57489SHong Zhang } 13742c57489SHong Zhang 1388cb82516SHong Zhang /* requires array of size pN, thus nonscalable version */ 13942c57489SHong Zhang #undef __FUNCT__ 140aa690a28SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ" 141aa690a28SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 1420d3441aeSHong Zhang { 1430d3441aeSHong Zhang PetscErrorCode ierr; 1440d3441aeSHong Zhang Mat_PtAPMPI *ptap; 145*681d504bSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 1462259aa2eSHong Zhang MPI_Comm comm; 1472259aa2eSHong Zhang PetscMPIInt size,rank; 14815a3b8e2SHong Zhang Mat Cmpi; 14915a3b8e2SHong Zhang PetscFreeSpaceList free_space=NULL,current_space=NULL; 15067a12041SHong Zhang PetscInt am=A->rmap->n,pm=P->rmap->n; 15167a12041SHong Zhang PetscInt *lnk,i,k,pnz,row,nnz; 152445158ffSHong Zhang PetscInt pN=P->cmap->N,pn=P->cmap->n; 15315a3b8e2SHong Zhang PetscBT lnkbt; 154f671be37SHong Zhang PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 15515a3b8e2SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 15615a3b8e2SHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 15767a12041SHong Zhang PetscInt nzi,nspacedouble; 15815a3b8e2SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 15915a3b8e2SHong Zhang MPI_Request *swaits,*rwaits; 16015a3b8e2SHong Zhang MPI_Status *sstatus,rstatus; 161445158ffSHong Zhang PetscLayout rowmap; 162445158ffSHong Zhang PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 163f671be37SHong Zhang PetscInt nsend; 164445158ffSHong Zhang PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 16567a12041SHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0; 16667a12041SHong Zhang PetscInt rmax,*aj,*ai,*pi; 16767a12041SHong Zhang Mat_SeqAIJ *p_loc,*p_oth,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*c_loc,*c_oth; 16867a12041SHong Zhang PetscScalar *apv; 169aa690a28SHong Zhang #if defined(PETSC_USE_INFO) 1708cb82516SHong Zhang PetscReal apfill; 171aa690a28SHong Zhang #endif 172*681d504bSHong Zhang #if defined(PTAP_PROFILE) 173*681d504bSHong Zhang PetscLogDouble t0,t1,t11,t12,t2,t3,t4; 174*681d504bSHong Zhang #endif 17567a12041SHong Zhang 17667a12041SHong Zhang PetscFunctionBegin; 17767a12041SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 17867a12041SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 17967a12041SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1808cb82516SHong Zhang #if defined(PTAP_PROFILE) 18180bb4639SHong Zhang ierr = PetscTime(&t0);CHKERRQ(ierr); 1828cb82516SHong Zhang #endif 1838cb82516SHong Zhang 18415a3b8e2SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 18515a3b8e2SHong Zhang ierr = PetscNew(&ptap);CHKERRQ(ierr); 18615a3b8e2SHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 18715a3b8e2SHong Zhang 18815a3b8e2SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 18915a3b8e2SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 19015a3b8e2SHong Zhang /* get P_loc by taking all local rows of P */ 19115a3b8e2SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 19215a3b8e2SHong Zhang 19367a12041SHong Zhang /* (0) compute Rd = Pd^T, Ro = Po^T */ 19467a12041SHong Zhang /* --------------------------------- */ 195de817e96SHong Zhang ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 196de817e96SHong Zhang ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 1978cb82516SHong Zhang #if defined(PTAP_PROFILE) 198445158ffSHong Zhang ierr = PetscTime(&t11);CHKERRQ(ierr); 1998cb82516SHong Zhang #endif 20015a3b8e2SHong Zhang 20167a12041SHong Zhang /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 20267a12041SHong Zhang /* ----------------------------------------------------------------- */ 20367a12041SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 20467a12041SHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 205445158ffSHong Zhang 2068cb82516SHong Zhang /* create and initialize a linked list - nonscalable version */ 207445158ffSHong Zhang ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 208445158ffSHong Zhang 2098cb82516SHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 21067a12041SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(ad->i[am]+ao->i[am]+p_loc->i[pm])),&free_space);CHKERRQ(ierr); 211445158ffSHong Zhang current_space = free_space; 21267a12041SHong Zhang nspacedouble = 0; 21367a12041SHong Zhang 21467a12041SHong Zhang ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 21567a12041SHong Zhang api[0] = 0; 216445158ffSHong Zhang for (i=0; i<am; i++) { 21767a12041SHong Zhang /* diagonal portion: Ad[i,:]*P */ 21867a12041SHong Zhang ai = ad->i; pi = p_loc->i; 21967a12041SHong Zhang nzi = ai[i+1] - ai[i]; 22067a12041SHong Zhang aj = ad->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_loc->j + pi[row]; 225445158ffSHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 226445158ffSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 227445158ffSHong Zhang } 22867a12041SHong Zhang /* off-diagonal portion: Ao[i,:]*P */ 22967a12041SHong Zhang ai = ao->i; pi = p_oth->i; 23067a12041SHong Zhang nzi = ai[i+1] - ai[i]; 23167a12041SHong Zhang aj = ao->j + ai[i]; 232445158ffSHong Zhang for (j=0; j<nzi; j++) { 233445158ffSHong Zhang row = aj[j]; 23467a12041SHong Zhang pnz = pi[row+1] - pi[row]; 23567a12041SHong Zhang Jptr = p_oth->j + pi[row]; 236445158ffSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 237445158ffSHong Zhang } 238445158ffSHong Zhang apnz = lnk[0]; 239445158ffSHong Zhang api[i+1] = api[i] + apnz; 240445158ffSHong Zhang if (ap_rmax < apnz) ap_rmax = apnz; 241445158ffSHong Zhang 242445158ffSHong Zhang /* if free space is not available, double the total space in the list */ 243445158ffSHong Zhang if (current_space->local_remaining<apnz) { 244445158ffSHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 245445158ffSHong Zhang nspacedouble++; 246445158ffSHong Zhang } 247445158ffSHong Zhang 248445158ffSHong Zhang /* Copy data into free space, then initialize lnk */ 249445158ffSHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 250445158ffSHong Zhang 251445158ffSHong Zhang current_space->array += apnz; 252445158ffSHong Zhang current_space->local_used += apnz; 253445158ffSHong Zhang current_space->local_remaining -= apnz; 254445158ffSHong Zhang } 255*681d504bSHong Zhang /* Allocate space for apj and apv, initialize apj, and */ 256445158ffSHong Zhang /* destroy list of free space and other temporary array(s) */ 257445158ffSHong Zhang ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 258445158ffSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 259445158ffSHong Zhang 260aa690a28SHong Zhang /* Create AP_loc for reuse */ 261445158ffSHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 262aa690a28SHong Zhang 263aa690a28SHong Zhang #if defined(PETSC_USE_INFO) 264aa690a28SHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 265aa690a28SHong Zhang ptap->AP_loc->info.mallocs = nspacedouble; 266aa690a28SHong Zhang ptap->AP_loc->info.fill_ratio_given = fill; 267aa690a28SHong Zhang ptap->AP_loc->info.fill_ratio_needed = apfill; 268aa690a28SHong Zhang 269aa690a28SHong Zhang if (api[am]) { 270aa690a28SHong Zhang ierr = PetscInfo3(ptap->AP_loc,"AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr); 271aa690a28SHong Zhang ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 272aa690a28SHong Zhang } else { 273aa690a28SHong Zhang ierr = PetscInfo(ptap->AP_loc,"AP_loc is empty \n");CHKERRQ(ierr); 274aa690a28SHong Zhang } 275aa690a28SHong Zhang #endif 276aa690a28SHong Zhang 2778cb82516SHong Zhang #if defined(PTAP_PROFILE) 278445158ffSHong Zhang ierr = PetscTime(&t12);CHKERRQ(ierr); 2798cb82516SHong Zhang #endif 28015a3b8e2SHong Zhang 281*681d504bSHong Zhang /* (2-1) compute symbolic Co = Ro*AP_loc */ 282*681d504bSHong Zhang /* ------------------------------------ */ 28367a12041SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 28467a12041SHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 2858cb82516SHong Zhang #if defined(PTAP_PROFILE) 28680bb4639SHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 2878cb82516SHong Zhang #endif 28815a3b8e2SHong Zhang 28967a12041SHong Zhang /* (3) send coj of C_oth to other processors */ 29067a12041SHong Zhang /* ------------------------------------------ */ 29115a3b8e2SHong Zhang /* determine row ownership */ 292445158ffSHong Zhang ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 293445158ffSHong Zhang rowmap->n = pn; 294445158ffSHong Zhang rowmap->bs = 1; 295445158ffSHong Zhang ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 296445158ffSHong Zhang owners = rowmap->range; 29715a3b8e2SHong Zhang 29815a3b8e2SHong Zhang /* determine the number of messages to send, their lengths */ 2998cb82516SHong Zhang ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 3008cb82516SHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 30115a3b8e2SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 30215a3b8e2SHong Zhang 30367a12041SHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 30467a12041SHong Zhang coi = c_oth->i; coj = c_oth->j; 30567a12041SHong Zhang con = ptap->C_oth->rmap->n; 30615a3b8e2SHong Zhang proc = 0; 30767a12041SHong Zhang for (i=0; i<con; i++) { 30815a3b8e2SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 30915a3b8e2SHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 31015a3b8e2SHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 31115a3b8e2SHong Zhang } 31215a3b8e2SHong Zhang 31315a3b8e2SHong Zhang len = 0; /* max length of buf_si[], see (4) */ 31415a3b8e2SHong Zhang owners_co[0] = 0; 31567a12041SHong Zhang nsend = 0; 31615a3b8e2SHong Zhang for (proc=0; proc<size; proc++) { 31715a3b8e2SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 31815a3b8e2SHong Zhang if (len_s[proc]) { 319445158ffSHong Zhang nsend++; 32015a3b8e2SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 32115a3b8e2SHong Zhang len += len_si[proc]; 32215a3b8e2SHong Zhang } 32315a3b8e2SHong Zhang } 32415a3b8e2SHong Zhang 32515a3b8e2SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 326445158ffSHong Zhang ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 327445158ffSHong Zhang ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 32815a3b8e2SHong Zhang 32915a3b8e2SHong Zhang /* post the Irecv and Isend of coj */ 33015a3b8e2SHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 331445158ffSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 332445158ffSHong Zhang ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 33315a3b8e2SHong Zhang for (proc=0, k=0; proc<size; proc++) { 33415a3b8e2SHong Zhang if (!len_s[proc]) continue; 33515a3b8e2SHong Zhang i = owners_co[proc]; 33615a3b8e2SHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 33715a3b8e2SHong Zhang k++; 33815a3b8e2SHong Zhang } 33915a3b8e2SHong Zhang 340*681d504bSHong Zhang /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 341*681d504bSHong Zhang /* ---------------------------------------- */ 342*681d504bSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 343*681d504bSHong Zhang c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 344*681d504bSHong Zhang 345*681d504bSHong Zhang /* receives 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 } 381*681d504bSHong 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 /* ------------------------------------------ */ 396*681d504bSHong Zhang /* set initial free space to be pN */ 3978cb82516SHong Zhang ierr = PetscFreeSpaceGet(pN,&free_space);CHKERRQ(ierr); /* non-scalable version */ 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 } 40715a3b8e2SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 408445158ffSHong Zhang 40915a3b8e2SHong Zhang rmax = 0; 41015a3b8e2SHong Zhang for (i=0; i<pn; i++) { 41167a12041SHong Zhang /* add C_loc into Cmpi */ 41267a12041SHong Zhang nzi = c_loc->i[i+1] - c_loc->i[i]; 41367a12041SHong Zhang Jptr = c_loc->j + c_loc->i[i]; 41467a12041SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 41515a3b8e2SHong Zhang 41615a3b8e2SHong Zhang /* add received col data into lnk */ 417445158ffSHong Zhang for (k=0; k<nrecv; k++) { /* k-th received message */ 41815a3b8e2SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 41915a3b8e2SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 42015a3b8e2SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 42115a3b8e2SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 42215a3b8e2SHong Zhang nextrow[k]++; nextci[k]++; 42315a3b8e2SHong Zhang } 42415a3b8e2SHong Zhang } 42515a3b8e2SHong Zhang nnz = lnk[0]; 4268cb82516SHong Zhang 42715a3b8e2SHong Zhang /* copy data into free space, then initialize lnk */ 42815a3b8e2SHong Zhang ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 42915a3b8e2SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 43015a3b8e2SHong Zhang if (nnz > rmax) rmax = nnz; 43115a3b8e2SHong Zhang } 43215a3b8e2SHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 43315a3b8e2SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 434445158ffSHong Zhang ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 4358cb82516SHong Zhang #if defined(PTAP_PROFILE) 43680bb4639SHong Zhang ierr = PetscTime(&t3);CHKERRQ(ierr); 4378cb82516SHong Zhang #endif 43880bb4639SHong Zhang 43915a3b8e2SHong Zhang /* (6) create symbolic parallel matrix Cmpi */ 44015a3b8e2SHong Zhang /*------------------------------------------*/ 44115a3b8e2SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 44215a3b8e2SHong Zhang ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 44315a3b8e2SHong Zhang ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 44415a3b8e2SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 44515a3b8e2SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 44615a3b8e2SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 44715a3b8e2SHong Zhang 448445158ffSHong Zhang /* members in merge */ 449445158ffSHong Zhang ierr = PetscFree(id_r);CHKERRQ(ierr); 450445158ffSHong Zhang ierr = PetscFree(len_r);CHKERRQ(ierr); 451445158ffSHong Zhang ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 452445158ffSHong Zhang ierr = PetscFree(buf_ri);CHKERRQ(ierr); 453445158ffSHong Zhang ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 454445158ffSHong Zhang ierr = PetscFree(buf_rj);CHKERRQ(ierr); 455445158ffSHong Zhang ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 45615a3b8e2SHong Zhang 45715a3b8e2SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 45815a3b8e2SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 45915a3b8e2SHong Zhang c->ptap = ptap; 46015a3b8e2SHong Zhang ptap->rmax = ap_rmax; 4611a47ec54SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 4621a47ec54SHong Zhang ptap->destroy = Cmpi->ops->destroy; 46315a3b8e2SHong Zhang 4648cb82516SHong Zhang /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ_new() */ 4658cb82516SHong Zhang ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 4662259aa2eSHong Zhang 4671a47ec54SHong Zhang 4681a47ec54SHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 4691a47ec54SHong Zhang Cmpi->assembled = PETSC_FALSE; 4701a47ec54SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 471aa690a28SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 472aa690a28SHong Zhang Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ; 4731a47ec54SHong Zhang *C = Cmpi; 4741a47ec54SHong Zhang 4758cb82516SHong Zhang #if defined(PTAP_PROFILE) 47680bb4639SHong Zhang ierr = PetscTime(&t4);CHKERRQ(ierr); 477e9c1f85fSHong Zhang if (rank == 1) { 478445158ffSHong 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); 479e9c1f85fSHong Zhang } 4808cb82516SHong Zhang #endif 4810d3441aeSHong Zhang PetscFunctionReturn(0); 4820d3441aeSHong Zhang } 4830d3441aeSHong Zhang 4840d3441aeSHong Zhang #undef __FUNCT__ 485aa690a28SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ" 486aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 4870d3441aeSHong Zhang { 4880d3441aeSHong Zhang PetscErrorCode ierr; 4890d3441aeSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 4900d3441aeSHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 4918cb82516SHong Zhang Mat_SeqAIJ *ap,*p_loc,*p_oth,*c_seq; 4920d3441aeSHong Zhang Mat_PtAPMPI *ptap = c->ptap; 4939ce11a7cSHong Zhang Mat AP_loc,C_loc,C_oth; 4940d3441aeSHong Zhang PetscInt i,rstart,rend,cm,ncols,row; 4958cb82516SHong Zhang PetscInt *api,*apj,am = A->rmap->n,j,col,apnz; 4968cb82516SHong Zhang PetscScalar *apa; 4970d3441aeSHong Zhang const PetscInt *cols; 4980d3441aeSHong Zhang const PetscScalar *vals; 4998cb82516SHong Zhang #if defined(PTAP_PROFILE) 5008cb82516SHong Zhang PetscMPIInt rank; 5018cb82516SHong Zhang MPI_Comm comm; 5020d3441aeSHong Zhang PetscLogDouble t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi; 5038cb82516SHong Zhang #endif 5040d3441aeSHong Zhang 5050d3441aeSHong Zhang PetscFunctionBegin; 5060d3441aeSHong Zhang ierr = MatZeroEntries(C);CHKERRQ(ierr); 5070d3441aeSHong Zhang 508e497d3c8SHong Zhang /* 1) get R = Pd^T,Ro = Po^T */ 5098cb82516SHong Zhang #if defined(PTAP_PROFILE) 5100d3441aeSHong Zhang ierr = PetscTime(&t0);CHKERRQ(ierr); 5118cb82516SHong Zhang #endif 512748c7196SHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 5130d3441aeSHong Zhang ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 5140d3441aeSHong Zhang ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 515748c7196SHong Zhang } 5168cb82516SHong Zhang #if defined(PTAP_PROFILE) 5170d3441aeSHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 5180d3441aeSHong Zhang eR = t1 - t0; 5198cb82516SHong Zhang #endif 5200d3441aeSHong Zhang 521e497d3c8SHong Zhang /* 2) get AP_loc */ 5220d3441aeSHong Zhang AP_loc = ptap->AP_loc; 5238cb82516SHong Zhang ap = (Mat_SeqAIJ*)AP_loc->data; 5240d3441aeSHong Zhang 525e497d3c8SHong Zhang /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 5260d3441aeSHong Zhang /*-----------------------------------------------------*/ 527748c7196SHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 528748c7196SHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 5290d3441aeSHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 5300d3441aeSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 531e497d3c8SHong Zhang } 5320d3441aeSHong Zhang 5338cb82516SHong Zhang /* 2-2) compute numeric A_loc*P - dominating part */ 5348cb82516SHong Zhang /* ---------------------------------------------- */ 5350d3441aeSHong Zhang /* get data from symbolic products */ 5368cb82516SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 5378cb82516SHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 5388cb82516SHong Zhang apa = ptap->apa; 539*681d504bSHong Zhang api = ap->i; 540*681d504bSHong Zhang apj = ap->j; 541e497d3c8SHong Zhang for (i=0; i<am; i++) { 542445158ffSHong Zhang /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 543e497d3c8SHong Zhang AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa); 544e497d3c8SHong Zhang apnz = api[i+1] - api[i]; 545e497d3c8SHong Zhang for (j=0; j<apnz; j++) { 546e497d3c8SHong Zhang col = apj[j+api[i]]; 547e497d3c8SHong Zhang ap->a[j+ap->i[i]] = apa[col]; 548e497d3c8SHong Zhang apa[col] = 0.0; 5490d3441aeSHong Zhang } 550aa690a28SHong Zhang ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 5510d3441aeSHong Zhang } 5528cb82516SHong Zhang #if defined(PTAP_PROFILE) 5530d3441aeSHong Zhang ierr = PetscTime(&t2);CHKERRQ(ierr); 5540d3441aeSHong Zhang eAP = t2 - t1; 5558cb82516SHong Zhang #endif 5560d3441aeSHong Zhang 5578cb82516SHong Zhang /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 558445158ffSHong Zhang ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 559445158ffSHong Zhang ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 5609ce11a7cSHong Zhang C_loc = ptap->C_loc; 5619ce11a7cSHong Zhang C_oth = ptap->C_oth; 5628cb82516SHong Zhang #if defined(PTAP_PROFILE) 5630d3441aeSHong Zhang ierr = PetscTime(&t3);CHKERRQ(ierr); 5640d3441aeSHong Zhang eCseq = t3 - t2; 5658cb82516SHong Zhang #endif 5660d3441aeSHong Zhang 5670d3441aeSHong Zhang /* add C_loc and Co to to C */ 5680d3441aeSHong Zhang ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 5690d3441aeSHong Zhang 5700d3441aeSHong Zhang /* C_loc -> C */ 5710d3441aeSHong Zhang cm = C_loc->rmap->N; 5729ce11a7cSHong Zhang c_seq = (Mat_SeqAIJ*)C_loc->data; 5738cb82516SHong Zhang cols = c_seq->j; 5748cb82516SHong Zhang vals = c_seq->a; 5750d3441aeSHong Zhang for (i=0; i<cm; i++) { 5769ce11a7cSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 5770d3441aeSHong Zhang row = rstart + i; 5780d3441aeSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 5798cb82516SHong Zhang cols += ncols; vals += ncols; 5800d3441aeSHong Zhang } 5810d3441aeSHong Zhang 5820d3441aeSHong Zhang /* Co -> C, off-processor part */ 5839ce11a7cSHong Zhang cm = C_oth->rmap->N; 5849ce11a7cSHong Zhang c_seq = (Mat_SeqAIJ*)C_oth->data; 5858cb82516SHong Zhang cols = c_seq->j; 5868cb82516SHong Zhang vals = c_seq->a; 5879ce11a7cSHong Zhang for (i=0; i<cm; i++) { 5889ce11a7cSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 5890d3441aeSHong Zhang row = p->garray[i]; 5900d3441aeSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 5918cb82516SHong Zhang cols += ncols; vals += ncols; 5920d3441aeSHong Zhang } 5930d3441aeSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5940d3441aeSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 5958cb82516SHong Zhang #if defined(PTAP_PROFILE) 5960d3441aeSHong Zhang ierr = PetscTime(&t4);CHKERRQ(ierr); 5970d3441aeSHong Zhang eCmpi = t4 - t3; 5980d3441aeSHong Zhang 5998cb82516SHong Zhang ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 6008cb82516SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 6010d3441aeSHong Zhang if (rank==1) { 6020d3441aeSHong 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); 6030d3441aeSHong Zhang } 6048cb82516SHong Zhang #endif 605748c7196SHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 6060d3441aeSHong Zhang PetscFunctionReturn(0); 6070d3441aeSHong Zhang } 6080d3441aeSHong Zhang 6090d3441aeSHong Zhang #undef __FUNCT__ 610aa690a28SHong Zhang #define __FUNCT__ "MatPtAPSymbolic_MPIAIJ_MPIAIJ_old" 611aa690a28SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_old(Mat A,Mat P,PetscReal fill,Mat *C) 61242c57489SHong Zhang { 61342c57489SHong Zhang PetscErrorCode ierr; 61477584fe6SHong Zhang Mat Cmpi; 615f8487c73SHong Zhang Mat_PtAPMPI *ptap; 6160298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 617f8487c73SHong Zhang Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 61842c57489SHong Zhang Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 61942c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 620ded4f561SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ; 62177584fe6SHong Zhang PetscInt *adi=ad->i,*aj,*aoi=ao->i,nnz; 622a914927fSHong Zhang PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 623d16ca5f0SHong Zhang PetscInt am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n; 62442c57489SHong Zhang PetscBT lnkbt; 625ce94432eSBarry Smith MPI_Comm comm; 626a914927fSHong Zhang PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0; 62742c57489SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 62842c57489SHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 62924ecddacSHong Zhang PetscInt nzi,*pti,*ptj; 63042c57489SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 631ded4f561SHong Zhang MPI_Request *swaits,*rwaits; 632ded4f561SHong Zhang MPI_Status *sstatus,rstatus; 63342c57489SHong Zhang Mat_Merge_SeqsToMPI *merge; 63477584fe6SHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0; 635d16ca5f0SHong Zhang PetscReal afill=1.0,afill_tmp; 636a914927fSHong Zhang PetscInt rmax; 63742c57489SHong Zhang 63842c57489SHong Zhang PetscFunctionBegin; 639ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 64083445d95SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 64183445d95SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 64283445d95SHong Zhang 643f8487c73SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 644b00a9115SJed Brown ierr = PetscNew(&ptap);CHKERRQ(ierr); 645b00a9115SJed Brown ierr = PetscNew(&merge);CHKERRQ(ierr); 6469f4155fbSHong Zhang ptap->merge = merge; 647f8487c73SHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 6486c6e00e0SHong Zhang 6496c6e00e0SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 650b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 651fe615159SHong Zhang 6526c6e00e0SHong Zhang /* get P_loc by taking all local rows of P */ 653a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 6546c6e00e0SHong Zhang 655a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 656a1a4e84aSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 6576c6e00e0SHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 6586c6e00e0SHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 65942c57489SHong Zhang 6602addb229SHong Zhang /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */ 6612addb229SHong Zhang /*--------------------------------------------------------------------------*/ 662854ce69bSBarry Smith ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 66382412ba7SHong Zhang api[0] = 0; 66483445d95SHong Zhang 66583445d95SHong Zhang /* create and initialize a linked list */ 666a914927fSHong Zhang ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 66783445d95SHong Zhang 668a914927fSHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */ 669d16ca5f0SHong Zhang ierr = PetscFreeSpaceGet((PetscInt)(fill*(adi[am]+aoi[am]+pi_loc[pm])),&free_space);CHKERRQ(ierr); 670f4a743e1SHong Zhang current_space = free_space; 671f4a743e1SHong Zhang 672f4a743e1SHong Zhang for (i=0; i<am; i++) { 673f4a743e1SHong Zhang /* diagonal portion of A */ 674ded4f561SHong Zhang nzi = adi[i+1] - adi[i]; 67577584fe6SHong Zhang aj = ad->j + adi[i]; 676ded4f561SHong Zhang for (j=0; j<nzi; j++) { 67777584fe6SHong Zhang row = aj[j]; 67882412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 679ded4f561SHong Zhang Jptr = pj_loc + pi_loc[row]; 68083445d95SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 681a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 682f4a743e1SHong Zhang } 683f4a743e1SHong Zhang /* off-diagonal portion of A */ 684ded4f561SHong Zhang nzi = aoi[i+1] - aoi[i]; 68577584fe6SHong Zhang aj = ao->j + aoi[i]; 686ded4f561SHong Zhang for (j=0; j<nzi; j++) { 68777584fe6SHong Zhang row = aj[j]; 68882412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 689ded4f561SHong Zhang Jptr = pj_oth + pi_oth[row]; 690a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 691f4a743e1SHong Zhang } 692a914927fSHong Zhang apnz = lnk[0]; 69382412ba7SHong Zhang api[i+1] = api[i] + apnz; 69477584fe6SHong Zhang if (ap_rmax < apnz) ap_rmax = apnz; 695f4a743e1SHong Zhang 69683445d95SHong Zhang /* if free space is not available, double the total space in the list */ 69782412ba7SHong Zhang if (current_space->local_remaining<apnz) { 6982ba1acd5SHong Zhang ierr = PetscFreeSpaceGet(apnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 699f2b054eeSHong Zhang nspacedouble++; 700f4a743e1SHong Zhang } 701f4a743e1SHong Zhang 702f4a743e1SHong Zhang /* Copy data into free space, then initialize lnk */ 703a914927fSHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 7042205254eSKarl Rupp 70582412ba7SHong Zhang current_space->array += apnz; 70682412ba7SHong Zhang current_space->local_used += apnz; 70782412ba7SHong Zhang current_space->local_remaining -= apnz; 708f4a743e1SHong Zhang } 709a914927fSHong Zhang 71082412ba7SHong Zhang /* Allocate space for apj, initialize apj, and */ 711f4a743e1SHong Zhang /* destroy list of free space and other temporary array(s) */ 712854ce69bSBarry Smith ierr = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr); 71377584fe6SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 714118727c9SMark F. Adams afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1); 715d16ca5f0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 716f4a743e1SHong Zhang 7172addb229SHong Zhang /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/ 7182addb229SHong Zhang /*-----------------------------------------------------------------*/ 719de0260b3SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 72042c57489SHong Zhang 721ded4f561SHong Zhang /* then, compute symbolic Co = (p->B)^T*AP */ 722d0f46423SBarry Smith pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 72383445d95SHong Zhang >= (num of nonzero rows of C_seq) - pn */ 724854ce69bSBarry Smith ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 725de0260b3SHong Zhang coi[0] = 0; 72642c57489SHong Zhang 727d16ca5f0SHong Zhang /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */ 728d16ca5f0SHong Zhang nnz = fill*(poti[pon] + api[am]); 72922d28d08SBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 73042c57489SHong Zhang current_space = free_space; 73142c57489SHong Zhang 732de0260b3SHong Zhang for (i=0; i<pon; i++) { 73382412ba7SHong Zhang pnz = poti[i+1] - poti[i]; 73477584fe6SHong Zhang ptJ = potj + poti[i]; 73577584fe6SHong Zhang for (j=0; j<pnz; j++) { 73677584fe6SHong Zhang row = ptJ[j]; /* row of AP == col of Pot */ 73782412ba7SHong Zhang apnz = api[row+1] - api[row]; 738ded4f561SHong Zhang Jptr = apj + api[row]; 73983445d95SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 740a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 74142c57489SHong Zhang } 742a914927fSHong Zhang nnz = lnk[0]; 74342c57489SHong Zhang 74483445d95SHong Zhang /* If free space is not available, double the total space in the list */ 745ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 7464238b7adSHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 747d16ca5f0SHong Zhang nspacedouble++; 74842c57489SHong Zhang } 74942c57489SHong Zhang 75042c57489SHong Zhang /* Copy data into free space, and zero out denserows */ 751a914927fSHong Zhang ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 7522205254eSKarl Rupp 753ded4f561SHong Zhang current_space->array += nnz; 754ded4f561SHong Zhang current_space->local_used += nnz; 755ded4f561SHong Zhang current_space->local_remaining -= nnz; 7562205254eSKarl Rupp 757ded4f561SHong Zhang coi[i+1] = coi[i] + nnz; 75842c57489SHong Zhang } 7596b911d16SHong Zhang 7606b911d16SHong Zhang ierr = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr); 761a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 762118727c9SMark F. Adams afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1); 763d16ca5f0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 764de0260b3SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 76542c57489SHong Zhang 7666b911d16SHong Zhang /* (3) send j-array (coj) of Co to other processors */ 7676b911d16SHong Zhang /*--------------------------------------------------*/ 7686b911d16SHong Zhang ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr); 7696b911d16SHong Zhang len_s = merge->len_s; 7706b911d16SHong Zhang merge->nsend = 0; 7716b911d16SHong Zhang 7726b911d16SHong Zhang 77342c57489SHong Zhang /* determine row ownership */ 77426283091SBarry Smith ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 7757a2fc3feSBarry Smith merge->rowmap->n = pn; 7767a2fc3feSBarry Smith merge->rowmap->bs = 1; 7772205254eSKarl Rupp 77826283091SBarry Smith ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 7797a2fc3feSBarry Smith owners = merge->rowmap->range; 78042c57489SHong Zhang 78142c57489SHong Zhang /* determine the number of messages to send, their lengths */ 782dcca6d9dSJed Brown ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr); 78383445d95SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 784854ce69bSBarry Smith ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 785de0260b3SHong Zhang 78683445d95SHong Zhang proc = 0; 787de0260b3SHong Zhang for (i=0; i<pon; i++) { 788de0260b3SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 7892addb229SHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 790ee6e4b5bSHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 79183445d95SHong Zhang } 792de0260b3SHong Zhang 7932addb229SHong Zhang len = 0; /* max length of buf_si[], see (4) */ 794de0260b3SHong Zhang owners_co[0] = 0; 79542c57489SHong Zhang for (proc=0; proc<size; proc++) { 796de0260b3SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 797ee6e4b5bSHong Zhang if (len_s[proc]) { 79842c57489SHong Zhang merge->nsend++; 7992addb229SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 80042c57489SHong Zhang len += len_si[proc]; 80142c57489SHong Zhang } 80242c57489SHong Zhang } 80342c57489SHong Zhang 804ded4f561SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 8050298fd71SBarry Smith ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 80642c57489SHong Zhang ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 80742c57489SHong Zhang 808ded4f561SHong Zhang /* post the Irecv and Isend of coj */ 809529bc5dcSHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 810529bc5dcSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 811854ce69bSBarry Smith ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 81242c57489SHong Zhang for (proc=0, k=0; proc<size; proc++) { 81342c57489SHong Zhang if (!len_s[proc]) continue; 814de0260b3SHong Zhang i = owners_co[proc]; 815529bc5dcSHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 81642c57489SHong Zhang k++; 81742c57489SHong Zhang } 81842c57489SHong Zhang 819ded4f561SHong Zhang /* receives and sends of coj are complete */ 82077584fe6SHong Zhang for (i=0; i<merge->nrecv; i++) { 821c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 822ded4f561SHong Zhang } 823ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 8240c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 82582412ba7SHong Zhang 8266b911d16SHong Zhang /* (4) send and recv coi */ 8276b911d16SHong Zhang /*-----------------------*/ 828529bc5dcSHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 829529bc5dcSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 830854ce69bSBarry Smith ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 83142c57489SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 83242c57489SHong Zhang for (proc=0,k=0; proc<size; proc++) { 83342c57489SHong Zhang if (!len_s[proc]) continue; 83442c57489SHong Zhang /* form outgoing message for i-structure: 83542c57489SHong Zhang buf_si[0]: nrows to be sent 83642c57489SHong Zhang [1:nrows]: row index (global) 83742c57489SHong Zhang [nrows+1:2*nrows+1]: i-structure index 83842c57489SHong Zhang */ 83942c57489SHong Zhang /*-------------------------------------------*/ 8402addb229SHong Zhang nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 84142c57489SHong Zhang buf_si_i = buf_si + nrows+1; 84242c57489SHong Zhang buf_si[0] = nrows; 84342c57489SHong Zhang buf_si_i[0] = 0; 84442c57489SHong Zhang nrows = 0; 845de0260b3SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 846ded4f561SHong Zhang nzi = coi[i+1] - coi[i]; 847ded4f561SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 848de0260b3SHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 84942c57489SHong Zhang nrows++; 85042c57489SHong Zhang } 851529bc5dcSHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 85242c57489SHong Zhang k++; 85342c57489SHong Zhang buf_si += len_si[proc]; 85442c57489SHong Zhang } 855ded4f561SHong Zhang i = merge->nrecv; 856ded4f561SHong Zhang while (i--) { 857c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 858ded4f561SHong Zhang } 859ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 8600c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 861a914927fSHong Zhang 86224ecddacSHong Zhang ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr); 86342c57489SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 864ded4f561SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 86542c57489SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 86642c57489SHong Zhang 8676b911d16SHong Zhang /* (5) compute the local portion of C (mpi mat) */ 8686b911d16SHong Zhang /*----------------------------------------------*/ 869ded4f561SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 870ded4f561SHong Zhang 87124ecddacSHong Zhang /* allocate pti array and free space for accumulating nonzero column info */ 872854ce69bSBarry Smith ierr = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr); 87324ecddacSHong Zhang pti[0] = 0; 87442c57489SHong Zhang 875d16ca5f0SHong Zhang /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 876d16ca5f0SHong Zhang nnz = fill*(pi_loc[pm] + api[am]); 87722d28d08SBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 87842c57489SHong Zhang current_space = free_space; 87942c57489SHong Zhang 880dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 88142c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { 88242c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 88342c57489SHong Zhang nrows = *buf_ri_k[k]; 88442c57489SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 88542c57489SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 88642c57489SHong Zhang } 88742c57489SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 88808cb4532SHong Zhang rmax = 0; 88942c57489SHong Zhang for (i=0; i<pn; i++) { 890ded4f561SHong Zhang /* add pdt[i,:]*AP into lnk */ 891ded4f561SHong Zhang pnz = pdti[i+1] - pdti[i]; 89277584fe6SHong Zhang ptJ = pdtj + pdti[i]; 89377584fe6SHong Zhang for (j=0; j<pnz; j++) { 89477584fe6SHong Zhang row = ptJ[j]; /* row of AP == col of Pt */ 895ded4f561SHong Zhang apnz = api[row+1] - api[row]; 896ded4f561SHong Zhang Jptr = apj + api[row]; 897ded4f561SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 898a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 899ded4f561SHong Zhang } 900a914927fSHong Zhang 90142c57489SHong Zhang /* add received col data into lnk */ 90242c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 90342c57489SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 904ded4f561SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 905ded4f561SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 906a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 90742c57489SHong Zhang nextrow[k]++; nextci[k]++; 90842c57489SHong Zhang } 90942c57489SHong Zhang } 910a914927fSHong Zhang nnz = lnk[0]; 91142c57489SHong Zhang 91242c57489SHong Zhang /* if free space is not available, make more free space */ 913ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 9144238b7adSHong Zhang ierr = PetscFreeSpaceGet(nnz+current_space->total_array_size,¤t_space);CHKERRQ(ierr); 915d16ca5f0SHong Zhang nspacedouble++; 91642c57489SHong Zhang } 91742c57489SHong Zhang /* copy data into free space, then initialize lnk */ 918a914927fSHong Zhang ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 919ded4f561SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 9202205254eSKarl Rupp 921ded4f561SHong Zhang current_space->array += nnz; 922ded4f561SHong Zhang current_space->local_used += nnz; 923ded4f561SHong Zhang current_space->local_remaining -= nnz; 9242205254eSKarl Rupp 92524ecddacSHong Zhang pti[i+1] = pti[i] + nnz; 92608cb4532SHong Zhang if (nnz > rmax) rmax = nnz; 92742c57489SHong Zhang } 928ded4f561SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 9290572522cSBarry Smith ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 93042c57489SHong Zhang 931854ce69bSBarry Smith ierr = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr); 93224ecddacSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr); 93324ecddacSHong Zhang afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1); 934d16ca5f0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 93542c57489SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 93642c57489SHong Zhang 9376b911d16SHong Zhang /* (6) create symbolic parallel matrix Cmpi */ 9386b911d16SHong Zhang /*------------------------------------------*/ 93977584fe6SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 94077584fe6SHong Zhang ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 94133d57670SJed Brown ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 94277584fe6SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 94377584fe6SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 94442c57489SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 945a2f3521dSMark F. Adams 946b091e509SHong Zhang merge->bi = pti; /* Cseq->i */ 947b091e509SHong Zhang merge->bj = ptj; /* Cseq->j */ 948b091e509SHong Zhang merge->coi = coi; /* Co->i */ 949b091e509SHong Zhang merge->coj = coj; /* Co->j */ 95042c57489SHong Zhang merge->buf_ri = buf_ri; 95142c57489SHong Zhang merge->buf_rj = buf_rj; 952de0260b3SHong Zhang merge->owners_co = owners_co; 95377584fe6SHong Zhang merge->destroy = Cmpi->ops->destroy; 95477584fe6SHong Zhang merge->duplicate = Cmpi->ops->duplicate; 955cc6cb787SHong Zhang 95677584fe6SHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 957a914927fSHong Zhang Cmpi->assembled = PETSC_FALSE; 95877584fe6SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 959aa690a28SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP_old; 960aa690a28SHong Zhang Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_old; 96142c57489SHong Zhang 96277584fe6SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 96377584fe6SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 964f8487c73SHong Zhang c->ptap = ptap; 96577584fe6SHong Zhang ptap->api = api; 96677584fe6SHong Zhang ptap->apj = apj; 967d6ab1dc0SHong Zhang ptap->rmax = ap_rmax; 96877584fe6SHong Zhang *C = Cmpi; 969414894bdSHong Zhang 970414894bdSHong Zhang /* flag 'scalable' determines which implementations to be used: 971414894bdSHong Zhang 0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa; 972414894bdSHong Zhang 1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */ 973414894bdSHong Zhang /* set default scalable */ 974da0a95b2SSatish Balay ptap->scalable = PETSC_FALSE; /* PETSC_TRUE; */ 9752205254eSKarl Rupp 9760298fd71SBarry Smith ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr); 977414894bdSHong Zhang if (!ptap->scalable) { /* Do dense axpy */ 9781795a4d1SJed Brown ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 979414894bdSHong Zhang } else { 9801795a4d1SJed Brown ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr); 981414894bdSHong Zhang } 982414894bdSHong Zhang 983f2b054eeSHong Zhang #if defined(PETSC_USE_INFO) 98424ecddacSHong Zhang if (pti[pn] != 0) { 98557622a8eSBarry Smith ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 98657622a8eSBarry Smith ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 987f2b054eeSHong Zhang } else { 98877584fe6SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 989f2b054eeSHong Zhang } 990f2b054eeSHong Zhang #endif 99142c57489SHong Zhang PetscFunctionReturn(0); 99242c57489SHong Zhang } 99342c57489SHong Zhang 99442c57489SHong Zhang #undef __FUNCT__ 995aa690a28SHong Zhang #define __FUNCT__ "MatPtAPNumeric_MPIAIJ_MPIAIJ_old" 996aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_old(Mat A,Mat P,Mat C) 99742c57489SHong Zhang { 99842c57489SHong Zhang PetscErrorCode ierr; 999f8487c73SHong Zhang Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 100042c57489SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 1001de0260b3SHong Zhang Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 100242c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 1003f8487c73SHong Zhang Mat_PtAPMPI *ptap; 10049f4155fbSHong Zhang Mat_Merge_SeqsToMPI *merge; 10051d633602SHong Zhang PetscInt *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp; 100682412ba7SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj; 100782412ba7SHong Zhang PetscInt i,j,k,anz,pnz,apnz,nextap,row,*cj; 1008e900a757SHong Zhang MatScalar *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp; 1009d0f46423SBarry Smith PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1010ce94432eSBarry Smith MPI_Comm comm; 101142c57489SHong Zhang PetscMPIInt size,rank,taga,*len_s; 101282412ba7SHong Zhang PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 101342c57489SHong Zhang PetscInt **buf_ri,**buf_rj; 101450f5bf86SHong Zhang PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 101542c57489SHong Zhang MPI_Request *s_waits,*r_waits; 101642c57489SHong Zhang MPI_Status *status; 101782412ba7SHong Zhang MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 101882412ba7SHong Zhang PetscInt *api,*apj,*coi,*coj; 1019d0f46423SBarry Smith PetscInt *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend; 10206ce94e8fSHong Zhang PetscBool scalable; 102138571addSHong Zhang #if defined(PTAP_PROFILE) 1022e497d3c8SHong Zhang PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2; 102338571addSHong Zhang #endif 102442c57489SHong Zhang 102542c57489SHong Zhang PetscFunctionBegin; 1026ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 102738571addSHong Zhang #if defined(PTAP_PROFILE) 10288563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 102938571addSHong Zhang #endif 103042c57489SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 103142c57489SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 103242c57489SHong Zhang 1033f8487c73SHong Zhang ptap = c->ptap; 1034ce94432eSBarry 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"); 1035f8487c73SHong Zhang merge = ptap->merge; 1036414894bdSHong Zhang apa = ptap->apa; 10376ce94e8fSHong Zhang scalable = ptap->scalable; 10386c6e00e0SHong Zhang 1039a1a4e84aSHong Zhang /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 10406b911d16SHong Zhang /*-----------------------------------------------------*/ 1041f8487c73SHong Zhang if (ptap->reuse == MAT_INITIAL_MATRIX) { 10429f4155fbSHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */ 1043f8487c73SHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 10446c6e00e0SHong Zhang } else { /* update numerical values of P_oth and P_loc */ 1045b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1046a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 10476c6e00e0SHong Zhang } 104838571addSHong Zhang #if defined(PTAP_PROFILE) 10498563dfccSBarry Smith ierr = PetscTime(&t1);CHKERRQ(ierr); 105005d62848SHong Zhang eP = t1-t0; 105138571addSHong Zhang #endif 105205d62848SHong Zhang /* 105305d62848SHong Zhang printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank, 105405d62848SHong Zhang a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N, 105505d62848SHong Zhang ptap->P_loc->rmap->N,ptap->P_loc->cmap->N, 105605d62848SHong Zhang ptap->P_oth->rmap->N,ptap->P_oth->cmap->N); 105705d62848SHong Zhang */ 1058f8487c73SHong Zhang 1059f9473cf7SHong Zhang /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */ 1060f9473cf7SHong Zhang /*--------------------------------------------------------------*/ 106142c57489SHong Zhang /* get data from symbolic products */ 1062a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1063a1a4e84aSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 1064a9d06231SHong Zhang pi_loc=p_loc->i; pj_loc=p_loc->j; pJ=pj_loc; pa_loc=p_loc->a; 106542c57489SHong Zhang pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 106642c57489SHong Zhang 1067de0260b3SHong Zhang coi = merge->coi; coj = merge->coj; 10681795a4d1SJed Brown ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 1069de0260b3SHong Zhang 1070de0260b3SHong Zhang bi = merge->bi; bj = merge->bj; 10717a2fc3feSBarry Smith owners = merge->rowmap->range; 10721795a4d1SJed Brown ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); /* ba: Cseq->a */ 107342c57489SHong Zhang 1074a1a4e84aSHong Zhang api = ptap->api; apj = ptap->apj; 1075d9824c15SHong Zhang 10769516a85cSHong Zhang if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */ 1077b38d1a2dSJed Brown ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr); 107805d62848SHong Zhang #if 0 10790d3441aeSHong Zhang /* ------ 10x slower -------------- */ 10800d3441aeSHong Zhang /*==================================*/ 10819516a85cSHong Zhang Mat R = ptap->R; 10829516a85cSHong Zhang Mat_SeqAIJ *r = (Mat_SeqAIJ*)R->data; 10839516a85cSHong Zhang PetscInt *ri=r->i,*rj=r->j,rnz,arow,l,prow,pcol,pN=P->cmap->N; 10849516a85cSHong Zhang PetscScalar *ra=r->a,tmp,cdense[pN]; 10859516a85cSHong Zhang 10869516a85cSHong Zhang ierr = PetscMemzero(cdense,pN*sizeof(PetscScalar));CHKERRQ(ierr); 10879516a85cSHong Zhang for (i=0; i<cm; i++) { /* each row of C or R */ 10889516a85cSHong Zhang rnz = ri[i+1] - ri[i]; 10899516a85cSHong Zhang 10909516a85cSHong Zhang for (j=0; j<rnz; j++) { /* each nz of R */ 10919516a85cSHong Zhang arow = rj[ri[i] + j]; 10929516a85cSHong Zhang 10939516a85cSHong Zhang /* diagonal portion of A */ 10949516a85cSHong Zhang anz = ad->i[arow+1] - ad->i[arow]; 10959516a85cSHong Zhang for (k=0; k<anz; k++) { /* each nz of Ad */ 10969516a85cSHong Zhang tmp = ra[ri[i] + j]*ad->a[ad->i[arow] + k]; 10979516a85cSHong Zhang prow = ad->j[ad->i[arow] + k]; 10989516a85cSHong Zhang pnz = pi_loc[prow+1] - pi_loc[prow]; 10999516a85cSHong Zhang 11009516a85cSHong Zhang for (l=0; l<pnz; l++) { /* each nz of P_loc */ 11019516a85cSHong Zhang pcol = pj_loc[pi_loc[prow] + l]; 11029516a85cSHong Zhang cdense[pcol] += tmp*pa_loc[pi_loc[prow] + l]; 11039516a85cSHong Zhang } 11049516a85cSHong Zhang } 11059516a85cSHong Zhang 11069516a85cSHong Zhang /* off-diagonal portion of A */ 11079516a85cSHong Zhang anz = ao->i[arow+1] - ao->i[arow]; 11089516a85cSHong Zhang for (k=0; k<anz; k++) { /* each nz of Ao */ 11099516a85cSHong Zhang tmp = ra[ri[i] + j]*ao->a[ao->i[arow] + k]; 11109516a85cSHong Zhang prow = ao->j[ao->i[arow] + k]; 11119516a85cSHong Zhang pnz = pi_oth[prow+1] - pi_oth[prow]; 11129516a85cSHong Zhang 11139516a85cSHong Zhang for (l=0; l<pnz; l++) { /* each nz of P_oth */ 11149516a85cSHong Zhang pcol = pj_oth[pi_oth[prow] + l]; 11159516a85cSHong Zhang cdense[pcol] += tmp*pa_oth[pi_oth[prow] + l]; 11169516a85cSHong Zhang } 11179516a85cSHong Zhang } 11189516a85cSHong Zhang 1119da0a95b2SSatish Balay } /* for (j=0; j<rnz; j++) */ 11209516a85cSHong Zhang 11219516a85cSHong Zhang /* copy cdense[] into ca; zero cdense[] */ 11229516a85cSHong Zhang cnz = bi[i+1] - bi[i]; 11239516a85cSHong Zhang cj = bj + bi[i]; 11249516a85cSHong Zhang ca = ba + bi[i]; 11259516a85cSHong Zhang for (j=0; j<cnz; j++) { 11269516a85cSHong Zhang ca[j] += cdense[cj[j]]; 11279516a85cSHong Zhang cdense[cj[j]] = 0.0; 11289516a85cSHong Zhang } 11299516a85cSHong Zhang #if 0 11309516a85cSHong Zhang if (rank == 0) { 11319516a85cSHong Zhang printf("[%d] row %d: ",rank,i); 11329516a85cSHong Zhang for (j=0; j<pN; j++) printf(" %g,",cdense[j]); 11339516a85cSHong Zhang printf("\n"); 11349516a85cSHong Zhang } 1135da0a95b2SSatish Balay for (j=0; j<pN; j++) cdense[j]=0.0; /* zero cdnese[] */ 11369516a85cSHong Zhang #endif 1137da0a95b2SSatish Balay } /* for (i=0; i<cm; i++) { */ 113805d62848SHong Zhang #endif 11399516a85cSHong Zhang 1140f671be37SHong Zhang /* ========================================== */ 11418cb82516SHong Zhang #if defined(PTAP_PROFILE) 114205d62848SHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 11438cb82516SHong Zhang #endif 1144a9d06231SHong Zhang for (i=0; i<am; i++) { 1145b091e509SHong Zhang #if defined(PTAP_PROFILE) 11468563dfccSBarry Smith ierr = PetscTime(&t2_0);CHKERRQ(ierr); 1147b091e509SHong Zhang #endif 1148a9d06231SHong Zhang /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */ 1149a9d06231SHong Zhang /*------------------------------------------------------------*/ 1150a9d06231SHong Zhang apJ = apj + api[i]; 1151a9d06231SHong Zhang 1152a9d06231SHong Zhang /* diagonal portion of A */ 1153a9d06231SHong Zhang anz = adi[i+1] - adi[i]; 1154a9d06231SHong Zhang adj = ad->j + adi[i]; 1155a9d06231SHong Zhang ada = ad->a + adi[i]; 1156a9d06231SHong Zhang for (j=0; j<anz; j++) { 1157a9d06231SHong Zhang row = adj[j]; 1158a9d06231SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 1159a9d06231SHong Zhang pj = pj_loc + pi_loc[row]; 1160a9d06231SHong Zhang pa = pa_loc + pi_loc[row]; 1161a9d06231SHong Zhang 1162a9d06231SHong Zhang /* perform dense axpy */ 1163e900a757SHong Zhang valtmp = ada[j]; 1164a9d06231SHong Zhang for (k=0; k<pnz; k++) { 1165e900a757SHong Zhang apa[pj[k]] += valtmp*pa[k]; 1166a9d06231SHong Zhang } 1167a9d06231SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1168a9d06231SHong Zhang } 1169a9d06231SHong Zhang 1170a9d06231SHong Zhang /* off-diagonal portion of A */ 1171a9d06231SHong Zhang anz = aoi[i+1] - aoi[i]; 1172a9d06231SHong Zhang aoj = ao->j + aoi[i]; 1173a9d06231SHong Zhang aoa = ao->a + aoi[i]; 1174a9d06231SHong Zhang for (j=0; j<anz; j++) { 1175a9d06231SHong Zhang row = aoj[j]; 1176a9d06231SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 1177a9d06231SHong Zhang pj = pj_oth + pi_oth[row]; 1178a9d06231SHong Zhang pa = pa_oth + pi_oth[row]; 1179a9d06231SHong Zhang 1180a9d06231SHong Zhang /* perform dense axpy */ 1181e900a757SHong Zhang valtmp = aoa[j]; 1182a9d06231SHong Zhang for (k=0; k<pnz; k++) { 1183e900a757SHong Zhang apa[pj[k]] += valtmp*pa[k]; 1184a9d06231SHong Zhang } 1185a9d06231SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1186a9d06231SHong Zhang } 1187b091e509SHong Zhang #if defined(PTAP_PROFILE) 11888563dfccSBarry Smith ierr = PetscTime(&t2_1);CHKERRQ(ierr); 1189b091e509SHong Zhang et2_AP += t2_1 - t2_0; 1190b091e509SHong Zhang #endif 1191a9d06231SHong Zhang 1192a9d06231SHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 1193a9d06231SHong Zhang /*--------------------------------------------------------------*/ 1194a9d06231SHong Zhang apnz = api[i+1] - api[i]; 1195a9d06231SHong Zhang /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */ 1196a9d06231SHong Zhang pnz = po->i[i+1] - po->i[i]; 1197e900a757SHong Zhang poJ = po->j + po->i[i]; 1198a9d06231SHong Zhang pA = po->a + po->i[i]; 1199a9d06231SHong Zhang for (j=0; j<pnz; j++) { 1200e900a757SHong Zhang row = poJ[j]; 1201e900a757SHong Zhang cnz = coi[row+1] - coi[row]; 1202e900a757SHong Zhang cj = coj + coi[row]; 1203e900a757SHong Zhang ca = coa + coi[row]; 1204a9d06231SHong Zhang /* perform dense axpy */ 1205e900a757SHong Zhang valtmp = pA[j]; 1206a9d06231SHong Zhang for (k=0; k<cnz; k++) { 1207e900a757SHong Zhang ca[k] += valtmp*apa[cj[k]]; 1208a9d06231SHong Zhang } 1209a9d06231SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1210a9d06231SHong Zhang } 121105d62848SHong Zhang #if 1 1212a9d06231SHong Zhang /* put the value into Cd (diagonal part) */ 1213a9d06231SHong Zhang pnz = pd->i[i+1] - pd->i[i]; 1214e900a757SHong Zhang pdJ = pd->j + pd->i[i]; 1215a9d06231SHong Zhang pA = pd->a + pd->i[i]; 1216a9d06231SHong Zhang for (j=0; j<pnz; j++) { 1217e900a757SHong Zhang row = pdJ[j]; 1218e900a757SHong Zhang cnz = bi[row+1] - bi[row]; 1219e900a757SHong Zhang cj = bj + bi[row]; 1220e900a757SHong Zhang ca = ba + bi[row]; 1221a9d06231SHong Zhang /* perform dense axpy */ 1222e900a757SHong Zhang valtmp = pA[j]; 1223a9d06231SHong Zhang for (k=0; k<cnz; k++) { 1224e900a757SHong Zhang ca[k] += valtmp*apa[cj[k]]; 1225a9d06231SHong Zhang } 1226a9d06231SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1227a9d06231SHong Zhang } 12289516a85cSHong Zhang #endif 1229a9d06231SHong Zhang /* zero the current row of A*P */ 1230a9d06231SHong Zhang for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0; 1231b091e509SHong Zhang #if defined(PTAP_PROFILE) 12328563dfccSBarry Smith ierr = PetscTime(&t2_2);CHKERRQ(ierr); 123305d62848SHong Zhang ePtAP += t2_2 - t2_1; 1234b091e509SHong Zhang #endif 1235a9d06231SHong Zhang } 123638571addSHong Zhang } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */ 1237b38d1a2dSJed Brown ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr); 123838571addSHong Zhang /*-----------------------------------------------------------------------------------------*/ 1239a9d06231SHong Zhang pA=pa_loc; 124042c57489SHong Zhang for (i=0; i<am; i++) { 1241b091e509SHong Zhang #if defined(PTAP_PROFILE) 12428563dfccSBarry Smith ierr = PetscTime(&t2_0);CHKERRQ(ierr); 1243b091e509SHong Zhang #endif 1244f4a743e1SHong Zhang /* form i-th sparse row of A*P */ 124582412ba7SHong Zhang apnz = api[i+1] - api[i]; 124682412ba7SHong Zhang apJ = apj + api[i]; 124742c57489SHong Zhang /* diagonal portion of A */ 124882412ba7SHong Zhang anz = adi[i+1] - adi[i]; 12491d633602SHong Zhang adj = ad->j + adi[i]; 12501d633602SHong Zhang ada = ad->a + adi[i]; 125182412ba7SHong Zhang for (j=0; j<anz; j++) { 12521d633602SHong Zhang row = adj[j]; 125382412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 125482412ba7SHong Zhang pj = pj_loc + pi_loc[row]; 125582412ba7SHong Zhang pa = pa_loc + pi_loc[row]; 1256e900a757SHong Zhang valtmp = ada[j]; 1257f4a743e1SHong Zhang nextp = 0; 125882412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 125982412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 1260e900a757SHong Zhang apa[k] += valtmp*pa[nextp++]; 126142c57489SHong Zhang } 126242c57489SHong Zhang } 1263dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 126442c57489SHong Zhang } 126542c57489SHong Zhang /* off-diagonal portion of A */ 126682412ba7SHong Zhang anz = aoi[i+1] - aoi[i]; 12671d633602SHong Zhang aoj = ao->j + aoi[i]; 12681d633602SHong Zhang aoa = ao->a + aoi[i]; 126982412ba7SHong Zhang for (j=0; j<anz; j++) { 12701d633602SHong Zhang row = aoj[j]; 127182412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 127282412ba7SHong Zhang pj = pj_oth + pi_oth[row]; 127382412ba7SHong Zhang pa = pa_oth + pi_oth[row]; 1274e900a757SHong Zhang valtmp = aoa[j]; 1275f4a743e1SHong Zhang nextp = 0; 127682412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 127782412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 1278e900a757SHong Zhang apa[k] += valtmp*pa[nextp++]; 127942c57489SHong Zhang } 128042c57489SHong Zhang } 1281dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 128242c57489SHong Zhang } 1283b091e509SHong Zhang #if defined(PTAP_PROFILE) 12848563dfccSBarry Smith ierr = PetscTime(&t2_1);CHKERRQ(ierr); 1285b091e509SHong Zhang et2_AP += t2_1 - t2_0; 1286b091e509SHong Zhang #endif 128742c57489SHong Zhang 1288a9d06231SHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 1289a9d06231SHong Zhang /*--------------------------------------------------------------*/ 129082412ba7SHong Zhang pnz = pi_loc[i+1] - pi_loc[i]; 1291f9473cf7SHong Zhang pJ = pj_loc + pi_loc[i]; 129282412ba7SHong Zhang for (j=0; j<pnz; j++) { 129342c57489SHong Zhang nextap = 0; 1294f9473cf7SHong Zhang row = pJ[j]; /* global index */ 129582412ba7SHong Zhang if (row < pcstart || row >=pcend) { /* put the value into Co */ 1296e900a757SHong Zhang row = *poJ; 1297e900a757SHong Zhang cj = coj + coi[row]; 1298e900a757SHong Zhang ca = coa + coi[row]; poJ++; 129982412ba7SHong Zhang } else { /* put the value into Cd */ 1300e900a757SHong Zhang row = *pdJ; 1301e900a757SHong Zhang cj = bj + bi[row]; 1302e900a757SHong Zhang ca = ba + bi[row]; pdJ++; 130342c57489SHong Zhang } 1304e900a757SHong Zhang valtmp = pA[j]; 130582412ba7SHong Zhang for (k=0; nextap<apnz; k++) { 1306e900a757SHong Zhang if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++]; 130742c57489SHong Zhang } 1308dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 1309de0260b3SHong Zhang } 13100496f32cSHong Zhang pA += pnz; 131142c57489SHong Zhang /* zero the current row info for A*P */ 131282412ba7SHong Zhang ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr); 1313b091e509SHong Zhang #if defined(PTAP_PROFILE) 13148563dfccSBarry Smith ierr = PetscTime(&t2_2);CHKERRQ(ierr); 131505d62848SHong Zhang ePtAP += t2_2 - t2_1; 1316b091e509SHong Zhang #endif 131742c57489SHong Zhang } 1318d9824c15SHong Zhang } 131938571addSHong Zhang #if defined(PTAP_PROFILE) 13208563dfccSBarry Smith ierr = PetscTime(&t2);CHKERRQ(ierr); 132138571addSHong Zhang #endif 132242c57489SHong Zhang 1323a9d06231SHong Zhang /* 3) send and recv matrix values coa */ 1324a9d06231SHong Zhang /*------------------------------------*/ 132542c57489SHong Zhang buf_ri = merge->buf_ri; 132642c57489SHong Zhang buf_rj = merge->buf_rj; 132742c57489SHong Zhang len_s = merge->len_s; 1328fc42d0c8SSatish Balay ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 132942c57489SHong Zhang ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 133042c57489SHong Zhang 1331dcca6d9dSJed Brown ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 133242c57489SHong Zhang for (proc=0,k=0; proc<size; proc++) { 133342c57489SHong Zhang if (!len_s[proc]) continue; 1334de0260b3SHong Zhang i = merge->owners_co[proc]; 1335de0260b3SHong Zhang ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 133642c57489SHong Zhang k++; 133742c57489SHong Zhang } 13380c468ba9SBarry Smith if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 13390c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 134042c57489SHong Zhang 1341a9d06231SHong Zhang ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 134242c57489SHong Zhang ierr = PetscFree(r_waits);CHKERRQ(ierr); 134382412ba7SHong Zhang ierr = PetscFree(coa);CHKERRQ(ierr); 134438571addSHong Zhang #if defined(PTAP_PROFILE) 13458563dfccSBarry Smith ierr = PetscTime(&t3);CHKERRQ(ierr); 134638571addSHong Zhang #endif 134742c57489SHong Zhang 1348a9d06231SHong Zhang /* 4) insert local Cseq and received values into Cmpi */ 1349a9d06231SHong Zhang /*------------------------------------------------------*/ 1350dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 135142c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { 135242c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 135342c57489SHong Zhang nrows = *(buf_ri_k[k]); 135442c57489SHong Zhang nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 135582412ba7SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 135642c57489SHong Zhang } 135742c57489SHong Zhang 1358de0260b3SHong Zhang for (i=0; i<cm; i++) { 135982412ba7SHong Zhang row = owners[rank] + i; /* global row index of C_seq */ 136042c57489SHong Zhang bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1361de0260b3SHong Zhang ba_i = ba + bi[i]; 136282412ba7SHong Zhang bnz = bi[i+1] - bi[i]; 136342c57489SHong Zhang /* add received vals into ba */ 136442c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 136542c57489SHong Zhang /* i-th row */ 136642c57489SHong Zhang if (i == *nextrow[k]) { 136782412ba7SHong Zhang cnz = *(nextci[k]+1) - *nextci[k]; 136882412ba7SHong Zhang cj = buf_rj[k] + *(nextci[k]); 136982412ba7SHong Zhang ca = abuf_r[k] + *(nextci[k]); 137082412ba7SHong Zhang nextcj = 0; 137182412ba7SHong Zhang for (j=0; nextcj<cnz; j++) { 137282412ba7SHong Zhang if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 137382412ba7SHong Zhang ba_i[j] += ca[nextcj++]; 137442c57489SHong Zhang } 137542c57489SHong Zhang } 137682412ba7SHong Zhang nextrow[k]++; nextci[k]++; 1377c5af039cSHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 137842c57489SHong Zhang } 137942c57489SHong Zhang } 138082412ba7SHong Zhang ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 138142c57489SHong Zhang } 138242c57489SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 138342c57489SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 138442c57489SHong Zhang 1385de0260b3SHong Zhang ierr = PetscFree(ba);CHKERRQ(ierr); 1386c05d87d6SBarry Smith ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 138742c57489SHong Zhang ierr = PetscFree(abuf_r);CHKERRQ(ierr); 13880572522cSBarry Smith ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 138938571addSHong Zhang #if defined(PTAP_PROFILE) 13908563dfccSBarry Smith ierr = PetscTime(&t4);CHKERRQ(ierr); 13919516a85cSHong Zhang if (rank==1) { 139205d62848SHong 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); 13939516a85cSHong Zhang } 139438571addSHong Zhang #endif 139542c57489SHong Zhang PetscFunctionReturn(0); 139642c57489SHong Zhang } 1397