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> 12694f88d4SFande Kong #include <petsc/private/hashmapiv.h> 134a56b808SFande Kong #include <petsc/private/hashseti.h> 144a56b808SFande Kong #include <petscsf.h> 1542c57489SHong Zhang 1624ecddacSHong Zhang 17cf97cf3cSHong Zhang PetscErrorCode MatView_MPIAIJ_PtAP(Mat A,PetscViewer viewer) 18cf97cf3cSHong Zhang { 19cf97cf3cSHong Zhang PetscErrorCode ierr; 20cf97cf3cSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 213cdca5ebSHong Zhang Mat_APMPI *ptap=a->ap; 22cf97cf3cSHong Zhang PetscBool iascii; 23cf97cf3cSHong Zhang PetscViewerFormat format; 24cf97cf3cSHong Zhang 25cf97cf3cSHong Zhang PetscFunctionBegin; 26baa1ba78SHong Zhang if (!ptap) { 27baa1ba78SHong Zhang /* hack: MatDuplicate() sets oldmat->ops->view to newmat which is a base mat class with null ptpa! */ 28baa1ba78SHong Zhang A->ops->view = MatView_MPIAIJ; 29baa1ba78SHong Zhang ierr = (A->ops->view)(A,viewer);CHKERRQ(ierr); 30baa1ba78SHong Zhang PetscFunctionReturn(0); 31baa1ba78SHong Zhang } 32baa1ba78SHong Zhang 33cf97cf3cSHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 34cf97cf3cSHong Zhang if (iascii) { 35cf97cf3cSHong Zhang ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 36cf97cf3cSHong Zhang if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 37cf97cf3cSHong Zhang if (ptap->algType == 0) { 38cf97cf3cSHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"using scalable MatPtAP() implementation\n");CHKERRQ(ierr); 39cf97cf3cSHong Zhang } else if (ptap->algType == 1) { 40cf97cf3cSHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"using nonscalable MatPtAP() implementation\n");CHKERRQ(ierr); 414a56b808SFande Kong } else if (ptap->algType == 2) { 424a56b808SFande Kong ierr = PetscViewerASCIIPrintf(viewer,"using allatonce MatPtAP() implementation\n");CHKERRQ(ierr); 434a56b808SFande Kong } else if (ptap->algType == 3) { 444a56b808SFande Kong ierr = PetscViewerASCIIPrintf(viewer,"using merged allatonce MatPtAP() implementation\n");CHKERRQ(ierr); 45cf97cf3cSHong Zhang } 46cf97cf3cSHong Zhang } 47cf97cf3cSHong Zhang } 48cf97cf3cSHong Zhang ierr = (ptap->view)(A,viewer);CHKERRQ(ierr); 49cf97cf3cSHong Zhang PetscFunctionReturn(0); 50cf97cf3cSHong Zhang } 51cf97cf3cSHong Zhang 524624976aSHong Zhang PetscErrorCode MatFreeIntermediateDataStructures_MPIAIJ_AP(Mat A) 53cc6cb787SHong Zhang { 54cc6cb787SHong Zhang PetscErrorCode ierr; 55f8487c73SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 563cdca5ebSHong Zhang Mat_APMPI *ptap=a->ap; 5760552ceaSHong Zhang Mat_Merge_SeqsToMPI *merge; 58cc6cb787SHong Zhang 59cc6cb787SHong Zhang PetscFunctionBegin; 60dd4011a9SHong Zhang if (!ptap) PetscFunctionReturn(0); 61dd4011a9SHong Zhang 62b7f45c76SHong Zhang ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr); 63f8487c73SHong Zhang ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 64a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 65a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 66c5af039cSHong Zhang ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */ 6705d62848SHong Zhang ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr); 6805d62848SHong Zhang ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr); 698403a639SHong Zhang if (ptap->AP_loc) { /* used by alg_rap */ 70681d504bSHong Zhang Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data; 71681d504bSHong Zhang ierr = PetscFree(ap->i);CHKERRQ(ierr); 72681d504bSHong Zhang ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr); 730d3441aeSHong Zhang ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr); 748403a639SHong Zhang } else { /* used by alg_ptap */ 758403a639SHong Zhang ierr = PetscFree(ptap->api);CHKERRQ(ierr); 768403a639SHong Zhang ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 77681d504bSHong Zhang } 782259aa2eSHong Zhang ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr); 792259aa2eSHong Zhang ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr); 80414894bdSHong Zhang if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);} 81a560ef98SHong Zhang 8260552ceaSHong Zhang ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr); 8360552ceaSHong Zhang 8460552ceaSHong Zhang merge=ptap->merge; 858403a639SHong Zhang if (merge) { /* used by alg_ptap */ 86cc6cb787SHong Zhang ierr = PetscFree(merge->id_r);CHKERRQ(ierr); 87cc6cb787SHong Zhang ierr = PetscFree(merge->len_s);CHKERRQ(ierr); 88cc6cb787SHong Zhang ierr = PetscFree(merge->len_r);CHKERRQ(ierr); 89cc6cb787SHong Zhang ierr = PetscFree(merge->bi);CHKERRQ(ierr); 90cc6cb787SHong Zhang ierr = PetscFree(merge->bj);CHKERRQ(ierr); 91c05d87d6SBarry Smith ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr); 92cc6cb787SHong Zhang ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr); 93c05d87d6SBarry Smith ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr); 94cc6cb787SHong Zhang ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr); 95445158ffSHong Zhang ierr = PetscFree(merge->coi);CHKERRQ(ierr); 96445158ffSHong Zhang ierr = PetscFree(merge->coj);CHKERRQ(ierr); 9705b42c5fSBarry Smith ierr = PetscFree(merge->owners_co);CHKERRQ(ierr); 986bf464f9SBarry Smith ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr); 99f8487c73SHong Zhang ierr = PetscFree(ptap->merge);CHKERRQ(ierr); 100bf0cc555SLisandro Dalcin } 101a3bb6f32SFande Kong ierr = ISLocalToGlobalMappingDestroy(&ptap->ltog);CHKERRQ(ierr); 1024a56b808SFande Kong 1034a56b808SFande Kong ierr = PetscSFDestroy(&ptap->sf);CHKERRQ(ierr); 1044a56b808SFande Kong ierr = PetscFree(ptap->c_othi);CHKERRQ(ierr); 1054a56b808SFande Kong ierr = PetscFree(ptap->c_rmti);CHKERRQ(ierr); 1063697aca6SHong Zhang PetscFunctionReturn(0); 107c8b0795cSMark F. Adams } 108dce485f0SBarry Smith 1093697aca6SHong Zhang PetscErrorCode MatDestroy_MPIAIJ_PtAP(Mat A) 1103697aca6SHong Zhang { 1113697aca6SHong Zhang PetscErrorCode ierr; 1129b3d8a68SHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data; 1139b3d8a68SHong Zhang Mat_APMPI *ptap=a->ap; 1143697aca6SHong Zhang 1153697aca6SHong Zhang PetscFunctionBegin; 116dd4011a9SHong Zhang ierr = (*A->ops->freeintermediatedatastructures)(A);CHKERRQ(ierr); 1179b3d8a68SHong Zhang ierr = (*ptap->destroy)(A);CHKERRQ(ierr); /* MatDestroy_MPIAIJ(A) */ 1189b3d8a68SHong Zhang ierr = PetscFree(ptap);CHKERRQ(ierr); 119cc6cb787SHong Zhang PetscFunctionReturn(0); 120cc6cb787SHong Zhang } 121cc6cb787SHong Zhang 122150d2497SBarry Smith PETSC_INTERN PetscErrorCode MatPtAP_MPIAIJ_MPIAIJ(Mat A,Mat P,MatReuse scall,PetscReal fill,Mat *C) 12342c57489SHong Zhang { 12442c57489SHong Zhang PetscErrorCode ierr; 125a4ffb656SHong Zhang PetscBool flg; 12667a12041SHong Zhang MPI_Comm comm; 127a4ffb656SHong Zhang #if !defined(PETSC_HAVE_HYPRE) 1284a56b808SFande Kong const char *algTypes[4] = {"scalable","nonscalable","allatonce","allatonce_merged"}; 1294a56b808SFande Kong PetscInt nalg=4; 130a4ffb656SHong Zhang #else 1314a56b808SFande Kong const char *algTypes[5] = {"scalable","nonscalable","allatonce","allatonce_merged","hypre"}; 1324a56b808SFande Kong PetscInt nalg=5; 133a4ffb656SHong Zhang #endif 134a4ffb656SHong Zhang PetscInt pN=P->cmap->N,alg=1; /* set default algorithm */ 13542c57489SHong Zhang 13642c57489SHong Zhang PetscFunctionBegin; 13767a12041SHong Zhang /* check if matrix local sizes are compatible */ 13867a12041SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1396c4ed002SBarry 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); 1406c4ed002SBarry 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); 14167a12041SHong Zhang 142cf3ca8ceSHong Zhang if (scall == MAT_INITIAL_MATRIX) { 143a4ffb656SHong Zhang /* pick an algorithm */ 144715a5346SHong Zhang ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatPtAP","Mat");CHKERRQ(ierr); 145a4ffb656SHong Zhang ierr = PetscOptionsEList("-matptap_via","Algorithmic approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 146a4ffb656SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 147a4ffb656SHong Zhang 148a4ffb656SHong Zhang if (!flg && pN > 100000) { /* may switch to scalable algorithm as default */ 149a4ffb656SHong Zhang MatInfo Ainfo,Pinfo; 150a4ffb656SHong Zhang PetscInt nz_local; 151a4ffb656SHong Zhang PetscBool alg_scalable_loc=PETSC_FALSE,alg_scalable; 152a4ffb656SHong Zhang 153a4ffb656SHong Zhang ierr = MatGetInfo(A,MAT_LOCAL,&Ainfo);CHKERRQ(ierr); 154a4ffb656SHong Zhang ierr = MatGetInfo(P,MAT_LOCAL,&Pinfo);CHKERRQ(ierr); 155a4ffb656SHong Zhang nz_local = (PetscInt)(Ainfo.nz_allocated + Pinfo.nz_allocated); 156a4ffb656SHong Zhang 157a4ffb656SHong Zhang if (pN > fill*nz_local) alg_scalable_loc = PETSC_TRUE; 158a4ffb656SHong Zhang ierr = MPIU_Allreduce(&alg_scalable_loc,&alg_scalable,1,MPIU_BOOL,MPI_LOR,comm);CHKERRQ(ierr); 159a4ffb656SHong Zhang 160a4ffb656SHong Zhang if (alg_scalable) { 161a4ffb656SHong Zhang alg = 0; /* scalable algorithm would 50% slower than nonscalable algorithm */ 1620d3441aeSHong Zhang } 163a4ffb656SHong Zhang } 164a4ffb656SHong Zhang 165a4ffb656SHong Zhang switch (alg) { 166a4ffb656SHong Zhang case 1: 16780da8df7SHong Zhang /* do R=P^T locally, then C=R*A*P -- nonscalable */ 168a4ffb656SHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 169a4ffb656SHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr); 1703ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 171a4ffb656SHong Zhang break; 172a4ffb656SHong Zhang case 2: 1734a56b808SFande Kong /* compute C=P^T*A*P allatonce */ 1744a56b808SFande Kong ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 1754a56b808SFande Kong ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce(A,P,fill,C);CHKERRQ(ierr); 1764a56b808SFande Kong ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 1774a56b808SFande Kong break; 1784a56b808SFande Kong case 3: 1794a56b808SFande Kong /* compute C=P^T*A*P allatonce */ 1804a56b808SFande Kong ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 1814a56b808SFande Kong ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce_merged(A,P,fill,C);CHKERRQ(ierr); 1824a56b808SFande Kong ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 1834a56b808SFande Kong break; 1844a56b808SFande Kong #if defined(PETSC_HAVE_HYPRE) 1854a56b808SFande Kong case 4: 186a4ffb656SHong Zhang /* Use boomerAMGBuildCoarseOperator */ 1870abfe76eSFande Kong ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 188a4ffb656SHong Zhang ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr); 1890abfe76eSFande Kong ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 190a4ffb656SHong Zhang break; 191a4ffb656SHong Zhang #endif 192a4ffb656SHong Zhang default: 193a4ffb656SHong Zhang /* do R=P^T locally, then C=R*A*P */ 194a4ffb656SHong Zhang ierr = PetscLogEventBegin(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 195a4ffb656SHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(A,P,fill,C);CHKERRQ(ierr); 196a4ffb656SHong Zhang ierr = PetscLogEventEnd(MAT_PtAPSymbolic,A,P,0,0);CHKERRQ(ierr); 197a4ffb656SHong Zhang break; 198a4ffb656SHong Zhang } 1997d0a89b7SHong Zhang 2004a56b808SFande Kong if (alg == 0 || alg == 1 || alg == 2 || alg == 3) { 2017d0a89b7SHong Zhang Mat_MPIAIJ *c = (Mat_MPIAIJ*)(*C)->data; 2027d0a89b7SHong Zhang Mat_APMPI *ap = c->ap; 2037d0a89b7SHong Zhang ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)(*C)),((PetscObject)(*C))->prefix,"MatFreeIntermediateDataStructures","Mat");CHKERRQ(ierr); 2047d0a89b7SHong Zhang ap->freestruct = PETSC_FALSE; 2057d0a89b7SHong Zhang ierr = PetscOptionsBool("-mat_freeintermediatedatastructures","Free intermediate data structures", "MatFreeIntermediateDataStructures",ap->freestruct,&ap->freestruct, NULL);CHKERRQ(ierr); 2067d0a89b7SHong Zhang ierr = PetscOptionsEnd();CHKERRQ(ierr); 2077a7894deSKris Buschelman } 2087d0a89b7SHong Zhang } 2097d0a89b7SHong Zhang 2103ff4c91cSHong Zhang ierr = PetscLogEventBegin(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 2118403a639SHong Zhang ierr = (*(*C)->ops->ptapnumeric)(A,P,*C);CHKERRQ(ierr); 2123ff4c91cSHong Zhang ierr = PetscLogEventEnd(MAT_PtAPNumeric,A,P,0,0);CHKERRQ(ierr); 21342c57489SHong Zhang PetscFunctionReturn(0); 21442c57489SHong Zhang } 21542c57489SHong Zhang 216cf742fccSHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,Mat C) 217cf742fccSHong Zhang { 218cf742fccSHong Zhang PetscErrorCode ierr; 219cf742fccSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 220cf742fccSHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 22192ec70a1SHong Zhang Mat_SeqAIJ *ap,*p_loc,*p_oth=NULL,*c_seq; 2223cdca5ebSHong Zhang Mat_APMPI *ptap = c->ap; 223cf742fccSHong Zhang Mat AP_loc,C_loc,C_oth; 224a3bb6f32SFande Kong PetscInt i,rstart,rend,cm,ncols,row,*api,*apj,am = A->rmap->n,apnz,nout; 225cf742fccSHong Zhang PetscScalar *apa; 226cf742fccSHong Zhang const PetscInt *cols; 227cf742fccSHong Zhang const PetscScalar *vals; 228cf742fccSHong Zhang 229cf742fccSHong Zhang PetscFunctionBegin; 230fa2a5dcfSHong Zhang if (!ptap->AP_loc) { 23180da8df7SHong Zhang MPI_Comm comm; 23280da8df7SHong Zhang ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 23380da8df7SHong Zhang SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 23480da8df7SHong Zhang } 23580da8df7SHong Zhang 236cf742fccSHong Zhang ierr = MatZeroEntries(C);CHKERRQ(ierr); 237cf742fccSHong Zhang 238cf742fccSHong Zhang /* 1) get R = Pd^T,Ro = Po^T */ 239cf742fccSHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 240419ecdd9Sandi selinger ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 241419ecdd9Sandi selinger ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 242cf742fccSHong Zhang } 243cf742fccSHong Zhang 244cf742fccSHong Zhang /* 2) get AP_loc */ 245cf742fccSHong Zhang AP_loc = ptap->AP_loc; 246cf742fccSHong Zhang ap = (Mat_SeqAIJ*)AP_loc->data; 247cf742fccSHong Zhang 248cf742fccSHong Zhang /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 249cf742fccSHong Zhang /*-----------------------------------------------------*/ 250cf742fccSHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 251cf742fccSHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 252cf742fccSHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 253cf742fccSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 254cf742fccSHong Zhang } 255cf742fccSHong Zhang 256cf742fccSHong Zhang /* 2-2) compute numeric A_loc*P - dominating part */ 257cf742fccSHong Zhang /* ---------------------------------------------- */ 258cf742fccSHong Zhang /* get data from symbolic products */ 259cf742fccSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 26052f7967eSHong Zhang if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 26152f7967eSHong Zhang 262cf742fccSHong Zhang api = ap->i; 263cf742fccSHong Zhang apj = ap->j; 264a3bb6f32SFande Kong ierr = ISLocalToGlobalMappingApply(ptap->ltog,api[AP_loc->rmap->n],apj,apj);CHKERRQ(ierr); 265cf742fccSHong Zhang for (i=0; i<am; i++) { 266cf742fccSHong Zhang /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 267cf742fccSHong Zhang apnz = api[i+1] - api[i]; 268b4e8d1b6SHong Zhang apa = ap->a + api[i]; 269b4e8d1b6SHong Zhang ierr = PetscMemzero(apa,sizeof(PetscScalar)*apnz);CHKERRQ(ierr); 270b4e8d1b6SHong Zhang AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa); 271cf742fccSHong Zhang ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 272cf742fccSHong Zhang } 273a3bb6f32SFande Kong ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,api[AP_loc->rmap->n],apj,&nout,apj);CHKERRQ(ierr); 274a3bb6f32SFande Kong if (api[AP_loc->rmap->n] != nout) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incorrect mapping %D != %D\n",api[AP_loc->rmap->n],nout); 275cf742fccSHong Zhang 276cf742fccSHong Zhang /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 277154d0d78SFande Kong /* Always use scalable version since we are in the MPI scalable version */ 278cf742fccSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 279cf742fccSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 280cf742fccSHong Zhang 281cf742fccSHong Zhang C_loc = ptap->C_loc; 282cf742fccSHong Zhang C_oth = ptap->C_oth; 283cf742fccSHong Zhang 284cf742fccSHong Zhang /* add C_loc and Co to to C */ 285cf742fccSHong Zhang ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 286cf742fccSHong Zhang 287cf742fccSHong Zhang /* C_loc -> C */ 288cf742fccSHong Zhang cm = C_loc->rmap->N; 289cf742fccSHong Zhang c_seq = (Mat_SeqAIJ*)C_loc->data; 290cf742fccSHong Zhang cols = c_seq->j; 291cf742fccSHong Zhang vals = c_seq->a; 292a3bb6f32SFande Kong ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_seq->i[C_loc->rmap->n],c_seq->j,c_seq->j);CHKERRQ(ierr); 293edeac6deSandi selinger 294e9ede7d0Sandi selinger /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */ 295edeac6deSandi selinger /* when there are no off-processor parts. */ 2961de21080Sandi selinger /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */ 2971de21080Sandi selinger /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */ 2981de21080Sandi selinger /* a table, and other, more complex stuff has to be done. */ 299edeac6deSandi selinger if (C->assembled) { 300edeac6deSandi selinger C->was_assembled = PETSC_TRUE; 301edeac6deSandi selinger C->assembled = PETSC_FALSE; 302edeac6deSandi selinger } 303edeac6deSandi selinger if (C->was_assembled) { 304cf742fccSHong Zhang for (i=0; i<cm; i++) { 305cf742fccSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 306cf742fccSHong Zhang row = rstart + i; 307edeac6deSandi selinger ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 308cf742fccSHong Zhang cols += ncols; vals += ncols; 309cf742fccSHong Zhang } 310edeac6deSandi selinger } else { 311e9ede7d0Sandi selinger ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr); 312edeac6deSandi selinger } 313a3bb6f32SFande Kong ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_seq->i[C_loc->rmap->n],c_seq->j,&nout,c_seq->j);CHKERRQ(ierr); 314a3bb6f32SFande Kong if (c_seq->i[C_loc->rmap->n] != nout) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incorrect mapping %D != %D\n",c_seq->i[C_loc->rmap->n],nout); 315cf742fccSHong Zhang 316cf742fccSHong Zhang /* Co -> C, off-processor part */ 317cf742fccSHong Zhang cm = C_oth->rmap->N; 318cf742fccSHong Zhang c_seq = (Mat_SeqAIJ*)C_oth->data; 319cf742fccSHong Zhang cols = c_seq->j; 320cf742fccSHong Zhang vals = c_seq->a; 321a3bb6f32SFande Kong ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_seq->i[C_oth->rmap->n],c_seq->j,c_seq->j);CHKERRQ(ierr); 322cf742fccSHong Zhang for (i=0; i<cm; i++) { 323cf742fccSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 324cf742fccSHong Zhang row = p->garray[i]; 325cf742fccSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 326cf742fccSHong Zhang cols += ncols; vals += ncols; 327cf742fccSHong Zhang } 328cf742fccSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 329cf742fccSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 330cf742fccSHong Zhang 331cf742fccSHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 33280da8df7SHong Zhang 333a3bb6f32SFande Kong ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_seq->i[C_oth->rmap->n],c_seq->j,&nout,c_seq->j);CHKERRQ(ierr); 334a3bb6f32SFande Kong if (c_seq->i[C_oth->rmap->n] != nout) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incorrect mapping %D != %D\n",c_seq->i[C_loc->rmap->n],nout); 335a3bb6f32SFande Kong 33680da8df7SHong 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 */ 3377d0a89b7SHong Zhang if (ptap->freestruct) { 33880da8df7SHong Zhang ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 33980da8df7SHong Zhang } 340cf742fccSHong Zhang PetscFunctionReturn(0); 341cf742fccSHong Zhang } 342cf742fccSHong Zhang 343e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat *C) 3440d3441aeSHong Zhang { 3450d3441aeSHong Zhang PetscErrorCode ierr; 3463cdca5ebSHong Zhang Mat_APMPI *ptap; 347681d504bSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 3482259aa2eSHong Zhang MPI_Comm comm; 3492259aa2eSHong Zhang PetscMPIInt size,rank; 35076db11f6SHong Zhang Mat Cmpi,P_loc,P_oth; 35115a3b8e2SHong Zhang PetscFreeSpaceList free_space=NULL,current_space=NULL; 352d0e9a2f7SHong Zhang PetscInt am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n; 353d0e9a2f7SHong Zhang PetscInt *lnk,i,k,pnz,row,nsend; 354f671be37SHong Zhang PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 35515a3b8e2SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 356d0e9a2f7SHong Zhang PetscInt len,proc,*dnz,*onz,*owners,nzi,nspacedouble; 35715a3b8e2SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 35815a3b8e2SHong Zhang MPI_Request *swaits,*rwaits; 35915a3b8e2SHong Zhang MPI_Status *sstatus,rstatus; 360445158ffSHong Zhang PetscLayout rowmap; 361445158ffSHong Zhang PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 362445158ffSHong Zhang PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 363a3bb6f32SFande Kong PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,Crmax,*aj,*ai,*pi,nout; 36452f7967eSHong Zhang Mat_SeqAIJ *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth; 36567a12041SHong Zhang PetscScalar *apv; 36678d30b94SHong Zhang PetscTable ta; 367b92f168fSBarry Smith MatType mtype; 368e83e5f86SFande Kong const char *prefix; 369aa690a28SHong Zhang #if defined(PETSC_USE_INFO) 3708cb82516SHong Zhang PetscReal apfill; 371aa690a28SHong Zhang #endif 37267a12041SHong Zhang 37367a12041SHong Zhang PetscFunctionBegin; 37467a12041SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 37567a12041SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 37667a12041SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 377ae5f0867Sstefano_zampini 37852f7967eSHong Zhang if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data; 37952f7967eSHong Zhang 380ae5f0867Sstefano_zampini /* create symbolic parallel matrix Cmpi */ 381ae5f0867Sstefano_zampini ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 382b92f168fSBarry Smith ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 383b92f168fSBarry Smith ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 384ae5f0867Sstefano_zampini 3853cdca5ebSHong Zhang /* create struct Mat_APMPI and attached it to C later */ 386e953e47cSHong Zhang ierr = PetscNew(&ptap);CHKERRQ(ierr); 387e953e47cSHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 388cf97cf3cSHong Zhang ptap->algType = 0; 389e953e47cSHong Zhang 390e953e47cSHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 39176db11f6SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&P_oth);CHKERRQ(ierr); 392e953e47cSHong Zhang /* get P_loc by taking all local rows of P */ 39376db11f6SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&P_loc);CHKERRQ(ierr); 39476db11f6SHong Zhang 39576db11f6SHong Zhang ptap->P_loc = P_loc; 39676db11f6SHong Zhang ptap->P_oth = P_oth; 397e953e47cSHong Zhang 398e953e47cSHong Zhang /* (0) compute Rd = Pd^T, Ro = Po^T */ 399e953e47cSHong Zhang /* --------------------------------- */ 400419ecdd9Sandi selinger ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 401419ecdd9Sandi selinger ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 402e953e47cSHong Zhang 403e953e47cSHong Zhang /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 404e953e47cSHong Zhang /* ----------------------------------------------------------------- */ 40576db11f6SHong Zhang p_loc = (Mat_SeqAIJ*)P_loc->data; 40652f7967eSHong Zhang if (P_oth) p_oth = (Mat_SeqAIJ*)P_oth->data; 407e953e47cSHong Zhang 408e953e47cSHong Zhang /* create and initialize a linked list */ 409e953e47cSHong Zhang ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */ 41076db11f6SHong Zhang MatRowMergeMax_SeqAIJ(p_loc,P_loc->rmap->N,ta); 41176db11f6SHong Zhang MatRowMergeMax_SeqAIJ(p_oth,P_oth->rmap->N,ta); 412e953e47cSHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */ 413e953e47cSHong Zhang 41476db11f6SHong Zhang ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 415e953e47cSHong Zhang 416e953e47cSHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 41752f7967eSHong Zhang if (ao) { 418e953e47cSHong Zhang ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr); 41952f7967eSHong Zhang } else { 42052f7967eSHong Zhang ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr); 42152f7967eSHong Zhang } 422e953e47cSHong Zhang current_space = free_space; 423e953e47cSHong Zhang nspacedouble = 0; 424e953e47cSHong Zhang 425e953e47cSHong Zhang ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 426e953e47cSHong Zhang api[0] = 0; 427e953e47cSHong Zhang for (i=0; i<am; i++) { 428e953e47cSHong Zhang /* diagonal portion: Ad[i,:]*P */ 429e953e47cSHong Zhang ai = ad->i; pi = p_loc->i; 430e953e47cSHong Zhang nzi = ai[i+1] - ai[i]; 431e953e47cSHong Zhang aj = ad->j + ai[i]; 432e953e47cSHong Zhang for (j=0; j<nzi; j++) { 433e953e47cSHong Zhang row = aj[j]; 434e953e47cSHong Zhang pnz = pi[row+1] - pi[row]; 435e953e47cSHong Zhang Jptr = p_loc->j + pi[row]; 436e953e47cSHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 43776db11f6SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 438e953e47cSHong Zhang } 439e953e47cSHong Zhang /* off-diagonal portion: Ao[i,:]*P */ 44052f7967eSHong Zhang if (ao) { 441e953e47cSHong Zhang ai = ao->i; pi = p_oth->i; 442e953e47cSHong Zhang nzi = ai[i+1] - ai[i]; 443e953e47cSHong Zhang aj = ao->j + ai[i]; 444e953e47cSHong Zhang for (j=0; j<nzi; j++) { 445e953e47cSHong Zhang row = aj[j]; 446e953e47cSHong Zhang pnz = pi[row+1] - pi[row]; 447e953e47cSHong Zhang Jptr = p_oth->j + pi[row]; 44876db11f6SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 449e953e47cSHong Zhang } 45052f7967eSHong Zhang } 451e953e47cSHong Zhang apnz = lnk[0]; 452e953e47cSHong Zhang api[i+1] = api[i] + apnz; 453e953e47cSHong Zhang 454e953e47cSHong Zhang /* if free space is not available, double the total space in the list */ 455e953e47cSHong Zhang if (current_space->local_remaining<apnz) { 456e953e47cSHong Zhang ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 457e953e47cSHong Zhang nspacedouble++; 458e953e47cSHong Zhang } 459e953e47cSHong Zhang 460e953e47cSHong Zhang /* Copy data into free space, then initialize lnk */ 46176db11f6SHong Zhang ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr); 462e953e47cSHong Zhang 463e953e47cSHong Zhang current_space->array += apnz; 464e953e47cSHong Zhang current_space->local_used += apnz; 465e953e47cSHong Zhang current_space->local_remaining -= apnz; 466e953e47cSHong Zhang } 467e953e47cSHong Zhang /* Allocate space for apj and apv, initialize apj, and */ 468e953e47cSHong Zhang /* destroy list of free space and other temporary array(s) */ 469a3bb6f32SFande Kong ierr = PetscCalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 470e953e47cSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 47176db11f6SHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 472e953e47cSHong Zhang 473e953e47cSHong Zhang /* Create AP_loc for reuse */ 474e953e47cSHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 475a3bb6f32SFande Kong ierr = MatSeqAIJCompactOutExtraColumns_SeqAIJ(ptap->AP_loc, &ptap->ltog);CHKERRQ(ierr); 476e953e47cSHong Zhang 477e953e47cSHong Zhang #if defined(PETSC_USE_INFO) 47852f7967eSHong Zhang if (ao) { 479e953e47cSHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 48052f7967eSHong Zhang } else { 48152f7967eSHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1); 48252f7967eSHong Zhang } 483e953e47cSHong Zhang ptap->AP_loc->info.mallocs = nspacedouble; 484e953e47cSHong Zhang ptap->AP_loc->info.fill_ratio_given = fill; 485e953e47cSHong Zhang ptap->AP_loc->info.fill_ratio_needed = apfill; 486e953e47cSHong Zhang 487e953e47cSHong Zhang if (api[am]) { 488b11c1ec8SHong 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); 489e953e47cSHong Zhang ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 490e953e47cSHong Zhang } else { 491b11c1ec8SHong Zhang ierr = PetscInfo(ptap->AP_loc,"Scalable algorithm, AP_loc is empty \n");CHKERRQ(ierr); 492e953e47cSHong Zhang } 493e953e47cSHong Zhang #endif 494e953e47cSHong Zhang 495e953e47cSHong Zhang /* (2-1) compute symbolic Co = Ro*AP_loc */ 496e953e47cSHong Zhang /* ------------------------------------ */ 497e83e5f86SFande Kong ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 498e83e5f86SFande Kong ierr = MatSetOptionsPrefix(ptap->Ro,prefix);CHKERRQ(ierr); 499e83e5f86SFande Kong ierr = MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");CHKERRQ(ierr); 500e953e47cSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 501e953e47cSHong Zhang 502e953e47cSHong Zhang /* (3) send coj of C_oth to other processors */ 503e953e47cSHong Zhang /* ------------------------------------------ */ 504e953e47cSHong Zhang /* determine row ownership */ 505e953e47cSHong Zhang ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 506e953e47cSHong Zhang rowmap->n = pn; 507e953e47cSHong Zhang rowmap->bs = 1; 508e953e47cSHong Zhang ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 509e953e47cSHong Zhang owners = rowmap->range; 510e953e47cSHong Zhang 511e953e47cSHong Zhang /* determine the number of messages to send, their lengths */ 512e953e47cSHong Zhang ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 513e953e47cSHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 514e953e47cSHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 515e953e47cSHong Zhang 516e953e47cSHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 517e953e47cSHong Zhang coi = c_oth->i; coj = c_oth->j; 518e953e47cSHong Zhang con = ptap->C_oth->rmap->n; 519e953e47cSHong Zhang proc = 0; 520a3bb6f32SFande Kong ierr = ISLocalToGlobalMappingApply(ptap->ltog,coi[con],coj,coj);CHKERRQ(ierr); 521e953e47cSHong Zhang for (i=0; i<con; i++) { 522e953e47cSHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 523e953e47cSHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 524e953e47cSHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 525e953e47cSHong Zhang } 526e953e47cSHong Zhang 527e953e47cSHong Zhang len = 0; /* max length of buf_si[], see (4) */ 528e953e47cSHong Zhang owners_co[0] = 0; 529e953e47cSHong Zhang nsend = 0; 530e953e47cSHong Zhang for (proc=0; proc<size; proc++) { 531e953e47cSHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 532e953e47cSHong Zhang if (len_s[proc]) { 533e953e47cSHong Zhang nsend++; 534e953e47cSHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 535e953e47cSHong Zhang len += len_si[proc]; 536e953e47cSHong Zhang } 537e953e47cSHong Zhang } 538e953e47cSHong Zhang 539e953e47cSHong Zhang /* determine the number and length of messages to receive for coi and coj */ 540e953e47cSHong Zhang ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 541e953e47cSHong Zhang ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 542e953e47cSHong Zhang 543e953e47cSHong Zhang /* post the Irecv and Isend of coj */ 544e953e47cSHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 545e953e47cSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 546e953e47cSHong Zhang ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 547e953e47cSHong Zhang for (proc=0, k=0; proc<size; proc++) { 548e953e47cSHong Zhang if (!len_s[proc]) continue; 549e953e47cSHong Zhang i = owners_co[proc]; 550e953e47cSHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 551e953e47cSHong Zhang k++; 552e953e47cSHong Zhang } 553e953e47cSHong Zhang 554e953e47cSHong Zhang /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 555e953e47cSHong Zhang /* ---------------------------------------- */ 556e83e5f86SFande Kong ierr = MatSetOptionsPrefix(ptap->Rd,prefix);CHKERRQ(ierr); 557e83e5f86SFande Kong ierr = MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");CHKERRQ(ierr); 558e953e47cSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 559e953e47cSHong Zhang c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 560a3bb6f32SFande Kong ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_loc->i[ptap->C_loc->rmap->n],c_loc->j,c_loc->j);CHKERRQ(ierr); 561e953e47cSHong Zhang 562e953e47cSHong Zhang /* receives coj are complete */ 563e953e47cSHong Zhang for (i=0; i<nrecv; i++) { 564e953e47cSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 565e953e47cSHong Zhang } 566e953e47cSHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 567e953e47cSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 568e953e47cSHong Zhang 569e953e47cSHong Zhang /* add received column indices into ta to update Crmax */ 570e953e47cSHong Zhang for (k=0; k<nrecv; k++) {/* k-th received message */ 571e953e47cSHong Zhang Jptr = buf_rj[k]; 572e953e47cSHong Zhang for (j=0; j<len_r[k]; j++) { 573e953e47cSHong Zhang ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 574e953e47cSHong Zhang } 575e953e47cSHong Zhang } 576e953e47cSHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 577e953e47cSHong Zhang ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 578e953e47cSHong Zhang 579e953e47cSHong Zhang /* (4) send and recv coi */ 580e953e47cSHong Zhang /*-----------------------*/ 581e953e47cSHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 582e953e47cSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 583e953e47cSHong Zhang ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 584e953e47cSHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 585e953e47cSHong Zhang for (proc=0,k=0; proc<size; proc++) { 586e953e47cSHong Zhang if (!len_s[proc]) continue; 587e953e47cSHong Zhang /* form outgoing message for i-structure: 588e953e47cSHong Zhang buf_si[0]: nrows to be sent 589e953e47cSHong Zhang [1:nrows]: row index (global) 590e953e47cSHong Zhang [nrows+1:2*nrows+1]: i-structure index 591e953e47cSHong Zhang */ 592e953e47cSHong Zhang /*-------------------------------------------*/ 593e953e47cSHong Zhang nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 594e953e47cSHong Zhang buf_si_i = buf_si + nrows+1; 595e953e47cSHong Zhang buf_si[0] = nrows; 596e953e47cSHong Zhang buf_si_i[0] = 0; 597e953e47cSHong Zhang nrows = 0; 598e953e47cSHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 599e953e47cSHong Zhang nzi = coi[i+1] - coi[i]; 600e953e47cSHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 601e953e47cSHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 602e953e47cSHong Zhang nrows++; 603e953e47cSHong Zhang } 604e953e47cSHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 605e953e47cSHong Zhang k++; 606e953e47cSHong Zhang buf_si += len_si[proc]; 607e953e47cSHong Zhang } 608e953e47cSHong Zhang for (i=0; i<nrecv; i++) { 609e953e47cSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 610e953e47cSHong Zhang } 611e953e47cSHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 612e953e47cSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 613e953e47cSHong Zhang 614e953e47cSHong Zhang ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 615e953e47cSHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 616e953e47cSHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 617e953e47cSHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 618b4e8d1b6SHong Zhang 619e953e47cSHong Zhang /* (5) compute the local portion of Cmpi */ 620e953e47cSHong Zhang /* ------------------------------------------ */ 621e953e47cSHong Zhang /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 622e953e47cSHong Zhang ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 623e953e47cSHong Zhang current_space = free_space; 624e953e47cSHong Zhang 625e953e47cSHong Zhang ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 626e953e47cSHong Zhang for (k=0; k<nrecv; k++) { 627e953e47cSHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 628e953e47cSHong Zhang nrows = *buf_ri_k[k]; 629e953e47cSHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 630e953e47cSHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 631e953e47cSHong Zhang } 632e953e47cSHong Zhang 633e953e47cSHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 634cf742fccSHong Zhang ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 635e953e47cSHong Zhang for (i=0; i<pn; i++) { 636e953e47cSHong Zhang /* add C_loc into Cmpi */ 637e953e47cSHong Zhang nzi = c_loc->i[i+1] - c_loc->i[i]; 638e953e47cSHong Zhang Jptr = c_loc->j + c_loc->i[i]; 639cf742fccSHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 640e953e47cSHong Zhang 641e953e47cSHong Zhang /* add received col data into lnk */ 642e953e47cSHong Zhang for (k=0; k<nrecv; k++) { /* k-th received message */ 643e953e47cSHong Zhang if (i == *nextrow[k]) { /* i-th row */ 644e953e47cSHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 645e953e47cSHong Zhang Jptr = buf_rj[k] + *nextci[k]; 646cf742fccSHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 647e953e47cSHong Zhang nextrow[k]++; nextci[k]++; 648e953e47cSHong Zhang } 649e953e47cSHong Zhang } 650e953e47cSHong Zhang nzi = lnk[0]; 651e953e47cSHong Zhang 652e953e47cSHong Zhang /* copy data into free space, then initialize lnk */ 653cf742fccSHong Zhang ierr = PetscLLCondensedClean_Scalable(nzi,current_space->array,lnk);CHKERRQ(ierr); 654e953e47cSHong Zhang ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 655e953e47cSHong Zhang } 656e953e47cSHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 657cf742fccSHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 658e953e47cSHong Zhang ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 659e953e47cSHong Zhang 660e953e47cSHong Zhang /* local sizes and preallocation */ 661e953e47cSHong Zhang ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 662ac94c67aSTristan Konolige if (P->cmap->bs > 0) { 663ac94c67aSTristan Konolige ierr = PetscLayoutSetBlockSize(Cmpi->rmap,P->cmap->bs);CHKERRQ(ierr); 664ac94c67aSTristan Konolige ierr = PetscLayoutSetBlockSize(Cmpi->cmap,P->cmap->bs);CHKERRQ(ierr); 665ac94c67aSTristan Konolige } 666e953e47cSHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 667e953e47cSHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 668e953e47cSHong Zhang 669e953e47cSHong Zhang /* members in merge */ 670e953e47cSHong Zhang ierr = PetscFree(id_r);CHKERRQ(ierr); 671e953e47cSHong Zhang ierr = PetscFree(len_r);CHKERRQ(ierr); 672e953e47cSHong Zhang ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 673e953e47cSHong Zhang ierr = PetscFree(buf_ri);CHKERRQ(ierr); 674e953e47cSHong Zhang ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 675e953e47cSHong Zhang ierr = PetscFree(buf_rj);CHKERRQ(ierr); 676e953e47cSHong Zhang ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 677e953e47cSHong Zhang 678e953e47cSHong Zhang /* attach the supporting struct to Cmpi for reuse */ 679e953e47cSHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 6803cdca5ebSHong Zhang c->ap = ptap; 681e953e47cSHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 682e953e47cSHong Zhang ptap->destroy = Cmpi->ops->destroy; 683cf97cf3cSHong Zhang ptap->view = Cmpi->ops->view; 684e953e47cSHong Zhang 685e953e47cSHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 686e953e47cSHong Zhang Cmpi->assembled = PETSC_FALSE; 687a4ffb656SHong Zhang Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable; 688e953e47cSHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 689cf97cf3cSHong Zhang Cmpi->ops->view = MatView_MPIAIJ_PtAP; 6904624976aSHong Zhang Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 691e953e47cSHong Zhang *C = Cmpi; 692a3bb6f32SFande Kong 693a3bb6f32SFande Kong nout = 0; 694a3bb6f32SFande Kong ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_oth->i[ptap->C_oth->rmap->n],c_oth->j,&nout,c_oth->j);CHKERRQ(ierr); 695a3bb6f32SFande Kong if (c_oth->i[ptap->C_oth->rmap->n] != nout) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incorrect mapping %D != %D\n",c_oth->i[ptap->C_oth->rmap->n],nout); 696a3bb6f32SFande Kong ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_loc->i[ptap->C_loc->rmap->n],c_loc->j,&nout,c_loc->j);CHKERRQ(ierr); 697a3bb6f32SFande Kong if (c_loc->i[ptap->C_loc->rmap->n] != nout) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Incorrect mapping %D != %D\n",c_loc->i[ptap->C_loc->rmap->n],nout); 698a3bb6f32SFande Kong 699e953e47cSHong Zhang PetscFunctionReturn(0); 700e953e47cSHong Zhang } 701e953e47cSHong Zhang 7024a56b808SFande Kong PETSC_STATIC_INLINE PetscErrorCode MatPtAPSymbolicComputeOneRowOfAP_private(Mat A,Mat P,Mat P_oth,PetscInt i,PetscHSetI dht,PetscHSetI oht) 7034a56b808SFande Kong { 7044a56b808SFande Kong Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 7054a56b808SFande Kong Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_oth=(Mat_SeqAIJ*)P_oth->data,*pd=(Mat_SeqAIJ*)p->A->data,*po=(Mat_SeqAIJ*)p->B->data; 7064a56b808SFande Kong PetscInt *ai,nzi,j,*aj,row,col,*pi,*pj,pnz,nzpi,*p_othcols,k; 7074a56b808SFande Kong PetscInt pcstart,pcend,column; 7084a56b808SFande Kong PetscErrorCode ierr; 7094a56b808SFande Kong 7104a56b808SFande Kong PetscFunctionBegin; 7114a56b808SFande Kong 7124a56b808SFande Kong pcstart = P->cmap->rstart; 7134a56b808SFande Kong pcend = P->cmap->rend; 7144a56b808SFande Kong /* diagonal portion: Ad[i,:]*P */ 7154a56b808SFande Kong ai = ad->i; 7164a56b808SFande Kong nzi = ai[i+1] - ai[i]; 7174a56b808SFande Kong aj = ad->j + ai[i]; 7184a56b808SFande Kong for (j=0; j<nzi; j++) { 7194a56b808SFande Kong row = aj[j]; 7204a56b808SFande Kong nzpi = pd->i[row+1] - pd->i[row]; 7214a56b808SFande Kong pj = pd->j + pd->i[row]; 7224a56b808SFande Kong for (k=0; k<nzpi; k++) { 7234a56b808SFande Kong ierr = PetscHSetIAdd(dht,pj[k]+pcstart);CHKERRQ(ierr); 7244a56b808SFande Kong } 7254a56b808SFande Kong } 7264a56b808SFande Kong for (j=0; j<nzi; j++) { 7274a56b808SFande Kong row = aj[j]; 7284a56b808SFande Kong nzpi = po->i[row+1] - po->i[row]; 7294a56b808SFande Kong pj = po->j + po->i[row]; 7304a56b808SFande Kong for (k=0; k<nzpi; k++) { 7314a56b808SFande Kong ierr = PetscHSetIAdd(oht,p->garray[pj[k]]);CHKERRQ(ierr); 7324a56b808SFande Kong } 7334a56b808SFande Kong } 7344a56b808SFande Kong 7354a56b808SFande Kong /* off diagonal part: Ao[i, :]*P_oth */ 7364a56b808SFande Kong if (ao) { 7374a56b808SFande Kong ai = ao->i; 7384a56b808SFande Kong pi = p_oth->i; 7394a56b808SFande Kong nzi = ai[i+1] - ai[i]; 7404a56b808SFande Kong aj = ao->j + ai[i]; 7414a56b808SFande Kong for (j=0; j<nzi; j++) { 7424a56b808SFande Kong row = aj[j]; 7434a56b808SFande Kong pnz = pi[row+1] - pi[row]; 7444a56b808SFande Kong p_othcols = p_oth->j + pi[row]; 7454a56b808SFande Kong for (col=0; col<pnz; col++) { 7464a56b808SFande Kong column = p_othcols[col]; 7474a56b808SFande Kong if (column>=pcstart && column<pcend) { 7484a56b808SFande Kong ierr = PetscHSetIAdd(dht,column);CHKERRQ(ierr); 7494a56b808SFande Kong } else { 7504a56b808SFande Kong ierr = PetscHSetIAdd(oht,column);CHKERRQ(ierr); 7514a56b808SFande Kong } 7524a56b808SFande Kong } 7534a56b808SFande Kong } 7544a56b808SFande Kong } /* end if (ao) */ 7554a56b808SFande Kong PetscFunctionReturn(0); 7564a56b808SFande Kong } 7574a56b808SFande Kong 7584a56b808SFande Kong PETSC_STATIC_INLINE PetscErrorCode MatPtAPNumericComputeOneRowOfAP_private(Mat A,Mat P,Mat P_oth,PetscInt i,PetscHMapIV hmap) 7594a56b808SFande Kong { 7604a56b808SFande Kong Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 7614a56b808SFande Kong Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data,*p_oth=(Mat_SeqAIJ*)P_oth->data,*pd=(Mat_SeqAIJ*)p->A->data,*po=(Mat_SeqAIJ*)p->B->data; 7624a56b808SFande Kong PetscInt *ai,nzi,j,*aj,row,col,*pi,pnz,*p_othcols,pcstart,*pj,k,nzpi; 7634a56b808SFande Kong PetscScalar ra,*aa,*pa; 764*bb5ddf68SFande Kong PetscLogDouble flops=0.0; 7654a56b808SFande Kong PetscErrorCode ierr; 7664a56b808SFande Kong 7674a56b808SFande Kong PetscFunctionBegin; 7684a56b808SFande Kong pcstart = P->cmap->rstart; 7694a56b808SFande Kong /* diagonal portion: Ad[i,:]*P */ 7704a56b808SFande Kong ai = ad->i; 7714a56b808SFande Kong nzi = ai[i+1] - ai[i]; 7724a56b808SFande Kong aj = ad->j + ai[i]; 7734a56b808SFande Kong aa = ad->a + ai[i]; 7744a56b808SFande Kong 7754a56b808SFande Kong for (j=0; j<nzi; j++) { 7764a56b808SFande Kong ra = aa[j]; 7774a56b808SFande Kong row = aj[j]; 7784a56b808SFande Kong nzpi = pd->i[row+1] - pd->i[row]; 7794a56b808SFande Kong pj = pd->j + pd->i[row]; 7804a56b808SFande Kong pa = pd->a + pd->i[row]; 7814a56b808SFande Kong for (k=0; k<nzpi; k++) { 782694f88d4SFande Kong ierr = PetscHMapIVAddValue(hmap,pj[k]+pcstart,ra*pa[k]);CHKERRQ(ierr); 7834a56b808SFande Kong } 784*bb5ddf68SFande Kong flops += nzpi; 7854a56b808SFande Kong } 7864a56b808SFande Kong for (j=0; j<nzi; j++) { 7874a56b808SFande Kong ra = aa[j]; 7884a56b808SFande Kong row = aj[j]; 7894a56b808SFande Kong nzpi = po->i[row+1] - po->i[row]; 7904a56b808SFande Kong pj = po->j + po->i[row]; 7914a56b808SFande Kong pa = po->a + po->i[row]; 7924a56b808SFande Kong for (k=0; k<nzpi; k++) { 793694f88d4SFande Kong ierr = PetscHMapIVAddValue(hmap,p->garray[pj[k]],ra*pa[k]);CHKERRQ(ierr); 7944a56b808SFande Kong } 795*bb5ddf68SFande Kong flops += nzpi; 7964a56b808SFande Kong } 7974a56b808SFande Kong 7984a56b808SFande Kong 7994a56b808SFande Kong /* off diagonal part: Ao[i, :]*P_oth */ 8004a56b808SFande Kong if (ao) { 8014a56b808SFande Kong ai = ao->i; 8024a56b808SFande Kong pi = p_oth->i; 8034a56b808SFande Kong nzi = ai[i+1] - ai[i]; 8044a56b808SFande Kong aj = ao->j + ai[i]; 8054a56b808SFande Kong aa = ao->a + ai[i]; 8064a56b808SFande Kong for (j=0; j<nzi; j++) { 8074a56b808SFande Kong row = aj[j]; 8084a56b808SFande Kong ra = aa[j]; 8094a56b808SFande Kong pnz = pi[row+1] - pi[row]; 8104a56b808SFande Kong p_othcols = p_oth->j + pi[row]; 8114a56b808SFande Kong pa = p_oth->a + pi[row]; 8124a56b808SFande Kong for (col=0; col<pnz; col++) { 813694f88d4SFande Kong ierr = PetscHMapIVAddValue(hmap,p_othcols[col],ra*pa[col]);CHKERRQ(ierr); 8144a56b808SFande Kong } 815*bb5ddf68SFande Kong flops += pnz; 8164a56b808SFande Kong } 8174a56b808SFande Kong } /* end if (ao) */ 818*bb5ddf68SFande Kong 819*bb5ddf68SFande Kong ierr = PetscLogFlops(flops);CHKERRQ(ierr); 8204a56b808SFande Kong PetscFunctionReturn(0); 8214a56b808SFande Kong } 8224a56b808SFande Kong 8234a56b808SFande Kong PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce(Mat A,Mat P,Mat C) 8244a56b808SFande Kong { 8254a56b808SFande Kong PetscErrorCode ierr; 8264a56b808SFande Kong Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 8274a56b808SFande Kong Mat_SeqAIJ *cd,*co,*po,*pd; 8284a56b808SFande Kong Mat_APMPI *ptap = c->ap; 8294a56b808SFande Kong PetscHMapIV hmap; 8304a56b808SFande Kong PetscInt i,j,jj,kk,nzi,*c_rmtj,voff,*c_othj,pn,pon,pcstart,pcend,ccstart,ccend,row,am,*poj,*pdj,*apindices,cmaxr,*c_rmtc,*c_rmtjj,*dcc,*occ,loc; 8314a56b808SFande Kong PetscScalar *c_rmta,*c_otha,*poa,*pda,*apvalues,*apvaluestmp,*c_rmtaa; 832*bb5ddf68SFande Kong PetscLogDouble flops=0.0; 8334a56b808SFande Kong MPI_Comm comm; 8344a56b808SFande Kong 8354a56b808SFande Kong PetscFunctionBegin; 8364a56b808SFande Kong ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 837bb9bda99SHong Zhang if (!ptap->P_oth) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 8384a56b808SFande Kong 8394a56b808SFande Kong ierr = MatZeroEntries(C);CHKERRQ(ierr); 8404a56b808SFande Kong 8414a56b808SFande Kong /* Get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 8424a56b808SFande Kong /*-----------------------------------------------------*/ 8434a56b808SFande Kong if (ptap->reuse == MAT_REUSE_MATRIX) { 8444a56b808SFande Kong /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 8454a56b808SFande Kong ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 8464a56b808SFande Kong } 8474a56b808SFande Kong 8484a56b808SFande Kong po = (Mat_SeqAIJ*) p->B->data; 8494a56b808SFande Kong pd = (Mat_SeqAIJ*) p->A->data; 8504a56b808SFande Kong 8514a56b808SFande Kong ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr); 8524a56b808SFande Kong 8534a56b808SFande Kong ierr = PetscCalloc2(ptap->c_rmti[pon],&c_rmtj,ptap->c_rmti[pon],&c_rmta);CHKERRQ(ierr); 8544a56b808SFande Kong ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr); 8554a56b808SFande Kong cmaxr = 0; 8564a56b808SFande Kong for (i=0; i<pon; i++) { 8574a56b808SFande Kong cmaxr = PetscMax(cmaxr,ptap->c_rmti[i+1]-ptap->c_rmti[i]); 8584a56b808SFande Kong } 8594a56b808SFande Kong ierr = PetscCalloc4(cmaxr,&apindices,cmaxr,&apvalues,cmaxr,&apvaluestmp,pon,&c_rmtc);CHKERRQ(ierr); 8604a56b808SFande Kong ierr = PetscHMapIVCreate(&hmap);CHKERRQ(ierr); 8614a56b808SFande Kong ierr = PetscHMapIVResize(hmap,cmaxr);CHKERRQ(ierr); 8624a56b808SFande Kong for (i=0; i<am && pon; i++) { 8634a56b808SFande Kong ierr = PetscHMapIVClear(hmap);CHKERRQ(ierr); 8644a56b808SFande Kong nzi = po->i[i+1] - po->i[i]; 8654a56b808SFande Kong if (!nzi) continue; 8664a56b808SFande Kong ierr = MatPtAPNumericComputeOneRowOfAP_private(A,P,ptap->P_oth,i,hmap);CHKERRQ(ierr); 8674a56b808SFande Kong voff = 0; 8684a56b808SFande Kong ierr = PetscHMapIVGetPairs(hmap,&voff,apindices,apvalues);CHKERRQ(ierr); 8694a56b808SFande Kong if (!voff) continue; 8704a56b808SFande Kong /*ierr = PetscMemzero(c_rmtc,sizeof(PetscInt)*pon);CHKERRQ(ierr);*/ 8714a56b808SFande Kong 8724a56b808SFande Kong /* Form C(ii, :) */ 8734a56b808SFande Kong poj = po->j + po->i[i]; 8744a56b808SFande Kong poa = po->a + po->i[i]; 8754a56b808SFande Kong for (j=0; j<nzi; j++) { 8764a56b808SFande Kong c_rmtjj = c_rmtj + ptap->c_rmti[poj[j]]; 8774a56b808SFande Kong c_rmtaa = c_rmta + ptap->c_rmti[poj[j]]; 8784a56b808SFande Kong for (jj=0; jj<voff; jj++) { 8794a56b808SFande Kong apvaluestmp[jj] = apvalues[jj]*poa[j]; 8804a56b808SFande Kong /*If the row is empty */ 8814a56b808SFande Kong if (!c_rmtc[poj[j]]) { 8824a56b808SFande Kong c_rmtjj[jj] = apindices[jj]; 8834a56b808SFande Kong c_rmtaa[jj] = apvaluestmp[jj]; 8844a56b808SFande Kong c_rmtc[poj[j]]++; 8854a56b808SFande Kong } else { 8864a56b808SFande Kong ierr = PetscFindInt(apindices[jj],c_rmtc[poj[j]],c_rmtjj,&loc);CHKERRQ(ierr); 8874a56b808SFande Kong if (loc>=0){ /* hit */ 8884a56b808SFande Kong c_rmtaa[loc] += apvaluestmp[jj]; 8894a56b808SFande Kong } else { /* new element */ 8904a56b808SFande Kong loc = -(loc+1); 8914a56b808SFande Kong /* Move data backward */ 8924a56b808SFande Kong for (kk=c_rmtc[poj[j]]; kk>loc; kk--) { 8934a56b808SFande Kong c_rmtjj[kk] = c_rmtjj[kk-1]; 8944a56b808SFande Kong c_rmtaa[kk] = c_rmtaa[kk-1]; 8954a56b808SFande Kong }/* End kk */ 8964a56b808SFande Kong c_rmtjj[loc] = apindices[jj]; 8974a56b808SFande Kong c_rmtaa[loc] = apvaluestmp[jj]; 8984a56b808SFande Kong c_rmtc[poj[j]]++; 8994a56b808SFande Kong } 9004a56b808SFande Kong } 901*bb5ddf68SFande Kong flops += voff; 9024a56b808SFande Kong } /* End jj */ 9034a56b808SFande Kong } /* End j */ 9044a56b808SFande Kong } /* End i */ 9054a56b808SFande Kong 9064a56b808SFande Kong ierr = PetscFree4(apindices,apvalues,apvaluestmp,c_rmtc);CHKERRQ(ierr); 9074a56b808SFande Kong 9084a56b808SFande Kong ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr); 9094a56b808SFande Kong ierr = PetscCalloc2(ptap->c_othi[pn],&c_othj,ptap->c_othi[pn],&c_otha);CHKERRQ(ierr); 9104a56b808SFande Kong 9114a56b808SFande Kong ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr); 9124a56b808SFande Kong ierr = PetscSFReduceBegin(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPIU_REPLACE);CHKERRQ(ierr); 9134a56b808SFande Kong ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr); 9144a56b808SFande Kong cd = (Mat_SeqAIJ*)(c->A)->data; 9154a56b808SFande Kong co = (Mat_SeqAIJ*)(c->B)->data; 9164a56b808SFande Kong 9174a56b808SFande Kong cmaxr = 0; 9184a56b808SFande Kong for (i=0; i<pn; i++) { 9194a56b808SFande Kong cmaxr = PetscMax(cmaxr,(cd->i[i+1]-cd->i[i])+(co->i[i+1]-co->i[i])); 9204a56b808SFande Kong } 9214a56b808SFande Kong ierr = PetscCalloc5(cmaxr,&apindices,cmaxr,&apvalues,cmaxr,&apvaluestmp,pn,&dcc,pn,&occ);CHKERRQ(ierr); 9224a56b808SFande Kong ierr = PetscHMapIVCreate(&hmap);CHKERRQ(ierr); 9234a56b808SFande Kong ierr = PetscHMapIVResize(hmap,cmaxr);CHKERRQ(ierr); 9244a56b808SFande Kong for (i=0; i<am && pn; i++) { 9254a56b808SFande Kong ierr = PetscHMapIVClear(hmap);CHKERRQ(ierr); 9264a56b808SFande Kong nzi = pd->i[i+1] - pd->i[i]; 9274a56b808SFande Kong if (!nzi) continue; 9284a56b808SFande Kong ierr = MatPtAPNumericComputeOneRowOfAP_private(A,P,ptap->P_oth,i,hmap);CHKERRQ(ierr); 9294a56b808SFande Kong voff = 0; 9304a56b808SFande Kong ierr = PetscHMapIVGetPairs(hmap,&voff,apindices,apvalues);CHKERRQ(ierr); 9314a56b808SFande Kong if (!voff) continue; 9324a56b808SFande Kong /* Form C(ii, :) */ 9334a56b808SFande Kong pdj = pd->j + pd->i[i]; 9344a56b808SFande Kong pda = pd->a + pd->i[i]; 9354a56b808SFande Kong for (j=0; j<nzi; j++) { 9364a56b808SFande Kong row = pcstart + pdj[j]; 9374a56b808SFande Kong for (jj=0; jj<voff; jj++) { 9384a56b808SFande Kong apvaluestmp[jj] = apvalues[jj]*pda[j]; 9394a56b808SFande Kong } 9404a56b808SFande Kong ierr = MatSetValues(C,1,&row,voff,apindices,apvaluestmp,ADD_VALUES);CHKERRQ(ierr); 941*bb5ddf68SFande Kong flops += voff; 9424a56b808SFande Kong } 9434a56b808SFande Kong } 9444a56b808SFande Kong 9454a56b808SFande Kong ierr = MatGetOwnershipRangeColumn(C,&ccstart,&ccend);CHKERRQ(ierr); 9464a56b808SFande Kong ierr = PetscFree5(apindices,apvalues,apvaluestmp,dcc,occ);CHKERRQ(ierr); 9474a56b808SFande Kong ierr = PetscHMapIVDestroy(&hmap);CHKERRQ(ierr); 9484a56b808SFande Kong ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPI_SUM);CHKERRQ(ierr); 9494a56b808SFande Kong ierr = PetscSFReduceEnd(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPIU_REPLACE);CHKERRQ(ierr); 9504a56b808SFande Kong ierr = PetscFree2(c_rmtj,c_rmta);CHKERRQ(ierr); 9514a56b808SFande Kong 9524a56b808SFande Kong /* Add contributions from remote */ 9534a56b808SFande Kong for (i = 0; i < pn; i++) { 9544a56b808SFande Kong row = i + pcstart; 9554a56b808SFande Kong ierr = MatSetValues(C,1,&row,ptap->c_othi[i+1]-ptap->c_othi[i],c_othj+ptap->c_othi[i],c_otha+ptap->c_othi[i],ADD_VALUES);CHKERRQ(ierr); 9564a56b808SFande Kong } 9574a56b808SFande Kong ierr = PetscFree2(c_othj,c_otha);CHKERRQ(ierr); 9584a56b808SFande Kong 9594a56b808SFande Kong ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9604a56b808SFande Kong ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 961*bb5ddf68SFande Kong ierr = PetscLogFlops(flops);CHKERRQ(ierr); 9624a56b808SFande Kong 9634a56b808SFande Kong ptap->reuse = MAT_REUSE_MATRIX; 9644a56b808SFande Kong 9654a56b808SFande Kong /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */ 9664a56b808SFande Kong if (ptap->freestruct) { 9674a56b808SFande Kong ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 9684a56b808SFande Kong } 9694a56b808SFande Kong PetscFunctionReturn(0); 9704a56b808SFande Kong } 9714a56b808SFande Kong 9724a56b808SFande Kong PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged(Mat A,Mat P,Mat C) 9734a56b808SFande Kong { 9744a56b808SFande Kong PetscErrorCode ierr; 9754a56b808SFande Kong Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 9764a56b808SFande Kong Mat_SeqAIJ *cd,*co,*po,*pd; 9774a56b808SFande Kong Mat_APMPI *ptap = c->ap; 9784a56b808SFande Kong PetscHMapIV hmap; 9794a56b808SFande Kong PetscInt i,j,jj,kk,nzi,dnzi,*c_rmtj,voff,*c_othj,pn,pon,pcstart,pcend,row,am,*poj,*pdj,*apindices,cmaxr,*c_rmtc,*c_rmtjj,loc; 9804a56b808SFande Kong PetscScalar *c_rmta,*c_otha,*poa,*pda,*apvalues,*apvaluestmp,*c_rmtaa; 981*bb5ddf68SFande Kong PetscLogDouble flops=0.0; 9824a56b808SFande Kong MPI_Comm comm; 9834a56b808SFande Kong 9844a56b808SFande Kong PetscFunctionBegin; 9854a56b808SFande Kong ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 986bb9bda99SHong Zhang if (!ptap->P_oth) SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 9874a56b808SFande Kong 9884a56b808SFande Kong ierr = MatZeroEntries(C);CHKERRQ(ierr); 9894a56b808SFande Kong 9904a56b808SFande Kong /* Get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 9914a56b808SFande Kong /*-----------------------------------------------------*/ 9924a56b808SFande Kong if (ptap->reuse == MAT_REUSE_MATRIX) { 9934a56b808SFande Kong /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 9944a56b808SFande Kong ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 9954a56b808SFande Kong } 9964a56b808SFande Kong 9974a56b808SFande Kong po = (Mat_SeqAIJ*) p->B->data; 9984a56b808SFande Kong pd = (Mat_SeqAIJ*) p->A->data; 9994a56b808SFande Kong 10004a56b808SFande Kong ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr); 10014a56b808SFande Kong ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr); 10024a56b808SFande Kong 10034a56b808SFande Kong ierr = PetscCalloc2(ptap->c_rmti[pon],&c_rmtj,ptap->c_rmti[pon],&c_rmta);CHKERRQ(ierr); 10044a56b808SFande Kong ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr); 10054a56b808SFande Kong ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr); 10064a56b808SFande Kong cmaxr = 0; 10074a56b808SFande Kong for (i=0; i<pon; i++) { 10084a56b808SFande Kong cmaxr = PetscMax(cmaxr,ptap->c_rmti[i+1]-ptap->c_rmti[i]); 10094a56b808SFande Kong } 10104a56b808SFande Kong cd = (Mat_SeqAIJ*)(c->A)->data; 10114a56b808SFande Kong co = (Mat_SeqAIJ*)(c->B)->data; 10124a56b808SFande Kong for (i=0; i<pn; i++) { 10134a56b808SFande Kong cmaxr = PetscMax(cmaxr,(cd->i[i+1]-cd->i[i])+(co->i[i+1]-co->i[i])); 10144a56b808SFande Kong } 10154a56b808SFande Kong ierr = PetscCalloc4(cmaxr,&apindices,cmaxr,&apvalues,cmaxr,&apvaluestmp,pon,&c_rmtc);CHKERRQ(ierr); 10164a56b808SFande Kong ierr = PetscHMapIVCreate(&hmap);CHKERRQ(ierr); 10174a56b808SFande Kong ierr = PetscHMapIVResize(hmap,cmaxr);CHKERRQ(ierr); 10184a56b808SFande Kong for (i=0; i<am && (pon || pn); i++) { 10194a56b808SFande Kong ierr = PetscHMapIVClear(hmap);CHKERRQ(ierr); 10204a56b808SFande Kong nzi = po->i[i+1] - po->i[i]; 10214a56b808SFande Kong dnzi = pd->i[i+1] - pd->i[i]; 10224a56b808SFande Kong if (!nzi && !dnzi) continue; 10234a56b808SFande Kong ierr = MatPtAPNumericComputeOneRowOfAP_private(A,P,ptap->P_oth,i,hmap);CHKERRQ(ierr); 10244a56b808SFande Kong voff = 0; 10254a56b808SFande Kong ierr = PetscHMapIVGetPairs(hmap,&voff,apindices,apvalues);CHKERRQ(ierr); 10264a56b808SFande Kong if (!voff) continue; 10274a56b808SFande Kong 10284a56b808SFande Kong /* Form remote C(ii, :) */ 10294a56b808SFande Kong poj = po->j + po->i[i]; 10304a56b808SFande Kong poa = po->a + po->i[i]; 10314a56b808SFande Kong for (j=0; j<nzi; j++) { 10324a56b808SFande Kong c_rmtjj = c_rmtj + ptap->c_rmti[poj[j]]; 10334a56b808SFande Kong c_rmtaa = c_rmta + ptap->c_rmti[poj[j]]; 10344a56b808SFande Kong for (jj=0; jj<voff; jj++) { 10354a56b808SFande Kong apvaluestmp[jj] = apvalues[jj]*poa[j]; 10364a56b808SFande Kong /*If the row is empty */ 10374a56b808SFande Kong if (!c_rmtc[poj[j]]) { 10384a56b808SFande Kong c_rmtjj[jj] = apindices[jj]; 10394a56b808SFande Kong c_rmtaa[jj] = apvaluestmp[jj]; 10404a56b808SFande Kong c_rmtc[poj[j]]++; 10414a56b808SFande Kong } else { 10424a56b808SFande Kong ierr = PetscFindInt(apindices[jj],c_rmtc[poj[j]],c_rmtjj,&loc);CHKERRQ(ierr); 10434a56b808SFande Kong if (loc>=0){ /* hit */ 10444a56b808SFande Kong c_rmtaa[loc] += apvaluestmp[jj]; 10454a56b808SFande Kong } else { /* new element */ 10464a56b808SFande Kong loc = -(loc+1); 10474a56b808SFande Kong /* Move data backward */ 10484a56b808SFande Kong for (kk=c_rmtc[poj[j]]; kk>loc; kk--) { 10494a56b808SFande Kong c_rmtjj[kk] = c_rmtjj[kk-1]; 10504a56b808SFande Kong c_rmtaa[kk] = c_rmtaa[kk-1]; 10514a56b808SFande Kong }/* End kk */ 10524a56b808SFande Kong c_rmtjj[loc] = apindices[jj]; 10534a56b808SFande Kong c_rmtaa[loc] = apvaluestmp[jj]; 10544a56b808SFande Kong c_rmtc[poj[j]]++; 10554a56b808SFande Kong } 10564a56b808SFande Kong } 10574a56b808SFande Kong } /* End jj */ 1058*bb5ddf68SFande Kong flops += voff; 10594a56b808SFande Kong } /* End j */ 10604a56b808SFande Kong 10614a56b808SFande Kong /* Form local C(ii, :) */ 10624a56b808SFande Kong pdj = pd->j + pd->i[i]; 10634a56b808SFande Kong pda = pd->a + pd->i[i]; 10644a56b808SFande Kong for (j=0; j<dnzi; j++) { 10654a56b808SFande Kong row = pcstart + pdj[j]; 10664a56b808SFande Kong for (jj=0; jj<voff; jj++) { 10674a56b808SFande Kong apvaluestmp[jj] = apvalues[jj]*pda[j]; 10684a56b808SFande Kong }/* End kk */ 10694a56b808SFande Kong ierr = MatSetValues(C,1,&row,voff,apindices,apvaluestmp,ADD_VALUES);CHKERRQ(ierr); 1070*bb5ddf68SFande Kong flops += voff; 10714a56b808SFande Kong }/* End j */ 10724a56b808SFande Kong } /* End i */ 10734a56b808SFande Kong 10744a56b808SFande Kong ierr = PetscFree4(apindices,apvalues,apvaluestmp,c_rmtc);CHKERRQ(ierr); 10754a56b808SFande Kong ierr = PetscHMapIVDestroy(&hmap);CHKERRQ(ierr); 10764a56b808SFande Kong ierr = PetscCalloc2(ptap->c_othi[pn],&c_othj,ptap->c_othi[pn],&c_otha);CHKERRQ(ierr); 10774a56b808SFande Kong 10784a56b808SFande Kong ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr); 10794a56b808SFande Kong ierr = PetscSFReduceBegin(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPIU_REPLACE);CHKERRQ(ierr); 10804a56b808SFande Kong ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPI_SUM);CHKERRQ(ierr); 10814a56b808SFande Kong ierr = PetscSFReduceEnd(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPIU_REPLACE);CHKERRQ(ierr); 10824a56b808SFande Kong ierr = PetscFree2(c_rmtj,c_rmta);CHKERRQ(ierr); 10834a56b808SFande Kong 10844a56b808SFande Kong /* Add contributions from remote */ 10854a56b808SFande Kong for (i = 0; i < pn; i++) { 10864a56b808SFande Kong row = i + pcstart; 10874a56b808SFande Kong ierr = MatSetValues(C,1,&row,ptap->c_othi[i+1]-ptap->c_othi[i],c_othj+ptap->c_othi[i],c_otha+ptap->c_othi[i],ADD_VALUES);CHKERRQ(ierr); 10884a56b808SFande Kong } 10894a56b808SFande Kong ierr = PetscFree2(c_othj,c_otha);CHKERRQ(ierr); 10904a56b808SFande Kong 10914a56b808SFande Kong ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10924a56b808SFande Kong ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1093*bb5ddf68SFande Kong ierr = PetscLogFlops(flops);CHKERRQ(ierr); 10944a56b808SFande Kong 10954a56b808SFande Kong ptap->reuse = MAT_REUSE_MATRIX; 10964a56b808SFande Kong 10974a56b808SFande Kong /* supporting struct ptap consumes almost same amount of memory as C=PtAP, release it if C will not be updated by A and P */ 10984a56b808SFande Kong if (ptap->freestruct) { 10994a56b808SFande Kong ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 11004a56b808SFande Kong } 11014a56b808SFande Kong PetscFunctionReturn(0); 11024a56b808SFande Kong } 11034a56b808SFande Kong 11044a56b808SFande Kong PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce(Mat A,Mat P,PetscReal fill,Mat *C) 11054a56b808SFande Kong { 11064a56b808SFande Kong Mat_APMPI *ptap; 11074a56b808SFande Kong Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c; 11084a56b808SFande Kong MPI_Comm comm; 11094a56b808SFande Kong Mat Cmpi; 11104a56b808SFande Kong Mat_SeqAIJ *pd,*po; 11114a56b808SFande Kong MatType mtype; 11124a56b808SFande Kong PetscSF sf; 11134a56b808SFande Kong PetscSFNode *iremote; 11144a56b808SFande Kong PetscInt rootspacesize,*rootspace,*rootspaceoffsets,nleaves; 11154a56b808SFande Kong const PetscInt *rootdegrees; 11164a56b808SFande Kong PetscHSetI ht,oht,*hta,*hto; 11174a56b808SFande Kong PetscInt pn,pon,*c_rmtc,i,j,nzi,htsize,htosize,*c_rmtj,off,*c_othj,rcvncols,sendncols,*c_rmtoffsets; 11184a56b808SFande Kong PetscInt owner,lidx,*rdj,col,pcstart,pcend,*dnz,*onz,am,arstart,arend,*poj,*pdj; 11194a56b808SFande Kong PetscInt nalg=2,alg=0; 11204a56b808SFande Kong PetscBool flg; 11214a56b808SFande Kong const char *algTypes[2] = {"overlapping","merged"}; 11224a56b808SFande Kong PetscErrorCode ierr; 11234a56b808SFande Kong 11244a56b808SFande Kong PetscFunctionBegin; 11254a56b808SFande Kong ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 11264a56b808SFande Kong 11274a56b808SFande Kong /* Create symbolic parallel matrix Cmpi */ 11284a56b808SFande Kong ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr); 11294a56b808SFande Kong ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 11304a56b808SFande Kong ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 11314a56b808SFande Kong ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 11324a56b808SFande Kong ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 11334a56b808SFande Kong 11344a56b808SFande Kong ierr = PetscNew(&ptap);CHKERRQ(ierr); 11354a56b808SFande Kong ptap->reuse = MAT_INITIAL_MATRIX; 11364a56b808SFande Kong ptap->algType = 2; 11374a56b808SFande Kong 11384a56b808SFande Kong /* Get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 11394a56b808SFande Kong ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 11404a56b808SFande Kong 11414a56b808SFande Kong po = (Mat_SeqAIJ*)p->B->data; 11424a56b808SFande Kong pd = (Mat_SeqAIJ*)p->A->data; 11434a56b808SFande Kong 11444a56b808SFande Kong /* This equals to the number of offdiag columns in P */ 11454a56b808SFande Kong ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr); 11464a56b808SFande Kong /* offsets */ 11474a56b808SFande Kong ierr = PetscMalloc1(pon+1,&ptap->c_rmti);CHKERRQ(ierr); 11484a56b808SFande Kong /* The number of columns we will send to remote ranks */ 11494a56b808SFande Kong ierr = PetscMalloc1(pon,&c_rmtc);CHKERRQ(ierr); 11504a56b808SFande Kong ierr = PetscMalloc1(pon,&hta);CHKERRQ(ierr); 11514a56b808SFande Kong for (i=0; i<pon; i++) { 11524a56b808SFande Kong ierr = PetscHSetICreate(&hta[i]);CHKERRQ(ierr); 11534a56b808SFande Kong } 11544a56b808SFande Kong ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr); 11554a56b808SFande Kong ierr = MatGetOwnershipRange(A,&arstart,&arend);CHKERRQ(ierr); 11564a56b808SFande Kong /* Create hash table to merge all columns for C(i, :) */ 11574a56b808SFande Kong ierr = PetscHSetICreate(&ht);CHKERRQ(ierr); 11584a56b808SFande Kong 11594a56b808SFande Kong ptap->c_rmti[0] = 0; 11604a56b808SFande Kong /* 2) Pass 1: calculate the size for C_rmt (a matrix need to be sent to other processors) */ 11614a56b808SFande Kong for (i=0; i<am && pon; i++) { 11624a56b808SFande Kong /* Form one row of AP */ 11634a56b808SFande Kong ierr = PetscHSetIClear(ht);CHKERRQ(ierr); 11644a56b808SFande Kong /* If the off diag is empty, we should not do any calculation */ 11654a56b808SFande Kong nzi = po->i[i+1] - po->i[i]; 11664a56b808SFande Kong if (!nzi) continue; 11674a56b808SFande Kong 11684a56b808SFande Kong ierr = MatPtAPSymbolicComputeOneRowOfAP_private(A,P,ptap->P_oth,i,ht,ht);CHKERRQ(ierr); 11694a56b808SFande Kong ierr = PetscHSetIGetSize(ht,&htsize);CHKERRQ(ierr); 11704a56b808SFande Kong /* If AP is empty, just continue */ 11714a56b808SFande Kong if (!htsize) continue; 11724a56b808SFande Kong /* Form C(ii, :) */ 11734a56b808SFande Kong poj = po->j + po->i[i]; 11744a56b808SFande Kong for (j=0; j<nzi; j++) { 1175694f88d4SFande Kong ierr = PetscHSetIUpdate(hta[poj[j]],ht);CHKERRQ(ierr); 11764a56b808SFande Kong } 11774a56b808SFande Kong } 11784a56b808SFande Kong 11794a56b808SFande Kong for (i=0; i<pon; i++) { 11804a56b808SFande Kong ierr = PetscHSetIGetSize(hta[i],&htsize);CHKERRQ(ierr); 11814a56b808SFande Kong ptap->c_rmti[i+1] = ptap->c_rmti[i] + htsize; 11824a56b808SFande Kong c_rmtc[i] = htsize; 11834a56b808SFande Kong } 11844a56b808SFande Kong 11854a56b808SFande Kong ierr = PetscMalloc1(ptap->c_rmti[pon],&c_rmtj);CHKERRQ(ierr); 11864a56b808SFande Kong 11874a56b808SFande Kong for (i=0; i<pon; i++) { 11884a56b808SFande Kong off = 0; 11894a56b808SFande Kong ierr = PetscHSetIGetElems(hta[i],&off,c_rmtj+ptap->c_rmti[i]);CHKERRQ(ierr); 11904a56b808SFande Kong ierr = PetscHSetIDestroy(&hta[i]);CHKERRQ(ierr); 11914a56b808SFande Kong } 11924a56b808SFande Kong ierr = PetscFree(hta);CHKERRQ(ierr); 11934a56b808SFande Kong 11944a56b808SFande Kong ierr = PetscCalloc1(pon,&iremote);CHKERRQ(ierr); 11954a56b808SFande Kong for (i=0; i<pon; i++) { 1196ef7a94f2SFande Kong owner = 0; lidx = 0; 11974a56b808SFande Kong ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[i],&owner,&lidx);CHKERRQ(ierr); 11984a56b808SFande Kong iremote[i].index = lidx; 11994a56b808SFande Kong iremote[i].rank = owner; 12004a56b808SFande Kong } 12014a56b808SFande Kong 12024a56b808SFande Kong ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr); 12034a56b808SFande Kong ierr = PetscSFSetGraph(sf,pn,pon,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 12044a56b808SFande Kong /* Reorder ranks properly so that the data handled by gather and scatter have the same order */ 12054a56b808SFande Kong ierr = PetscSFSetRankOrder(sf,PETSC_TRUE);CHKERRQ(ierr); 12064a56b808SFande Kong ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 12074a56b808SFande Kong ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 12084a56b808SFande Kong /* How many neighbors have contributions to my rows? */ 12094a56b808SFande Kong ierr = PetscSFComputeDegreeBegin(sf,&rootdegrees);CHKERRQ(ierr); 12104a56b808SFande Kong ierr = PetscSFComputeDegreeEnd(sf,&rootdegrees);CHKERRQ(ierr); 12114a56b808SFande Kong rootspacesize = 0; 12124a56b808SFande Kong for (i = 0; i < pn; i++) { 12134a56b808SFande Kong rootspacesize += rootdegrees[i]; 12144a56b808SFande Kong } 12154a56b808SFande Kong ierr = PetscCalloc1(rootspacesize,&rootspace);CHKERRQ(ierr); 12164a56b808SFande Kong ierr = PetscCalloc1(rootspacesize+1,&rootspaceoffsets);CHKERRQ(ierr); 12174a56b808SFande Kong /* Get information from leaves 12184a56b808SFande Kong * Number of columns other people contribute to my rows 12194a56b808SFande Kong * */ 12204a56b808SFande Kong ierr = PetscSFGatherBegin(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr); 12214a56b808SFande Kong ierr = PetscSFGatherEnd(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr); 12224a56b808SFande Kong ierr = PetscFree(c_rmtc);CHKERRQ(ierr); 12234a56b808SFande Kong ierr = PetscCalloc1(pn+1,&ptap->c_othi);CHKERRQ(ierr); 12244a56b808SFande Kong /* The number of columns is received for each row */ 12254a56b808SFande Kong ptap->c_othi[0] = 0; 12264a56b808SFande Kong rootspacesize = 0; 12274a56b808SFande Kong rootspaceoffsets[0] = 0; 12284a56b808SFande Kong for (i = 0; i < pn; i++) { 12294a56b808SFande Kong rcvncols = 0; 12304a56b808SFande Kong for (j = 0; j<rootdegrees[i]; j++) { 12314a56b808SFande Kong rcvncols += rootspace[rootspacesize]; 12324a56b808SFande Kong rootspaceoffsets[rootspacesize+1] = rootspaceoffsets[rootspacesize] + rootspace[rootspacesize]; 12334a56b808SFande Kong rootspacesize++; 12344a56b808SFande Kong } 12354a56b808SFande Kong ptap->c_othi[i+1] = ptap->c_othi[i] + rcvncols; 12364a56b808SFande Kong } 12374a56b808SFande Kong ierr = PetscFree(rootspace);CHKERRQ(ierr); 12384a56b808SFande Kong 12394a56b808SFande Kong ierr = PetscCalloc1(pon,&c_rmtoffsets);CHKERRQ(ierr); 12404a56b808SFande Kong ierr = PetscSFScatterBegin(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr); 12414a56b808SFande Kong ierr = PetscSFScatterEnd(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr); 12424a56b808SFande Kong ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 12434a56b808SFande Kong ierr = PetscFree(rootspaceoffsets);CHKERRQ(ierr); 12444a56b808SFande Kong 12454a56b808SFande Kong ierr = PetscCalloc1(ptap->c_rmti[pon],&iremote);CHKERRQ(ierr); 12464a56b808SFande Kong nleaves = 0; 12474a56b808SFande Kong for (i = 0; i<pon; i++) { 1248ef7a94f2SFande Kong owner = 0; lidx = 0; 12494a56b808SFande Kong ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[i],&owner,&lidx);CHKERRQ(ierr); 12504a56b808SFande Kong sendncols = ptap->c_rmti[i+1] - ptap->c_rmti[i]; 12514a56b808SFande Kong for (j=0; j<sendncols; j++) { 12524a56b808SFande Kong iremote[nleaves].rank = owner; 12534a56b808SFande Kong iremote[nleaves++].index = c_rmtoffsets[i] + j; 12544a56b808SFande Kong } 12554a56b808SFande Kong } 12564a56b808SFande Kong ierr = PetscFree(c_rmtoffsets);CHKERRQ(ierr); 12574a56b808SFande Kong ierr = PetscCalloc1(ptap->c_othi[pn],&c_othj);CHKERRQ(ierr); 12584a56b808SFande Kong 12594a56b808SFande Kong ierr = PetscSFCreate(comm,&ptap->sf);CHKERRQ(ierr); 12604a56b808SFande Kong ierr = PetscSFSetGraph(ptap->sf,ptap->c_othi[pn],nleaves,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 12614a56b808SFande Kong ierr = PetscSFSetFromOptions(ptap->sf);CHKERRQ(ierr); 12624a56b808SFande Kong /* One to one map */ 12634a56b808SFande Kong ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr); 12644a56b808SFande Kong 12654a56b808SFande Kong ierr = PetscCalloc2(pn,&dnz,pn,&onz);CHKERRQ(ierr); 12664a56b808SFande Kong ierr = PetscHSetICreate(&oht);CHKERRQ(ierr); 12674a56b808SFande Kong ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr); 12684a56b808SFande Kong ierr = PetscMalloc2(pn,&hta,pn,&hto);CHKERRQ(ierr); 12694a56b808SFande Kong for (i=0; i<pn; i++) { 12704a56b808SFande Kong ierr = PetscHSetICreate(&hta[i]);CHKERRQ(ierr); 12714a56b808SFande Kong ierr = PetscHSetICreate(&hto[i]);CHKERRQ(ierr); 12724a56b808SFande Kong } 12734a56b808SFande Kong /* Work on local part */ 12744a56b808SFande Kong /* 4) Pass 1: Estimate memory for C_loc */ 12754a56b808SFande Kong for (i=0; i<am && pn; i++) { 12764a56b808SFande Kong ierr = PetscHSetIClear(ht);CHKERRQ(ierr); 12774a56b808SFande Kong ierr = PetscHSetIClear(oht);CHKERRQ(ierr); 12784a56b808SFande Kong nzi = pd->i[i+1] - pd->i[i]; 12794a56b808SFande Kong if (!nzi) continue; 12804a56b808SFande Kong 12814a56b808SFande Kong ierr = MatPtAPSymbolicComputeOneRowOfAP_private(A,P,ptap->P_oth,i,ht,oht);CHKERRQ(ierr); 12824a56b808SFande Kong ierr = PetscHSetIGetSize(ht,&htsize);CHKERRQ(ierr); 12834a56b808SFande Kong ierr = PetscHSetIGetSize(oht,&htosize);CHKERRQ(ierr); 12844a56b808SFande Kong if (!(htsize+htosize)) continue; 12854a56b808SFande Kong /* Form C(ii, :) */ 12864a56b808SFande Kong pdj = pd->j + pd->i[i]; 12874a56b808SFande Kong for (j=0; j<nzi; j++) { 1288694f88d4SFande Kong ierr = PetscHSetIUpdate(hta[pdj[j]],ht);CHKERRQ(ierr); 1289694f88d4SFande Kong ierr = PetscHSetIUpdate(hto[pdj[j]],oht);CHKERRQ(ierr); 12904a56b808SFande Kong } 12914a56b808SFande Kong } 12924a56b808SFande Kong 12934a56b808SFande Kong ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr); 12944a56b808SFande Kong ierr = PetscHSetIDestroy(&oht);CHKERRQ(ierr); 12954a56b808SFande Kong 12964a56b808SFande Kong /* Get remote data */ 12974a56b808SFande Kong ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr); 12984a56b808SFande Kong ierr = PetscFree(c_rmtj);CHKERRQ(ierr); 12994a56b808SFande Kong 13004a56b808SFande Kong for (i = 0; i < pn; i++) { 13014a56b808SFande Kong nzi = ptap->c_othi[i+1] - ptap->c_othi[i]; 13024a56b808SFande Kong rdj = c_othj + ptap->c_othi[i]; 13034a56b808SFande Kong for (j = 0; j < nzi; j++) { 13044a56b808SFande Kong col = rdj[j]; 13054a56b808SFande Kong /* diag part */ 13064a56b808SFande Kong if (col>=pcstart && col<pcend) { 13074a56b808SFande Kong ierr = PetscHSetIAdd(hta[i],col);CHKERRQ(ierr); 13084a56b808SFande Kong } else { /* off diag */ 13094a56b808SFande Kong ierr = PetscHSetIAdd(hto[i],col);CHKERRQ(ierr); 13104a56b808SFande Kong } 13114a56b808SFande Kong } 13124a56b808SFande Kong ierr = PetscHSetIGetSize(hta[i],&htsize);CHKERRQ(ierr); 13134a56b808SFande Kong dnz[i] = htsize; 13144a56b808SFande Kong ierr = PetscHSetIDestroy(&hta[i]);CHKERRQ(ierr); 13154a56b808SFande Kong ierr = PetscHSetIGetSize(hto[i],&htsize);CHKERRQ(ierr); 13164a56b808SFande Kong onz[i] = htsize; 13174a56b808SFande Kong ierr = PetscHSetIDestroy(&hto[i]);CHKERRQ(ierr); 13184a56b808SFande Kong } 13194a56b808SFande Kong 13204a56b808SFande Kong ierr = PetscFree2(hta,hto);CHKERRQ(ierr); 13214a56b808SFande Kong ierr = PetscFree(c_othj);CHKERRQ(ierr); 13224a56b808SFande Kong 13234a56b808SFande Kong /* local sizes and preallocation */ 13244a56b808SFande Kong ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 13254a56b808SFande Kong ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 13264a56b808SFande Kong ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 13274a56b808SFande Kong ierr = PetscFree2(dnz,onz);CHKERRQ(ierr); 13284a56b808SFande Kong 13294a56b808SFande Kong /* attach the supporting struct to Cmpi for reuse */ 13304a56b808SFande Kong c = (Mat_MPIAIJ*)Cmpi->data; 13314a56b808SFande Kong c->ap = ptap; 13324a56b808SFande Kong ptap->duplicate = Cmpi->ops->duplicate; 13334a56b808SFande Kong ptap->destroy = Cmpi->ops->destroy; 13344a56b808SFande Kong ptap->view = Cmpi->ops->view; 13354a56b808SFande Kong 13364a56b808SFande Kong /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 13374a56b808SFande Kong Cmpi->assembled = PETSC_FALSE; 13384a56b808SFande Kong /* pick an algorithm */ 13394a56b808SFande Kong ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatPtAP","Mat");CHKERRQ(ierr); 13404a56b808SFande Kong alg = 0; 13414a56b808SFande Kong ierr = PetscOptionsEList("-matptap_allatonce_via","PtAP allatonce numeric approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 13424a56b808SFande Kong ierr = PetscOptionsEnd();CHKERRQ(ierr); 13434a56b808SFande Kong switch (alg) { 13444a56b808SFande Kong case 0: 13454a56b808SFande Kong Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce; 13464a56b808SFande Kong break; 13474a56b808SFande Kong case 1: 13484a56b808SFande Kong Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged; 13494a56b808SFande Kong break; 13504a56b808SFande Kong default: 13514a56b808SFande Kong SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG," Unsupported allatonce numerical algorithm \n"); 13524a56b808SFande Kong } 13534a56b808SFande Kong Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 13544a56b808SFande Kong Cmpi->ops->view = MatView_MPIAIJ_PtAP; 13554a56b808SFande Kong Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 13564a56b808SFande Kong *C = Cmpi; 13574a56b808SFande Kong PetscFunctionReturn(0); 13584a56b808SFande Kong } 13594a56b808SFande Kong 13604a56b808SFande Kong PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce_merged(Mat A,Mat P,PetscReal fill,Mat *C) 13614a56b808SFande Kong { 13624a56b808SFande Kong Mat_APMPI *ptap; 13634a56b808SFande Kong Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c; 13644a56b808SFande Kong MPI_Comm comm; 13654a56b808SFande Kong Mat Cmpi; 13664a56b808SFande Kong Mat_SeqAIJ *pd,*po; 13674a56b808SFande Kong MatType mtype; 13684a56b808SFande Kong PetscSF sf; 13694a56b808SFande Kong PetscSFNode *iremote; 13704a56b808SFande Kong PetscInt rootspacesize,*rootspace,*rootspaceoffsets,nleaves; 13714a56b808SFande Kong const PetscInt *rootdegrees; 13724a56b808SFande Kong PetscHSetI ht,oht,*hta,*hto,*htd; 13734a56b808SFande Kong PetscInt pn,pon,*c_rmtc,i,j,nzi,dnzi,htsize,htosize,*c_rmtj,off,*c_othj,rcvncols,sendncols,*c_rmtoffsets; 13744a56b808SFande Kong PetscInt owner,lidx,*rdj,col,pcstart,pcend,*dnz,*onz,am,arstart,arend,*poj,*pdj; 13754a56b808SFande Kong PetscInt nalg=2,alg=0; 13764a56b808SFande Kong PetscBool flg; 13774a56b808SFande Kong const char *algTypes[2] = {"merged","overlapping"}; 13784a56b808SFande Kong PetscErrorCode ierr; 13794a56b808SFande Kong 13804a56b808SFande Kong PetscFunctionBegin; 13814a56b808SFande Kong ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 13824a56b808SFande Kong 13834a56b808SFande Kong /* Create symbolic parallel matrix Cmpi */ 13844a56b808SFande Kong ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr); 13854a56b808SFande Kong ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 13864a56b808SFande Kong ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 13874a56b808SFande Kong ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 13884a56b808SFande Kong ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 13894a56b808SFande Kong 13904a56b808SFande Kong ierr = PetscNew(&ptap);CHKERRQ(ierr); 13914a56b808SFande Kong ptap->reuse = MAT_INITIAL_MATRIX; 13924a56b808SFande Kong ptap->algType = 3; 13934a56b808SFande Kong 13944a56b808SFande Kong /* 0) Get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 13954a56b808SFande Kong ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 13964a56b808SFande Kong 13974a56b808SFande Kong po = (Mat_SeqAIJ*)p->B->data; 13984a56b808SFande Kong pd = (Mat_SeqAIJ*)p->A->data; 13994a56b808SFande Kong 14004a56b808SFande Kong /* This equals to the number of offdiag columns in P */ 14014a56b808SFande Kong ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr); 14024a56b808SFande Kong /* offsets */ 14034a56b808SFande Kong ierr = PetscMalloc1(pon+1,&ptap->c_rmti);CHKERRQ(ierr); 14044a56b808SFande Kong /* The number of columns we will send to remote ranks */ 14054a56b808SFande Kong ierr = PetscMalloc1(pon,&c_rmtc);CHKERRQ(ierr); 14064a56b808SFande Kong ierr = PetscMalloc1(pon,&hta);CHKERRQ(ierr); 14074a56b808SFande Kong for (i=0; i<pon; i++) { 14084a56b808SFande Kong ierr = PetscHSetICreate(&hta[i]);CHKERRQ(ierr); 14094a56b808SFande Kong } 14104a56b808SFande Kong ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr); 14114a56b808SFande Kong ierr = MatGetOwnershipRange(A,&arstart,&arend);CHKERRQ(ierr); 14124a56b808SFande Kong /* Create hash table to merge all columns for C(i, :) */ 14134a56b808SFande Kong ierr = PetscHSetICreate(&ht);CHKERRQ(ierr); 14144a56b808SFande Kong ierr = PetscHSetICreate(&oht);CHKERRQ(ierr); 14154a56b808SFande Kong ierr = PetscMalloc2(pn,&htd,pn,&hto);CHKERRQ(ierr); 14164a56b808SFande Kong for (i=0; i<pn; i++) { 14174a56b808SFande Kong ierr = PetscHSetICreate(&htd[i]);CHKERRQ(ierr); 14184a56b808SFande Kong ierr = PetscHSetICreate(&hto[i]);CHKERRQ(ierr); 14194a56b808SFande Kong } 14204a56b808SFande Kong ptap->c_rmti[0] = 0; 14214a56b808SFande Kong /* 2) Pass 1: calculate the size for C_rmt (a matrix need to be sent to other processors) */ 14224a56b808SFande Kong for (i=0; i<am && (pon || pn); i++) { 14234a56b808SFande Kong /* Form one row of AP */ 14244a56b808SFande Kong ierr = PetscHSetIClear(ht);CHKERRQ(ierr); 14254a56b808SFande Kong ierr = PetscHSetIClear(oht);CHKERRQ(ierr); 14264a56b808SFande Kong /* If the off diag is empty, we should not do any calculation */ 14274a56b808SFande Kong nzi = po->i[i+1] - po->i[i]; 14284a56b808SFande Kong dnzi = pd->i[i+1] - pd->i[i]; 14294a56b808SFande Kong if (!nzi && !dnzi) continue; 14304a56b808SFande Kong 14314a56b808SFande Kong ierr = MatPtAPSymbolicComputeOneRowOfAP_private(A,P,ptap->P_oth,i,ht,oht);CHKERRQ(ierr); 14324a56b808SFande Kong ierr = PetscHSetIGetSize(ht,&htsize);CHKERRQ(ierr); 14334a56b808SFande Kong ierr = PetscHSetIGetSize(oht,&htosize);CHKERRQ(ierr); 14344a56b808SFande Kong /* If AP is empty, just continue */ 14354a56b808SFande Kong if (!(htsize+htosize)) continue; 14364a56b808SFande Kong 14374a56b808SFande Kong /* Form remote C(ii, :) */ 14384a56b808SFande Kong poj = po->j + po->i[i]; 14394a56b808SFande Kong for (j=0; j<nzi; j++) { 1440694f88d4SFande Kong ierr = PetscHSetIUpdate(hta[poj[j]],ht);CHKERRQ(ierr); 1441694f88d4SFande Kong ierr = PetscHSetIUpdate(hta[poj[j]],oht);CHKERRQ(ierr); 14424a56b808SFande Kong } 14434a56b808SFande Kong 14444a56b808SFande Kong /* Form local C(ii, :) */ 14454a56b808SFande Kong pdj = pd->j + pd->i[i]; 14464a56b808SFande Kong for (j=0; j<dnzi; j++) { 1447694f88d4SFande Kong ierr = PetscHSetIUpdate(htd[pdj[j]],ht);CHKERRQ(ierr); 1448694f88d4SFande Kong ierr = PetscHSetIUpdate(hto[pdj[j]],oht);CHKERRQ(ierr); 14494a56b808SFande Kong } 14504a56b808SFande Kong } 14514a56b808SFande Kong 14524a56b808SFande Kong ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr); 14534a56b808SFande Kong ierr = PetscHSetIDestroy(&oht);CHKERRQ(ierr); 14544a56b808SFande Kong 14554a56b808SFande Kong for (i=0; i<pon; i++) { 14564a56b808SFande Kong ierr = PetscHSetIGetSize(hta[i],&htsize);CHKERRQ(ierr); 14574a56b808SFande Kong ptap->c_rmti[i+1] = ptap->c_rmti[i] + htsize; 14584a56b808SFande Kong c_rmtc[i] = htsize; 14594a56b808SFande Kong } 14604a56b808SFande Kong 14614a56b808SFande Kong ierr = PetscMalloc1(ptap->c_rmti[pon],&c_rmtj);CHKERRQ(ierr); 14624a56b808SFande Kong 14634a56b808SFande Kong for (i=0; i<pon; i++) { 14644a56b808SFande Kong off = 0; 14654a56b808SFande Kong ierr = PetscHSetIGetElems(hta[i],&off,c_rmtj+ptap->c_rmti[i]);CHKERRQ(ierr); 14664a56b808SFande Kong ierr = PetscHSetIDestroy(&hta[i]);CHKERRQ(ierr); 14674a56b808SFande Kong } 14684a56b808SFande Kong ierr = PetscFree(hta);CHKERRQ(ierr); 14694a56b808SFande Kong 14704a56b808SFande Kong ierr = PetscCalloc1(pon,&iremote);CHKERRQ(ierr); 14714a56b808SFande Kong for (i=0; i<pon; i++) { 1472ef7a94f2SFande Kong owner = 0; lidx = 0; 14734a56b808SFande Kong ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[i],&owner,&lidx);CHKERRQ(ierr); 14744a56b808SFande Kong iremote[i].index = lidx; 14754a56b808SFande Kong iremote[i].rank = owner; 14764a56b808SFande Kong } 14774a56b808SFande Kong 14784a56b808SFande Kong ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr); 14794a56b808SFande Kong ierr = PetscSFSetGraph(sf,pn,pon,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 14804a56b808SFande Kong /* Reorder ranks properly so that the data handled by gather and scatter have the same order */ 14814a56b808SFande Kong ierr = PetscSFSetRankOrder(sf,PETSC_TRUE);CHKERRQ(ierr); 14824a56b808SFande Kong ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 14834a56b808SFande Kong ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 14844a56b808SFande Kong /* How many neighbors have contributions to my rows? */ 14854a56b808SFande Kong ierr = PetscSFComputeDegreeBegin(sf,&rootdegrees);CHKERRQ(ierr); 14864a56b808SFande Kong ierr = PetscSFComputeDegreeEnd(sf,&rootdegrees);CHKERRQ(ierr); 14874a56b808SFande Kong rootspacesize = 0; 14884a56b808SFande Kong for (i = 0; i < pn; i++) { 14894a56b808SFande Kong rootspacesize += rootdegrees[i]; 14904a56b808SFande Kong } 14914a56b808SFande Kong ierr = PetscCalloc1(rootspacesize,&rootspace);CHKERRQ(ierr); 14924a56b808SFande Kong ierr = PetscCalloc1(rootspacesize+1,&rootspaceoffsets);CHKERRQ(ierr); 14934a56b808SFande Kong /* Get information from leaves 14944a56b808SFande Kong * Number of columns other people contribute to my rows 14954a56b808SFande Kong * */ 14964a56b808SFande Kong ierr = PetscSFGatherBegin(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr); 14974a56b808SFande Kong ierr = PetscSFGatherEnd(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr); 14984a56b808SFande Kong ierr = PetscFree(c_rmtc);CHKERRQ(ierr); 14994a56b808SFande Kong ierr = PetscCalloc1(pn+1,&ptap->c_othi);CHKERRQ(ierr); 15004a56b808SFande Kong /* The number of columns is received for each row */ 15014a56b808SFande Kong ptap->c_othi[0] = 0; 15024a56b808SFande Kong rootspacesize = 0; 15034a56b808SFande Kong rootspaceoffsets[0] = 0; 15044a56b808SFande Kong for (i = 0; i < pn; i++) { 15054a56b808SFande Kong rcvncols = 0; 15064a56b808SFande Kong for (j = 0; j<rootdegrees[i]; j++) { 15074a56b808SFande Kong rcvncols += rootspace[rootspacesize]; 15084a56b808SFande Kong rootspaceoffsets[rootspacesize+1] = rootspaceoffsets[rootspacesize] + rootspace[rootspacesize]; 15094a56b808SFande Kong rootspacesize++; 15104a56b808SFande Kong } 15114a56b808SFande Kong ptap->c_othi[i+1] = ptap->c_othi[i] + rcvncols; 15124a56b808SFande Kong } 15134a56b808SFande Kong ierr = PetscFree(rootspace);CHKERRQ(ierr); 15144a56b808SFande Kong 15154a56b808SFande Kong ierr = PetscCalloc1(pon,&c_rmtoffsets);CHKERRQ(ierr); 15164a56b808SFande Kong ierr = PetscSFScatterBegin(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr); 15174a56b808SFande Kong ierr = PetscSFScatterEnd(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr); 15184a56b808SFande Kong ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 15194a56b808SFande Kong ierr = PetscFree(rootspaceoffsets);CHKERRQ(ierr); 15204a56b808SFande Kong 15214a56b808SFande Kong ierr = PetscCalloc1(ptap->c_rmti[pon],&iremote);CHKERRQ(ierr); 15224a56b808SFande Kong nleaves = 0; 15234a56b808SFande Kong for (i = 0; i<pon; i++) { 1524ef7a94f2SFande Kong owner = 0; lidx = 0; 15254a56b808SFande Kong ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[i],&owner,&lidx);CHKERRQ(ierr); 15264a56b808SFande Kong sendncols = ptap->c_rmti[i+1] - ptap->c_rmti[i]; 15274a56b808SFande Kong for (j=0; j<sendncols; j++) { 15284a56b808SFande Kong iremote[nleaves].rank = owner; 15294a56b808SFande Kong iremote[nleaves++].index = c_rmtoffsets[i] + j; 15304a56b808SFande Kong } 15314a56b808SFande Kong } 15324a56b808SFande Kong ierr = PetscFree(c_rmtoffsets);CHKERRQ(ierr); 15334a56b808SFande Kong ierr = PetscCalloc1(ptap->c_othi[pn],&c_othj);CHKERRQ(ierr); 15344a56b808SFande Kong 15354a56b808SFande Kong ierr = PetscSFCreate(comm,&ptap->sf);CHKERRQ(ierr); 15364a56b808SFande Kong ierr = PetscSFSetGraph(ptap->sf,ptap->c_othi[pn],nleaves,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 15374a56b808SFande Kong ierr = PetscSFSetFromOptions(ptap->sf);CHKERRQ(ierr); 15384a56b808SFande Kong /* One to one map */ 15394a56b808SFande Kong ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr); 15404a56b808SFande Kong /* Get remote data */ 15414a56b808SFande Kong ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPIU_REPLACE);CHKERRQ(ierr); 15424a56b808SFande Kong ierr = PetscFree(c_rmtj);CHKERRQ(ierr); 15434a56b808SFande Kong ierr = PetscCalloc2(pn,&dnz,pn,&onz);CHKERRQ(ierr); 15444a56b808SFande Kong ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr); 15454a56b808SFande Kong 15464a56b808SFande Kong for (i = 0; i < pn; i++) { 15474a56b808SFande Kong nzi = ptap->c_othi[i+1] - ptap->c_othi[i]; 15484a56b808SFande Kong rdj = c_othj + ptap->c_othi[i]; 15494a56b808SFande Kong for (j = 0; j < nzi; j++) { 15504a56b808SFande Kong col = rdj[j]; 15514a56b808SFande Kong /* diag part */ 15524a56b808SFande Kong if (col>=pcstart && col<pcend) { 15534a56b808SFande Kong ierr = PetscHSetIAdd(htd[i],col);CHKERRQ(ierr); 15544a56b808SFande Kong } else { /* off diag */ 15554a56b808SFande Kong ierr = PetscHSetIAdd(hto[i],col);CHKERRQ(ierr); 15564a56b808SFande Kong } 15574a56b808SFande Kong } 15584a56b808SFande Kong ierr = PetscHSetIGetSize(htd[i],&htsize);CHKERRQ(ierr); 15594a56b808SFande Kong dnz[i] = htsize; 15604a56b808SFande Kong ierr = PetscHSetIDestroy(&htd[i]);CHKERRQ(ierr); 15614a56b808SFande Kong ierr = PetscHSetIGetSize(hto[i],&htsize);CHKERRQ(ierr); 15624a56b808SFande Kong onz[i] = htsize; 15634a56b808SFande Kong ierr = PetscHSetIDestroy(&hto[i]);CHKERRQ(ierr); 15644a56b808SFande Kong } 15654a56b808SFande Kong 15664a56b808SFande Kong ierr = PetscFree2(htd,hto);CHKERRQ(ierr); 15674a56b808SFande Kong ierr = PetscFree(c_othj);CHKERRQ(ierr); 15684a56b808SFande Kong 15694a56b808SFande Kong /* local sizes and preallocation */ 15704a56b808SFande Kong ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 15714a56b808SFande Kong ierr = MatSetBlockSizes(Cmpi,PetscAbs(P->cmap->bs),PetscAbs(P->cmap->bs));CHKERRQ(ierr); 15724a56b808SFande Kong ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 15734a56b808SFande Kong ierr = PetscFree2(dnz,onz);CHKERRQ(ierr); 15744a56b808SFande Kong 15754a56b808SFande Kong /* attach the supporting struct to Cmpi for reuse */ 15764a56b808SFande Kong c = (Mat_MPIAIJ*)Cmpi->data; 15774a56b808SFande Kong c->ap = ptap; 15784a56b808SFande Kong ptap->duplicate = Cmpi->ops->duplicate; 15794a56b808SFande Kong ptap->destroy = Cmpi->ops->destroy; 15804a56b808SFande Kong ptap->view = Cmpi->ops->view; 15814a56b808SFande Kong 15824a56b808SFande Kong /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 15834a56b808SFande Kong Cmpi->assembled = PETSC_FALSE; 15844a56b808SFande Kong /* pick an algorithm */ 15854a56b808SFande Kong ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatPtAP","Mat");CHKERRQ(ierr); 15864a56b808SFande Kong alg = 0; 15874a56b808SFande Kong ierr = PetscOptionsEList("-matptap_allatonce_via","PtAP allatonce numeric approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 15884a56b808SFande Kong ierr = PetscOptionsEnd();CHKERRQ(ierr); 15894a56b808SFande Kong switch (alg) { 15904a56b808SFande Kong case 0: 15914a56b808SFande Kong Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged; 15924a56b808SFande Kong break; 15934a56b808SFande Kong case 1: 15944a56b808SFande Kong Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce; 15954a56b808SFande Kong break; 15964a56b808SFande Kong default: 15974a56b808SFande Kong SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG," Unsupported allatonce numerical algorithm \n"); 15984a56b808SFande Kong } 15994a56b808SFande Kong Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 16004a56b808SFande Kong Cmpi->ops->view = MatView_MPIAIJ_PtAP; 16014a56b808SFande Kong Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 16024a56b808SFande Kong *C = Cmpi; 16034a56b808SFande Kong PetscFunctionReturn(0); 16044a56b808SFande Kong } 16054a56b808SFande Kong 1606e953e47cSHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat *C) 1607e953e47cSHong Zhang { 1608e953e47cSHong Zhang PetscErrorCode ierr; 16093cdca5ebSHong Zhang Mat_APMPI *ptap; 1610e953e47cSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c; 1611e953e47cSHong Zhang MPI_Comm comm; 1612e953e47cSHong Zhang PetscMPIInt size,rank; 1613e953e47cSHong Zhang Mat Cmpi; 1614e953e47cSHong Zhang PetscFreeSpaceList free_space=NULL,current_space=NULL; 1615e953e47cSHong Zhang PetscInt am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n; 1616e953e47cSHong Zhang PetscInt *lnk,i,k,pnz,row,nsend; 1617e953e47cSHong Zhang PetscBT lnkbt; 1618e953e47cSHong Zhang PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,icompleted=0,nrecv; 1619e953e47cSHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1620e953e47cSHong Zhang PetscInt len,proc,*dnz,*onz,*owners,nzi,nspacedouble; 1621e953e47cSHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1622e953e47cSHong Zhang MPI_Request *swaits,*rwaits; 1623e953e47cSHong Zhang MPI_Status *sstatus,rstatus; 1624e953e47cSHong Zhang PetscLayout rowmap; 1625e953e47cSHong Zhang PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 1626e953e47cSHong Zhang PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 1627e953e47cSHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi; 1628ec07b8f8SHong Zhang Mat_SeqAIJ *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth; 1629e953e47cSHong Zhang PetscScalar *apv; 1630e953e47cSHong Zhang PetscTable ta; 1631b92f168fSBarry Smith MatType mtype; 1632e83e5f86SFande Kong const char *prefix; 1633e953e47cSHong Zhang #if defined(PETSC_USE_INFO) 1634e953e47cSHong Zhang PetscReal apfill; 1635e953e47cSHong Zhang #endif 1636e953e47cSHong Zhang 1637e953e47cSHong Zhang PetscFunctionBegin; 1638b4e8d1b6SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1639e953e47cSHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1640e953e47cSHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 1641e953e47cSHong Zhang 164252f7967eSHong Zhang if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data; 1643ec07b8f8SHong Zhang 1644e953e47cSHong Zhang /* create symbolic parallel matrix Cmpi */ 1645e953e47cSHong Zhang ierr = MatCreate(comm,&Cmpi);CHKERRQ(ierr); 1646b92f168fSBarry Smith ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 1647b92f168fSBarry Smith ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 1648e953e47cSHong Zhang 1649e953e47cSHong Zhang /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */ 1650e953e47cSHong Zhang Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ; 1651e953e47cSHong Zhang 16523cdca5ebSHong Zhang /* create struct Mat_APMPI and attached it to C later */ 165315a3b8e2SHong Zhang ierr = PetscNew(&ptap);CHKERRQ(ierr); 165415a3b8e2SHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 1655a4ffb656SHong Zhang ptap->algType = 1; 165615a3b8e2SHong Zhang 165715a3b8e2SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 165815a3b8e2SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 165915a3b8e2SHong Zhang /* get P_loc by taking all local rows of P */ 166015a3b8e2SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 166115a3b8e2SHong Zhang 166267a12041SHong Zhang /* (0) compute Rd = Pd^T, Ro = Po^T */ 166367a12041SHong Zhang /* --------------------------------- */ 1664419ecdd9Sandi selinger ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 1665419ecdd9Sandi selinger ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 166615a3b8e2SHong Zhang 166767a12041SHong Zhang /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 166867a12041SHong Zhang /* ----------------------------------------------------------------- */ 166967a12041SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 167052f7967eSHong Zhang if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 1671445158ffSHong Zhang 16729a6dcab7SHong Zhang /* create and initialize a linked list */ 167345d00d1dSHong Zhang ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */ 16744b38b95cSHong Zhang MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta); 16754b38b95cSHong Zhang MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta); 167678d30b94SHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */ 1677d0e9a2f7SHong 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); */ 167878d30b94SHong Zhang 167978d30b94SHong Zhang ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 1680445158ffSHong Zhang 16818cb82516SHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 1682ec07b8f8SHong Zhang if (ao) { 1683f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr); 1684ec07b8f8SHong Zhang } else { 1685ec07b8f8SHong Zhang ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr); 1686ec07b8f8SHong Zhang } 1687445158ffSHong Zhang current_space = free_space; 168867a12041SHong Zhang nspacedouble = 0; 168967a12041SHong Zhang 169067a12041SHong Zhang ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 169167a12041SHong Zhang api[0] = 0; 1692445158ffSHong Zhang for (i=0; i<am; i++) { 169367a12041SHong Zhang /* diagonal portion: Ad[i,:]*P */ 169467a12041SHong Zhang ai = ad->i; pi = p_loc->i; 169567a12041SHong Zhang nzi = ai[i+1] - ai[i]; 169667a12041SHong Zhang aj = ad->j + ai[i]; 1697445158ffSHong Zhang for (j=0; j<nzi; j++) { 1698445158ffSHong Zhang row = aj[j]; 169967a12041SHong Zhang pnz = pi[row+1] - pi[row]; 170067a12041SHong Zhang Jptr = p_loc->j + pi[row]; 1701445158ffSHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 1702445158ffSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1703445158ffSHong Zhang } 170467a12041SHong Zhang /* off-diagonal portion: Ao[i,:]*P */ 1705ec07b8f8SHong Zhang if (ao) { 170667a12041SHong Zhang ai = ao->i; pi = p_oth->i; 170767a12041SHong Zhang nzi = ai[i+1] - ai[i]; 170867a12041SHong Zhang aj = ao->j + ai[i]; 1709445158ffSHong Zhang for (j=0; j<nzi; j++) { 1710445158ffSHong Zhang row = aj[j]; 171167a12041SHong Zhang pnz = pi[row+1] - pi[row]; 171267a12041SHong Zhang Jptr = p_oth->j + pi[row]; 1713445158ffSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1714445158ffSHong Zhang } 1715ec07b8f8SHong Zhang } 1716445158ffSHong Zhang apnz = lnk[0]; 1717445158ffSHong Zhang api[i+1] = api[i] + apnz; 1718445158ffSHong Zhang if (ap_rmax < apnz) ap_rmax = apnz; 1719445158ffSHong Zhang 1720445158ffSHong Zhang /* if free space is not available, double the total space in the list */ 1721445158ffSHong Zhang if (current_space->local_remaining<apnz) { 1722f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1723445158ffSHong Zhang nspacedouble++; 1724445158ffSHong Zhang } 1725445158ffSHong Zhang 1726445158ffSHong Zhang /* Copy data into free space, then initialize lnk */ 1727445158ffSHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1728445158ffSHong Zhang 1729445158ffSHong Zhang current_space->array += apnz; 1730445158ffSHong Zhang current_space->local_used += apnz; 1731445158ffSHong Zhang current_space->local_remaining -= apnz; 1732445158ffSHong Zhang } 1733681d504bSHong Zhang /* Allocate space for apj and apv, initialize apj, and */ 1734445158ffSHong Zhang /* destroy list of free space and other temporary array(s) */ 1735445158ffSHong Zhang ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 1736445158ffSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 17379a6dcab7SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 17389a6dcab7SHong Zhang 1739aa690a28SHong Zhang /* Create AP_loc for reuse */ 1740445158ffSHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 1741aa690a28SHong Zhang 1742aa690a28SHong Zhang #if defined(PETSC_USE_INFO) 1743ec07b8f8SHong Zhang if (ao) { 1744aa690a28SHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 1745ec07b8f8SHong Zhang } else { 1746ec07b8f8SHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1); 1747ec07b8f8SHong Zhang } 1748aa690a28SHong Zhang ptap->AP_loc->info.mallocs = nspacedouble; 1749aa690a28SHong Zhang ptap->AP_loc->info.fill_ratio_given = fill; 1750aa690a28SHong Zhang ptap->AP_loc->info.fill_ratio_needed = apfill; 1751aa690a28SHong Zhang 1752aa690a28SHong Zhang if (api[am]) { 1753b11c1ec8SHong 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); 1754aa690a28SHong Zhang ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 1755aa690a28SHong Zhang } else { 1756b11c1ec8SHong Zhang ierr = PetscInfo(ptap->AP_loc,"Nonscalable algorithm, AP_loc is empty \n");CHKERRQ(ierr); 1757aa690a28SHong Zhang } 1758aa690a28SHong Zhang #endif 1759aa690a28SHong Zhang 1760681d504bSHong Zhang /* (2-1) compute symbolic Co = Ro*AP_loc */ 1761681d504bSHong Zhang /* ------------------------------------ */ 1762e83e5f86SFande Kong ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1763e83e5f86SFande Kong ierr = MatSetOptionsPrefix(ptap->Ro,prefix);CHKERRQ(ierr); 1764e83e5f86SFande Kong ierr = MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");CHKERRQ(ierr); 176567a12041SHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Ro,ptap->AP_loc,fill,&ptap->C_oth);CHKERRQ(ierr); 176615a3b8e2SHong Zhang 176767a12041SHong Zhang /* (3) send coj of C_oth to other processors */ 176867a12041SHong Zhang /* ------------------------------------------ */ 176915a3b8e2SHong Zhang /* determine row ownership */ 1770445158ffSHong Zhang ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 1771445158ffSHong Zhang rowmap->n = pn; 1772445158ffSHong Zhang rowmap->bs = 1; 1773445158ffSHong Zhang ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 1774445158ffSHong Zhang owners = rowmap->range; 177515a3b8e2SHong Zhang 177615a3b8e2SHong Zhang /* determine the number of messages to send, their lengths */ 17778cb82516SHong Zhang ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 17788cb82516SHong Zhang ierr = PetscMemzero(len_s,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 177915a3b8e2SHong Zhang ierr = PetscMemzero(len_si,size*sizeof(PetscMPIInt));CHKERRQ(ierr); 178015a3b8e2SHong Zhang 178167a12041SHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 178267a12041SHong Zhang coi = c_oth->i; coj = c_oth->j; 178367a12041SHong Zhang con = ptap->C_oth->rmap->n; 178415a3b8e2SHong Zhang proc = 0; 178567a12041SHong Zhang for (i=0; i<con; i++) { 178615a3b8e2SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 178715a3b8e2SHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 178815a3b8e2SHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 178915a3b8e2SHong Zhang } 179015a3b8e2SHong Zhang 179115a3b8e2SHong Zhang len = 0; /* max length of buf_si[], see (4) */ 179215a3b8e2SHong Zhang owners_co[0] = 0; 179367a12041SHong Zhang nsend = 0; 179415a3b8e2SHong Zhang for (proc=0; proc<size; proc++) { 179515a3b8e2SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 179615a3b8e2SHong Zhang if (len_s[proc]) { 1797445158ffSHong Zhang nsend++; 179815a3b8e2SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 179915a3b8e2SHong Zhang len += len_si[proc]; 180015a3b8e2SHong Zhang } 180115a3b8e2SHong Zhang } 180215a3b8e2SHong Zhang 180315a3b8e2SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 1804445158ffSHong Zhang ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 1805445158ffSHong Zhang ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 180615a3b8e2SHong Zhang 180715a3b8e2SHong Zhang /* post the Irecv and Isend of coj */ 180815a3b8e2SHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1809445158ffSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1810445158ffSHong Zhang ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 181115a3b8e2SHong Zhang for (proc=0, k=0; proc<size; proc++) { 181215a3b8e2SHong Zhang if (!len_s[proc]) continue; 181315a3b8e2SHong Zhang i = owners_co[proc]; 181415a3b8e2SHong Zhang ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRQ(ierr); 181515a3b8e2SHong Zhang k++; 181615a3b8e2SHong Zhang } 181715a3b8e2SHong Zhang 1818681d504bSHong Zhang /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 1819681d504bSHong Zhang /* ---------------------------------------- */ 1820e83e5f86SFande Kong ierr = MatSetOptionsPrefix(ptap->Rd,prefix);CHKERRQ(ierr); 1821e83e5f86SFande Kong ierr = MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");CHKERRQ(ierr); 1822681d504bSHong Zhang ierr = MatMatMultSymbolic_SeqAIJ_SeqAIJ(ptap->Rd,ptap->AP_loc,fill,&ptap->C_loc);CHKERRQ(ierr); 1823681d504bSHong Zhang c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 1824681d504bSHong Zhang 1825681d504bSHong Zhang /* receives coj are complete */ 1826445158ffSHong Zhang for (i=0; i<nrecv; i++) { 1827445158ffSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 182815a3b8e2SHong Zhang } 182915a3b8e2SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1830445158ffSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 183115a3b8e2SHong Zhang 183278d30b94SHong Zhang /* add received column indices into ta to update Crmax */ 183378d30b94SHong Zhang for (k=0; k<nrecv; k++) {/* k-th received message */ 183478d30b94SHong Zhang Jptr = buf_rj[k]; 183578d30b94SHong Zhang for (j=0; j<len_r[k]; j++) { 183678d30b94SHong Zhang ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 183778d30b94SHong Zhang } 183878d30b94SHong Zhang } 183978d30b94SHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 184078d30b94SHong Zhang ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 18419a6dcab7SHong Zhang 184215a3b8e2SHong Zhang /* (4) send and recv coi */ 184315a3b8e2SHong Zhang /*-----------------------*/ 184415a3b8e2SHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1845445158ffSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 184615a3b8e2SHong Zhang ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 184715a3b8e2SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 184815a3b8e2SHong Zhang for (proc=0,k=0; proc<size; proc++) { 184915a3b8e2SHong Zhang if (!len_s[proc]) continue; 185015a3b8e2SHong Zhang /* form outgoing message for i-structure: 185115a3b8e2SHong Zhang buf_si[0]: nrows to be sent 185215a3b8e2SHong Zhang [1:nrows]: row index (global) 185315a3b8e2SHong Zhang [nrows+1:2*nrows+1]: i-structure index 185415a3b8e2SHong Zhang */ 185515a3b8e2SHong Zhang /*-------------------------------------------*/ 185615a3b8e2SHong Zhang nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 185715a3b8e2SHong Zhang buf_si_i = buf_si + nrows+1; 185815a3b8e2SHong Zhang buf_si[0] = nrows; 185915a3b8e2SHong Zhang buf_si_i[0] = 0; 186015a3b8e2SHong Zhang nrows = 0; 186115a3b8e2SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 186215a3b8e2SHong Zhang nzi = coi[i+1] - coi[i]; 186315a3b8e2SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 186415a3b8e2SHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 186515a3b8e2SHong Zhang nrows++; 186615a3b8e2SHong Zhang } 186715a3b8e2SHong Zhang ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRQ(ierr); 186815a3b8e2SHong Zhang k++; 186915a3b8e2SHong Zhang buf_si += len_si[proc]; 187015a3b8e2SHong Zhang } 1871681d504bSHong Zhang for (i=0; i<nrecv; i++) { 1872445158ffSHong Zhang ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRQ(ierr); 187315a3b8e2SHong Zhang } 187415a3b8e2SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1875445158ffSHong Zhang if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRQ(ierr);} 187615a3b8e2SHong Zhang 18778cb82516SHong Zhang ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 187815a3b8e2SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 187915a3b8e2SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 188015a3b8e2SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 1881a4ffb656SHong Zhang 188267a12041SHong Zhang /* (5) compute the local portion of Cmpi */ 188367a12041SHong Zhang /* ------------------------------------------ */ 188478d30b94SHong Zhang /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 188578d30b94SHong Zhang ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 188615a3b8e2SHong Zhang current_space = free_space; 188715a3b8e2SHong Zhang 1888445158ffSHong Zhang ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 1889445158ffSHong Zhang for (k=0; k<nrecv; k++) { 189015a3b8e2SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 189115a3b8e2SHong Zhang nrows = *buf_ri_k[k]; 189215a3b8e2SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 189315a3b8e2SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 189415a3b8e2SHong Zhang } 1895445158ffSHong Zhang 189678d30b94SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 189778d30b94SHong Zhang ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 189815a3b8e2SHong Zhang for (i=0; i<pn; i++) { 189967a12041SHong Zhang /* add C_loc into Cmpi */ 190067a12041SHong Zhang nzi = c_loc->i[i+1] - c_loc->i[i]; 190167a12041SHong Zhang Jptr = c_loc->j + c_loc->i[i]; 190267a12041SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 190315a3b8e2SHong Zhang 190415a3b8e2SHong Zhang /* add received col data into lnk */ 1905445158ffSHong Zhang for (k=0; k<nrecv; k++) { /* k-th received message */ 190615a3b8e2SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 190715a3b8e2SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 190815a3b8e2SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 190915a3b8e2SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 191015a3b8e2SHong Zhang nextrow[k]++; nextci[k]++; 191115a3b8e2SHong Zhang } 191215a3b8e2SHong Zhang } 1913d0e9a2f7SHong Zhang nzi = lnk[0]; 19148cb82516SHong Zhang 191515a3b8e2SHong Zhang /* copy data into free space, then initialize lnk */ 1916d0e9a2f7SHong Zhang ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1917d0e9a2f7SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 191815a3b8e2SHong Zhang } 191915a3b8e2SHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 192015a3b8e2SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1921445158ffSHong Zhang ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 192280bb4639SHong Zhang 1923ae5f0867Sstefano_zampini /* local sizes and preallocation */ 192415a3b8e2SHong Zhang ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1925ac94c67aSTristan Konolige if (P->cmap->bs > 0) { 1926ac94c67aSTristan Konolige ierr = PetscLayoutSetBlockSize(Cmpi->rmap,P->cmap->bs);CHKERRQ(ierr); 1927ac94c67aSTristan Konolige ierr = PetscLayoutSetBlockSize(Cmpi->cmap,P->cmap->bs);CHKERRQ(ierr); 1928ac94c67aSTristan Konolige } 192915a3b8e2SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 193015a3b8e2SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 193115a3b8e2SHong Zhang 1932445158ffSHong Zhang /* members in merge */ 1933445158ffSHong Zhang ierr = PetscFree(id_r);CHKERRQ(ierr); 1934445158ffSHong Zhang ierr = PetscFree(len_r);CHKERRQ(ierr); 1935445158ffSHong Zhang ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 1936445158ffSHong Zhang ierr = PetscFree(buf_ri);CHKERRQ(ierr); 1937445158ffSHong Zhang ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 1938445158ffSHong Zhang ierr = PetscFree(buf_rj);CHKERRQ(ierr); 1939445158ffSHong Zhang ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 194015a3b8e2SHong Zhang 194115a3b8e2SHong Zhang /* attach the supporting struct to Cmpi for reuse */ 194215a3b8e2SHong Zhang c = (Mat_MPIAIJ*)Cmpi->data; 19433cdca5ebSHong Zhang c->ap = ptap; 19441a47ec54SHong Zhang ptap->duplicate = Cmpi->ops->duplicate; 19451a47ec54SHong Zhang ptap->destroy = Cmpi->ops->destroy; 1946cf97cf3cSHong Zhang ptap->view = Cmpi->ops->view; 19478cb82516SHong Zhang ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 19482259aa2eSHong Zhang 19491a47ec54SHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 19501a47ec54SHong Zhang Cmpi->assembled = PETSC_FALSE; 19511a47ec54SHong Zhang Cmpi->ops->destroy = MatDestroy_MPIAIJ_PtAP; 1952cf97cf3cSHong Zhang Cmpi->ops->view = MatView_MPIAIJ_PtAP; 19534624976aSHong Zhang Cmpi->ops->freeintermediatedatastructures = MatFreeIntermediateDataStructures_MPIAIJ_AP; 19541a47ec54SHong Zhang *C = Cmpi; 19550d3441aeSHong Zhang PetscFunctionReturn(0); 19560d3441aeSHong Zhang } 19570d3441aeSHong Zhang 1958aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 19590d3441aeSHong Zhang { 19600d3441aeSHong Zhang PetscErrorCode ierr; 19610d3441aeSHong Zhang Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 19620d3441aeSHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 19632dd9e643SHong Zhang Mat_SeqAIJ *ap,*p_loc,*p_oth=NULL,*c_seq; 19643cdca5ebSHong Zhang Mat_APMPI *ptap = c->ap; 19659ce11a7cSHong Zhang Mat AP_loc,C_loc,C_oth; 19660d3441aeSHong Zhang PetscInt i,rstart,rend,cm,ncols,row; 19678cb82516SHong Zhang PetscInt *api,*apj,am = A->rmap->n,j,col,apnz; 19688cb82516SHong Zhang PetscScalar *apa; 19690d3441aeSHong Zhang const PetscInt *cols; 19700d3441aeSHong Zhang const PetscScalar *vals; 19710d3441aeSHong Zhang 19720d3441aeSHong Zhang PetscFunctionBegin; 1973bb9bda99SHong Zhang if (!ptap->AP_loc) { 1974a9899c97SHong Zhang MPI_Comm comm; 1975a9899c97SHong Zhang ierr = PetscObjectGetComm((PetscObject)C,&comm);CHKERRQ(ierr); 1976dd4011a9SHong Zhang SETERRQ(comm,PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatFreeIntermediateDataStructures() or use '-mat_freeintermediatedatastructures'"); 1977a9899c97SHong Zhang } 1978a9899c97SHong Zhang 19790d3441aeSHong Zhang ierr = MatZeroEntries(C);CHKERRQ(ierr); 1980e497d3c8SHong Zhang /* 1) get R = Pd^T,Ro = Po^T */ 1981748c7196SHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 1982419ecdd9Sandi selinger ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 1983419ecdd9Sandi selinger ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 1984748c7196SHong Zhang } 19850d3441aeSHong Zhang 1986e497d3c8SHong Zhang /* 2) get AP_loc */ 19870d3441aeSHong Zhang AP_loc = ptap->AP_loc; 19888cb82516SHong Zhang ap = (Mat_SeqAIJ*)AP_loc->data; 19890d3441aeSHong Zhang 1990e497d3c8SHong Zhang /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 19910d3441aeSHong Zhang /*-----------------------------------------------------*/ 1992748c7196SHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 1993748c7196SHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 19940d3441aeSHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 19950d3441aeSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 1996e497d3c8SHong Zhang } 19970d3441aeSHong Zhang 19988cb82516SHong Zhang /* 2-2) compute numeric A_loc*P - dominating part */ 19998cb82516SHong Zhang /* ---------------------------------------------- */ 20000d3441aeSHong Zhang /* get data from symbolic products */ 20018cb82516SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 20022dd9e643SHong Zhang if (ptap->P_oth) { 20038cb82516SHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 20042dd9e643SHong Zhang } 20058cb82516SHong Zhang apa = ptap->apa; 2006681d504bSHong Zhang api = ap->i; 2007681d504bSHong Zhang apj = ap->j; 2008e497d3c8SHong Zhang for (i=0; i<am; i++) { 2009445158ffSHong Zhang /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 2010e497d3c8SHong Zhang AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa); 2011e497d3c8SHong Zhang apnz = api[i+1] - api[i]; 2012e497d3c8SHong Zhang for (j=0; j<apnz; j++) { 2013e497d3c8SHong Zhang col = apj[j+api[i]]; 2014e497d3c8SHong Zhang ap->a[j+ap->i[i]] = apa[col]; 2015e497d3c8SHong Zhang apa[col] = 0.0; 20160d3441aeSHong Zhang } 2017aa690a28SHong Zhang ierr = PetscLogFlops(2.0*apnz);CHKERRQ(ierr); 20180d3441aeSHong Zhang } 20190d3441aeSHong Zhang 20208cb82516SHong Zhang /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 2021445158ffSHong Zhang ierr = ((ptap->C_loc)->ops->matmultnumeric)(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 2022445158ffSHong Zhang ierr = ((ptap->C_oth)->ops->matmultnumeric)(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 20239ce11a7cSHong Zhang C_loc = ptap->C_loc; 20249ce11a7cSHong Zhang C_oth = ptap->C_oth; 20250d3441aeSHong Zhang 20260d3441aeSHong Zhang /* add C_loc and Co to to C */ 20270d3441aeSHong Zhang ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 20280d3441aeSHong Zhang 20290d3441aeSHong Zhang /* C_loc -> C */ 20300d3441aeSHong Zhang cm = C_loc->rmap->N; 20319ce11a7cSHong Zhang c_seq = (Mat_SeqAIJ*)C_loc->data; 20328cb82516SHong Zhang cols = c_seq->j; 20338cb82516SHong Zhang vals = c_seq->a; 2034904d1e70Sandi selinger 2035904d1e70Sandi selinger 2036e9ede7d0Sandi selinger /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */ 2037904d1e70Sandi selinger /* when there are no off-processor parts. */ 20381de21080Sandi selinger /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */ 20391de21080Sandi selinger /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */ 20401de21080Sandi selinger /* a table, and other, more complex stuff has to be done. */ 2041904d1e70Sandi selinger if (C->assembled) { 2042904d1e70Sandi selinger C->was_assembled = PETSC_TRUE; 2043904d1e70Sandi selinger C->assembled = PETSC_FALSE; 2044904d1e70Sandi selinger } 2045904d1e70Sandi selinger if (C->was_assembled) { 20460d3441aeSHong Zhang for (i=0; i<cm; i++) { 20479ce11a7cSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 20480d3441aeSHong Zhang row = rstart + i; 2049904d1e70Sandi selinger ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 20508cb82516SHong Zhang cols += ncols; vals += ncols; 20510d3441aeSHong Zhang } 2052904d1e70Sandi selinger } else { 2053e9ede7d0Sandi selinger ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr); 2054904d1e70Sandi selinger } 20550d3441aeSHong Zhang 20560d3441aeSHong Zhang /* Co -> C, off-processor part */ 20579ce11a7cSHong Zhang cm = C_oth->rmap->N; 20589ce11a7cSHong Zhang c_seq = (Mat_SeqAIJ*)C_oth->data; 20598cb82516SHong Zhang cols = c_seq->j; 20608cb82516SHong Zhang vals = c_seq->a; 20619ce11a7cSHong Zhang for (i=0; i<cm; i++) { 20629ce11a7cSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 20630d3441aeSHong Zhang row = p->garray[i]; 20640d3441aeSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 20658cb82516SHong Zhang cols += ncols; vals += ncols; 20660d3441aeSHong Zhang } 2067904d1e70Sandi selinger 20680d3441aeSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20690d3441aeSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20700d3441aeSHong Zhang 2071748c7196SHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 20723697aca6SHong Zhang 2073dd4011a9SHong 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 */ 20747d0a89b7SHong Zhang if (ptap->freestruct) { 2075dd4011a9SHong Zhang ierr = MatFreeIntermediateDataStructures(C);CHKERRQ(ierr); 20763697aca6SHong Zhang } 20770d3441aeSHong Zhang PetscFunctionReturn(0); 20780d3441aeSHong Zhang } 2079