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 15cf97cf3cSHong Zhang PetscErrorCode MatView_MPIAIJ_PtAP(Mat A,PetscViewer viewer) 16cf97cf3cSHong Zhang { 17cf97cf3cSHong Zhang PetscErrorCode ierr; 18cf97cf3cSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 19cf97cf3cSHong Zhang Mat_PtAPMPI *ptap=a->ptap; 20cf97cf3cSHong Zhang PetscBool iascii; 21cf97cf3cSHong Zhang PetscViewerFormat format; 22cf97cf3cSHong Zhang 23cf97cf3cSHong Zhang PetscFunctionBegin; 24cf97cf3cSHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 25cf97cf3cSHong Zhang if (iascii) { 26cf97cf3cSHong Zhang ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 27cf97cf3cSHong Zhang if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 28cf97cf3cSHong Zhang if (ptap->algType == 0) { 29cf97cf3cSHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"using scalable MatPtAP() implementation\n");CHKERRQ(ierr); 30cf97cf3cSHong Zhang } else if (ptap->algType == 1) { 31cf97cf3cSHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"using nonscalable MatPtAP() implementation\n");CHKERRQ(ierr); 32cf97cf3cSHong Zhang } 33cf97cf3cSHong Zhang } 34cf97cf3cSHong Zhang } 35cf97cf3cSHong Zhang ierr = (ptap->view)(A,viewer);CHKERRQ(ierr); 36cf97cf3cSHong Zhang PetscFunctionReturn(0); 37cf97cf3cSHong Zhang } 38cf97cf3cSHong Zhang 39f8487c73SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A) 40cc6cb787SHong Zhang { 41cc6cb787SHong Zhang PetscErrorCode ierr; 42f8487c73SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 43f8487c73SHong Zhang Mat_PtAPMPI *ptap=a->ptap; 44cc6cb787SHong Zhang 45cc6cb787SHong Zhang PetscFunctionBegin; 46f8487c73SHong Zhang if (ptap) { 47c8b0795cSMark F. Adams Mat_Merge_SeqsToMPI *merge=ptap->merge; 48b7f45c76SHong Zhang ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr); 49f8487c73SHong Zhang ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 50a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 51a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 52c5af039cSHong Zhang ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */ 5305d62848SHong Zhang ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr); 5405d62848SHong Zhang ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr); 558403a639SHong Zhang if (ptap->AP_loc) { /* used by alg_rap */ 56681d504bSHong Zhang Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data; 57681d504bSHong Zhang ierr = PetscFree(ap->i);CHKERRQ(ierr); 58681d504bSHong Zhang ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr); 590d3441aeSHong Zhang ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr); 608403a639SHong Zhang } else { /* used by alg_ptap */ 618403a639SHong Zhang ierr = PetscFree(ptap->api);CHKERRQ(ierr); 628403a639SHong Zhang ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 63681d504bSHong Zhang } 642259aa2eSHong Zhang ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr); 652259aa2eSHong Zhang ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr); 66414894bdSHong Zhang if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);} 67a560ef98SHong Zhang 688403a639SHong Zhang if (merge) { /* used by alg_ptap */ 69cc6cb787SHong Zhang ierr = PetscFree(merge->id_r);CHKERRQ(ierr); 70cc6cb787SHong Zhang ierr = PetscFree(merge->len_s);CHKERRQ(ierr); 71cc6cb787SHong Zhang ierr = PetscFree(merge->len_r);CHKERRQ(ierr); 72cc6cb787SHong Zhang ierr = PetscFree(merge->bi);CHKERRQ(ierr); 73cc6cb787SHong Zhang ierr = PetscFree(merge->bj);CHKERRQ(ierr); 74c05d87d6SBarry Smith ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr); 75cc6cb787SHong Zhang ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr); 76c05d87d6SBarry Smith ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr); 77cc6cb787SHong Zhang ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr); 78445158ffSHong Zhang ierr = PetscFree(merge->coi);CHKERRQ(ierr); 79445158ffSHong Zhang ierr = PetscFree(merge->coj);CHKERRQ(ierr); 8005b42c5fSBarry Smith ierr = PetscFree(merge->owners_co);CHKERRQ(ierr); 816bf464f9SBarry Smith ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr); 82f8487c73SHong Zhang ierr = PetscFree(ptap->merge);CHKERRQ(ierr); 83bf0cc555SLisandro Dalcin } 84a560ef98SHong Zhang 85a560ef98SHong Zhang ierr = ptap->destroy(A);CHKERRQ(ierr); 86f8487c73SHong Zhang ierr = PetscFree(ptap);CHKERRQ(ierr); 87c8b0795cSMark F. Adams } 88cc6cb787SHong Zhang PetscFunctionReturn(0); 89cc6cb787SHong Zhang } 90cc6cb787SHong Zhang 91aa690a28SHong Zhang PetscErrorCode MatDuplicate_MPIAIJ_MatPtAP(Mat A, MatDuplicateOption op, Mat *M) 921a47ec54SHong Zhang { 931a47ec54SHong Zhang PetscErrorCode ierr; 941a47ec54SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 951a47ec54SHong Zhang Mat_PtAPMPI *ptap = a->ptap; 961a47ec54SHong Zhang 971a47ec54SHong Zhang PetscFunctionBegin; 981a47ec54SHong Zhang ierr = (*ptap->duplicate)(A,op,M);CHKERRQ(ierr); 991a47ec54SHong Zhang (*M)->ops->destroy = ptap->destroy; 1001a47ec54SHong Zhang (*M)->ops->duplicate = ptap->duplicate; 101*cff31925SHong Zhang (*M)->ops->view = ptap->view; 1021a47ec54SHong Zhang PetscFunctionReturn(0); 1031a47ec54SHong Zhang } 1041a47ec54SHong Zhang 105150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 10642c57489SHong Zhang { 10742c57489SHong Zhang PetscErrorCode ierr; 1088403a639SHong Zhang PetscBool rap=PETSC_TRUE; /* do R=P^T locally, then C=R*A*P */ 10967a12041SHong Zhang MPI_Comm comm; 11042c57489SHong Zhang 11142c57489SHong Zhang PetscFunctionBegin; 11267a12041SHong Zhang /* check if matrix local sizes are compatible */ 11367a12041SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1146c4ed002SBarry Smith if (A->rmap->rstart != P->rmap->rstart || A->rmap->rend != P->rmap->rend) 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); 1156c4ed002SBarry Smith if (A->cmap->rstart != P->rmap->rstart || A->cmap->rend != P->rmap->rend) 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 117c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,NULL,"-matrap",&rap,NULL);CHKERRQ(ierr); 118cf3ca8ceSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 1193ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 1208403a639SHong Zhang if (rap) { /* do R=P^T locally, then C=R*A*P */ 121cf3ca8ceSHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr); 1228403a639SHong Zhang } else { /* do P^T*A*P */ 1238403a639SHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(A,P,fill,C);CHKERRQ(ierr); 1240d3441aeSHong Zhang } 1253ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 1267a7894deSKris Buschelman } 1273ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 1288403a639SHong Zhang ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr); 1293ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 13042c57489SHong Zhang PetscFunctionReturn(0); 13142c57489SHong Zhang } 13242c57489SHong Zhang 133ae5f0867Sstefano_zampini #if defined(PETSC_HAVE_HYPRE) 134ae5f0867Sstefano_zampini PETSC_INTERN PetscErrorCode MatPtAPSymbolic_AIJ_AIJ_wHYPRE(Mat,Mat,PetscReal,Mat*); 135ae5f0867Sstefano_zampini #endif 136ae5f0867Sstefano_zampini 137cf742fccSHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,Mat C) 138cf742fccSHong Zhang { 139cf742fccSHong Zhang PetscErrorCode ierr; 140cf742fccSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 141cf742fccSHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 142cf742fccSHong Zhang Mat_SeqAIJ *ap,*p_loc,*p_oth,*c_seq; 143cf742fccSHong Zhang Mat_PtAPMPI *ptap = c->ptap; 144cf742fccSHong Zhang Mat AP_loc,C_loc,C_oth; 1455ca39647SHong Zhang PetscInt i,rstart,rend,cm,ncols,row,*api,*apj,am = A->rmap->n,apnz; 146cf742fccSHong Zhang PetscScalar *apa; 147cf742fccSHong Zhang const PetscInt *cols; 148cf742fccSHong Zhang const PetscScalar *vals; 149cf742fccSHong Zhang 150cf742fccSHong Zhang PetscFunctionBegin; 151cf742fccSHong Zhang ierr = MatZeroEntries(C);CHKERRQ(ierr); 152cf742fccSHong Zhang 153cf742fccSHong Zhang /* 1) get R = Pd^T,Ro = Po^T */ 154cf742fccSHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 155cf742fccSHong Zhang ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 156cf742fccSHong Zhang ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 157cf742fccSHong Zhang } 158cf742fccSHong Zhang 159cf742fccSHong Zhang /* 2) get AP_loc */ 160cf742fccSHong Zhang AP_loc = ptap->AP_loc; 161cf742fccSHong Zhang ap = (Mat_SeqAIJ*)AP_loc->data; 162cf742fccSHong Zhang 163cf742fccSHong Zhang /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 164cf742fccSHong Zhang /*-----------------------------------------------------*/ 165cf742fccSHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 166cf742fccSHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 167cf742fccSHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 168cf742fccSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 169cf742fccSHong Zhang } 170cf742fccSHong Zhang 171cf742fccSHong Zhang /* 2-2) compute numeric A_loc*P - dominating part */ 172cf742fccSHong Zhang /* ---------------------------------------------- */ 173cf742fccSHong Zhang /* get data from symbolic products */ 174cf742fccSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 17552f7967eSHong Zhang if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 176c072d3e2SSatish Balay else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ptap->P_oth is NULL. Cannot proceed!"); 17752f7967eSHong Zhang 178cf742fccSHong Zhang api = ap->i; 179cf742fccSHong Zhang apj = ap->j; 180cf742fccSHong Zhang for (i=0; i<am; i++) { 181cf742fccSHong Zhang /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 182cf742fccSHong Zhang apnz = api[i+1] - api[i]; 183b4e8d1b6SHong Zhang apa = ap->a + api[i]; 184b4e8d1b6SHong Zhang ierr = PetscMemzero(apa,sizeof(PetscScalar)*apnz);CHKERRQ(ierr); 185b4e8d1b6SHong Zhang AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa); 186cf742fccSHong Zhang ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 187cf742fccSHong Zhang } 188cf742fccSHong Zhang 189cf742fccSHong Zhang /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 190cf742fccSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 191cf742fccSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 192cf742fccSHong Zhang 193cf742fccSHong Zhang C_loc = ptap->C_loc; 194cf742fccSHong Zhang C_oth = ptap->C_oth; 195cf742fccSHong Zhang 196cf742fccSHong Zhang /* add C_loc and Co to to C */ 197cf742fccSHong Zhang ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 198cf742fccSHong Zhang 199cf742fccSHong Zhang /* C_loc -> C */ 200cf742fccSHong Zhang cm = C_loc->rmap->N; 201cf742fccSHong Zhang c_seq = (Mat_SeqAIJ*)C_loc->data; 202cf742fccSHong Zhang cols = c_seq->j; 203cf742fccSHong Zhang vals = c_seq->a; 204cf742fccSHong Zhang for (i=0; i<cm; i++) { 205cf742fccSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 206cf742fccSHong Zhang row = rstart + i; 207cf742fccSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 208cf742fccSHong Zhang cols += ncols; vals += ncols; 209cf742fccSHong Zhang } 210cf742fccSHong Zhang 211cf742fccSHong Zhang /* Co -> C, off-processor part */ 212cf742fccSHong Zhang cm = C_oth->rmap->N; 213cf742fccSHong Zhang c_seq = (Mat_SeqAIJ*)C_oth->data; 214cf742fccSHong Zhang cols = c_seq->j; 215cf742fccSHong Zhang vals = c_seq->a; 216cf742fccSHong Zhang for (i=0; i<cm; i++) { 217cf742fccSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 218cf742fccSHong Zhang row = p->garray[i]; 219cf742fccSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 220cf742fccSHong Zhang cols += ncols; vals += ncols; 221cf742fccSHong Zhang } 222cf742fccSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 223cf742fccSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 224cf742fccSHong Zhang 225cf742fccSHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 226cf742fccSHong Zhang PetscFunctionReturn(0); 227cf742fccSHong Zhang } 228cf742fccSHong Zhang 229e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat *C) 2300d3441aeSHong Zhang { 2310d3441aeSHong Zhang PetscErrorCode ierr; 2320d3441aeSHong Zhang Mat_PtAPMPI *ptap; 233681d504bSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 2342259aa2eSHong Zhang MPI_Comm comm; 2352259aa2eSHong Zhang PetscMPIInt size,rank; 23676db11f6SHong Zhang Mat Cmpi,P_loc,P_oth; 23715a3b8e2SHong Zhang PetscFreeSpaceList free_space=NULL,current_space=NULL; 238d0e9a2f7SHong Zhang PetscInt am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n; 239d0e9a2f7SHong Zhang PetscInt *lnk,i,k,pnz,row,nsend; 240f671be37SHong Zhang PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 24115a3b8e2SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 242d0e9a2f7SHong Zhang PetscInt len,proc,*dnz,*onz,*owners,nzi,nspacedouble; 24315a3b8e2SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 24415a3b8e2SHong Zhang MPI_Request *swaits,*rwaits; 24515a3b8e2SHong Zhang MPI_Status *sstatus,rstatus; 246445158ffSHong Zhang PetscLayout rowmap; 247445158ffSHong Zhang PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 248445158ffSHong Zhang PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 249b4e8d1b6SHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,Crmax,*aj,*ai,*pi; 25052f7967eSHong Zhang Mat_SeqAIJ *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth; 25167a12041SHong Zhang PetscScalar *apv; 25278d30b94SHong Zhang PetscTable ta; 253aa690a28SHong Zhang #if defined(PETSC_USE_INFO) 2548cb82516SHong Zhang PetscReal apfill; 255aa690a28SHong Zhang #endif 25667a12041SHong Zhang 25767a12041SHong Zhang PetscFunctionBegin; 25867a12041SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 25967a12041SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 26067a12041SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 261ae5f0867Sstefano_zampini 26252f7967eSHong Zhang if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data; 26352f7967eSHong Zhang 264ae5f0867Sstefano_zampini /* create symbolic parallel matrix Cmpi */ 265ae5f0867Sstefano_zampini ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 266ae5f0867Sstefano_zampini ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 267ae5f0867Sstefano_zampini 268e953e47cSHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 269e953e47cSHong Zhang ierr = PetscNew(&ptap);CHKERRQ(ierr); 270e953e47cSHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 271cf97cf3cSHong Zhang ptap->algType = 0; 272e953e47cSHong Zhang 273e953e47cSHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 27476db11f6SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&P_oth);CHKERRQ(ierr); 275e953e47cSHong Zhang /* get P_loc by taking all local rows of P */ 27676db11f6SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&P_loc);CHKERRQ(ierr); 27776db11f6SHong Zhang 27876db11f6SHong Zhang ptap->P_loc = P_loc; 27976db11f6SHong Zhang ptap->P_oth = P_oth; 280e953e47cSHong Zhang 281e953e47cSHong Zhang /* (0) compute Rd = Pd^T, Ro = Po^T */ 282e953e47cSHong Zhang /* --------------------------------- */ 283e953e47cSHong Zhang ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 284e953e47cSHong Zhang ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 285e953e47cSHong Zhang 286e953e47cSHong Zhang /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 287e953e47cSHong Zhang /* ----------------------------------------------------------------- */ 28876db11f6SHong Zhang p_loc = (Mat_SeqAIJ*)P_loc->data; 28952f7967eSHong Zhang if (P_oth) p_oth = (Mat_SeqAIJ*)P_oth->data; 290e953e47cSHong Zhang 291e953e47cSHong Zhang /* create and initialize a linked list */ 292e953e47cSHong Zhang ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */ 29376db11f6SHong Zhang MatRowMergeMax_SeqAIJ(p_loc,P_loc->rmap->N,ta); 29476db11f6SHong Zhang MatRowMergeMax_SeqAIJ(p_oth,P_oth->rmap->N,ta); 295e953e47cSHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */ 296e953e47cSHong Zhang 29776db11f6SHong Zhang ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 298e953e47cSHong Zhang 299e953e47cSHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 30052f7967eSHong Zhang if (ao) { 301e953e47cSHong Zhang ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr); 30252f7967eSHong Zhang } else { 30352f7967eSHong Zhang ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr); 30452f7967eSHong Zhang } 305e953e47cSHong Zhang current_space = free_space; 306e953e47cSHong Zhang nspacedouble = 0; 307e953e47cSHong Zhang 308e953e47cSHong Zhang ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 309e953e47cSHong Zhang api[0] = 0; 310e953e47cSHong Zhang for (i=0; i<am; i++) { 311e953e47cSHong Zhang /* diagonal portion: Ad[i,:]*P */ 312e953e47cSHong Zhang ai = ad->i; pi = p_loc->i; 313e953e47cSHong Zhang nzi = ai[i+1] - ai[i]; 314e953e47cSHong Zhang aj = ad->j + ai[i]; 315e953e47cSHong Zhang for (j=0; j<nzi; j++) { 316e953e47cSHong Zhang row = aj[j]; 317e953e47cSHong Zhang pnz = pi[row+1] - pi[row]; 318e953e47cSHong Zhang Jptr = p_loc->j + pi[row]; 319e953e47cSHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 32076db11f6SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 321e953e47cSHong Zhang } 322e953e47cSHong Zhang /* off-diagonal portion: Ao[i,:]*P */ 32352f7967eSHong Zhang if (ao) { 324e953e47cSHong Zhang ai = ao->i; pi = p_oth->i; 325e953e47cSHong Zhang nzi = ai[i+1] - ai[i]; 326e953e47cSHong Zhang aj = ao->j + ai[i]; 327e953e47cSHong Zhang for (j=0; j<nzi; j++) { 328e953e47cSHong Zhang row = aj[j]; 329e953e47cSHong Zhang pnz = pi[row+1] - pi[row]; 330e953e47cSHong Zhang Jptr = p_oth->j + pi[row]; 33176db11f6SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 332e953e47cSHong Zhang } 33352f7967eSHong Zhang } 334e953e47cSHong Zhang apnz = lnk[0]; 335e953e47cSHong Zhang api[i+1] = api[i] + apnz; 336e953e47cSHong Zhang 337e953e47cSHong Zhang /* if free space is not available, double the total space in the list */ 338e953e47cSHong Zhang if (current_space->local_remaining<apnz) { 339e953e47cSHong Zhang ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 340e953e47cSHong Zhang nspacedouble++; 341e953e47cSHong Zhang } 342e953e47cSHong Zhang 343e953e47cSHong Zhang /* Copy data into free space, then initialize lnk */ 34476db11f6SHong Zhang ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr); 345e953e47cSHong Zhang 346e953e47cSHong Zhang current_space->array += apnz; 347e953e47cSHong Zhang current_space->local_used += apnz; 348e953e47cSHong Zhang current_space->local_remaining -= apnz; 349e953e47cSHong Zhang } 350e953e47cSHong Zhang /* Allocate space for apj and apv, initialize apj, and */ 351e953e47cSHong Zhang /* destroy list of free space and other temporary array(s) */ 352e953e47cSHong Zhang ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 353e953e47cSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 35476db11f6SHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 355e953e47cSHong Zhang 356e953e47cSHong Zhang /* Create AP_loc for reuse */ 357e953e47cSHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 358e953e47cSHong Zhang 359e953e47cSHong Zhang #if defined(PETSC_USE_INFO) 36052f7967eSHong Zhang if (ao) { 361e953e47cSHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 36252f7967eSHong Zhang } else { 36352f7967eSHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1); 36452f7967eSHong Zhang } 365e953e47cSHong Zhang ptap->AP_loc->info.mallocs = nspacedouble; 366e953e47cSHong Zhang ptap->AP_loc->info.fill_ratio_given = fill; 367e953e47cSHong Zhang ptap->AP_loc->info.fill_ratio_needed = apfill; 368e953e47cSHong Zhang 369e953e47cSHong Zhang if (api[am]) { 370b11c1ec8SHong Zhang ierr = PetscInfo3(ptap->AP_loc,"Scalable algorithm, AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr); 371e953e47cSHong Zhang ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 372e953e47cSHong Zhang } else { 373b11c1ec8SHong Zhang ierr = PetscInfo(ptap->AP_loc,"Scalable algorithm, AP_loc is empty \n");CHKERRQ(ierr); 374e953e47cSHong Zhang } 375e953e47cSHong Zhang #endif 376e953e47cSHong Zhang 377e953e47cSHong Zhang /* (2-1) compute symbolic Co = Ro*AP_loc */ 378e953e47cSHong Zhang /* ------------------------------------ */ 379e953e47cSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 380e953e47cSHong Zhang 381e953e47cSHong Zhang /* (3) send coj of C_oth to other processors */ 382e953e47cSHong Zhang /* ------------------------------------------ */ 383e953e47cSHong Zhang /* determine row ownership */ 384e953e47cSHong Zhang ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 385e953e47cSHong Zhang rowmap->n = pn; 386e953e47cSHong Zhang rowmap->bs = 1; 387e953e47cSHong Zhang ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 388e953e47cSHong Zhang owners = rowmap->range; 389e953e47cSHong Zhang 390e953e47cSHong Zhang /* determine the number of messages to send, their lengths */ 391e953e47cSHong Zhang ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 392e953e47cSHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 393e953e47cSHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 394e953e47cSHong Zhang 395e953e47cSHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 396e953e47cSHong Zhang coi = c_oth->i; coj = c_oth->j; 397e953e47cSHong Zhang con = ptap->C_oth->rmap->n; 398e953e47cSHong Zhang proc = 0; 399e953e47cSHong Zhang for (i=0; i<con; i++) { 400e953e47cSHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 401e953e47cSHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 402e953e47cSHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 403e953e47cSHong Zhang } 404e953e47cSHong Zhang 405e953e47cSHong Zhang len = 0; /* max length of buf_si[], see (4) */ 406e953e47cSHong Zhang owners_co[0] = 0; 407e953e47cSHong Zhang nsend = 0; 408e953e47cSHong Zhang for (proc=0; proc<size; proc++) { 409e953e47cSHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 410e953e47cSHong Zhang if (len_s[proc]) { 411e953e47cSHong Zhang nsend++; 412e953e47cSHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 413e953e47cSHong Zhang len += len_si[proc]; 414e953e47cSHong Zhang } 415e953e47cSHong Zhang } 416e953e47cSHong Zhang 417e953e47cSHong Zhang /* determine the number and length of messages to receive for coi and coj */ 418e953e47cSHong Zhang ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 419e953e47cSHong Zhang ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 420e953e47cSHong Zhang 421e953e47cSHong Zhang /* post the Irecv and Isend of coj */ 422e953e47cSHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 423e953e47cSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 424e953e47cSHong Zhang ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 425e953e47cSHong Zhang for (proc=0, k=0; proc<size; proc++) { 426e953e47cSHong Zhang if (!len_s[proc]) continue; 427e953e47cSHong Zhang i = owners_co[proc]; 428e953e47cSHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 429e953e47cSHong Zhang k++; 430e953e47cSHong Zhang } 431e953e47cSHong Zhang 432e953e47cSHong Zhang /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 433e953e47cSHong Zhang /* ---------------------------------------- */ 434e953e47cSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 435e953e47cSHong Zhang c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 436e953e47cSHong Zhang 437e953e47cSHong Zhang /* receives coj are complete */ 438e953e47cSHong Zhang for (i=0; i<nrecv; i++) { 439e953e47cSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 440e953e47cSHong Zhang } 441e953e47cSHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 442e953e47cSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 443e953e47cSHong Zhang 444e953e47cSHong Zhang /* add received column indices into ta to update Crmax */ 445e953e47cSHong Zhang for (k=0; k<nrecv; k++) {/* k-th received message */ 446e953e47cSHong Zhang Jptr = buf_rj[k]; 447e953e47cSHong Zhang for (j=0; j<len_r[k]; j++) { 448e953e47cSHong Zhang ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 449e953e47cSHong Zhang } 450e953e47cSHong Zhang } 451e953e47cSHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 452e953e47cSHong Zhang ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 453e953e47cSHong Zhang 454e953e47cSHong Zhang /* (4) send and recv coi */ 455e953e47cSHong Zhang /*-----------------------*/ 456e953e47cSHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 457e953e47cSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 458e953e47cSHong Zhang ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 459e953e47cSHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 460e953e47cSHong Zhang for (proc=0,k=0; proc<size; proc++) { 461e953e47cSHong Zhang if (!len_s[proc]) continue; 462e953e47cSHong Zhang /* form outgoing message for i-structure: 463e953e47cSHong Zhang buf_si[0]: nrows to be sent 464e953e47cSHong Zhang [1:nrows]: row index (global) 465e953e47cSHong Zhang [nrows+1:2*nrows+1]: i-structure index 466e953e47cSHong Zhang */ 467e953e47cSHong Zhang /*-------------------------------------------*/ 468e953e47cSHong Zhang nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 469e953e47cSHong Zhang buf_si_i = buf_si + nrows+1; 470e953e47cSHong Zhang buf_si[0] = nrows; 471e953e47cSHong Zhang buf_si_i[0] = 0; 472e953e47cSHong Zhang nrows = 0; 473e953e47cSHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 474e953e47cSHong Zhang nzi = coi[i+1] - coi[i]; 475e953e47cSHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 476e953e47cSHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 477e953e47cSHong Zhang nrows++; 478e953e47cSHong Zhang } 479e953e47cSHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 480e953e47cSHong Zhang k++; 481e953e47cSHong Zhang buf_si += len_si[proc]; 482e953e47cSHong Zhang } 483e953e47cSHong Zhang for (i=0; i<nrecv; i++) { 484e953e47cSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 485e953e47cSHong Zhang } 486e953e47cSHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 487e953e47cSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 488e953e47cSHong Zhang 489e953e47cSHong Zhang ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 490e953e47cSHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 491e953e47cSHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 492e953e47cSHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 493b4e8d1b6SHong Zhang 494e953e47cSHong Zhang /* (5) compute the local portion of Cmpi */ 495e953e47cSHong Zhang /* ------------------------------------------ */ 496e953e47cSHong Zhang /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 497e953e47cSHong Zhang ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 498e953e47cSHong Zhang current_space = free_space; 499e953e47cSHong Zhang 500e953e47cSHong Zhang ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 501e953e47cSHong Zhang for (k=0; k<nrecv; k++) { 502e953e47cSHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 503e953e47cSHong Zhang nrows = *buf_ri_k[k]; 504e953e47cSHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 505e953e47cSHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 506e953e47cSHong Zhang } 507e953e47cSHong Zhang 508e953e47cSHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 509cf742fccSHong Zhang ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 510e953e47cSHong Zhang for (i=0; i<pn; i++) { 511e953e47cSHong Zhang /* add C_loc into Cmpi */ 512e953e47cSHong Zhang nzi = c_loc->i[i+1] - c_loc->i[i]; 513e953e47cSHong Zhang Jptr = c_loc->j + c_loc->i[i]; 514cf742fccSHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 515e953e47cSHong Zhang 516e953e47cSHong Zhang /* add received col data into lnk */ 517e953e47cSHong Zhang for (k=0; k<nrecv; k++) { /* k-th received message */ 518e953e47cSHong Zhang if (i == *nextrow[k]) { /* i-th row */ 519e953e47cSHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 520e953e47cSHong Zhang Jptr = buf_rj[k] + *nextci[k]; 521cf742fccSHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 522e953e47cSHong Zhang nextrow[k]++; nextci[k]++; 523e953e47cSHong Zhang } 524e953e47cSHong Zhang } 525e953e47cSHong Zhang nzi = lnk[0]; 526e953e47cSHong Zhang 527e953e47cSHong Zhang /* copy data into free space, then initialize lnk */ 528cf742fccSHong Zhang ierr = PetscLLCondensedClean_Scalable(nzi,current_space->array,lnk);CHKERRQ(ierr); 529e953e47cSHong Zhang ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 530e953e47cSHong Zhang } 531e953e47cSHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 532cf742fccSHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 533e953e47cSHong Zhang ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 534e953e47cSHong Zhang 535e953e47cSHong Zhang /* local sizes and preallocation */ 536e953e47cSHong Zhang ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 537e953e47cSHong Zhang ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 538e953e47cSHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 539e953e47cSHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 540e953e47cSHong Zhang 541e953e47cSHong Zhang /* members in merge */ 542e953e47cSHong Zhang ierr = PetscFree(id_r);CHKERRQ(ierr); 543e953e47cSHong Zhang ierr = PetscFree(len_r);CHKERRQ(ierr); 544e953e47cSHong Zhang ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 545e953e47cSHong Zhang ierr = PetscFree(buf_ri);CHKERRQ(ierr); 546e953e47cSHong Zhang ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 547e953e47cSHong Zhang ierr = PetscFree(buf_rj);CHKERRQ(ierr); 548e953e47cSHong Zhang ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 549e953e47cSHong Zhang 550e953e47cSHong Zhang /* attach the supporting struct to Cmpi for reuse */ 551e953e47cSHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 552e953e47cSHong Zhang c->ptap = ptap; 553e953e47cSHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 554e953e47cSHong Zhang ptap->destroy = Cmpi->ops->destroy; 555cf97cf3cSHong Zhang ptap->view = Cmpi->ops->view; 556e953e47cSHong Zhang 557e953e47cSHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 558e953e47cSHong Zhang Cmpi->assembled = PETSC_FALSE; 559e953e47cSHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 560e953e47cSHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 561cf97cf3cSHong Zhang Cmpi->ops->view = MatView_MPIAIJ_PtAP; 562e953e47cSHong Zhang *C = Cmpi; 563e953e47cSHong Zhang PetscFunctionReturn(0); 564e953e47cSHong Zhang } 565e953e47cSHong Zhang 566e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 567e953e47cSHong Zhang { 568e953e47cSHong Zhang PetscErrorCode ierr; 569e953e47cSHong Zhang Mat_PtAPMPI *ptap; 570e953e47cSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 571e953e47cSHong Zhang MPI_Comm comm; 572e953e47cSHong Zhang PetscMPIInt size,rank; 573e953e47cSHong Zhang Mat Cmpi; 574e953e47cSHong Zhang PetscFreeSpaceList free_space=NULL,current_space=NULL; 575e953e47cSHong Zhang PetscInt am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n; 576e953e47cSHong Zhang PetscInt *lnk,i,k,pnz,row,nsend; 577e953e47cSHong Zhang PetscBT lnkbt; 578e953e47cSHong Zhang PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 579e953e47cSHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 580e953e47cSHong Zhang PetscInt len,proc,*dnz,*onz,*owners,nzi,nspacedouble; 581e953e47cSHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 582e953e47cSHong Zhang MPI_Request *swaits,*rwaits; 583e953e47cSHong Zhang MPI_Status *sstatus,rstatus; 584e953e47cSHong Zhang PetscLayout rowmap; 585e953e47cSHong Zhang PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 586e953e47cSHong Zhang PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 587e953e47cSHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi; 588ec07b8f8SHong Zhang Mat_SeqAIJ *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth; 589e953e47cSHong Zhang PetscScalar *apv; 590e953e47cSHong Zhang PetscTable ta; 591e953e47cSHong Zhang #if defined(PETSC_HAVE_HYPRE) 592e953e47cSHong Zhang const char *algTypes[3] = {"scalable","nonscalable","hypre"}; 593e953e47cSHong Zhang PetscInt nalg = 3; 594e953e47cSHong Zhang #else 595e953e47cSHong Zhang const char *algTypes[2] = {"scalable","nonscalable"}; 596e953e47cSHong Zhang PetscInt nalg = 2; 597e953e47cSHong Zhang #endif 598b4e8d1b6SHong Zhang PetscInt alg = 1; /* set default algorithm */ 599e953e47cSHong Zhang #if defined(PETSC_USE_INFO) 600e953e47cSHong Zhang PetscReal apfill; 601e953e47cSHong Zhang #endif 602e953e47cSHong Zhang #if defined(PTAP_PROFILE) 603e953e47cSHong Zhang PetscLogDouble t0,t1,t11,t12,t2,t3,t4; 604e953e47cSHong Zhang #endif 605b4e8d1b6SHong Zhang PetscBool flg; 606e953e47cSHong Zhang 607e953e47cSHong Zhang PetscFunctionBegin; 608b4e8d1b6SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 609b4e8d1b6SHong Zhang 610e953e47cSHong Zhang /* pick an algorithm */ 611ae5f0867Sstefano_zampini ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 61268eacb73SHong Zhang PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */ 613b4e8d1b6SHong Zhang ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[1],&alg,&flg);CHKERRQ(ierr); 614ae5f0867Sstefano_zampini ierr = PetscOptionsEnd();CHKERRQ(ierr); 615ae5f0867Sstefano_zampini 616b4e8d1b6SHong Zhang if (!flg && pN > 100000) { /* may switch to scalable algorithm as default */ 617b4e8d1b6SHong Zhang MatInfo Ainfo,Pinfo; 618b4e8d1b6SHong Zhang PetscInt nz_local; 619b4e8d1b6SHong Zhang PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable; 620b4e8d1b6SHong Zhang 621b4e8d1b6SHong Zhang ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr); 622b4e8d1b6SHong Zhang ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr); 623b4e8d1b6SHong Zhang nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated); 624b4e8d1b6SHong Zhang 625b4e8d1b6SHong Zhang if (pN > fill*nz_local) alg_scalable_loc = PETSC_TRUE; 626b4e8d1b6SHong Zhang ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr); 627b4e8d1b6SHong Zhang 628b4e8d1b6SHong Zhang if (alg_scalable) { 629b4e8d1b6SHong Zhang alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */ 630b4e8d1b6SHong Zhang } 631b4e8d1b6SHong Zhang } 632b4e8d1b6SHong Zhang 633e953e47cSHong Zhang if (alg == 0) { 634e953e47cSHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(A,P,fill,C);CHKERRQ(ierr); 635cf742fccSHong Zhang (*C)->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable; 636e953e47cSHong Zhang PetscFunctionReturn(0); 637e953e47cSHong Zhang 638ae5f0867Sstefano_zampini #if defined(PETSC_HAVE_HYPRE) 639ae5f0867Sstefano_zampini } else if (alg == 2) { 640ae5f0867Sstefano_zampini /* Use boomerAMGBuildCoarseOperator */ 641ae5f0867Sstefano_zampini ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr); 642ae5f0867Sstefano_zampini PetscFunctionReturn(0); 643ae5f0867Sstefano_zampini #endif 644ae5f0867Sstefano_zampini } 645ae5f0867Sstefano_zampini 6468cb82516SHong Zhang #if defined(PTAP_PROFILE) 64780bb4639SHong Zhang ierr = PetscTime(&t0);CHKERRQ(ierr); 6488cb82516SHong Zhang #endif 6498cb82516SHong Zhang 650e953e47cSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 651e953e47cSHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 652e953e47cSHong Zhang 65352f7967eSHong Zhang if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data; 654ec07b8f8SHong Zhang 655e953e47cSHong Zhang /* create symbolic parallel matrix Cmpi */ 656e953e47cSHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 657e953e47cSHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 658e953e47cSHong Zhang 659e953e47cSHong Zhang /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */ 660e953e47cSHong Zhang Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ; 661e953e47cSHong Zhang 66215a3b8e2SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 66315a3b8e2SHong Zhang ierr = PetscNew(&ptap);CHKERRQ(ierr); 66415a3b8e2SHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 665cf97cf3cSHong Zhang ptap->algType = alg; 66615a3b8e2SHong Zhang 66715a3b8e2SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 66815a3b8e2SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 66915a3b8e2SHong Zhang /* get P_loc by taking all local rows of P */ 67015a3b8e2SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 67115a3b8e2SHong Zhang 67267a12041SHong Zhang /* (0) compute Rd = Pd^T, Ro = Po^T */ 67367a12041SHong Zhang /* --------------------------------- */ 674de817e96SHong Zhang ierr = MatTranspose_SeqAIJ(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 675de817e96SHong Zhang ierr = MatTranspose_SeqAIJ(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 6768cb82516SHong Zhang #if defined(PTAP_PROFILE) 677445158ffSHong Zhang ierr = PetscTime(&t11);CHKERRQ(ierr); 6788cb82516SHong Zhang #endif 67915a3b8e2SHong Zhang 68067a12041SHong Zhang /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 68167a12041SHong Zhang /* ----------------------------------------------------------------- */ 68267a12041SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 68352f7967eSHong Zhang if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 684445158ffSHong Zhang 6859a6dcab7SHong Zhang /* create and initialize a linked list */ 68645d00d1dSHong Zhang ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */ 6874b38b95cSHong Zhang MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta); 6884b38b95cSHong Zhang MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta); 68978d30b94SHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */ 690d0e9a2f7SHong Zhang /* printf("[%d] est %d, Crmax %d; pN %d\n",rank,5*(p_loc->rmax+p_oth->rmax + (PetscInt)(1.e-2*pN)),Crmax,pN); */ 69178d30b94SHong Zhang 69278d30b94SHong Zhang ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 693445158ffSHong Zhang 6948cb82516SHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 695ec07b8f8SHong Zhang if (ao) { 696f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr); 697ec07b8f8SHong Zhang } else { 698ec07b8f8SHong Zhang ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr); 699ec07b8f8SHong Zhang } 700445158ffSHong Zhang current_space = free_space; 70167a12041SHong Zhang nspacedouble = 0; 70267a12041SHong Zhang 70367a12041SHong Zhang ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 70467a12041SHong Zhang api[0] = 0; 705445158ffSHong Zhang for (i=0; i<am; i++) { 70667a12041SHong Zhang /* diagonal portion: Ad[i,:]*P */ 70767a12041SHong Zhang ai = ad->i; pi = p_loc->i; 70867a12041SHong Zhang nzi = ai[i+1] - ai[i]; 70967a12041SHong Zhang aj = ad->j + ai[i]; 710445158ffSHong Zhang for (j=0; j<nzi; j++) { 711445158ffSHong Zhang row = aj[j]; 71267a12041SHong Zhang pnz = pi[row+1] - pi[row]; 71367a12041SHong Zhang Jptr = p_loc->j + pi[row]; 714445158ffSHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 715445158ffSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 716445158ffSHong Zhang } 71767a12041SHong Zhang /* off-diagonal portion: Ao[i,:]*P */ 718ec07b8f8SHong Zhang if (ao) { 71967a12041SHong Zhang ai = ao->i; pi = p_oth->i; 72067a12041SHong Zhang nzi = ai[i+1] - ai[i]; 72167a12041SHong Zhang aj = ao->j + ai[i]; 722445158ffSHong Zhang for (j=0; j<nzi; j++) { 723445158ffSHong Zhang row = aj[j]; 72467a12041SHong Zhang pnz = pi[row+1] - pi[row]; 72567a12041SHong Zhang Jptr = p_oth->j + pi[row]; 726445158ffSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 727445158ffSHong Zhang } 728ec07b8f8SHong Zhang } 729445158ffSHong Zhang apnz = lnk[0]; 730445158ffSHong Zhang api[i+1] = api[i] + apnz; 731445158ffSHong Zhang if (ap_rmax < apnz) ap_rmax = apnz; 732445158ffSHong Zhang 733445158ffSHong Zhang /* if free space is not available, double the total space in the list */ 734445158ffSHong Zhang if (current_space->local_remaining<apnz) { 735f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 736445158ffSHong Zhang nspacedouble++; 737445158ffSHong Zhang } 738445158ffSHong Zhang 739445158ffSHong Zhang /* Copy data into free space, then initialize lnk */ 740445158ffSHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 741445158ffSHong Zhang 742445158ffSHong Zhang current_space->array += apnz; 743445158ffSHong Zhang current_space->local_used += apnz; 744445158ffSHong Zhang current_space->local_remaining -= apnz; 745445158ffSHong Zhang } 746681d504bSHong Zhang /* Allocate space for apj and apv, initialize apj, and */ 747445158ffSHong Zhang /* destroy list of free space and other temporary array(s) */ 748445158ffSHong Zhang ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 749445158ffSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 7509a6dcab7SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 7519a6dcab7SHong Zhang 752aa690a28SHong Zhang /* Create AP_loc for reuse */ 753445158ffSHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 754aa690a28SHong Zhang 755aa690a28SHong Zhang #if defined(PETSC_USE_INFO) 756ec07b8f8SHong Zhang if (ao) { 757aa690a28SHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 758ec07b8f8SHong Zhang } else { 759ec07b8f8SHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1); 760ec07b8f8SHong Zhang } 761aa690a28SHong Zhang ptap->AP_loc->info.mallocs = nspacedouble; 762aa690a28SHong Zhang ptap->AP_loc->info.fill_ratio_given = fill; 763aa690a28SHong Zhang ptap->AP_loc->info.fill_ratio_needed = apfill; 764aa690a28SHong Zhang 765aa690a28SHong Zhang if (api[am]) { 766b11c1ec8SHong Zhang ierr = PetscInfo3(ptap->AP_loc,"Nonscalable algorithm, AP_loc reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)apfill);CHKERRQ(ierr); 767aa690a28SHong Zhang ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 768aa690a28SHong Zhang } else { 769b11c1ec8SHong Zhang ierr = PetscInfo(ptap->AP_loc,"Nonscalable algorithm, AP_loc is empty \n");CHKERRQ(ierr); 770aa690a28SHong Zhang } 771aa690a28SHong Zhang #endif 772aa690a28SHong Zhang 7738cb82516SHong Zhang #if defined(PTAP_PROFILE) 774445158ffSHong Zhang ierr = PetscTime(&t12);CHKERRQ(ierr); 7758cb82516SHong Zhang #endif 77615a3b8e2SHong Zhang 777681d504bSHong Zhang /* (2-1) compute symbolic Co = Ro*AP_loc */ 778681d504bSHong Zhang /* ------------------------------------ */ 77967a12041SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 7808cb82516SHong Zhang #if defined(PTAP_PROFILE) 78180bb4639SHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 7828cb82516SHong Zhang #endif 78315a3b8e2SHong Zhang 78467a12041SHong Zhang /* (3) send coj of C_oth to other processors */ 78567a12041SHong Zhang /* ------------------------------------------ */ 78615a3b8e2SHong Zhang /* determine row ownership */ 787445158ffSHong Zhang ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 788445158ffSHong Zhang rowmap->n = pn; 789445158ffSHong Zhang rowmap->bs = 1; 790445158ffSHong Zhang ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 791445158ffSHong Zhang owners = rowmap->range; 79215a3b8e2SHong Zhang 79315a3b8e2SHong Zhang /* determine the number of messages to send, their lengths */ 7948cb82516SHong Zhang ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 7958cb82516SHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 79615a3b8e2SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 79715a3b8e2SHong Zhang 79867a12041SHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 79967a12041SHong Zhang coi = c_oth->i; coj = c_oth->j; 80067a12041SHong Zhang con = ptap->C_oth->rmap->n; 80115a3b8e2SHong Zhang proc = 0; 80267a12041SHong Zhang for (i=0; i<con; i++) { 80315a3b8e2SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 80415a3b8e2SHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 80515a3b8e2SHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 80615a3b8e2SHong Zhang } 80715a3b8e2SHong Zhang 80815a3b8e2SHong Zhang len = 0; /* max length of buf_si[], see (4) */ 80915a3b8e2SHong Zhang owners_co[0] = 0; 81067a12041SHong Zhang nsend = 0; 81115a3b8e2SHong Zhang for (proc=0; proc<size; proc++) { 81215a3b8e2SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 81315a3b8e2SHong Zhang if (len_s[proc]) { 814445158ffSHong Zhang nsend++; 81515a3b8e2SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 81615a3b8e2SHong Zhang len += len_si[proc]; 81715a3b8e2SHong Zhang } 81815a3b8e2SHong Zhang } 81915a3b8e2SHong Zhang 82015a3b8e2SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 821445158ffSHong Zhang ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 822445158ffSHong Zhang ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 82315a3b8e2SHong Zhang 82415a3b8e2SHong Zhang /* post the Irecv and Isend of coj */ 82515a3b8e2SHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 826445158ffSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 827445158ffSHong Zhang ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 82815a3b8e2SHong Zhang for (proc=0, k=0; proc<size; proc++) { 82915a3b8e2SHong Zhang if (!len_s[proc]) continue; 83015a3b8e2SHong Zhang i = owners_co[proc]; 83115a3b8e2SHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 83215a3b8e2SHong Zhang k++; 83315a3b8e2SHong Zhang } 83415a3b8e2SHong Zhang 835681d504bSHong Zhang /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 836681d504bSHong Zhang /* ---------------------------------------- */ 837681d504bSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 838681d504bSHong Zhang c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 839681d504bSHong Zhang 840681d504bSHong Zhang /* receives coj are complete */ 841445158ffSHong Zhang for (i=0; i<nrecv; i++) { 842445158ffSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 84315a3b8e2SHong Zhang } 84415a3b8e2SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 845445158ffSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 84615a3b8e2SHong Zhang 84778d30b94SHong Zhang /* add received column indices into ta to update Crmax */ 84878d30b94SHong Zhang for (k=0; k<nrecv; k++) {/* k-th received message */ 84978d30b94SHong Zhang Jptr = buf_rj[k]; 85078d30b94SHong Zhang for (j=0; j<len_r[k]; j++) { 85178d30b94SHong Zhang ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 85278d30b94SHong Zhang } 85378d30b94SHong Zhang } 85478d30b94SHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 85578d30b94SHong Zhang ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 8569a6dcab7SHong Zhang 85715a3b8e2SHong Zhang /* (4) send and recv coi */ 85815a3b8e2SHong Zhang /*-----------------------*/ 85915a3b8e2SHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 860445158ffSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 86115a3b8e2SHong Zhang ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 86215a3b8e2SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 86315a3b8e2SHong Zhang for (proc=0,k=0; proc<size; proc++) { 86415a3b8e2SHong Zhang if (!len_s[proc]) continue; 86515a3b8e2SHong Zhang /* form outgoing message for i-structure: 86615a3b8e2SHong Zhang buf_si[0]: nrows to be sent 86715a3b8e2SHong Zhang [1:nrows]: row index (global) 86815a3b8e2SHong Zhang [nrows+1:2*nrows+1]: i-structure index 86915a3b8e2SHong Zhang */ 87015a3b8e2SHong Zhang /*-------------------------------------------*/ 87115a3b8e2SHong Zhang nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 87215a3b8e2SHong Zhang buf_si_i = buf_si + nrows+1; 87315a3b8e2SHong Zhang buf_si[0] = nrows; 87415a3b8e2SHong Zhang buf_si_i[0] = 0; 87515a3b8e2SHong Zhang nrows = 0; 87615a3b8e2SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 87715a3b8e2SHong Zhang nzi = coi[i+1] - coi[i]; 87815a3b8e2SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 87915a3b8e2SHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 88015a3b8e2SHong Zhang nrows++; 88115a3b8e2SHong Zhang } 88215a3b8e2SHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 88315a3b8e2SHong Zhang k++; 88415a3b8e2SHong Zhang buf_si += len_si[proc]; 88515a3b8e2SHong Zhang } 886681d504bSHong Zhang for (i=0; i<nrecv; i++) { 887445158ffSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 88815a3b8e2SHong Zhang } 88915a3b8e2SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 890445158ffSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 89115a3b8e2SHong Zhang 8928cb82516SHong Zhang ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 89315a3b8e2SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 89415a3b8e2SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 89515a3b8e2SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 8968cb82516SHong Zhang #if defined(PTAP_PROFILE) 89780bb4639SHong Zhang ierr = PetscTime(&t2);CHKERRQ(ierr); 8988cb82516SHong Zhang #endif 89967a12041SHong Zhang /* (5) compute the local portion of Cmpi */ 90067a12041SHong Zhang /* ------------------------------------------ */ 90178d30b94SHong Zhang /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 90278d30b94SHong Zhang ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 90315a3b8e2SHong Zhang current_space = free_space; 90415a3b8e2SHong Zhang 905445158ffSHong Zhang ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 906445158ffSHong Zhang for (k=0; k<nrecv; k++) { 90715a3b8e2SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 90815a3b8e2SHong Zhang nrows = *buf_ri_k[k]; 90915a3b8e2SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 91015a3b8e2SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 91115a3b8e2SHong Zhang } 912445158ffSHong Zhang 91378d30b94SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 91478d30b94SHong Zhang ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 91515a3b8e2SHong Zhang for (i=0; i<pn; i++) { 91667a12041SHong Zhang /* add C_loc into Cmpi */ 91767a12041SHong Zhang nzi = c_loc->i[i+1] - c_loc->i[i]; 91867a12041SHong Zhang Jptr = c_loc->j + c_loc->i[i]; 91967a12041SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 92015a3b8e2SHong Zhang 92115a3b8e2SHong Zhang /* add received col data into lnk */ 922445158ffSHong Zhang for (k=0; k<nrecv; k++) { /* k-th received message */ 92315a3b8e2SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 92415a3b8e2SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 92515a3b8e2SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 92615a3b8e2SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 92715a3b8e2SHong Zhang nextrow[k]++; nextci[k]++; 92815a3b8e2SHong Zhang } 92915a3b8e2SHong Zhang } 930d0e9a2f7SHong Zhang nzi = lnk[0]; 9318cb82516SHong Zhang 93215a3b8e2SHong Zhang /* copy data into free space, then initialize lnk */ 933d0e9a2f7SHong Zhang ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 934d0e9a2f7SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 93515a3b8e2SHong Zhang } 93615a3b8e2SHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 93715a3b8e2SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 938445158ffSHong Zhang ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 9398cb82516SHong Zhang #if defined(PTAP_PROFILE) 94080bb4639SHong Zhang ierr = PetscTime(&t3);CHKERRQ(ierr); 9418cb82516SHong Zhang #endif 94280bb4639SHong Zhang 943ae5f0867Sstefano_zampini /* local sizes and preallocation */ 94415a3b8e2SHong Zhang ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 94515a3b8e2SHong Zhang ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 94615a3b8e2SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 94715a3b8e2SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 94815a3b8e2SHong Zhang 949445158ffSHong Zhang /* members in merge */ 950445158ffSHong Zhang ierr = PetscFree(id_r);CHKERRQ(ierr); 951445158ffSHong Zhang ierr = PetscFree(len_r);CHKERRQ(ierr); 952445158ffSHong Zhang ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 953445158ffSHong Zhang ierr = PetscFree(buf_ri);CHKERRQ(ierr); 954445158ffSHong Zhang ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 955445158ffSHong Zhang ierr = PetscFree(buf_rj);CHKERRQ(ierr); 956445158ffSHong Zhang ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 95715a3b8e2SHong Zhang 95815a3b8e2SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 95915a3b8e2SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 96015a3b8e2SHong Zhang c->ptap = ptap; 9611a47ec54SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 9621a47ec54SHong Zhang ptap->destroy = Cmpi->ops->destroy; 963cf97cf3cSHong Zhang ptap->view = Cmpi->ops->view; 96415a3b8e2SHong Zhang 96578d30b94SHong Zhang if (alg == 1) { 9668cb82516SHong Zhang ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 96778d30b94SHong Zhang } 9682259aa2eSHong Zhang 9691a47ec54SHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 9701a47ec54SHong Zhang Cmpi->assembled = PETSC_FALSE; 9711a47ec54SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 972aa690a28SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 973cf97cf3cSHong Zhang Cmpi->ops->view = MatView_MPIAIJ_PtAP; 9741a47ec54SHong Zhang *C = Cmpi; 9751a47ec54SHong Zhang 9768cb82516SHong Zhang #if defined(PTAP_PROFILE) 97780bb4639SHong Zhang ierr = PetscTime(&t4);CHKERRQ(ierr); 978e9c1f85fSHong Zhang if (rank == 1) { 979445158ffSHong 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); 980e9c1f85fSHong Zhang } 9818cb82516SHong Zhang #endif 9820d3441aeSHong Zhang PetscFunctionReturn(0); 9830d3441aeSHong Zhang } 9840d3441aeSHong Zhang 985aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 9860d3441aeSHong Zhang { 9870d3441aeSHong Zhang PetscErrorCode ierr; 9880d3441aeSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 9890d3441aeSHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 9902dd9e643SHong Zhang Mat_SeqAIJ *ap,*p_loc,*p_oth=NULL,*c_seq; 9910d3441aeSHong Zhang Mat_PtAPMPI *ptap = c->ptap; 9929ce11a7cSHong Zhang Mat AP_loc,C_loc,C_oth; 9930d3441aeSHong Zhang PetscInt i,rstart,rend,cm,ncols,row; 9948cb82516SHong Zhang PetscInt *api,*apj,am = A->rmap->n,j,col,apnz; 9958cb82516SHong Zhang PetscScalar *apa; 9960d3441aeSHong Zhang const PetscInt *cols; 9970d3441aeSHong Zhang const PetscScalar *vals; 9988cb82516SHong Zhang #if defined(PTAP_PROFILE) 9998cb82516SHong Zhang PetscMPIInt rank; 10008cb82516SHong Zhang MPI_Comm comm; 10010d3441aeSHong Zhang PetscLogDouble t0,t1,t2,t3,t4,eR,eAP,eCseq,eCmpi; 10028cb82516SHong Zhang #endif 10030d3441aeSHong Zhang 10040d3441aeSHong Zhang PetscFunctionBegin; 10050d3441aeSHong Zhang ierr = MatZeroEntries(C);CHKERRQ(ierr); 10060d3441aeSHong Zhang 1007e497d3c8SHong Zhang /* 1) get R = Pd^T,Ro = Po^T */ 10088cb82516SHong Zhang #if defined(PTAP_PROFILE) 10090d3441aeSHong Zhang ierr = PetscTime(&t0);CHKERRQ(ierr); 10108cb82516SHong Zhang #endif 1011748c7196SHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 10120d3441aeSHong Zhang ierr = MatTranspose_SeqAIJ(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 10130d3441aeSHong Zhang ierr = MatTranspose_SeqAIJ(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 1014748c7196SHong Zhang } 10158cb82516SHong Zhang #if defined(PTAP_PROFILE) 10160d3441aeSHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 10170d3441aeSHong Zhang eR = t1 - t0; 10188cb82516SHong Zhang #endif 10190d3441aeSHong Zhang 1020e497d3c8SHong Zhang /* 2) get AP_loc */ 10210d3441aeSHong Zhang AP_loc = ptap->AP_loc; 10228cb82516SHong Zhang ap = (Mat_SeqAIJ*)AP_loc->data; 10230d3441aeSHong Zhang 1024e497d3c8SHong Zhang /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 10250d3441aeSHong Zhang /*-----------------------------------------------------*/ 1026748c7196SHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 1027748c7196SHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 10280d3441aeSHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 10290d3441aeSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 1030e497d3c8SHong Zhang } 10310d3441aeSHong Zhang 10328cb82516SHong Zhang /* 2-2) compute numeric A_loc*P - dominating part */ 10338cb82516SHong Zhang /* ---------------------------------------------- */ 10340d3441aeSHong Zhang /* get data from symbolic products */ 10358cb82516SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 10362dd9e643SHong Zhang if (ptap->P_oth) { 10378cb82516SHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 10382dd9e643SHong Zhang } 10398cb82516SHong Zhang apa = ptap->apa; 1040681d504bSHong Zhang api = ap->i; 1041681d504bSHong Zhang apj = ap->j; 1042e497d3c8SHong Zhang for (i=0; i<am; i++) { 1043445158ffSHong Zhang /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 1044e497d3c8SHong Zhang AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa); 1045e497d3c8SHong Zhang apnz = api[i+1] - api[i]; 1046e497d3c8SHong Zhang for (j=0; j<apnz; j++) { 1047e497d3c8SHong Zhang col = apj[j+api[i]]; 1048e497d3c8SHong Zhang ap->a[j+ap->i[i]] = apa[col]; 1049e497d3c8SHong Zhang apa[col] = 0.0; 10500d3441aeSHong Zhang } 1051aa690a28SHong Zhang ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 10520d3441aeSHong Zhang } 10538cb82516SHong Zhang #if defined(PTAP_PROFILE) 10540d3441aeSHong Zhang ierr = PetscTime(&t2);CHKERRQ(ierr); 10550d3441aeSHong Zhang eAP = t2 - t1; 10568cb82516SHong Zhang #endif 10570d3441aeSHong Zhang 10588cb82516SHong Zhang /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 1059445158ffSHong Zhang ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 1060445158ffSHong Zhang ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 10619ce11a7cSHong Zhang C_loc = ptap->C_loc; 10629ce11a7cSHong Zhang C_oth = ptap->C_oth; 10638cb82516SHong Zhang #if defined(PTAP_PROFILE) 10640d3441aeSHong Zhang ierr = PetscTime(&t3);CHKERRQ(ierr); 10650d3441aeSHong Zhang eCseq = t3 - t2; 10668cb82516SHong Zhang #endif 10670d3441aeSHong Zhang 10680d3441aeSHong Zhang /* add C_loc and Co to to C */ 10690d3441aeSHong Zhang ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 10700d3441aeSHong Zhang 10710d3441aeSHong Zhang /* C_loc -> C */ 10720d3441aeSHong Zhang cm = C_loc->rmap->N; 10739ce11a7cSHong Zhang c_seq = (Mat_SeqAIJ*)C_loc->data; 10748cb82516SHong Zhang cols = c_seq->j; 10758cb82516SHong Zhang vals = c_seq->a; 10760d3441aeSHong Zhang for (i=0; i<cm; i++) { 10779ce11a7cSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 10780d3441aeSHong Zhang row = rstart + i; 10790d3441aeSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 10808cb82516SHong Zhang cols += ncols; vals += ncols; 10810d3441aeSHong Zhang } 10820d3441aeSHong Zhang 10830d3441aeSHong Zhang /* Co -> C, off-processor part */ 10849ce11a7cSHong Zhang cm = C_oth->rmap->N; 10859ce11a7cSHong Zhang c_seq = (Mat_SeqAIJ*)C_oth->data; 10868cb82516SHong Zhang cols = c_seq->j; 10878cb82516SHong Zhang vals = c_seq->a; 10889ce11a7cSHong Zhang for (i=0; i<cm; i++) { 10899ce11a7cSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 10900d3441aeSHong Zhang row = p->garray[i]; 10910d3441aeSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 10928cb82516SHong Zhang cols += ncols; vals += ncols; 10930d3441aeSHong Zhang } 10940d3441aeSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10950d3441aeSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10968cb82516SHong Zhang #if defined(PTAP_PROFILE) 10970d3441aeSHong Zhang ierr = PetscTime(&t4);CHKERRQ(ierr); 10980d3441aeSHong Zhang eCmpi = t4 - t3; 10990d3441aeSHong Zhang 11008cb82516SHong Zhang ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 11018cb82516SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 11020d3441aeSHong Zhang if (rank==1) { 11030d3441aeSHong 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); 11040d3441aeSHong Zhang } 11058cb82516SHong Zhang #endif 1106748c7196SHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 11070d3441aeSHong Zhang PetscFunctionReturn(0); 11080d3441aeSHong Zhang } 11090d3441aeSHong Zhang 11108403a639SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,PetscReal fill,Mat *C) 111142c57489SHong Zhang { 111242c57489SHong Zhang PetscErrorCode ierr; 111377584fe6SHong Zhang Mat Cmpi; 1114f8487c73SHong Zhang Mat_PtAPMPI *ptap; 11150298fd71SBarry Smith PetscFreeSpaceList free_space=NULL,current_space=NULL; 1116f8487c73SHong Zhang Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 111742c57489SHong Zhang Mat_SeqAIJ *ad =(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 111842c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 1119ded4f561SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pdti,*pdtj,*poti,*potj,*ptJ; 112077584fe6SHong Zhang PetscInt *adi=ad->i,*aj,*aoi=ao->i,nnz; 1121a914927fSHong Zhang PetscInt *lnk,*owners_co,*coi,*coj,i,k,pnz,row; 1122d16ca5f0SHong Zhang PetscInt am=A->rmap->n,pN=P->cmap->N,pm=P->rmap->n,pn=P->cmap->n; 112342c57489SHong Zhang PetscBT lnkbt; 1124ce94432eSBarry Smith MPI_Comm comm; 1125a914927fSHong Zhang PetscMPIInt size,rank,tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0; 112642c57489SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 112742c57489SHong Zhang PetscInt len,proc,*dnz,*onz,*owners; 112824ecddacSHong Zhang PetscInt nzi,*pti,*ptj; 112942c57489SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1130ded4f561SHong Zhang MPI_Request *swaits,*rwaits; 1131ded4f561SHong Zhang MPI_Status *sstatus,rstatus; 113242c57489SHong Zhang Mat_Merge_SeqsToMPI *merge; 113377584fe6SHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,pon,nspacedouble=0,j,ap_rmax=0; 1134d16ca5f0SHong Zhang PetscReal afill=1.0,afill_tmp; 113542c57489SHong Zhang 113642c57489SHong Zhang PetscFunctionBegin; 1137ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 113883445d95SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 113983445d95SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 114083445d95SHong Zhang 1141f8487c73SHong Zhang /* create struct Mat_PtAPMPI and attached it to C later */ 1142b00a9115SJed Brown ierr = PetscNew(&ptap);CHKERRQ(ierr); 1143b00a9115SJed Brown ierr = PetscNew(&merge);CHKERRQ(ierr); 11449f4155fbSHong Zhang ptap->merge = merge; 1145f8487c73SHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 11466c6e00e0SHong Zhang 11476c6e00e0SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 1148b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1149fe615159SHong Zhang 11506c6e00e0SHong Zhang /* get P_loc by taking all local rows of P */ 1151a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 11526c6e00e0SHong Zhang 1153a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1154a1a4e84aSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 11556c6e00e0SHong Zhang pi_loc = p_loc->i; pj_loc = p_loc->j; 11566c6e00e0SHong Zhang pi_oth = p_oth->i; pj_oth = p_oth->j; 115742c57489SHong Zhang 11582addb229SHong Zhang /* (1) compute symbolic AP = A_loc*P = A_diag*P_loc + A_off*P_oth (api,apj) */ 11592addb229SHong Zhang /*--------------------------------------------------------------------------*/ 1160854ce69bSBarry Smith ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 116182412ba7SHong Zhang api[0] = 0; 116283445d95SHong Zhang 116383445d95SHong Zhang /* create and initialize a linked list */ 1164a914927fSHong Zhang ierr = PetscLLCondensedCreate(pN,pN,&lnk,&lnkbt);CHKERRQ(ierr); 116583445d95SHong Zhang 1166a914927fSHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) -OOM for ex56, np=8k on Intrepid! */ 1167f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(adi[am],PetscIntSumTruncate(aoi[am],pi_loc[pm]))),&free_space);CHKERRQ(ierr); 1168f4a743e1SHong Zhang current_space = free_space; 1169f4a743e1SHong Zhang 1170f4a743e1SHong Zhang for (i=0; i<am; i++) { 1171f4a743e1SHong Zhang /* diagonal portion of A */ 1172ded4f561SHong Zhang nzi = adi[i+1] - adi[i]; 117377584fe6SHong Zhang aj = ad->j + adi[i]; 1174ded4f561SHong Zhang for (j=0; j<nzi; j++) { 117577584fe6SHong Zhang row = aj[j]; 117682412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 1177ded4f561SHong Zhang Jptr = pj_loc + pi_loc[row]; 117883445d95SHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 1179a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1180f4a743e1SHong Zhang } 1181f4a743e1SHong Zhang /* off-diagonal portion of A */ 1182ded4f561SHong Zhang nzi = aoi[i+1] - aoi[i]; 118377584fe6SHong Zhang aj = ao->j + aoi[i]; 1184ded4f561SHong Zhang for (j=0; j<nzi; j++) { 118577584fe6SHong Zhang row = aj[j]; 118682412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 1187ded4f561SHong Zhang Jptr = pj_oth + pi_oth[row]; 1188a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1189f4a743e1SHong Zhang } 1190a914927fSHong Zhang apnz = lnk[0]; 119182412ba7SHong Zhang api[i+1] = api[i] + apnz; 119277584fe6SHong Zhang if (ap_rmax < apnz) ap_rmax = apnz; 1193f4a743e1SHong Zhang 119483445d95SHong Zhang /* if free space is not available, double the total space in the list */ 119582412ba7SHong Zhang if (current_space->local_remaining<apnz) { 1196f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1197f2b054eeSHong Zhang nspacedouble++; 1198f4a743e1SHong Zhang } 1199f4a743e1SHong Zhang 1200f4a743e1SHong Zhang /* Copy data into free space, then initialize lnk */ 1201a914927fSHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 12022205254eSKarl Rupp 120382412ba7SHong Zhang current_space->array += apnz; 120482412ba7SHong Zhang current_space->local_used += apnz; 120582412ba7SHong Zhang current_space->local_remaining -= apnz; 1206f4a743e1SHong Zhang } 1207a914927fSHong Zhang 120882412ba7SHong Zhang /* Allocate space for apj, initialize apj, and */ 1209f4a743e1SHong Zhang /* destroy list of free space and other temporary array(s) */ 1210854ce69bSBarry Smith ierr = PetscMalloc1(api[am]+1,&apj);CHKERRQ(ierr); 121177584fe6SHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 1212118727c9SMark F. Adams afill_tmp = (PetscReal)api[am]/(adi[am]+aoi[am]+pi_loc[pm]+1); 1213d16ca5f0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 1214f4a743e1SHong Zhang 12152addb229SHong Zhang /* (2) determine symbolic Co=(p->B)^T*AP - send to others (coi,coj)*/ 12162addb229SHong Zhang /*-----------------------------------------------------------------*/ 1217de0260b3SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 121842c57489SHong Zhang 1219ded4f561SHong Zhang /* then, compute symbolic Co = (p->B)^T*AP */ 1220d0f46423SBarry Smith pon = (p->B)->cmap->n; /* total num of rows to be sent to other processors 122183445d95SHong Zhang >= (num of nonzero rows of C_seq) - pn */ 1222854ce69bSBarry Smith ierr = PetscMalloc1(pon+1,&coi);CHKERRQ(ierr); 1223de0260b3SHong Zhang coi[0] = 0; 122442c57489SHong Zhang 1225d16ca5f0SHong Zhang /* set initial free space to be fill*(nnz(p->B) + nnz(AP)) */ 1226f91af8c7SBarry Smith nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(poti[pon],api[am])); 122722d28d08SBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 122842c57489SHong Zhang current_space = free_space; 122942c57489SHong Zhang 1230de0260b3SHong Zhang for (i=0; i<pon; i++) { 123182412ba7SHong Zhang pnz = poti[i+1] - poti[i]; 123277584fe6SHong Zhang ptJ = potj + poti[i]; 123377584fe6SHong Zhang for (j=0; j<pnz; j++) { 123477584fe6SHong Zhang row = ptJ[j]; /* row of AP == col of Pot */ 123582412ba7SHong Zhang apnz = api[row+1] - api[row]; 1236ded4f561SHong Zhang Jptr = apj + api[row]; 123783445d95SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 1238a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 123942c57489SHong Zhang } 1240a914927fSHong Zhang nnz = lnk[0]; 124142c57489SHong Zhang 124283445d95SHong Zhang /* If free space is not available, double the total space in the list */ 1243ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 1244f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1245d16ca5f0SHong Zhang nspacedouble++; 124642c57489SHong Zhang } 124742c57489SHong Zhang 124842c57489SHong Zhang /* Copy data into free space, and zero out denserows */ 1249a914927fSHong Zhang ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 12502205254eSKarl Rupp 1251ded4f561SHong Zhang current_space->array += nnz; 1252ded4f561SHong Zhang current_space->local_used += nnz; 1253ded4f561SHong Zhang current_space->local_remaining -= nnz; 12542205254eSKarl Rupp 1255ded4f561SHong Zhang coi[i+1] = coi[i] + nnz; 125642c57489SHong Zhang } 12576b911d16SHong Zhang 12586b911d16SHong Zhang ierr = PetscMalloc1(coi[pon],&coj);CHKERRQ(ierr); 1259a1a86e44SBarry Smith ierr = PetscFreeSpaceContiguous(&free_space,coj);CHKERRQ(ierr); 1260118727c9SMark F. Adams afill_tmp = (PetscReal)coi[pon]/(poti[pon] + api[am]+1); 1261d16ca5f0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 1262de0260b3SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->B,&poti,&potj);CHKERRQ(ierr); 126342c57489SHong Zhang 12646b911d16SHong Zhang /* (3) send j-array (coj) of Co to other processors */ 12656b911d16SHong Zhang /*--------------------------------------------------*/ 12666b911d16SHong Zhang ierr = PetscCalloc1(size,&merge->len_s);CHKERRQ(ierr); 12676b911d16SHong Zhang len_s = merge->len_s; 12686b911d16SHong Zhang merge->nsend = 0; 12696b911d16SHong Zhang 12706b911d16SHong Zhang 127142c57489SHong Zhang /* determine row ownership */ 127226283091SBarry Smith ierr = PetscLayoutCreate(comm,&merge->rowmap);CHKERRQ(ierr); 12737a2fc3feSBarry Smith merge->rowmap->n = pn; 12747a2fc3feSBarry Smith merge->rowmap->bs = 1; 12752205254eSKarl Rupp 127626283091SBarry Smith ierr = PetscLayoutSetUp(merge->rowmap);CHKERRQ(ierr); 12777a2fc3feSBarry Smith owners = merge->rowmap->range; 127842c57489SHong Zhang 127942c57489SHong Zhang /* determine the number of messages to send, their lengths */ 1280dcca6d9dSJed Brown ierr = PetscMalloc2(size,&len_si,size,&sstatus);CHKERRQ(ierr); 128183445d95SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 1282854ce69bSBarry Smith ierr = PetscMalloc1(size+2,&owners_co);CHKERRQ(ierr); 1283de0260b3SHong Zhang 128483445d95SHong Zhang proc = 0; 1285de0260b3SHong Zhang for (i=0; i<pon; i++) { 1286de0260b3SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 12872addb229SHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 1288ee6e4b5bSHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 128983445d95SHong Zhang } 1290de0260b3SHong Zhang 12912addb229SHong Zhang len = 0; /* max length of buf_si[], see (4) */ 1292de0260b3SHong Zhang owners_co[0] = 0; 129342c57489SHong Zhang for (proc=0; proc<size; proc++) { 1294de0260b3SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 1295ee6e4b5bSHong Zhang if (len_s[proc]) { 129642c57489SHong Zhang merge->nsend++; 12972addb229SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 129842c57489SHong Zhang len += len_si[proc]; 129942c57489SHong Zhang } 130042c57489SHong Zhang } 130142c57489SHong Zhang 1302ded4f561SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 13030298fd71SBarry Smith ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&merge->nrecv);CHKERRQ(ierr); 130442c57489SHong Zhang ierr = PetscGatherMessageLengths2(comm,merge->nsend,merge->nrecv,len_s,len_si,&merge->id_r,&merge->len_r,&len_ri);CHKERRQ(ierr); 130542c57489SHong Zhang 1306ded4f561SHong Zhang /* post the Irecv and Isend of coj */ 1307529bc5dcSHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1308529bc5dcSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,merge->nrecv,merge->id_r,merge->len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1309854ce69bSBarry Smith ierr = PetscMalloc1(merge->nsend+1,&swaits);CHKERRQ(ierr); 131042c57489SHong Zhang for (proc=0, k=0; proc<size; proc++) { 131142c57489SHong Zhang if (!len_s[proc]) continue; 1312de0260b3SHong Zhang i = owners_co[proc]; 1313529bc5dcSHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 131442c57489SHong Zhang k++; 131542c57489SHong Zhang } 131642c57489SHong Zhang 1317ded4f561SHong Zhang /* receives and sends of coj are complete */ 131877584fe6SHong Zhang for (i=0; i<merge->nrecv; i++) { 1319c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1320ded4f561SHong Zhang } 1321ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 13220c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 132382412ba7SHong Zhang 13246b911d16SHong Zhang /* (4) send and recv coi */ 13256b911d16SHong Zhang /*-----------------------*/ 1326529bc5dcSHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1327529bc5dcSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,merge->nrecv,merge->id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 1328854ce69bSBarry Smith ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 132942c57489SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 133042c57489SHong Zhang for (proc=0,k=0; proc<size; proc++) { 133142c57489SHong Zhang if (!len_s[proc]) continue; 133242c57489SHong Zhang /* form outgoing message for i-structure: 133342c57489SHong Zhang buf_si[0]: nrows to be sent 133442c57489SHong Zhang [1:nrows]: row index (global) 133542c57489SHong Zhang [nrows+1:2*nrows+1]: i-structure index 133642c57489SHong Zhang */ 133742c57489SHong Zhang /*-------------------------------------------*/ 13382addb229SHong Zhang nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 133942c57489SHong Zhang buf_si_i = buf_si + nrows+1; 134042c57489SHong Zhang buf_si[0] = nrows; 134142c57489SHong Zhang buf_si_i[0] = 0; 134242c57489SHong Zhang nrows = 0; 1343de0260b3SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 1344ded4f561SHong Zhang nzi = coi[i+1] - coi[i]; 1345ded4f561SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 1346de0260b3SHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 134742c57489SHong Zhang nrows++; 134842c57489SHong Zhang } 1349529bc5dcSHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 135042c57489SHong Zhang k++; 135142c57489SHong Zhang buf_si += len_si[proc]; 135242c57489SHong Zhang } 1353ded4f561SHong Zhang i = merge->nrecv; 1354ded4f561SHong Zhang while (i--) { 1355c280ed6aSJed Brown ierr = MPI_Waitany(merge->nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 1356ded4f561SHong Zhang } 1357ded4f561SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 13580c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,swaits,sstatus);CHKERRQ(ierr);} 1359a914927fSHong Zhang 136024ecddacSHong Zhang ierr = PetscFree2(len_si,sstatus);CHKERRQ(ierr); 136142c57489SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 1362ded4f561SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 136342c57489SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 136442c57489SHong Zhang 13656b911d16SHong Zhang /* (5) compute the local portion of C (mpi mat) */ 13666b911d16SHong Zhang /*----------------------------------------------*/ 1367ded4f561SHong Zhang ierr = MatGetSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 1368ded4f561SHong Zhang 136924ecddacSHong Zhang /* allocate pti array and free space for accumulating nonzero column info */ 1370854ce69bSBarry Smith ierr = PetscMalloc1(pn+1,&pti);CHKERRQ(ierr); 137124ecddacSHong Zhang pti[0] = 0; 137242c57489SHong Zhang 1373d16ca5f0SHong Zhang /* set initial free space to be fill*(nnz(P) + nnz(AP)) */ 1374f91af8c7SBarry Smith nnz = PetscRealIntMultTruncate(fill,PetscIntSumTruncate(pi_loc[pm],api[am])); 137522d28d08SBarry Smith ierr = PetscFreeSpaceGet(nnz,&free_space);CHKERRQ(ierr); 137642c57489SHong Zhang current_space = free_space; 137742c57489SHong Zhang 1378dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 137942c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { 138042c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 138142c57489SHong Zhang nrows = *buf_ri_k[k]; 138242c57489SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 138342c57489SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 138442c57489SHong Zhang } 138542c57489SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 1386936ad1eaSHong Zhang 138742c57489SHong Zhang for (i=0; i<pn; i++) { 1388ded4f561SHong Zhang /* add pdt[i,:]*AP into lnk */ 1389ded4f561SHong Zhang pnz = pdti[i+1] - pdti[i]; 139077584fe6SHong Zhang ptJ = pdtj + pdti[i]; 139177584fe6SHong Zhang for (j=0; j<pnz; j++) { 139277584fe6SHong Zhang row = ptJ[j]; /* row of AP == col of Pt */ 1393ded4f561SHong Zhang apnz = api[row+1] - api[row]; 1394ded4f561SHong Zhang Jptr = apj + api[row]; 1395ded4f561SHong Zhang /* add non-zero cols of AP into the sorted linked list lnk */ 1396a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(apnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1397ded4f561SHong Zhang } 1398a914927fSHong Zhang 139942c57489SHong Zhang /* add received col data into lnk */ 140042c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 140142c57489SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 1402ded4f561SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 1403ded4f561SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 1404a914927fSHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 140542c57489SHong Zhang nextrow[k]++; nextci[k]++; 140642c57489SHong Zhang } 140742c57489SHong Zhang } 1408a914927fSHong Zhang nnz = lnk[0]; 140942c57489SHong Zhang 141042c57489SHong Zhang /* if free space is not available, make more free space */ 1411ded4f561SHong Zhang if (current_space->local_remaining<nnz) { 1412f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscIntSumTruncate(nnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1413d16ca5f0SHong Zhang nspacedouble++; 141442c57489SHong Zhang } 141542c57489SHong Zhang /* copy data into free space, then initialize lnk */ 1416a914927fSHong Zhang ierr = PetscLLCondensedClean(pN,nnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1417ded4f561SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nnz,current_space->array,dnz,onz);CHKERRQ(ierr); 14182205254eSKarl Rupp 1419ded4f561SHong Zhang current_space->array += nnz; 1420ded4f561SHong Zhang current_space->local_used += nnz; 1421ded4f561SHong Zhang current_space->local_remaining -= nnz; 14222205254eSKarl Rupp 142324ecddacSHong Zhang pti[i+1] = pti[i] + nnz; 142442c57489SHong Zhang } 1425ded4f561SHong Zhang ierr = MatRestoreSymbolicTranspose_SeqAIJ(p->A,&pdti,&pdtj);CHKERRQ(ierr); 14260572522cSBarry Smith ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 142742c57489SHong Zhang 1428854ce69bSBarry Smith ierr = PetscMalloc1(pti[pn]+1,&ptj);CHKERRQ(ierr); 142924ecddacSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,ptj);CHKERRQ(ierr); 143024ecddacSHong Zhang afill_tmp = (PetscReal)pti[pn]/(pi_loc[pm] + api[am]+1); 1431d16ca5f0SHong Zhang if (afill_tmp > afill) afill = afill_tmp; 143242c57489SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 143342c57489SHong Zhang 14346b911d16SHong Zhang /* (6) create symbolic parallel matrix Cmpi */ 14356b911d16SHong Zhang /*------------------------------------------*/ 143677584fe6SHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 143777584fe6SHong Zhang ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 143833d57670SJed Brown ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 143977584fe6SHong Zhang ierr = MatSetType(Cmpi,MATMPIAIJ);CHKERRQ(ierr); 144077584fe6SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 144142c57489SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 1442a2f3521dSMark F. Adams 1443b091e509SHong Zhang merge->bi = pti; /* Cseq->i */ 1444b091e509SHong Zhang merge->bj = ptj; /* Cseq->j */ 1445b091e509SHong Zhang merge->coi = coi; /* Co->i */ 1446b091e509SHong Zhang merge->coj = coj; /* Co->j */ 144742c57489SHong Zhang merge->buf_ri = buf_ri; 144842c57489SHong Zhang merge->buf_rj = buf_rj; 1449de0260b3SHong Zhang merge->owners_co = owners_co; 145077584fe6SHong Zhang merge->destroy = Cmpi->ops->destroy; 145142c57489SHong Zhang 145277584fe6SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 145377584fe6SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 1454f8487c73SHong Zhang c->ptap = ptap; 145577584fe6SHong Zhang ptap->api = api; 145677584fe6SHong Zhang ptap->apj = apj; 14578403a639SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 14588403a639SHong Zhang ptap->destroy = Cmpi->ops->destroy; 14598403a639SHong Zhang 14608403a639SHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 14618403a639SHong Zhang Cmpi->assembled = PETSC_FALSE; 14628403a639SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 14638403a639SHong Zhang Cmpi->ops->duplicate = MatDuplicate_MPIAIJ_MatPtAP; 14648403a639SHong Zhang Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap; 146577584fe6SHong Zhang *C = Cmpi; 1466414894bdSHong Zhang 1467414894bdSHong Zhang /* flag 'scalable' determines which implementations to be used: 1468414894bdSHong Zhang 0: do dense axpy in MatPtAPNumeric() - fast, but requires storage of a nonscalable dense array apa; 1469414894bdSHong Zhang 1: do sparse axpy in MatPtAPNumeric() - might slow, uses a sparse array apa */ 1470414894bdSHong Zhang /* set default scalable */ 1471da0a95b2SSatish Balay ptap->scalable = PETSC_FALSE; /* PETSC_TRUE; */ 14722205254eSKarl Rupp 1473c5929fdfSBarry Smith ierr = PetscOptionsGetBool(((PetscObject)Cmpi)->options,((PetscObject)Cmpi)->prefix,"-matptap_scalable",&ptap->scalable,NULL);CHKERRQ(ierr); 1474414894bdSHong Zhang if (!ptap->scalable) { /* Do dense axpy */ 14751795a4d1SJed Brown ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 1476414894bdSHong Zhang } else { 14771795a4d1SJed Brown ierr = PetscCalloc1(ap_rmax+1,&ptap->apa);CHKERRQ(ierr); 1478414894bdSHong Zhang } 1479414894bdSHong Zhang 1480f2b054eeSHong Zhang #if defined(PETSC_USE_INFO) 148124ecddacSHong Zhang if (pti[pn] != 0) { 148257622a8eSBarry Smith ierr = PetscInfo3(Cmpi,"Reallocs %D; Fill ratio: given %g needed %g.\n",nspacedouble,(double)fill,(double)afill);CHKERRQ(ierr); 148357622a8eSBarry Smith ierr = PetscInfo1(Cmpi,"Use MatPtAP(A,P,MatReuse,%g,&C) for best performance.\n",(double)afill);CHKERRQ(ierr); 1484f2b054eeSHong Zhang } else { 148577584fe6SHong Zhang ierr = PetscInfo(Cmpi,"Empty matrix product\n");CHKERRQ(ierr); 1486f2b054eeSHong Zhang } 1487f2b054eeSHong Zhang #endif 148842c57489SHong Zhang PetscFunctionReturn(0); 148942c57489SHong Zhang } 149042c57489SHong Zhang 14918403a639SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_ptap(Mat A,Mat P,Mat C) 149242c57489SHong Zhang { 149342c57489SHong Zhang PetscErrorCode ierr; 1494f8487c73SHong Zhang Mat_MPIAIJ *a =(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 149542c57489SHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 1496de0260b3SHong Zhang Mat_SeqAIJ *pd=(Mat_SeqAIJ*)(p->A)->data,*po=(Mat_SeqAIJ*)(p->B)->data; 149742c57489SHong Zhang Mat_SeqAIJ *p_loc,*p_oth; 1498f8487c73SHong Zhang Mat_PtAPMPI *ptap; 14999f4155fbSHong Zhang Mat_Merge_SeqsToMPI *merge; 15001d633602SHong Zhang PetscInt *adi=ad->i,*aoi=ao->i,*adj,*aoj,*apJ,nextp; 150182412ba7SHong Zhang PetscInt *pi_loc,*pj_loc,*pi_oth,*pj_oth,*pJ,*pj; 150282412ba7SHong Zhang PetscInt i,j,k,anz,pnz,apnz,nextap,row,*cj; 1503e900a757SHong Zhang MatScalar *ada,*aoa,*apa,*pa,*ca,*pa_loc,*pa_oth,valtmp; 1504d0f46423SBarry Smith PetscInt am =A->rmap->n,cm=C->rmap->n,pon=(p->B)->cmap->n; 1505ce94432eSBarry Smith MPI_Comm comm; 150642c57489SHong Zhang PetscMPIInt size,rank,taga,*len_s; 150782412ba7SHong Zhang PetscInt *owners,proc,nrows,**buf_ri_k,**nextrow,**nextci; 150842c57489SHong Zhang PetscInt **buf_ri,**buf_rj; 150950f5bf86SHong Zhang PetscInt cnz=0,*bj_i,*bi,*bj,bnz,nextcj; /* bi,bj,ba: local array of C(mpi mat) */ 151042c57489SHong Zhang MPI_Request *s_waits,*r_waits; 151142c57489SHong Zhang MPI_Status *status; 151282412ba7SHong Zhang MatScalar **abuf_r,*ba_i,*pA,*coa,*ba; 151382412ba7SHong Zhang PetscInt *api,*apj,*coi,*coj; 1514d0f46423SBarry Smith PetscInt *poJ=po->j,*pdJ=pd->j,pcstart=P->cmap->rstart,pcend=P->cmap->rend; 15156ce94e8fSHong Zhang PetscBool scalable; 151638571addSHong Zhang #if defined(PTAP_PROFILE) 1517e497d3c8SHong Zhang PetscLogDouble t0,t1,t2,eP,t3,t4,et2_AP=0.0,ePtAP=0.0,t2_0,t2_1,t2_2; 151838571addSHong Zhang #endif 151942c57489SHong Zhang 152042c57489SHong Zhang PetscFunctionBegin; 1521ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 152238571addSHong Zhang #if defined(PTAP_PROFILE) 15238563dfccSBarry Smith ierr = PetscTime(&t0);CHKERRQ(ierr); 152438571addSHong Zhang #endif 152542c57489SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 152642c57489SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 152742c57489SHong Zhang 1528f8487c73SHong Zhang ptap = c->ptap; 1529ce94432eSBarry 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"); 1530f8487c73SHong Zhang merge = ptap->merge; 1531414894bdSHong Zhang apa = ptap->apa; 15326ce94e8fSHong Zhang scalable = ptap->scalable; 15336c6e00e0SHong Zhang 1534a1a4e84aSHong Zhang /* 1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 15356b911d16SHong Zhang /*-----------------------------------------------------*/ 1536f8487c73SHong Zhang if (ptap->reuse == MAT_INITIAL_MATRIX) { 15379f4155fbSHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic(), skip calling MatGetBrowsOfAoCols() and MatMPIAIJGetLocalMat() */ 1538f8487c73SHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 15396c6e00e0SHong Zhang } else { /* update numerical values of P_oth and P_loc */ 1540b7f45c76SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 1541a1a4e84aSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 15426c6e00e0SHong Zhang } 154338571addSHong Zhang #if defined(PTAP_PROFILE) 15448563dfccSBarry Smith ierr = PetscTime(&t1);CHKERRQ(ierr); 154505d62848SHong Zhang eP = t1-t0; 154638571addSHong Zhang #endif 154705d62848SHong Zhang /* 154805d62848SHong Zhang printf("[%d] Ad: %d, %d; Ao: %d, %d; P_loc: %d, %d; P_oth %d, %d;\n",rank, 154905d62848SHong Zhang a->A->rmap->N,a->A->cmap->N,a->B->rmap->N,a->B->cmap->N, 155005d62848SHong Zhang ptap->P_loc->rmap->N,ptap->P_loc->cmap->N, 155105d62848SHong Zhang ptap->P_oth->rmap->N,ptap->P_oth->cmap->N); 155205d62848SHong Zhang */ 1553f8487c73SHong Zhang 1554f9473cf7SHong Zhang /* 2) compute numeric C_seq = P_loc^T*A_loc*P - dominating part */ 1555f9473cf7SHong Zhang /*--------------------------------------------------------------*/ 155642c57489SHong Zhang /* get data from symbolic products */ 1557a1a4e84aSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 1558a1a4e84aSHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 155989ae1891SBarry Smith pi_loc=p_loc->i; pj_loc=p_loc->j; pa_loc=p_loc->a; 156042c57489SHong Zhang pi_oth=p_oth->i; pj_oth=p_oth->j; pa_oth=p_oth->a; 156142c57489SHong Zhang 1562de0260b3SHong Zhang coi = merge->coi; coj = merge->coj; 15631795a4d1SJed Brown ierr = PetscCalloc1(coi[pon]+1,&coa);CHKERRQ(ierr); 1564de0260b3SHong Zhang 1565de0260b3SHong Zhang bi = merge->bi; bj = merge->bj; 15667a2fc3feSBarry Smith owners = merge->rowmap->range; 15671795a4d1SJed Brown ierr = PetscCalloc1(bi[cm]+1,&ba);CHKERRQ(ierr); /* ba: Cseq->a */ 156842c57489SHong Zhang 1569a1a4e84aSHong Zhang api = ptap->api; apj = ptap->apj; 1570d9824c15SHong Zhang 15719516a85cSHong Zhang if (!scalable) { /* Do dense axpy on apa (length of pN, stores A[i,:]*P) - nonscalable, but faster (could take 1/3 scalable time) */ 1572b38d1a2dSJed Brown ierr = PetscInfo(C,"Using non-scalable dense axpy\n");CHKERRQ(ierr); 15738cb82516SHong Zhang #if defined(PTAP_PROFILE) 157405d62848SHong Zhang ierr = PetscTime(&t1);CHKERRQ(ierr); 15758cb82516SHong Zhang #endif 1576a9d06231SHong Zhang for (i=0; i<am; i++) { 1577b091e509SHong Zhang #if defined(PTAP_PROFILE) 15788563dfccSBarry Smith ierr = PetscTime(&t2_0);CHKERRQ(ierr); 1579b091e509SHong Zhang #endif 1580a9d06231SHong Zhang /* 2-a) form i-th sparse row of A_loc*P = Ad*P_loc + Ao*P_oth */ 1581a9d06231SHong Zhang /*------------------------------------------------------------*/ 1582a9d06231SHong Zhang apJ = apj + api[i]; 1583a9d06231SHong Zhang 1584a9d06231SHong Zhang /* diagonal portion of A */ 1585a9d06231SHong Zhang anz = adi[i+1] - adi[i]; 1586a9d06231SHong Zhang adj = ad->j + adi[i]; 1587a9d06231SHong Zhang ada = ad->a + adi[i]; 1588a9d06231SHong Zhang for (j=0; j<anz; j++) { 1589a9d06231SHong Zhang row = adj[j]; 1590a9d06231SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 1591a9d06231SHong Zhang pj = pj_loc + pi_loc[row]; 1592a9d06231SHong Zhang pa = pa_loc + pi_loc[row]; 1593a9d06231SHong Zhang 1594a9d06231SHong Zhang /* perform dense axpy */ 1595e900a757SHong Zhang valtmp = ada[j]; 1596a9d06231SHong Zhang for (k=0; k<pnz; k++) { 1597e900a757SHong Zhang apa[pj[k]] += valtmp*pa[k]; 1598a9d06231SHong Zhang } 1599a9d06231SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1600a9d06231SHong Zhang } 1601a9d06231SHong Zhang 1602a9d06231SHong Zhang /* off-diagonal portion of A */ 1603a9d06231SHong Zhang anz = aoi[i+1] - aoi[i]; 1604a9d06231SHong Zhang aoj = ao->j + aoi[i]; 1605a9d06231SHong Zhang aoa = ao->a + aoi[i]; 1606a9d06231SHong Zhang for (j=0; j<anz; j++) { 1607a9d06231SHong Zhang row = aoj[j]; 1608a9d06231SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 1609a9d06231SHong Zhang pj = pj_oth + pi_oth[row]; 1610a9d06231SHong Zhang pa = pa_oth + pi_oth[row]; 1611a9d06231SHong Zhang 1612a9d06231SHong Zhang /* perform dense axpy */ 1613e900a757SHong Zhang valtmp = aoa[j]; 1614a9d06231SHong Zhang for (k=0; k<pnz; k++) { 1615e900a757SHong Zhang apa[pj[k]] += valtmp*pa[k]; 1616a9d06231SHong Zhang } 1617a9d06231SHong Zhang ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 1618a9d06231SHong Zhang } 1619b091e509SHong Zhang #if defined(PTAP_PROFILE) 16208563dfccSBarry Smith ierr = PetscTime(&t2_1);CHKERRQ(ierr); 1621b091e509SHong Zhang et2_AP += t2_1 - t2_0; 1622b091e509SHong Zhang #endif 1623a9d06231SHong Zhang 1624a9d06231SHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 1625a9d06231SHong Zhang /*--------------------------------------------------------------*/ 1626a9d06231SHong Zhang apnz = api[i+1] - api[i]; 1627a9d06231SHong Zhang /* put the value into Co=(p->B)^T*AP (off-diagonal part, send to others) */ 1628a9d06231SHong Zhang pnz = po->i[i+1] - po->i[i]; 1629e900a757SHong Zhang poJ = po->j + po->i[i]; 1630a9d06231SHong Zhang pA = po->a + po->i[i]; 1631a9d06231SHong Zhang for (j=0; j<pnz; j++) { 1632e900a757SHong Zhang row = poJ[j]; 1633e900a757SHong Zhang cnz = coi[row+1] - coi[row]; 1634e900a757SHong Zhang cj = coj + coi[row]; 1635e900a757SHong Zhang ca = coa + coi[row]; 1636a9d06231SHong Zhang /* perform dense axpy */ 1637e900a757SHong Zhang valtmp = pA[j]; 1638a9d06231SHong Zhang for (k=0; k<cnz; k++) { 1639e900a757SHong Zhang ca[k] += valtmp*apa[cj[k]]; 1640a9d06231SHong Zhang } 1641a9d06231SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1642a9d06231SHong Zhang } 1643a9d06231SHong Zhang /* put the value into Cd (diagonal part) */ 1644a9d06231SHong Zhang pnz = pd->i[i+1] - pd->i[i]; 1645e900a757SHong Zhang pdJ = pd->j + pd->i[i]; 1646a9d06231SHong Zhang pA = pd->a + pd->i[i]; 1647a9d06231SHong Zhang for (j=0; j<pnz; j++) { 1648e900a757SHong Zhang row = pdJ[j]; 1649e900a757SHong Zhang cnz = bi[row+1] - bi[row]; 1650e900a757SHong Zhang cj = bj + bi[row]; 1651e900a757SHong Zhang ca = ba + bi[row]; 1652a9d06231SHong Zhang /* perform dense axpy */ 1653e900a757SHong Zhang valtmp = pA[j]; 1654a9d06231SHong Zhang for (k=0; k<cnz; k++) { 1655e900a757SHong Zhang ca[k] += valtmp*apa[cj[k]]; 1656a9d06231SHong Zhang } 1657a9d06231SHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 1658a9d06231SHong Zhang } 16598403a639SHong Zhang 1660a9d06231SHong Zhang /* zero the current row of A*P */ 1661a9d06231SHong Zhang for (k=0; k<apnz; k++) apa[apJ[k]] = 0.0; 1662b091e509SHong Zhang #if defined(PTAP_PROFILE) 16638563dfccSBarry Smith ierr = PetscTime(&t2_2);CHKERRQ(ierr); 166405d62848SHong Zhang ePtAP += t2_2 - t2_1; 1665b091e509SHong Zhang #endif 1666a9d06231SHong Zhang } 166738571addSHong Zhang } else { /* Do sparse axpy on apa (length of ap_rmax, stores A[i,:]*P) - scalable, but slower */ 1668b38d1a2dSJed Brown ierr = PetscInfo(C,"Using scalable sparse axpy\n");CHKERRQ(ierr); 166938571addSHong Zhang /*-----------------------------------------------------------------------------------------*/ 1670a9d06231SHong Zhang pA=pa_loc; 167142c57489SHong Zhang for (i=0; i<am; i++) { 1672b091e509SHong Zhang #if defined(PTAP_PROFILE) 16738563dfccSBarry Smith ierr = PetscTime(&t2_0);CHKERRQ(ierr); 1674b091e509SHong Zhang #endif 1675f4a743e1SHong Zhang /* form i-th sparse row of A*P */ 167682412ba7SHong Zhang apnz = api[i+1] - api[i]; 167782412ba7SHong Zhang apJ = apj + api[i]; 167842c57489SHong Zhang /* diagonal portion of A */ 167982412ba7SHong Zhang anz = adi[i+1] - adi[i]; 16801d633602SHong Zhang adj = ad->j + adi[i]; 16811d633602SHong Zhang ada = ad->a + adi[i]; 168282412ba7SHong Zhang for (j=0; j<anz; j++) { 16831d633602SHong Zhang row = adj[j]; 168482412ba7SHong Zhang pnz = pi_loc[row+1] - pi_loc[row]; 168582412ba7SHong Zhang pj = pj_loc + pi_loc[row]; 168682412ba7SHong Zhang pa = pa_loc + pi_loc[row]; 1687e900a757SHong Zhang valtmp = ada[j]; 1688f4a743e1SHong Zhang nextp = 0; 168982412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 169082412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 1691e900a757SHong Zhang apa[k] += valtmp*pa[nextp++]; 169242c57489SHong Zhang } 169342c57489SHong Zhang } 1694dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 169542c57489SHong Zhang } 169642c57489SHong Zhang /* off-diagonal portion of A */ 169782412ba7SHong Zhang anz = aoi[i+1] - aoi[i]; 16981d633602SHong Zhang aoj = ao->j + aoi[i]; 16991d633602SHong Zhang aoa = ao->a + aoi[i]; 170082412ba7SHong Zhang for (j=0; j<anz; j++) { 17011d633602SHong Zhang row = aoj[j]; 170282412ba7SHong Zhang pnz = pi_oth[row+1] - pi_oth[row]; 170382412ba7SHong Zhang pj = pj_oth + pi_oth[row]; 170482412ba7SHong Zhang pa = pa_oth + pi_oth[row]; 1705e900a757SHong Zhang valtmp = aoa[j]; 1706f4a743e1SHong Zhang nextp = 0; 170782412ba7SHong Zhang for (k=0; nextp<pnz; k++) { 170882412ba7SHong Zhang if (apJ[k] == pj[nextp]) { /* col of AP == col of P */ 1709e900a757SHong Zhang apa[k] += valtmp*pa[nextp++]; 171042c57489SHong Zhang } 171142c57489SHong Zhang } 1712dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 171342c57489SHong Zhang } 1714b091e509SHong Zhang #if defined(PTAP_PROFILE) 17158563dfccSBarry Smith ierr = PetscTime(&t2_1);CHKERRQ(ierr); 1716b091e509SHong Zhang et2_AP += t2_1 - t2_0; 1717b091e509SHong Zhang #endif 171842c57489SHong Zhang 1719a9d06231SHong Zhang /* 2-b) Compute Cseq = P_loc[i,:]^T*AP[i,:] using outer product */ 1720a9d06231SHong Zhang /*--------------------------------------------------------------*/ 172182412ba7SHong Zhang pnz = pi_loc[i+1] - pi_loc[i]; 1722f9473cf7SHong Zhang pJ = pj_loc + pi_loc[i]; 172382412ba7SHong Zhang for (j=0; j<pnz; j++) { 172442c57489SHong Zhang nextap = 0; 1725f9473cf7SHong Zhang row = pJ[j]; /* global index */ 172682412ba7SHong Zhang if (row < pcstart || row >=pcend) { /* put the value into Co */ 1727e900a757SHong Zhang row = *poJ; 1728e900a757SHong Zhang cj = coj + coi[row]; 1729e900a757SHong Zhang ca = coa + coi[row]; poJ++; 173082412ba7SHong Zhang } else { /* put the value into Cd */ 1731e900a757SHong Zhang row = *pdJ; 1732e900a757SHong Zhang cj = bj + bi[row]; 1733e900a757SHong Zhang ca = ba + bi[row]; pdJ++; 173442c57489SHong Zhang } 1735e900a757SHong Zhang valtmp = pA[j]; 173682412ba7SHong Zhang for (k=0; nextap<apnz; k++) { 1737e900a757SHong Zhang if (cj[k]==apJ[nextap]) ca[k] += valtmp*apa[nextap++]; 173842c57489SHong Zhang } 1739dc0b31edSSatish Balay ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 1740de0260b3SHong Zhang } 17410496f32cSHong Zhang pA += pnz; 174242c57489SHong Zhang /* zero the current row info for A*P */ 174382412ba7SHong Zhang ierr = PetscMemzero(apa,apnz*sizeof(MatScalar));CHKERRQ(ierr); 1744b091e509SHong Zhang #if defined(PTAP_PROFILE) 17458563dfccSBarry Smith ierr = PetscTime(&t2_2);CHKERRQ(ierr); 174605d62848SHong Zhang ePtAP += t2_2 - t2_1; 1747b091e509SHong Zhang #endif 174842c57489SHong Zhang } 1749d9824c15SHong Zhang } 175038571addSHong Zhang #if defined(PTAP_PROFILE) 17518563dfccSBarry Smith ierr = PetscTime(&t2);CHKERRQ(ierr); 175238571addSHong Zhang #endif 175342c57489SHong Zhang 1754a9d06231SHong Zhang /* 3) send and recv matrix values coa */ 1755a9d06231SHong Zhang /*------------------------------------*/ 175642c57489SHong Zhang buf_ri = merge->buf_ri; 175742c57489SHong Zhang buf_rj = merge->buf_rj; 175842c57489SHong Zhang len_s = merge->len_s; 1759fc42d0c8SSatish Balay ierr = PetscCommGetNewTag(comm,&taga);CHKERRQ(ierr); 176042c57489SHong Zhang ierr = PetscPostIrecvScalar(comm,taga,merge->nrecv,merge->id_r,merge->len_r,&abuf_r,&r_waits);CHKERRQ(ierr); 176142c57489SHong Zhang 1762dcca6d9dSJed Brown ierr = PetscMalloc2(merge->nsend+1,&s_waits,size,&status);CHKERRQ(ierr); 176342c57489SHong Zhang for (proc=0,k=0; proc<size; proc++) { 176442c57489SHong Zhang if (!len_s[proc]) continue; 1765de0260b3SHong Zhang i = merge->owners_co[proc]; 1766de0260b3SHong Zhang ierr = MPI_Isend(coa+coi[i],len_s[proc],MPIU_MATSCALAR,proc,taga,comm,s_waits+k);CHKERRQ(ierr); 176742c57489SHong Zhang k++; 176842c57489SHong Zhang } 17690c468ba9SBarry Smith if (merge->nrecv) {ierr = MPI_Waitall(merge->nrecv,r_waits,status);CHKERRQ(ierr);} 17700c468ba9SBarry Smith if (merge->nsend) {ierr = MPI_Waitall(merge->nsend,s_waits,status);CHKERRQ(ierr);} 177142c57489SHong Zhang 1772a9d06231SHong Zhang ierr = PetscFree2(s_waits,status);CHKERRQ(ierr); 177342c57489SHong Zhang ierr = PetscFree(r_waits);CHKERRQ(ierr); 177482412ba7SHong Zhang ierr = PetscFree(coa);CHKERRQ(ierr); 177538571addSHong Zhang #if defined(PTAP_PROFILE) 17768563dfccSBarry Smith ierr = PetscTime(&t3);CHKERRQ(ierr); 177738571addSHong Zhang #endif 177842c57489SHong Zhang 1779a9d06231SHong Zhang /* 4) insert local Cseq and received values into Cmpi */ 1780a9d06231SHong Zhang /*------------------------------------------------------*/ 1781dcca6d9dSJed Brown ierr = PetscMalloc3(merge->nrecv,&buf_ri_k,merge->nrecv,&nextrow,merge->nrecv,&nextci);CHKERRQ(ierr); 178242c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { 178342c57489SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 178442c57489SHong Zhang nrows = *(buf_ri_k[k]); 178542c57489SHong Zhang nextrow[k] = buf_ri_k[k]+1; /* next row number of k-th recved i-structure */ 178682412ba7SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 178742c57489SHong Zhang } 178842c57489SHong Zhang 1789de0260b3SHong Zhang for (i=0; i<cm; i++) { 179082412ba7SHong Zhang row = owners[rank] + i; /* global row index of C_seq */ 179142c57489SHong Zhang bj_i = bj + bi[i]; /* col indices of the i-th row of C */ 1792de0260b3SHong Zhang ba_i = ba + bi[i]; 179382412ba7SHong Zhang bnz = bi[i+1] - bi[i]; 179442c57489SHong Zhang /* add received vals into ba */ 179542c57489SHong Zhang for (k=0; k<merge->nrecv; k++) { /* k-th received message */ 179642c57489SHong Zhang /* i-th row */ 179742c57489SHong Zhang if (i == *nextrow[k]) { 179882412ba7SHong Zhang cnz = *(nextci[k]+1) - *nextci[k]; 179982412ba7SHong Zhang cj = buf_rj[k] + *(nextci[k]); 180082412ba7SHong Zhang ca = abuf_r[k] + *(nextci[k]); 180182412ba7SHong Zhang nextcj = 0; 180282412ba7SHong Zhang for (j=0; nextcj<cnz; j++) { 180382412ba7SHong Zhang if (bj_i[j] == cj[nextcj]) { /* bcol == ccol */ 180482412ba7SHong Zhang ba_i[j] += ca[nextcj++]; 180542c57489SHong Zhang } 180642c57489SHong Zhang } 180782412ba7SHong Zhang nextrow[k]++; nextci[k]++; 1808c5af039cSHong Zhang ierr = PetscLogFlops(2.0*cnz);CHKERRQ(ierr); 180942c57489SHong Zhang } 181042c57489SHong Zhang } 181182412ba7SHong Zhang ierr = MatSetValues(C,1,&row,bnz,bj_i,ba_i,INSERT_VALUES);CHKERRQ(ierr); 181242c57489SHong Zhang } 181342c57489SHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 181442c57489SHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 181542c57489SHong Zhang 1816de0260b3SHong Zhang ierr = PetscFree(ba);CHKERRQ(ierr); 1817c05d87d6SBarry Smith ierr = PetscFree(abuf_r[0]);CHKERRQ(ierr); 181842c57489SHong Zhang ierr = PetscFree(abuf_r);CHKERRQ(ierr); 18190572522cSBarry Smith ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 182038571addSHong Zhang #if defined(PTAP_PROFILE) 18218563dfccSBarry Smith ierr = PetscTime(&t4);CHKERRQ(ierr); 18229516a85cSHong Zhang if (rank==1) { 182305d62848SHong 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); 18249516a85cSHong Zhang } 182538571addSHong Zhang #endif 182642c57489SHong Zhang PetscFunctionReturn(0); 182742c57489SHong Zhang } 1828