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; 19*3cdca5ebSHong Zhang Mat_APMPI *ptap=a->ap; 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 39dd4011a9SHong Zhang PetscErrorCode MatFreeIntermediateDataStructures_MPIAIJ_PtAP(Mat A) 40cc6cb787SHong Zhang { 41cc6cb787SHong Zhang PetscErrorCode ierr; 42dd4011a9SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 43*3cdca5ebSHong Zhang Mat_APMPI *ptap=a->ap; 443697aca6SHong Zhang Mat_Merge_SeqsToMPI *merge=ptap->merge; 45cc6cb787SHong Zhang 46cc6cb787SHong Zhang PetscFunctionBegin; 47dd4011a9SHong Zhang if (!ptap) PetscFunctionReturn(0); 48dd4011a9SHong Zhang 49b7f45c76SHong Zhang ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr); 50f8487c73SHong Zhang ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 51a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 52a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 53c5af039cSHong Zhang ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */ 5405d62848SHong Zhang ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr); 5505d62848SHong Zhang ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr); 568403a639SHong Zhang if (ptap->AP_loc) { /* used by alg_rap */ 57681d504bSHong Zhang Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data; 58681d504bSHong Zhang ierr = PetscFree(ap->i);CHKERRQ(ierr); 59681d504bSHong Zhang ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr); 600d3441aeSHong Zhang ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr); 618403a639SHong Zhang } else { /* used by alg_ptap */ 628403a639SHong Zhang ierr = PetscFree(ptap->api);CHKERRQ(ierr); 638403a639SHong Zhang ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 64681d504bSHong Zhang } 652259aa2eSHong Zhang ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr); 662259aa2eSHong Zhang ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr); 67414894bdSHong Zhang if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);} 68a560ef98SHong Zhang 698403a639SHong Zhang if (merge) { /* used by alg_ptap */ 70cc6cb787SHong Zhang ierr = PetscFree(merge->id_r);CHKERRQ(ierr); 71cc6cb787SHong Zhang ierr = PetscFree(merge->len_s);CHKERRQ(ierr); 72cc6cb787SHong Zhang ierr = PetscFree(merge->len_r);CHKERRQ(ierr); 73cc6cb787SHong Zhang ierr = PetscFree(merge->bi);CHKERRQ(ierr); 74cc6cb787SHong Zhang ierr = PetscFree(merge->bj);CHKERRQ(ierr); 75c05d87d6SBarry Smith ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr); 76cc6cb787SHong Zhang ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr); 77c05d87d6SBarry Smith ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr); 78cc6cb787SHong Zhang ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr); 79445158ffSHong Zhang ierr = PetscFree(merge->coi);CHKERRQ(ierr); 80445158ffSHong Zhang ierr = PetscFree(merge->coj);CHKERRQ(ierr); 8105b42c5fSBarry Smith ierr = PetscFree(merge->owners_co);CHKERRQ(ierr); 826bf464f9SBarry Smith ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr); 83f8487c73SHong Zhang ierr = PetscFree(ptap->merge);CHKERRQ(ierr); 84bf0cc555SLisandro Dalcin } 85dd4011a9SHong Zhang 86dd4011a9SHong Zhang A->ops->destroy = ptap->destroy; 87*3cdca5ebSHong Zhang ierr = PetscFree(a->ap);CHKERRQ(ierr); 883697aca6SHong Zhang PetscFunctionReturn(0); 893697aca6SHong Zhang } 90a560ef98SHong Zhang 913697aca6SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A) 923697aca6SHong Zhang { 933697aca6SHong Zhang PetscErrorCode ierr; 943697aca6SHong Zhang 953697aca6SHong Zhang PetscFunctionBegin; 96dd4011a9SHong Zhang ierr = (*A->ops->freeintermediatedatastructures)(A);CHKERRQ(ierr); 97dd4011a9SHong Zhang ierr = (*A->ops->destroy)(A);CHKERRQ(ierr); 98cc6cb787SHong Zhang PetscFunctionReturn(0); 99cc6cb787SHong Zhang } 100cc6cb787SHong Zhang 101150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 10242c57489SHong Zhang { 10342c57489SHong Zhang PetscErrorCode ierr; 104a4ffb656SHong Zhang PetscBool flg; 10567a12041SHong Zhang MPI_Comm comm; 106a4ffb656SHong Zhang #if !defined(PETSC_HAVE_HYPRE) 107a4ffb656SHong Zhang const char *algTypes[2] = {"scalable","nonscalable"}; 108a4ffb656SHong Zhang PetscInt nalg=2; 109a4ffb656SHong Zhang #else 110a4ffb656SHong Zhang const char *algTypes[3] = {"scalable","nonscalable","hypre"}; 111a4ffb656SHong Zhang PetscInt nalg=3; 112a4ffb656SHong Zhang #endif 113a4ffb656SHong Zhang PetscInt pN=P->cmap->N,alg=1; /* set default algorithm */ 11442c57489SHong Zhang 11542c57489SHong Zhang PetscFunctionBegin; 11667a12041SHong Zhang /* check if matrix local sizes are compatible */ 11767a12041SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1186c4ed002SBarry 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); 1196c4ed002SBarry 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); 12067a12041SHong Zhang 121cf3ca8ceSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 122a4ffb656SHong Zhang /* pick an algorithm */ 123a4ffb656SHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)A);CHKERRQ(ierr); 124a4ffb656SHong Zhang PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */ 125a4ffb656SHong Zhang ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 126a4ffb656SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 127a4ffb656SHong Zhang 128a4ffb656SHong Zhang if (!flg && pN > 100000) { /* may switch to scalable algorithm as default */ 129a4ffb656SHong Zhang MatInfo Ainfo,Pinfo; 130a4ffb656SHong Zhang PetscInt nz_local; 131a4ffb656SHong Zhang PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable; 132a4ffb656SHong Zhang 133a4ffb656SHong Zhang ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr); 134a4ffb656SHong Zhang ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr); 135a4ffb656SHong Zhang nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated); 136a4ffb656SHong Zhang 137a4ffb656SHong Zhang if (pN > fill*nz_local) alg_scalable_loc = PETSC_TRUE; 138a4ffb656SHong Zhang ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr); 139a4ffb656SHong Zhang 140a4ffb656SHong Zhang if (alg_scalable) { 141a4ffb656SHong Zhang alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */ 1420d3441aeSHong Zhang } 143a4ffb656SHong Zhang } 144a4ffb656SHong Zhang 145a4ffb656SHong Zhang switch (alg) { 146a4ffb656SHong Zhang case 1: 14780da8df7SHong Zhang /* do R=P^T locally, then C=R*A*P -- nonscalable */ 148a4ffb656SHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 149a4ffb656SHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr); 1503ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 151a4ffb656SHong Zhang break; 152a4ffb656SHong Zhang #if defined(PETSC_HAVE_HYPRE) 153a4ffb656SHong Zhang case 2: 154a4ffb656SHong Zhang /* Use boomerAMGBuildCoarseOperator */ 155a4ffb656SHong Zhang ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr); 156a4ffb656SHong Zhang PetscFunctionReturn(0); 157a4ffb656SHong Zhang break; 158a4ffb656SHong Zhang #endif 159a4ffb656SHong Zhang default: 160a4ffb656SHong Zhang /* do R=P^T locally, then C=R*A*P */ 161a4ffb656SHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 162a4ffb656SHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(A,P,fill,C);CHKERRQ(ierr); 163a4ffb656SHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 164a4ffb656SHong Zhang break; 165a4ffb656SHong Zhang } 1667a7894deSKris Buschelman } 1673ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 1688403a639SHong Zhang ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr); 1693ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 17042c57489SHong Zhang PetscFunctionReturn(0); 17142c57489SHong Zhang } 17242c57489SHong Zhang 173cf742fccSHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,Mat C) 174cf742fccSHong Zhang { 175cf742fccSHong Zhang PetscErrorCode ierr; 176cf742fccSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 177cf742fccSHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 178cf742fccSHong Zhang Mat_SeqAIJ *ap,*p_loc,*p_oth,*c_seq; 179*3cdca5ebSHong Zhang Mat_APMPI *ptap = c->ap; 180cf742fccSHong Zhang Mat AP_loc,C_loc,C_oth; 1815ca39647SHong Zhang PetscInt i,rstart,rend,cm,ncols,row,*api,*apj,am = A->rmap->n,apnz; 182cf742fccSHong Zhang PetscScalar *apa; 183cf742fccSHong Zhang const PetscInt *cols; 184cf742fccSHong Zhang const PetscScalar *vals; 18580da8df7SHong Zhang PetscBool freestruct; 186cf742fccSHong Zhang 187cf742fccSHong Zhang PetscFunctionBegin; 18880da8df7SHong Zhang if (!ptap) { 18980da8df7SHong Zhang MPI_Comm comm; 19080da8df7SHong Zhang ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 19180da8df7SHong Zhang SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 19280da8df7SHong Zhang } 19380da8df7SHong Zhang 194cf742fccSHong Zhang ierr = MatZeroEntries(C);CHKERRQ(ierr); 195cf742fccSHong Zhang 196cf742fccSHong Zhang /* 1) get R = Pd^T,Ro = Po^T */ 197cf742fccSHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 198419ecdd9Sandi selinger ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 199419ecdd9Sandi selinger ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 200cf742fccSHong Zhang } 201cf742fccSHong Zhang 202cf742fccSHong Zhang /* 2) get AP_loc */ 203cf742fccSHong Zhang AP_loc = ptap->AP_loc; 204cf742fccSHong Zhang ap = (Mat_SeqAIJ*)AP_loc->data; 205cf742fccSHong Zhang 206cf742fccSHong Zhang /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 207cf742fccSHong Zhang /*-----------------------------------------------------*/ 208cf742fccSHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 209cf742fccSHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 210cf742fccSHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 211cf742fccSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 212cf742fccSHong Zhang } 213cf742fccSHong Zhang 214cf742fccSHong Zhang /* 2-2) compute numeric A_loc*P - dominating part */ 215cf742fccSHong Zhang /* ---------------------------------------------- */ 216cf742fccSHong Zhang /* get data from symbolic products */ 217cf742fccSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 21852f7967eSHong Zhang if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 219c072d3e2SSatish Balay else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"ptap->P_oth is NULL. Cannot proceed!"); 22052f7967eSHong Zhang 221cf742fccSHong Zhang api = ap->i; 222cf742fccSHong Zhang apj = ap->j; 223cf742fccSHong Zhang for (i=0; i<am; i++) { 224cf742fccSHong Zhang /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 225cf742fccSHong Zhang apnz = api[i+1] - api[i]; 226b4e8d1b6SHong Zhang apa = ap->a + api[i]; 227b4e8d1b6SHong Zhang ierr = PetscMemzero(apa,sizeof(PetscScalar)*apnz);CHKERRQ(ierr); 228b4e8d1b6SHong Zhang AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa); 229cf742fccSHong Zhang ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 230cf742fccSHong Zhang } 231cf742fccSHong Zhang 232cf742fccSHong Zhang /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 233cf742fccSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 234cf742fccSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 235cf742fccSHong Zhang 236cf742fccSHong Zhang C_loc = ptap->C_loc; 237cf742fccSHong Zhang C_oth = ptap->C_oth; 238cf742fccSHong Zhang 239cf742fccSHong Zhang /* add C_loc and Co to to C */ 240cf742fccSHong Zhang ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 241cf742fccSHong Zhang 242cf742fccSHong Zhang /* C_loc -> C */ 243cf742fccSHong Zhang cm = C_loc->rmap->N; 244cf742fccSHong Zhang c_seq = (Mat_SeqAIJ*)C_loc->data; 245cf742fccSHong Zhang cols = c_seq->j; 246cf742fccSHong Zhang vals = c_seq->a; 247edeac6deSandi selinger 248e9ede7d0Sandi selinger /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */ 249edeac6deSandi selinger /* when there are no off-processor parts. */ 2501de21080Sandi selinger /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */ 2511de21080Sandi selinger /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */ 2521de21080Sandi selinger /* a table, and other, more complex stuff has to be done. */ 253edeac6deSandi selinger if (C->assembled) { 254edeac6deSandi selinger C->was_assembled = PETSC_TRUE; 255edeac6deSandi selinger C->assembled = PETSC_FALSE; 256edeac6deSandi selinger } 257edeac6deSandi selinger if (C->was_assembled) { 258cf742fccSHong Zhang for (i=0; i<cm; i++) { 259cf742fccSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 260cf742fccSHong Zhang row = rstart + i; 261edeac6deSandi selinger ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 262cf742fccSHong Zhang cols += ncols; vals += ncols; 263cf742fccSHong Zhang } 264edeac6deSandi selinger } else { 265e9ede7d0Sandi selinger ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr); 266edeac6deSandi selinger } 267cf742fccSHong Zhang 268cf742fccSHong Zhang /* Co -> C, off-processor part */ 269cf742fccSHong Zhang cm = C_oth->rmap->N; 270cf742fccSHong Zhang c_seq = (Mat_SeqAIJ*)C_oth->data; 271cf742fccSHong Zhang cols = c_seq->j; 272cf742fccSHong Zhang vals = c_seq->a; 273cf742fccSHong Zhang for (i=0; i<cm; i++) { 274cf742fccSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 275cf742fccSHong Zhang row = p->garray[i]; 276cf742fccSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 277cf742fccSHong Zhang cols += ncols; vals += ncols; 278cf742fccSHong Zhang } 279cf742fccSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 280cf742fccSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 281cf742fccSHong Zhang 282cf742fccSHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 28380da8df7SHong Zhang 28480da8df7SHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)C);CHKERRQ(ierr); 28580da8df7SHong Zhang PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */ 28680da8df7SHong Zhang /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */ 28780da8df7SHong Zhang freestruct = PETSC_FALSE; 28880da8df7SHong Zhang ierr = PetscOptionsBool("-mat_freeintermediatedatastructures","Free intermediate data structures", "MatFreeIntermediateDataStructures",freestruct, &freestruct, NULL);CHKERRQ(ierr); 28980da8df7SHong Zhang if (freestruct) { 29080da8df7SHong Zhang ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 29180da8df7SHong Zhang } 29280da8df7SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 293cf742fccSHong Zhang PetscFunctionReturn(0); 294cf742fccSHong Zhang } 295cf742fccSHong Zhang 296e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat *C) 2970d3441aeSHong Zhang { 2980d3441aeSHong Zhang PetscErrorCode ierr; 299*3cdca5ebSHong Zhang Mat_APMPI *ptap; 300681d504bSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 3012259aa2eSHong Zhang MPI_Comm comm; 3022259aa2eSHong Zhang PetscMPIInt size,rank; 30376db11f6SHong Zhang Mat Cmpi,P_loc,P_oth; 30415a3b8e2SHong Zhang PetscFreeSpaceList free_space=NULL,current_space=NULL; 305d0e9a2f7SHong Zhang PetscInt am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n; 306d0e9a2f7SHong Zhang PetscInt *lnk,i,k,pnz,row,nsend; 307f671be37SHong Zhang PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 30815a3b8e2SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 309d0e9a2f7SHong Zhang PetscInt len,proc,*dnz,*onz,*owners,nzi,nspacedouble; 31015a3b8e2SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 31115a3b8e2SHong Zhang MPI_Request *swaits,*rwaits; 31215a3b8e2SHong Zhang MPI_Status *sstatus,rstatus; 313445158ffSHong Zhang PetscLayout rowmap; 314445158ffSHong Zhang PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 315445158ffSHong Zhang PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 316b4e8d1b6SHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,Crmax,*aj,*ai,*pi; 31752f7967eSHong Zhang Mat_SeqAIJ *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth; 31867a12041SHong Zhang PetscScalar *apv; 31978d30b94SHong Zhang PetscTable ta; 320b92f168fSBarry Smith MatType mtype; 321aa690a28SHong Zhang #if defined(PETSC_USE_INFO) 3228cb82516SHong Zhang PetscReal apfill; 323aa690a28SHong Zhang #endif 32467a12041SHong Zhang 32567a12041SHong Zhang PetscFunctionBegin; 32667a12041SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 32767a12041SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 32867a12041SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 329ae5f0867Sstefano_zampini 33052f7967eSHong Zhang if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data; 33152f7967eSHong Zhang 332ae5f0867Sstefano_zampini /* create symbolic parallel matrix Cmpi */ 333ae5f0867Sstefano_zampini ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 334b92f168fSBarry Smith ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 335b92f168fSBarry Smith ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 336ae5f0867Sstefano_zampini 337*3cdca5ebSHong Zhang /* create struct Mat_APMPI and attached it to C later */ 338e953e47cSHong Zhang ierr = PetscNew(&ptap);CHKERRQ(ierr); 339e953e47cSHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 340cf97cf3cSHong Zhang ptap->algType = 0; 341e953e47cSHong Zhang 342e953e47cSHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 34376db11f6SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&P_oth);CHKERRQ(ierr); 344e953e47cSHong Zhang /* get P_loc by taking all local rows of P */ 34576db11f6SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&P_loc);CHKERRQ(ierr); 34676db11f6SHong Zhang 34776db11f6SHong Zhang ptap->P_loc = P_loc; 34876db11f6SHong Zhang ptap->P_oth = P_oth; 349e953e47cSHong Zhang 350e953e47cSHong Zhang /* (0) compute Rd = Pd^T, Ro = Po^T */ 351e953e47cSHong Zhang /* --------------------------------- */ 352419ecdd9Sandi selinger ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 353419ecdd9Sandi selinger ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 354e953e47cSHong Zhang 355e953e47cSHong Zhang /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 356e953e47cSHong Zhang /* ----------------------------------------------------------------- */ 35776db11f6SHong Zhang p_loc = (Mat_SeqAIJ*)P_loc->data; 35852f7967eSHong Zhang if (P_oth) p_oth = (Mat_SeqAIJ*)P_oth->data; 359e953e47cSHong Zhang 360e953e47cSHong Zhang /* create and initialize a linked list */ 361e953e47cSHong Zhang ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */ 36276db11f6SHong Zhang MatRowMergeMax_SeqAIJ(p_loc,P_loc->rmap->N,ta); 36376db11f6SHong Zhang MatRowMergeMax_SeqAIJ(p_oth,P_oth->rmap->N,ta); 364e953e47cSHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */ 365e953e47cSHong Zhang 36676db11f6SHong Zhang ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 367e953e47cSHong Zhang 368e953e47cSHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 36952f7967eSHong Zhang if (ao) { 370e953e47cSHong Zhang ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr); 37152f7967eSHong Zhang } else { 37252f7967eSHong Zhang ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr); 37352f7967eSHong Zhang } 374e953e47cSHong Zhang current_space = free_space; 375e953e47cSHong Zhang nspacedouble = 0; 376e953e47cSHong Zhang 377e953e47cSHong Zhang ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 378e953e47cSHong Zhang api[0] = 0; 379e953e47cSHong Zhang for (i=0; i<am; i++) { 380e953e47cSHong Zhang /* diagonal portion: Ad[i,:]*P */ 381e953e47cSHong Zhang ai = ad->i; pi = p_loc->i; 382e953e47cSHong Zhang nzi = ai[i+1] - ai[i]; 383e953e47cSHong Zhang aj = ad->j + ai[i]; 384e953e47cSHong Zhang for (j=0; j<nzi; j++) { 385e953e47cSHong Zhang row = aj[j]; 386e953e47cSHong Zhang pnz = pi[row+1] - pi[row]; 387e953e47cSHong Zhang Jptr = p_loc->j + pi[row]; 388e953e47cSHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 38976db11f6SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 390e953e47cSHong Zhang } 391e953e47cSHong Zhang /* off-diagonal portion: Ao[i,:]*P */ 39252f7967eSHong Zhang if (ao) { 393e953e47cSHong Zhang ai = ao->i; pi = p_oth->i; 394e953e47cSHong Zhang nzi = ai[i+1] - ai[i]; 395e953e47cSHong Zhang aj = ao->j + ai[i]; 396e953e47cSHong Zhang for (j=0; j<nzi; j++) { 397e953e47cSHong Zhang row = aj[j]; 398e953e47cSHong Zhang pnz = pi[row+1] - pi[row]; 399e953e47cSHong Zhang Jptr = p_oth->j + pi[row]; 40076db11f6SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 401e953e47cSHong Zhang } 40252f7967eSHong Zhang } 403e953e47cSHong Zhang apnz = lnk[0]; 404e953e47cSHong Zhang api[i+1] = api[i] + apnz; 405e953e47cSHong Zhang 406e953e47cSHong Zhang /* if free space is not available, double the total space in the list */ 407e953e47cSHong Zhang if (current_space->local_remaining<apnz) { 408e953e47cSHong Zhang ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 409e953e47cSHong Zhang nspacedouble++; 410e953e47cSHong Zhang } 411e953e47cSHong Zhang 412e953e47cSHong Zhang /* Copy data into free space, then initialize lnk */ 41376db11f6SHong Zhang ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr); 414e953e47cSHong Zhang 415e953e47cSHong Zhang current_space->array += apnz; 416e953e47cSHong Zhang current_space->local_used += apnz; 417e953e47cSHong Zhang current_space->local_remaining -= apnz; 418e953e47cSHong Zhang } 419e953e47cSHong Zhang /* Allocate space for apj and apv, initialize apj, and */ 420e953e47cSHong Zhang /* destroy list of free space and other temporary array(s) */ 421e953e47cSHong Zhang ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 422e953e47cSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 42376db11f6SHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 424e953e47cSHong Zhang 425e953e47cSHong Zhang /* Create AP_loc for reuse */ 426e953e47cSHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 427e953e47cSHong Zhang 428e953e47cSHong Zhang #if defined(PETSC_USE_INFO) 42952f7967eSHong Zhang if (ao) { 430e953e47cSHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 43152f7967eSHong Zhang } else { 43252f7967eSHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1); 43352f7967eSHong Zhang } 434e953e47cSHong Zhang ptap->AP_loc->info.mallocs = nspacedouble; 435e953e47cSHong Zhang ptap->AP_loc->info.fill_ratio_given = fill; 436e953e47cSHong Zhang ptap->AP_loc->info.fill_ratio_needed = apfill; 437e953e47cSHong Zhang 438e953e47cSHong Zhang if (api[am]) { 439b11c1ec8SHong 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); 440e953e47cSHong Zhang ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 441e953e47cSHong Zhang } else { 442b11c1ec8SHong Zhang ierr = PetscInfo(ptap->AP_loc,"Scalable algorithm, AP_loc is empty \n");CHKERRQ(ierr); 443e953e47cSHong Zhang } 444e953e47cSHong Zhang #endif 445e953e47cSHong Zhang 446e953e47cSHong Zhang /* (2-1) compute symbolic Co = Ro*AP_loc */ 447e953e47cSHong Zhang /* ------------------------------------ */ 448e953e47cSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 449e953e47cSHong Zhang 450e953e47cSHong Zhang /* (3) send coj of C_oth to other processors */ 451e953e47cSHong Zhang /* ------------------------------------------ */ 452e953e47cSHong Zhang /* determine row ownership */ 453e953e47cSHong Zhang ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 454e953e47cSHong Zhang rowmap->n = pn; 455e953e47cSHong Zhang rowmap->bs = 1; 456e953e47cSHong Zhang ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 457e953e47cSHong Zhang owners = rowmap->range; 458e953e47cSHong Zhang 459e953e47cSHong Zhang /* determine the number of messages to send, their lengths */ 460e953e47cSHong Zhang ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 461e953e47cSHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 462e953e47cSHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 463e953e47cSHong Zhang 464e953e47cSHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 465e953e47cSHong Zhang coi = c_oth->i; coj = c_oth->j; 466e953e47cSHong Zhang con = ptap->C_oth->rmap->n; 467e953e47cSHong Zhang proc = 0; 468e953e47cSHong Zhang for (i=0; i<con; i++) { 469e953e47cSHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 470e953e47cSHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 471e953e47cSHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 472e953e47cSHong Zhang } 473e953e47cSHong Zhang 474e953e47cSHong Zhang len = 0; /* max length of buf_si[], see (4) */ 475e953e47cSHong Zhang owners_co[0] = 0; 476e953e47cSHong Zhang nsend = 0; 477e953e47cSHong Zhang for (proc=0; proc<size; proc++) { 478e953e47cSHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 479e953e47cSHong Zhang if (len_s[proc]) { 480e953e47cSHong Zhang nsend++; 481e953e47cSHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 482e953e47cSHong Zhang len += len_si[proc]; 483e953e47cSHong Zhang } 484e953e47cSHong Zhang } 485e953e47cSHong Zhang 486e953e47cSHong Zhang /* determine the number and length of messages to receive for coi and coj */ 487e953e47cSHong Zhang ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 488e953e47cSHong Zhang ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 489e953e47cSHong Zhang 490e953e47cSHong Zhang /* post the Irecv and Isend of coj */ 491e953e47cSHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 492e953e47cSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 493e953e47cSHong Zhang ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 494e953e47cSHong Zhang for (proc=0, k=0; proc<size; proc++) { 495e953e47cSHong Zhang if (!len_s[proc]) continue; 496e953e47cSHong Zhang i = owners_co[proc]; 497e953e47cSHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 498e953e47cSHong Zhang k++; 499e953e47cSHong Zhang } 500e953e47cSHong Zhang 501e953e47cSHong Zhang /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 502e953e47cSHong Zhang /* ---------------------------------------- */ 503e953e47cSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 504e953e47cSHong Zhang c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 505e953e47cSHong Zhang 506e953e47cSHong Zhang /* receives coj are complete */ 507e953e47cSHong Zhang for (i=0; i<nrecv; i++) { 508e953e47cSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 509e953e47cSHong Zhang } 510e953e47cSHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 511e953e47cSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 512e953e47cSHong Zhang 513e953e47cSHong Zhang /* add received column indices into ta to update Crmax */ 514e953e47cSHong Zhang for (k=0; k<nrecv; k++) {/* k-th received message */ 515e953e47cSHong Zhang Jptr = buf_rj[k]; 516e953e47cSHong Zhang for (j=0; j<len_r[k]; j++) { 517e953e47cSHong Zhang ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 518e953e47cSHong Zhang } 519e953e47cSHong Zhang } 520e953e47cSHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 521e953e47cSHong Zhang ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 522e953e47cSHong Zhang 523e953e47cSHong Zhang /* (4) send and recv coi */ 524e953e47cSHong Zhang /*-----------------------*/ 525e953e47cSHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 526e953e47cSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 527e953e47cSHong Zhang ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 528e953e47cSHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 529e953e47cSHong Zhang for (proc=0,k=0; proc<size; proc++) { 530e953e47cSHong Zhang if (!len_s[proc]) continue; 531e953e47cSHong Zhang /* form outgoing message for i-structure: 532e953e47cSHong Zhang buf_si[0]: nrows to be sent 533e953e47cSHong Zhang [1:nrows]: row index (global) 534e953e47cSHong Zhang [nrows+1:2*nrows+1]: i-structure index 535e953e47cSHong Zhang */ 536e953e47cSHong Zhang /*-------------------------------------------*/ 537e953e47cSHong Zhang nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 538e953e47cSHong Zhang buf_si_i = buf_si + nrows+1; 539e953e47cSHong Zhang buf_si[0] = nrows; 540e953e47cSHong Zhang buf_si_i[0] = 0; 541e953e47cSHong Zhang nrows = 0; 542e953e47cSHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 543e953e47cSHong Zhang nzi = coi[i+1] - coi[i]; 544e953e47cSHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 545e953e47cSHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 546e953e47cSHong Zhang nrows++; 547e953e47cSHong Zhang } 548e953e47cSHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 549e953e47cSHong Zhang k++; 550e953e47cSHong Zhang buf_si += len_si[proc]; 551e953e47cSHong Zhang } 552e953e47cSHong Zhang for (i=0; i<nrecv; i++) { 553e953e47cSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 554e953e47cSHong Zhang } 555e953e47cSHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 556e953e47cSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 557e953e47cSHong Zhang 558e953e47cSHong Zhang ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 559e953e47cSHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 560e953e47cSHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 561e953e47cSHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 562b4e8d1b6SHong Zhang 563e953e47cSHong Zhang /* (5) compute the local portion of Cmpi */ 564e953e47cSHong Zhang /* ------------------------------------------ */ 565e953e47cSHong Zhang /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 566e953e47cSHong Zhang ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 567e953e47cSHong Zhang current_space = free_space; 568e953e47cSHong Zhang 569e953e47cSHong Zhang ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 570e953e47cSHong Zhang for (k=0; k<nrecv; k++) { 571e953e47cSHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 572e953e47cSHong Zhang nrows = *buf_ri_k[k]; 573e953e47cSHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 574e953e47cSHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 575e953e47cSHong Zhang } 576e953e47cSHong Zhang 577e953e47cSHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 578cf742fccSHong Zhang ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 579e953e47cSHong Zhang for (i=0; i<pn; i++) { 580e953e47cSHong Zhang /* add C_loc into Cmpi */ 581e953e47cSHong Zhang nzi = c_loc->i[i+1] - c_loc->i[i]; 582e953e47cSHong Zhang Jptr = c_loc->j + c_loc->i[i]; 583cf742fccSHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 584e953e47cSHong Zhang 585e953e47cSHong Zhang /* add received col data into lnk */ 586e953e47cSHong Zhang for (k=0; k<nrecv; k++) { /* k-th received message */ 587e953e47cSHong Zhang if (i == *nextrow[k]) { /* i-th row */ 588e953e47cSHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 589e953e47cSHong Zhang Jptr = buf_rj[k] + *nextci[k]; 590cf742fccSHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 591e953e47cSHong Zhang nextrow[k]++; nextci[k]++; 592e953e47cSHong Zhang } 593e953e47cSHong Zhang } 594e953e47cSHong Zhang nzi = lnk[0]; 595e953e47cSHong Zhang 596e953e47cSHong Zhang /* copy data into free space, then initialize lnk */ 597cf742fccSHong Zhang ierr = PetscLLCondensedClean_Scalable(nzi,current_space->array,lnk);CHKERRQ(ierr); 598e953e47cSHong Zhang ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 599e953e47cSHong Zhang } 600e953e47cSHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 601cf742fccSHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 602e953e47cSHong Zhang ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 603e953e47cSHong Zhang 604e953e47cSHong Zhang /* local sizes and preallocation */ 605e953e47cSHong Zhang ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 606e953e47cSHong Zhang ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 607e953e47cSHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 608e953e47cSHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 609e953e47cSHong Zhang 610e953e47cSHong Zhang /* members in merge */ 611e953e47cSHong Zhang ierr = PetscFree(id_r);CHKERRQ(ierr); 612e953e47cSHong Zhang ierr = PetscFree(len_r);CHKERRQ(ierr); 613e953e47cSHong Zhang ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 614e953e47cSHong Zhang ierr = PetscFree(buf_ri);CHKERRQ(ierr); 615e953e47cSHong Zhang ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 616e953e47cSHong Zhang ierr = PetscFree(buf_rj);CHKERRQ(ierr); 617e953e47cSHong Zhang ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 618e953e47cSHong Zhang 619e953e47cSHong Zhang /* attach the supporting struct to Cmpi for reuse */ 620e953e47cSHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 621*3cdca5ebSHong Zhang c->ap = ptap; 622e953e47cSHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 623e953e47cSHong Zhang ptap->destroy = Cmpi->ops->destroy; 624cf97cf3cSHong Zhang ptap->view = Cmpi->ops->view; 625e953e47cSHong Zhang 626e953e47cSHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 627e953e47cSHong Zhang Cmpi->assembled = PETSC_FALSE; 628a4ffb656SHong Zhang Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable; 629e953e47cSHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 630cf97cf3cSHong Zhang Cmpi->ops->view = MatView_MPIAIJ_PtAP; 63180da8df7SHong Zhang Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_PtAP; 632e953e47cSHong Zhang *C = Cmpi; 633e953e47cSHong Zhang PetscFunctionReturn(0); 634e953e47cSHong Zhang } 635e953e47cSHong Zhang 636e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 637e953e47cSHong Zhang { 638e953e47cSHong Zhang PetscErrorCode ierr; 639*3cdca5ebSHong Zhang Mat_APMPI *ptap; 640e953e47cSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 641e953e47cSHong Zhang MPI_Comm comm; 642e953e47cSHong Zhang PetscMPIInt size,rank; 643e953e47cSHong Zhang Mat Cmpi; 644e953e47cSHong Zhang PetscFreeSpaceList free_space=NULL,current_space=NULL; 645e953e47cSHong Zhang PetscInt am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n; 646e953e47cSHong Zhang PetscInt *lnk,i,k,pnz,row,nsend; 647e953e47cSHong Zhang PetscBT lnkbt; 648e953e47cSHong Zhang PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 649e953e47cSHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 650e953e47cSHong Zhang PetscInt len,proc,*dnz,*onz,*owners,nzi,nspacedouble; 651e953e47cSHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 652e953e47cSHong Zhang MPI_Request *swaits,*rwaits; 653e953e47cSHong Zhang MPI_Status *sstatus,rstatus; 654e953e47cSHong Zhang PetscLayout rowmap; 655e953e47cSHong Zhang PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 656e953e47cSHong Zhang PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 657e953e47cSHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi; 658ec07b8f8SHong Zhang Mat_SeqAIJ *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth; 659e953e47cSHong Zhang PetscScalar *apv; 660e953e47cSHong Zhang PetscTable ta; 661b92f168fSBarry Smith MatType mtype; 662e953e47cSHong Zhang #if defined(PETSC_USE_INFO) 663e953e47cSHong Zhang PetscReal apfill; 664e953e47cSHong Zhang #endif 665e953e47cSHong Zhang 666e953e47cSHong Zhang PetscFunctionBegin; 667b4e8d1b6SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 668e953e47cSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 669e953e47cSHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 670e953e47cSHong Zhang 67152f7967eSHong Zhang if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data; 672ec07b8f8SHong Zhang 673e953e47cSHong Zhang /* create symbolic parallel matrix Cmpi */ 674e953e47cSHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 675b92f168fSBarry Smith ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 676b92f168fSBarry Smith ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 6773697aca6SHong Zhang 678e953e47cSHong Zhang /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */ 679e953e47cSHong Zhang Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ; 680e953e47cSHong Zhang 681*3cdca5ebSHong Zhang /* create struct Mat_APMPI and attached it to C later */ 68215a3b8e2SHong Zhang ierr = PetscNew(&ptap);CHKERRQ(ierr); 68315a3b8e2SHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 684a4ffb656SHong Zhang ptap->algType = 1; 68515a3b8e2SHong Zhang 68615a3b8e2SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 68715a3b8e2SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 68815a3b8e2SHong Zhang /* get P_loc by taking all local rows of P */ 68915a3b8e2SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 69015a3b8e2SHong Zhang 69167a12041SHong Zhang /* (0) compute Rd = Pd^T, Ro = Po^T */ 69267a12041SHong Zhang /* --------------------------------- */ 693419ecdd9Sandi selinger ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 694419ecdd9Sandi selinger ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 69515a3b8e2SHong Zhang 69667a12041SHong Zhang /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 69767a12041SHong Zhang /* ----------------------------------------------------------------- */ 69867a12041SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 69952f7967eSHong Zhang if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 700445158ffSHong Zhang 7019a6dcab7SHong Zhang /* create and initialize a linked list */ 70245d00d1dSHong Zhang ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */ 7034b38b95cSHong Zhang MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta); 7044b38b95cSHong Zhang MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta); 70578d30b94SHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */ 706d0e9a2f7SHong 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); */ 70778d30b94SHong Zhang 70878d30b94SHong Zhang ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 709445158ffSHong Zhang 7108cb82516SHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 711ec07b8f8SHong Zhang if (ao) { 712f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr); 713ec07b8f8SHong Zhang } else { 714ec07b8f8SHong Zhang ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr); 715ec07b8f8SHong Zhang } 716445158ffSHong Zhang current_space = free_space; 71767a12041SHong Zhang nspacedouble = 0; 71867a12041SHong Zhang 71967a12041SHong Zhang ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 72067a12041SHong Zhang api[0] = 0; 721445158ffSHong Zhang for (i=0; i<am; i++) { 72267a12041SHong Zhang /* diagonal portion: Ad[i,:]*P */ 72367a12041SHong Zhang ai = ad->i; pi = p_loc->i; 72467a12041SHong Zhang nzi = ai[i+1] - ai[i]; 72567a12041SHong Zhang aj = ad->j + ai[i]; 726445158ffSHong Zhang for (j=0; j<nzi; j++) { 727445158ffSHong Zhang row = aj[j]; 72867a12041SHong Zhang pnz = pi[row+1] - pi[row]; 72967a12041SHong Zhang Jptr = p_loc->j + pi[row]; 730445158ffSHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 731445158ffSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 732445158ffSHong Zhang } 73367a12041SHong Zhang /* off-diagonal portion: Ao[i,:]*P */ 734ec07b8f8SHong Zhang if (ao) { 73567a12041SHong Zhang ai = ao->i; pi = p_oth->i; 73667a12041SHong Zhang nzi = ai[i+1] - ai[i]; 73767a12041SHong Zhang aj = ao->j + ai[i]; 738445158ffSHong Zhang for (j=0; j<nzi; j++) { 739445158ffSHong Zhang row = aj[j]; 74067a12041SHong Zhang pnz = pi[row+1] - pi[row]; 74167a12041SHong Zhang Jptr = p_oth->j + pi[row]; 742445158ffSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 743445158ffSHong Zhang } 744ec07b8f8SHong Zhang } 745445158ffSHong Zhang apnz = lnk[0]; 746445158ffSHong Zhang api[i+1] = api[i] + apnz; 747445158ffSHong Zhang if (ap_rmax < apnz) ap_rmax = apnz; 748445158ffSHong Zhang 749445158ffSHong Zhang /* if free space is not available, double the total space in the list */ 750445158ffSHong Zhang if (current_space->local_remaining<apnz) { 751f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 752445158ffSHong Zhang nspacedouble++; 753445158ffSHong Zhang } 754445158ffSHong Zhang 755445158ffSHong Zhang /* Copy data into free space, then initialize lnk */ 756445158ffSHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 757445158ffSHong Zhang 758445158ffSHong Zhang current_space->array += apnz; 759445158ffSHong Zhang current_space->local_used += apnz; 760445158ffSHong Zhang current_space->local_remaining -= apnz; 761445158ffSHong Zhang } 762681d504bSHong Zhang /* Allocate space for apj and apv, initialize apj, and */ 763445158ffSHong Zhang /* destroy list of free space and other temporary array(s) */ 764445158ffSHong Zhang ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 765445158ffSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 7669a6dcab7SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 7679a6dcab7SHong Zhang 768aa690a28SHong Zhang /* Create AP_loc for reuse */ 769445158ffSHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 770aa690a28SHong Zhang 771aa690a28SHong Zhang #if defined(PETSC_USE_INFO) 772ec07b8f8SHong Zhang if (ao) { 773aa690a28SHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 774ec07b8f8SHong Zhang } else { 775ec07b8f8SHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1); 776ec07b8f8SHong Zhang } 777aa690a28SHong Zhang ptap->AP_loc->info.mallocs = nspacedouble; 778aa690a28SHong Zhang ptap->AP_loc->info.fill_ratio_given = fill; 779aa690a28SHong Zhang ptap->AP_loc->info.fill_ratio_needed = apfill; 780aa690a28SHong Zhang 781aa690a28SHong Zhang if (api[am]) { 782b11c1ec8SHong 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); 783aa690a28SHong Zhang ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 784aa690a28SHong Zhang } else { 785b11c1ec8SHong Zhang ierr = PetscInfo(ptap->AP_loc,"Nonscalable algorithm, AP_loc is empty \n");CHKERRQ(ierr); 786aa690a28SHong Zhang } 787aa690a28SHong Zhang #endif 788aa690a28SHong Zhang 789681d504bSHong Zhang /* (2-1) compute symbolic Co = Ro*AP_loc */ 790681d504bSHong Zhang /* ------------------------------------ */ 79167a12041SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 79215a3b8e2SHong Zhang 79367a12041SHong Zhang /* (3) send coj of C_oth to other processors */ 79467a12041SHong Zhang /* ------------------------------------------ */ 79515a3b8e2SHong Zhang /* determine row ownership */ 796445158ffSHong Zhang ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 797445158ffSHong Zhang rowmap->n = pn; 798445158ffSHong Zhang rowmap->bs = 1; 799445158ffSHong Zhang ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 800445158ffSHong Zhang owners = rowmap->range; 80115a3b8e2SHong Zhang 80215a3b8e2SHong Zhang /* determine the number of messages to send, their lengths */ 8038cb82516SHong Zhang ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 8048cb82516SHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 80515a3b8e2SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 80615a3b8e2SHong Zhang 80767a12041SHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 80867a12041SHong Zhang coi = c_oth->i; coj = c_oth->j; 80967a12041SHong Zhang con = ptap->C_oth->rmap->n; 81015a3b8e2SHong Zhang proc = 0; 81167a12041SHong Zhang for (i=0; i<con; i++) { 81215a3b8e2SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 81315a3b8e2SHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 81415a3b8e2SHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 81515a3b8e2SHong Zhang } 81615a3b8e2SHong Zhang 81715a3b8e2SHong Zhang len = 0; /* max length of buf_si[], see (4) */ 81815a3b8e2SHong Zhang owners_co[0] = 0; 81967a12041SHong Zhang nsend = 0; 82015a3b8e2SHong Zhang for (proc=0; proc<size; proc++) { 82115a3b8e2SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 82215a3b8e2SHong Zhang if (len_s[proc]) { 823445158ffSHong Zhang nsend++; 82415a3b8e2SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 82515a3b8e2SHong Zhang len += len_si[proc]; 82615a3b8e2SHong Zhang } 82715a3b8e2SHong Zhang } 82815a3b8e2SHong Zhang 82915a3b8e2SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 830445158ffSHong Zhang ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 831445158ffSHong Zhang ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 83215a3b8e2SHong Zhang 83315a3b8e2SHong Zhang /* post the Irecv and Isend of coj */ 83415a3b8e2SHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 835445158ffSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 836445158ffSHong Zhang ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 83715a3b8e2SHong Zhang for (proc=0, k=0; proc<size; proc++) { 83815a3b8e2SHong Zhang if (!len_s[proc]) continue; 83915a3b8e2SHong Zhang i = owners_co[proc]; 84015a3b8e2SHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 84115a3b8e2SHong Zhang k++; 84215a3b8e2SHong Zhang } 84315a3b8e2SHong Zhang 844681d504bSHong Zhang /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 845681d504bSHong Zhang /* ---------------------------------------- */ 846681d504bSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 847681d504bSHong Zhang c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 848681d504bSHong Zhang 849681d504bSHong Zhang /* receives coj are complete */ 850445158ffSHong Zhang for (i=0; i<nrecv; i++) { 851445158ffSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 85215a3b8e2SHong Zhang } 85315a3b8e2SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 854445158ffSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 85515a3b8e2SHong Zhang 85678d30b94SHong Zhang /* add received column indices into ta to update Crmax */ 85778d30b94SHong Zhang for (k=0; k<nrecv; k++) {/* k-th received message */ 85878d30b94SHong Zhang Jptr = buf_rj[k]; 85978d30b94SHong Zhang for (j=0; j<len_r[k]; j++) { 86078d30b94SHong Zhang ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 86178d30b94SHong Zhang } 86278d30b94SHong Zhang } 86378d30b94SHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 86478d30b94SHong Zhang ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 8659a6dcab7SHong Zhang 86615a3b8e2SHong Zhang /* (4) send and recv coi */ 86715a3b8e2SHong Zhang /*-----------------------*/ 86815a3b8e2SHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 869445158ffSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 87015a3b8e2SHong Zhang ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 87115a3b8e2SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 87215a3b8e2SHong Zhang for (proc=0,k=0; proc<size; proc++) { 87315a3b8e2SHong Zhang if (!len_s[proc]) continue; 87415a3b8e2SHong Zhang /* form outgoing message for i-structure: 87515a3b8e2SHong Zhang buf_si[0]: nrows to be sent 87615a3b8e2SHong Zhang [1:nrows]: row index (global) 87715a3b8e2SHong Zhang [nrows+1:2*nrows+1]: i-structure index 87815a3b8e2SHong Zhang */ 87915a3b8e2SHong Zhang /*-------------------------------------------*/ 88015a3b8e2SHong Zhang nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 88115a3b8e2SHong Zhang buf_si_i = buf_si + nrows+1; 88215a3b8e2SHong Zhang buf_si[0] = nrows; 88315a3b8e2SHong Zhang buf_si_i[0] = 0; 88415a3b8e2SHong Zhang nrows = 0; 88515a3b8e2SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 88615a3b8e2SHong Zhang nzi = coi[i+1] - coi[i]; 88715a3b8e2SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 88815a3b8e2SHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 88915a3b8e2SHong Zhang nrows++; 89015a3b8e2SHong Zhang } 89115a3b8e2SHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 89215a3b8e2SHong Zhang k++; 89315a3b8e2SHong Zhang buf_si += len_si[proc]; 89415a3b8e2SHong Zhang } 895681d504bSHong Zhang for (i=0; i<nrecv; i++) { 896445158ffSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 89715a3b8e2SHong Zhang } 89815a3b8e2SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 899445158ffSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 90015a3b8e2SHong Zhang 9018cb82516SHong Zhang ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 90215a3b8e2SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 90315a3b8e2SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 90415a3b8e2SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 905a4ffb656SHong Zhang 90667a12041SHong Zhang /* (5) compute the local portion of Cmpi */ 90767a12041SHong Zhang /* ------------------------------------------ */ 90878d30b94SHong Zhang /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 90978d30b94SHong Zhang ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 91015a3b8e2SHong Zhang current_space = free_space; 91115a3b8e2SHong Zhang 912445158ffSHong Zhang ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 913445158ffSHong Zhang for (k=0; k<nrecv; k++) { 91415a3b8e2SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 91515a3b8e2SHong Zhang nrows = *buf_ri_k[k]; 91615a3b8e2SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 91715a3b8e2SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 91815a3b8e2SHong Zhang } 919445158ffSHong Zhang 92078d30b94SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 92178d30b94SHong Zhang ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 92215a3b8e2SHong Zhang for (i=0; i<pn; i++) { 92367a12041SHong Zhang /* add C_loc into Cmpi */ 92467a12041SHong Zhang nzi = c_loc->i[i+1] - c_loc->i[i]; 92567a12041SHong Zhang Jptr = c_loc->j + c_loc->i[i]; 92667a12041SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 92715a3b8e2SHong Zhang 92815a3b8e2SHong Zhang /* add received col data into lnk */ 929445158ffSHong Zhang for (k=0; k<nrecv; k++) { /* k-th received message */ 93015a3b8e2SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 93115a3b8e2SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 93215a3b8e2SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 93315a3b8e2SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 93415a3b8e2SHong Zhang nextrow[k]++; nextci[k]++; 93515a3b8e2SHong Zhang } 93615a3b8e2SHong Zhang } 937d0e9a2f7SHong Zhang nzi = lnk[0]; 9388cb82516SHong Zhang 93915a3b8e2SHong Zhang /* copy data into free space, then initialize lnk */ 940d0e9a2f7SHong Zhang ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 941d0e9a2f7SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 94215a3b8e2SHong Zhang } 94315a3b8e2SHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 94415a3b8e2SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 945445158ffSHong Zhang ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 94680bb4639SHong Zhang 947ae5f0867Sstefano_zampini /* local sizes and preallocation */ 94815a3b8e2SHong Zhang ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 94915a3b8e2SHong Zhang ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 95015a3b8e2SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 95115a3b8e2SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 95215a3b8e2SHong Zhang 953445158ffSHong Zhang /* members in merge */ 954445158ffSHong Zhang ierr = PetscFree(id_r);CHKERRQ(ierr); 955445158ffSHong Zhang ierr = PetscFree(len_r);CHKERRQ(ierr); 956445158ffSHong Zhang ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 957445158ffSHong Zhang ierr = PetscFree(buf_ri);CHKERRQ(ierr); 958445158ffSHong Zhang ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 959445158ffSHong Zhang ierr = PetscFree(buf_rj);CHKERRQ(ierr); 960445158ffSHong Zhang ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 96115a3b8e2SHong Zhang 96215a3b8e2SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 96315a3b8e2SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 964*3cdca5ebSHong Zhang c->ap = ptap; 9651a47ec54SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 9661a47ec54SHong Zhang ptap->destroy = Cmpi->ops->destroy; 967cf97cf3cSHong Zhang ptap->view = Cmpi->ops->view; 9688cb82516SHong Zhang ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 9692259aa2eSHong Zhang 9701a47ec54SHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 9711a47ec54SHong Zhang Cmpi->assembled = PETSC_FALSE; 9721a47ec54SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 973cf97cf3cSHong Zhang Cmpi->ops->view = MatView_MPIAIJ_PtAP; 974dd4011a9SHong Zhang Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_PtAP; 9751a47ec54SHong Zhang *C = Cmpi; 9760d3441aeSHong Zhang PetscFunctionReturn(0); 9770d3441aeSHong Zhang } 9780d3441aeSHong Zhang 979aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 9800d3441aeSHong Zhang { 9810d3441aeSHong Zhang PetscErrorCode ierr; 9820d3441aeSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 9830d3441aeSHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 9842dd9e643SHong Zhang Mat_SeqAIJ *ap,*p_loc,*p_oth=NULL,*c_seq; 985*3cdca5ebSHong Zhang Mat_APMPI *ptap = c->ap; 9869ce11a7cSHong Zhang Mat AP_loc,C_loc,C_oth; 9870d3441aeSHong Zhang PetscInt i,rstart,rend,cm,ncols,row; 9888cb82516SHong Zhang PetscInt *api,*apj,am = A->rmap->n,j,col,apnz; 9898cb82516SHong Zhang PetscScalar *apa; 9900d3441aeSHong Zhang const PetscInt *cols; 9910d3441aeSHong Zhang const PetscScalar *vals; 992dd4011a9SHong Zhang PetscBool freestruct; 9930d3441aeSHong Zhang 9940d3441aeSHong Zhang PetscFunctionBegin; 995dd4011a9SHong Zhang if (!ptap) { 996a9899c97SHong Zhang MPI_Comm comm; 997a9899c97SHong Zhang ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 998dd4011a9SHong Zhang SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 999a9899c97SHong Zhang } 1000a9899c97SHong Zhang 10010d3441aeSHong Zhang ierr = MatZeroEntries(C);CHKERRQ(ierr); 1002e497d3c8SHong Zhang /* 1) get R = Pd^T,Ro = Po^T */ 1003748c7196SHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 1004419ecdd9Sandi selinger ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 1005419ecdd9Sandi selinger ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 1006748c7196SHong Zhang } 10070d3441aeSHong Zhang 1008e497d3c8SHong Zhang /* 2) get AP_loc */ 10090d3441aeSHong Zhang AP_loc = ptap->AP_loc; 10108cb82516SHong Zhang ap = (Mat_SeqAIJ*)AP_loc->data; 10110d3441aeSHong Zhang 1012e497d3c8SHong Zhang /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 10130d3441aeSHong Zhang /*-----------------------------------------------------*/ 1014748c7196SHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 1015748c7196SHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 10160d3441aeSHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 10170d3441aeSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 1018e497d3c8SHong Zhang } 10190d3441aeSHong Zhang 10208cb82516SHong Zhang /* 2-2) compute numeric A_loc*P - dominating part */ 10218cb82516SHong Zhang /* ---------------------------------------------- */ 10220d3441aeSHong Zhang /* get data from symbolic products */ 10238cb82516SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 10242dd9e643SHong Zhang if (ptap->P_oth) { 10258cb82516SHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 10262dd9e643SHong Zhang } 10278cb82516SHong Zhang apa = ptap->apa; 1028681d504bSHong Zhang api = ap->i; 1029681d504bSHong Zhang apj = ap->j; 1030e497d3c8SHong Zhang for (i=0; i<am; i++) { 1031445158ffSHong Zhang /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 1032e497d3c8SHong Zhang AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa); 1033e497d3c8SHong Zhang apnz = api[i+1] - api[i]; 1034e497d3c8SHong Zhang for (j=0; j<apnz; j++) { 1035e497d3c8SHong Zhang col = apj[j+api[i]]; 1036e497d3c8SHong Zhang ap->a[j+ap->i[i]] = apa[col]; 1037e497d3c8SHong Zhang apa[col] = 0.0; 10380d3441aeSHong Zhang } 1039aa690a28SHong Zhang ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 10400d3441aeSHong Zhang } 10410d3441aeSHong Zhang 10428cb82516SHong Zhang /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 1043445158ffSHong Zhang ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 1044445158ffSHong Zhang ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 10459ce11a7cSHong Zhang C_loc = ptap->C_loc; 10469ce11a7cSHong Zhang C_oth = ptap->C_oth; 10470d3441aeSHong Zhang 10480d3441aeSHong Zhang /* add C_loc and Co to to C */ 10490d3441aeSHong Zhang ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 10500d3441aeSHong Zhang 10510d3441aeSHong Zhang /* C_loc -> C */ 10520d3441aeSHong Zhang cm = C_loc->rmap->N; 10539ce11a7cSHong Zhang c_seq = (Mat_SeqAIJ*)C_loc->data; 10548cb82516SHong Zhang cols = c_seq->j; 10558cb82516SHong Zhang vals = c_seq->a; 1056904d1e70Sandi selinger 1057904d1e70Sandi selinger 1058e9ede7d0Sandi selinger /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */ 1059904d1e70Sandi selinger /* when there are no off-processor parts. */ 10601de21080Sandi selinger /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */ 10611de21080Sandi selinger /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */ 10621de21080Sandi selinger /* a table, and other, more complex stuff has to be done. */ 1063904d1e70Sandi selinger if (C->assembled) { 1064904d1e70Sandi selinger C->was_assembled = PETSC_TRUE; 1065904d1e70Sandi selinger C->assembled = PETSC_FALSE; 1066904d1e70Sandi selinger } 1067904d1e70Sandi selinger if (C->was_assembled) { 10680d3441aeSHong Zhang for (i=0; i<cm; i++) { 10699ce11a7cSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 10700d3441aeSHong Zhang row = rstart + i; 1071904d1e70Sandi selinger ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 10728cb82516SHong Zhang cols += ncols; vals += ncols; 10730d3441aeSHong Zhang } 1074904d1e70Sandi selinger } else { 1075e9ede7d0Sandi selinger ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr); 1076904d1e70Sandi selinger } 10770d3441aeSHong Zhang 10780d3441aeSHong Zhang /* Co -> C, off-processor part */ 10799ce11a7cSHong Zhang cm = C_oth->rmap->N; 10809ce11a7cSHong Zhang c_seq = (Mat_SeqAIJ*)C_oth->data; 10818cb82516SHong Zhang cols = c_seq->j; 10828cb82516SHong Zhang vals = c_seq->a; 10839ce11a7cSHong Zhang for (i=0; i<cm; i++) { 10849ce11a7cSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 10850d3441aeSHong Zhang row = p->garray[i]; 10860d3441aeSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 10878cb82516SHong Zhang cols += ncols; vals += ncols; 10880d3441aeSHong Zhang } 1089904d1e70Sandi selinger 10900d3441aeSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10910d3441aeSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10920d3441aeSHong Zhang 1093748c7196SHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 10943697aca6SHong Zhang 1095dd4011a9SHong Zhang ierr = PetscObjectOptionsBegin((PetscObject)C);CHKERRQ(ierr); 1096dd4011a9SHong Zhang PetscOptionsObject->alreadyprinted = PETSC_FALSE; /* a hack to ensure the option shows in '-help' */ 1097dd4011a9SHong Zhang /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */ 1098dd4011a9SHong Zhang freestruct = PETSC_FALSE; 1099dd4011a9SHong Zhang ierr = PetscOptionsBool("-mat_freeintermediatedatastructures","Free intermediate data structures", "MatFreeIntermediateDataStructures",freestruct, &freestruct, NULL);CHKERRQ(ierr); 1100dd4011a9SHong Zhang if (freestruct) { 1101dd4011a9SHong Zhang ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 11023697aca6SHong Zhang } 110380da8df7SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 11040d3441aeSHong Zhang PetscFunctionReturn(0); 11050d3441aeSHong Zhang } 1106