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 16cf97cf3cSHong Zhang PetscErrorCode MatView_MPIAIJ_PtAP(Mat A,PetscViewer viewer) 17cf97cf3cSHong Zhang { 18cf97cf3cSHong Zhang PetscErrorCode ierr; 19cf97cf3cSHong Zhang PetscBool iascii; 20cf97cf3cSHong Zhang PetscViewerFormat format; 216718818eSStefano Zampini Mat_APMPI *ptap; 22cf97cf3cSHong Zhang 23cf97cf3cSHong Zhang PetscFunctionBegin; 246718818eSStefano Zampini MatCheckProduct(A,1); 256718818eSStefano Zampini ptap = (Mat_APMPI*)A->product->data; 26cf97cf3cSHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 27cf97cf3cSHong Zhang if (iascii) { 28cf97cf3cSHong Zhang ierr = PetscViewerGetFormat(viewer,&format);CHKERRQ(ierr); 29cf97cf3cSHong Zhang if (format == PETSC_VIEWER_ASCII_INFO || format == PETSC_VIEWER_ASCII_INFO_DETAIL) { 30cf97cf3cSHong Zhang if (ptap->algType == 0) { 31cf97cf3cSHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"using scalable MatPtAP() implementation\n");CHKERRQ(ierr); 32cf97cf3cSHong Zhang } else if (ptap->algType == 1) { 33cf97cf3cSHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"using nonscalable MatPtAP() implementation\n");CHKERRQ(ierr); 344a56b808SFande Kong } else if (ptap->algType == 2) { 354a56b808SFande Kong ierr = PetscViewerASCIIPrintf(viewer,"using allatonce MatPtAP() implementation\n");CHKERRQ(ierr); 364a56b808SFande Kong } else if (ptap->algType == 3) { 374a56b808SFande Kong ierr = PetscViewerASCIIPrintf(viewer,"using merged allatonce MatPtAP() implementation\n");CHKERRQ(ierr); 38cf97cf3cSHong Zhang } 39cf97cf3cSHong Zhang } 40cf97cf3cSHong Zhang } 41cf97cf3cSHong Zhang PetscFunctionReturn(0); 42cf97cf3cSHong Zhang } 43cf97cf3cSHong Zhang 446718818eSStefano Zampini PetscErrorCode MatDestroy_MPIAIJ_PtAP(void *data) 45cc6cb787SHong Zhang { 46cc6cb787SHong Zhang PetscErrorCode ierr; 476718818eSStefano Zampini Mat_APMPI *ptap = (Mat_APMPI*)data; 4860552ceaSHong Zhang Mat_Merge_SeqsToMPI *merge; 49cc6cb787SHong Zhang 50cc6cb787SHong Zhang PetscFunctionBegin; 51b7f45c76SHong Zhang ierr = PetscFree2(ptap->startsj_s,ptap->startsj_r);CHKERRQ(ierr); 52f8487c73SHong Zhang ierr = PetscFree(ptap->bufa);CHKERRQ(ierr); 53a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_loc);CHKERRQ(ierr); 54a1a4e84aSHong Zhang ierr = MatDestroy(&ptap->P_oth);CHKERRQ(ierr); 55c5af039cSHong Zhang ierr = MatDestroy(&ptap->A_loc);CHKERRQ(ierr); /* used by MatTransposeMatMult() */ 5605d62848SHong Zhang ierr = MatDestroy(&ptap->Rd);CHKERRQ(ierr); 5705d62848SHong Zhang ierr = MatDestroy(&ptap->Ro);CHKERRQ(ierr); 588403a639SHong Zhang if (ptap->AP_loc) { /* used by alg_rap */ 59681d504bSHong Zhang Mat_SeqAIJ *ap = (Mat_SeqAIJ*)(ptap->AP_loc)->data; 60681d504bSHong Zhang ierr = PetscFree(ap->i);CHKERRQ(ierr); 61681d504bSHong Zhang ierr = PetscFree2(ap->j,ap->a);CHKERRQ(ierr); 620d3441aeSHong Zhang ierr = MatDestroy(&ptap->AP_loc);CHKERRQ(ierr); 638403a639SHong Zhang } else { /* used by alg_ptap */ 648403a639SHong Zhang ierr = PetscFree(ptap->api);CHKERRQ(ierr); 658403a639SHong Zhang ierr = PetscFree(ptap->apj);CHKERRQ(ierr); 66681d504bSHong Zhang } 672259aa2eSHong Zhang ierr = MatDestroy(&ptap->C_loc);CHKERRQ(ierr); 682259aa2eSHong Zhang ierr = MatDestroy(&ptap->C_oth);CHKERRQ(ierr); 69414894bdSHong Zhang if (ptap->apa) {ierr = PetscFree(ptap->apa);CHKERRQ(ierr);} 70a560ef98SHong Zhang 7160552ceaSHong Zhang ierr = MatDestroy(&ptap->Pt);CHKERRQ(ierr); 7260552ceaSHong Zhang 7360552ceaSHong Zhang merge = ptap->merge; 748403a639SHong Zhang if (merge) { /* used by alg_ptap */ 75cc6cb787SHong Zhang ierr = PetscFree(merge->id_r);CHKERRQ(ierr); 76cc6cb787SHong Zhang ierr = PetscFree(merge->len_s);CHKERRQ(ierr); 77cc6cb787SHong Zhang ierr = PetscFree(merge->len_r);CHKERRQ(ierr); 78cc6cb787SHong Zhang ierr = PetscFree(merge->bi);CHKERRQ(ierr); 79cc6cb787SHong Zhang ierr = PetscFree(merge->bj);CHKERRQ(ierr); 80c05d87d6SBarry Smith ierr = PetscFree(merge->buf_ri[0]);CHKERRQ(ierr); 81cc6cb787SHong Zhang ierr = PetscFree(merge->buf_ri);CHKERRQ(ierr); 82c05d87d6SBarry Smith ierr = PetscFree(merge->buf_rj[0]);CHKERRQ(ierr); 83cc6cb787SHong Zhang ierr = PetscFree(merge->buf_rj);CHKERRQ(ierr); 84445158ffSHong Zhang ierr = PetscFree(merge->coi);CHKERRQ(ierr); 85445158ffSHong Zhang ierr = PetscFree(merge->coj);CHKERRQ(ierr); 8605b42c5fSBarry Smith ierr = PetscFree(merge->owners_co);CHKERRQ(ierr); 876bf464f9SBarry Smith ierr = PetscLayoutDestroy(&merge->rowmap);CHKERRQ(ierr); 88f8487c73SHong Zhang ierr = PetscFree(ptap->merge);CHKERRQ(ierr); 89bf0cc555SLisandro Dalcin } 90a3bb6f32SFande Kong ierr = ISLocalToGlobalMappingDestroy(&ptap->ltog);CHKERRQ(ierr); 914a56b808SFande Kong 924a56b808SFande Kong ierr = PetscSFDestroy(&ptap->sf);CHKERRQ(ierr); 934a56b808SFande Kong ierr = PetscFree(ptap->c_othi);CHKERRQ(ierr); 944a56b808SFande Kong ierr = PetscFree(ptap->c_rmti);CHKERRQ(ierr); 959b3d8a68SHong Zhang ierr = PetscFree(ptap);CHKERRQ(ierr); 96cc6cb787SHong Zhang PetscFunctionReturn(0); 97cc6cb787SHong Zhang } 98cc6cb787SHong Zhang 99cf742fccSHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,Mat C) 100cf742fccSHong Zhang { 101cf742fccSHong Zhang PetscErrorCode ierr; 1026718818eSStefano Zampini Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 103cf742fccSHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 10492ec70a1SHong Zhang Mat_SeqAIJ *ap,*p_loc,*p_oth=NULL,*c_seq; 1056718818eSStefano Zampini Mat_APMPI *ptap; 106cf742fccSHong Zhang Mat AP_loc,C_loc,C_oth; 107a3bb6f32SFande Kong PetscInt i,rstart,rend,cm,ncols,row,*api,*apj,am = A->rmap->n,apnz,nout; 108cf742fccSHong Zhang PetscScalar *apa; 109cf742fccSHong Zhang const PetscInt *cols; 110cf742fccSHong Zhang const PetscScalar *vals; 111cf742fccSHong Zhang 112cf742fccSHong Zhang PetscFunctionBegin; 1136718818eSStefano Zampini MatCheckProduct(C,3); 1146718818eSStefano Zampini ptap = (Mat_APMPI*)C->product->data; 1156718818eSStefano Zampini if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data"); 1166718818eSStefano Zampini if (!ptap->AP_loc) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatProductClear()"); 11780da8df7SHong Zhang 118cf742fccSHong Zhang ierr = MatZeroEntries(C);CHKERRQ(ierr); 119cf742fccSHong Zhang 120cf742fccSHong Zhang /* 1) get R = Pd^T,Ro = Po^T */ 121cf742fccSHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 122419ecdd9Sandi selinger ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 123419ecdd9Sandi selinger ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 124cf742fccSHong Zhang } 125cf742fccSHong Zhang 126cf742fccSHong Zhang /* 2) get AP_loc */ 127cf742fccSHong Zhang AP_loc = ptap->AP_loc; 128cf742fccSHong Zhang ap = (Mat_SeqAIJ*)AP_loc->data; 129cf742fccSHong Zhang 130cf742fccSHong Zhang /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 131cf742fccSHong Zhang /*-----------------------------------------------------*/ 132cf742fccSHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 133cf742fccSHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 134cf742fccSHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 135cf742fccSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 136cf742fccSHong Zhang } 137cf742fccSHong Zhang 138cf742fccSHong Zhang /* 2-2) compute numeric A_loc*P - dominating part */ 139cf742fccSHong Zhang /* ---------------------------------------------- */ 140cf742fccSHong Zhang /* get data from symbolic products */ 141cf742fccSHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 14252f7967eSHong Zhang if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 14352f7967eSHong Zhang 144cf742fccSHong Zhang api = ap->i; 145cf742fccSHong Zhang apj = ap->j; 146a3bb6f32SFande Kong ierr = ISLocalToGlobalMappingApply(ptap->ltog,api[AP_loc->rmap->n],apj,apj);CHKERRQ(ierr); 147cf742fccSHong Zhang for (i=0; i<am; i++) { 148cf742fccSHong Zhang /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 149cf742fccSHong Zhang apnz = api[i+1] - api[i]; 150b4e8d1b6SHong Zhang apa = ap->a + api[i]; 151580bdb30SBarry Smith ierr = PetscArrayzero(apa,apnz);CHKERRQ(ierr); 152b4e8d1b6SHong Zhang AProw_scalable(i,ad,ao,p_loc,p_oth,api,apj,apa); 153cf742fccSHong Zhang } 154a3bb6f32SFande Kong ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,api[AP_loc->rmap->n],apj,&nout,apj);CHKERRQ(ierr); 155a3bb6f32SFande 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); 156cf742fccSHong Zhang 157cf742fccSHong Zhang /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 158154d0d78SFande Kong /* Always use scalable version since we are in the MPI scalable version */ 159cf742fccSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Rd,AP_loc,ptap->C_loc);CHKERRQ(ierr); 160cf742fccSHong Zhang ierr = MatMatMultNumeric_SeqAIJ_SeqAIJ_Scalable(ptap->Ro,AP_loc,ptap->C_oth);CHKERRQ(ierr); 161cf742fccSHong Zhang 162cf742fccSHong Zhang C_loc = ptap->C_loc; 163cf742fccSHong Zhang C_oth = ptap->C_oth; 164cf742fccSHong Zhang 165cf742fccSHong Zhang /* add C_loc and Co to to C */ 166cf742fccSHong Zhang ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 167cf742fccSHong Zhang 168cf742fccSHong Zhang /* C_loc -> C */ 169cf742fccSHong Zhang cm = C_loc->rmap->N; 170cf742fccSHong Zhang c_seq = (Mat_SeqAIJ*)C_loc->data; 171cf742fccSHong Zhang cols = c_seq->j; 172cf742fccSHong Zhang vals = c_seq->a; 173a3bb6f32SFande Kong ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_seq->i[C_loc->rmap->n],c_seq->j,c_seq->j);CHKERRQ(ierr); 174edeac6deSandi selinger 175e9ede7d0Sandi selinger /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */ 176edeac6deSandi selinger /* when there are no off-processor parts. */ 1771de21080Sandi selinger /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */ 1781de21080Sandi selinger /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */ 1791de21080Sandi selinger /* a table, and other, more complex stuff has to be done. */ 180edeac6deSandi selinger if (C->assembled) { 181edeac6deSandi selinger C->was_assembled = PETSC_TRUE; 182edeac6deSandi selinger C->assembled = PETSC_FALSE; 183edeac6deSandi selinger } 184edeac6deSandi selinger if (C->was_assembled) { 185cf742fccSHong Zhang for (i=0; i<cm; i++) { 186cf742fccSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 187cf742fccSHong Zhang row = rstart + i; 188edeac6deSandi selinger ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 189cf742fccSHong Zhang cols += ncols; vals += ncols; 190cf742fccSHong Zhang } 191edeac6deSandi selinger } else { 192e9ede7d0Sandi selinger ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr); 193edeac6deSandi selinger } 194a3bb6f32SFande Kong ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_seq->i[C_loc->rmap->n],c_seq->j,&nout,c_seq->j);CHKERRQ(ierr); 195a3bb6f32SFande 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); 196cf742fccSHong Zhang 197cf742fccSHong Zhang /* Co -> C, off-processor part */ 198cf742fccSHong Zhang cm = C_oth->rmap->N; 199cf742fccSHong Zhang c_seq = (Mat_SeqAIJ*)C_oth->data; 200cf742fccSHong Zhang cols = c_seq->j; 201cf742fccSHong Zhang vals = c_seq->a; 202a3bb6f32SFande Kong ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_seq->i[C_oth->rmap->n],c_seq->j,c_seq->j);CHKERRQ(ierr); 203cf742fccSHong Zhang for (i=0; i<cm; i++) { 204cf742fccSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 205cf742fccSHong Zhang row = p->garray[i]; 206cf742fccSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 207cf742fccSHong Zhang cols += ncols; vals += ncols; 208cf742fccSHong Zhang } 209cf742fccSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 210cf742fccSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 211cf742fccSHong Zhang 212cf742fccSHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 21380da8df7SHong Zhang 214a3bb6f32SFande Kong ierr = ISGlobalToLocalMappingApply(ptap->ltog,IS_GTOLM_DROP,c_seq->i[C_oth->rmap->n],c_seq->j,&nout,c_seq->j);CHKERRQ(ierr); 215a3bb6f32SFande 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); 216cf742fccSHong Zhang PetscFunctionReturn(0); 217cf742fccSHong Zhang } 218cf742fccSHong Zhang 2194222ddf1SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(Mat A,Mat P,PetscReal fill,Mat Cmpi) 2200d3441aeSHong Zhang { 2210d3441aeSHong Zhang PetscErrorCode ierr; 2223cdca5ebSHong Zhang Mat_APMPI *ptap; 2236718818eSStefano Zampini Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 2242259aa2eSHong Zhang MPI_Comm comm; 2252259aa2eSHong Zhang PetscMPIInt size,rank; 2264222ddf1SHong Zhang Mat P_loc,P_oth; 22715a3b8e2SHong Zhang PetscFreeSpaceList free_space=NULL,current_space=NULL; 228d0e9a2f7SHong Zhang PetscInt am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n; 229d0e9a2f7SHong Zhang PetscInt *lnk,i,k,pnz,row,nsend; 230ec4bef21SJose E. Roman PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,nrecv; 231ec4bef21SJose E. Roman PETSC_UNUSED PetscMPIInt icompleted=0; 23215a3b8e2SHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 233d0e9a2f7SHong Zhang PetscInt len,proc,*dnz,*onz,*owners,nzi,nspacedouble; 23415a3b8e2SHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 23515a3b8e2SHong Zhang MPI_Request *swaits,*rwaits; 23615a3b8e2SHong Zhang MPI_Status *sstatus,rstatus; 237445158ffSHong Zhang PetscLayout rowmap; 238445158ffSHong Zhang PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 239445158ffSHong Zhang PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 240a3bb6f32SFande Kong PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,Crmax,*aj,*ai,*pi,nout; 24152f7967eSHong Zhang Mat_SeqAIJ *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth; 24267a12041SHong Zhang PetscScalar *apv; 24378d30b94SHong Zhang PetscTable ta; 244b92f168fSBarry Smith MatType mtype; 245e83e5f86SFande Kong const char *prefix; 246aa690a28SHong Zhang #if defined(PETSC_USE_INFO) 2478cb82516SHong Zhang PetscReal apfill; 248aa690a28SHong Zhang #endif 24967a12041SHong Zhang 25067a12041SHong Zhang PetscFunctionBegin; 2516718818eSStefano Zampini MatCheckProduct(Cmpi,4); 2526718818eSStefano Zampini if (Cmpi->product->data) SETERRQ(PetscObjectComm((PetscObject)Cmpi),PETSC_ERR_PLIB,"Product data not empty"); 25367a12041SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 254ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 255ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 256ae5f0867Sstefano_zampini 25752f7967eSHong Zhang if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data; 25852f7967eSHong Zhang 259ae5f0867Sstefano_zampini /* create symbolic parallel matrix Cmpi */ 260b92f168fSBarry Smith ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 261b92f168fSBarry Smith ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 262ae5f0867Sstefano_zampini 2633cdca5ebSHong Zhang /* create struct Mat_APMPI and attached it to C later */ 264e953e47cSHong Zhang ierr = PetscNew(&ptap);CHKERRQ(ierr); 265e953e47cSHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 266cf97cf3cSHong Zhang ptap->algType = 0; 267e953e47cSHong Zhang 268e953e47cSHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 26976db11f6SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&P_oth);CHKERRQ(ierr); 270e953e47cSHong Zhang /* get P_loc by taking all local rows of P */ 27176db11f6SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&P_loc);CHKERRQ(ierr); 27276db11f6SHong Zhang 27376db11f6SHong Zhang ptap->P_loc = P_loc; 27476db11f6SHong Zhang ptap->P_oth = P_oth; 275e953e47cSHong Zhang 276e953e47cSHong Zhang /* (0) compute Rd = Pd^T, Ro = Po^T */ 277e953e47cSHong Zhang /* --------------------------------- */ 278419ecdd9Sandi selinger ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 279419ecdd9Sandi selinger ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 280e953e47cSHong Zhang 281e953e47cSHong Zhang /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 282e953e47cSHong Zhang /* ----------------------------------------------------------------- */ 28376db11f6SHong Zhang p_loc = (Mat_SeqAIJ*)P_loc->data; 28452f7967eSHong Zhang if (P_oth) p_oth = (Mat_SeqAIJ*)P_oth->data; 285e953e47cSHong Zhang 286e953e47cSHong Zhang /* create and initialize a linked list */ 287e953e47cSHong Zhang ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */ 28876db11f6SHong Zhang MatRowMergeMax_SeqAIJ(p_loc,P_loc->rmap->N,ta); 28976db11f6SHong Zhang MatRowMergeMax_SeqAIJ(p_oth,P_oth->rmap->N,ta); 290e953e47cSHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */ 291e953e47cSHong Zhang 29276db11f6SHong Zhang ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 293e953e47cSHong Zhang 294e953e47cSHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 29552f7967eSHong Zhang if (ao) { 296e953e47cSHong Zhang ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr); 29752f7967eSHong Zhang } else { 29852f7967eSHong Zhang ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr); 29952f7967eSHong Zhang } 300e953e47cSHong Zhang current_space = free_space; 301e953e47cSHong Zhang nspacedouble = 0; 302e953e47cSHong Zhang 303e953e47cSHong Zhang ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 304e953e47cSHong Zhang api[0] = 0; 305e953e47cSHong Zhang for (i=0; i<am; i++) { 306e953e47cSHong Zhang /* diagonal portion: Ad[i,:]*P */ 307e953e47cSHong Zhang ai = ad->i; pi = p_loc->i; 308e953e47cSHong Zhang nzi = ai[i+1] - ai[i]; 309e953e47cSHong Zhang aj = ad->j + ai[i]; 310e953e47cSHong Zhang for (j=0; j<nzi; j++) { 311e953e47cSHong Zhang row = aj[j]; 312e953e47cSHong Zhang pnz = pi[row+1] - pi[row]; 313e953e47cSHong Zhang Jptr = p_loc->j + pi[row]; 314e953e47cSHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 31576db11f6SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 316e953e47cSHong Zhang } 317e953e47cSHong Zhang /* off-diagonal portion: Ao[i,:]*P */ 31852f7967eSHong Zhang if (ao) { 319e953e47cSHong Zhang ai = ao->i; pi = p_oth->i; 320e953e47cSHong Zhang nzi = ai[i+1] - ai[i]; 321e953e47cSHong Zhang aj = ao->j + ai[i]; 322e953e47cSHong Zhang for (j=0; j<nzi; j++) { 323e953e47cSHong Zhang row = aj[j]; 324e953e47cSHong Zhang pnz = pi[row+1] - pi[row]; 325e953e47cSHong Zhang Jptr = p_oth->j + pi[row]; 32676db11f6SHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(pnz,Jptr,lnk);CHKERRQ(ierr); 327e953e47cSHong Zhang } 32852f7967eSHong Zhang } 329e953e47cSHong Zhang apnz = lnk[0]; 330e953e47cSHong Zhang api[i+1] = api[i] + apnz; 331e953e47cSHong Zhang 332e953e47cSHong Zhang /* if free space is not available, double the total space in the list */ 333e953e47cSHong Zhang if (current_space->local_remaining<apnz) { 334e953e47cSHong Zhang ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 335e953e47cSHong Zhang nspacedouble++; 336e953e47cSHong Zhang } 337e953e47cSHong Zhang 338e953e47cSHong Zhang /* Copy data into free space, then initialize lnk */ 33976db11f6SHong Zhang ierr = PetscLLCondensedClean_Scalable(apnz,current_space->array,lnk);CHKERRQ(ierr); 340e953e47cSHong Zhang 341e953e47cSHong Zhang current_space->array += apnz; 342e953e47cSHong Zhang current_space->local_used += apnz; 343e953e47cSHong Zhang current_space->local_remaining -= apnz; 344e953e47cSHong Zhang } 345e953e47cSHong Zhang /* Allocate space for apj and apv, initialize apj, and */ 346e953e47cSHong Zhang /* destroy list of free space and other temporary array(s) */ 347a3bb6f32SFande Kong ierr = PetscCalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 348e953e47cSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 34976db11f6SHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 350e953e47cSHong Zhang 351e953e47cSHong Zhang /* Create AP_loc for reuse */ 352e953e47cSHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 353a3bb6f32SFande Kong ierr = MatSeqAIJCompactOutExtraColumns_SeqAIJ(ptap->AP_loc, &ptap->ltog);CHKERRQ(ierr); 354e953e47cSHong Zhang 355e953e47cSHong Zhang #if defined(PETSC_USE_INFO) 35652f7967eSHong Zhang if (ao) { 357e953e47cSHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 35852f7967eSHong Zhang } else { 35952f7967eSHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1); 36052f7967eSHong Zhang } 361e953e47cSHong Zhang ptap->AP_loc->info.mallocs = nspacedouble; 362e953e47cSHong Zhang ptap->AP_loc->info.fill_ratio_given = fill; 363e953e47cSHong Zhang ptap->AP_loc->info.fill_ratio_needed = apfill; 364e953e47cSHong Zhang 365e953e47cSHong Zhang if (api[am]) { 366b11c1ec8SHong 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); 367e953e47cSHong Zhang ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 368e953e47cSHong Zhang } else { 369b11c1ec8SHong Zhang ierr = PetscInfo(ptap->AP_loc,"Scalable algorithm, AP_loc is empty \n");CHKERRQ(ierr); 370e953e47cSHong Zhang } 371e953e47cSHong Zhang #endif 372e953e47cSHong Zhang 373e953e47cSHong Zhang /* (2-1) compute symbolic Co = Ro*AP_loc */ 3744222ddf1SHong Zhang /* -------------------------------------- */ 3754222ddf1SHong Zhang ierr = MatProductCreate(ptap->Ro,ptap->AP_loc,NULL,&ptap->C_oth);CHKERRQ(ierr); 376e83e5f86SFande Kong ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 3774222ddf1SHong Zhang ierr = MatSetOptionsPrefix(ptap->C_oth,prefix);CHKERRQ(ierr); 3784222ddf1SHong Zhang ierr = MatAppendOptionsPrefix(ptap->C_oth,"inner_offdiag_");CHKERRQ(ierr); 3794222ddf1SHong Zhang 3804222ddf1SHong Zhang ierr = MatProductSetType(ptap->C_oth,MATPRODUCT_AB);CHKERRQ(ierr); 3814222ddf1SHong Zhang ierr = MatProductSetAlgorithm(ptap->C_oth,"sorted");CHKERRQ(ierr); 3824222ddf1SHong Zhang ierr = MatProductSetFill(ptap->C_oth,fill);CHKERRQ(ierr); 3834222ddf1SHong Zhang ierr = MatProductSetFromOptions(ptap->C_oth);CHKERRQ(ierr); 3844222ddf1SHong Zhang ierr = MatProductSymbolic(ptap->C_oth);CHKERRQ(ierr); 385e953e47cSHong Zhang 386e953e47cSHong Zhang /* (3) send coj of C_oth to other processors */ 387e953e47cSHong Zhang /* ------------------------------------------ */ 388e953e47cSHong Zhang /* determine row ownership */ 389e953e47cSHong Zhang ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 390e953e47cSHong Zhang rowmap->n = pn; 391e953e47cSHong Zhang rowmap->bs = 1; 392e953e47cSHong Zhang ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 393e953e47cSHong Zhang owners = rowmap->range; 394e953e47cSHong Zhang 395e953e47cSHong Zhang /* determine the number of messages to send, their lengths */ 396e953e47cSHong Zhang ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 397580bdb30SBarry Smith ierr = PetscArrayzero(len_s,size);CHKERRQ(ierr); 398580bdb30SBarry Smith ierr = PetscArrayzero(len_si,size);CHKERRQ(ierr); 399e953e47cSHong Zhang 400e953e47cSHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 401e953e47cSHong Zhang coi = c_oth->i; coj = c_oth->j; 402e953e47cSHong Zhang con = ptap->C_oth->rmap->n; 403e953e47cSHong Zhang proc = 0; 404a3bb6f32SFande Kong ierr = ISLocalToGlobalMappingApply(ptap->ltog,coi[con],coj,coj);CHKERRQ(ierr); 405e953e47cSHong Zhang for (i=0; i<con; i++) { 406e953e47cSHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 407e953e47cSHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 408e953e47cSHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 409e953e47cSHong Zhang } 410e953e47cSHong Zhang 411e953e47cSHong Zhang len = 0; /* max length of buf_si[], see (4) */ 412e953e47cSHong Zhang owners_co[0] = 0; 413e953e47cSHong Zhang nsend = 0; 414e953e47cSHong Zhang for (proc=0; proc<size; proc++) { 415e953e47cSHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 416e953e47cSHong Zhang if (len_s[proc]) { 417e953e47cSHong Zhang nsend++; 418e953e47cSHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 419e953e47cSHong Zhang len += len_si[proc]; 420e953e47cSHong Zhang } 421e953e47cSHong Zhang } 422e953e47cSHong Zhang 423e953e47cSHong Zhang /* determine the number and length of messages to receive for coi and coj */ 424e953e47cSHong Zhang ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 425e953e47cSHong Zhang ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 426e953e47cSHong Zhang 427e953e47cSHong Zhang /* post the Irecv and Isend of coj */ 428e953e47cSHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 429e953e47cSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 430e953e47cSHong Zhang ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 431e953e47cSHong Zhang for (proc=0, k=0; proc<size; proc++) { 432e953e47cSHong Zhang if (!len_s[proc]) continue; 433e953e47cSHong Zhang i = owners_co[proc]; 434ffc4695bSBarry Smith ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRMPI(ierr); 435e953e47cSHong Zhang k++; 436e953e47cSHong Zhang } 437e953e47cSHong Zhang 438e953e47cSHong Zhang /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 439e953e47cSHong Zhang /* ---------------------------------------- */ 4404222ddf1SHong Zhang ierr = MatProductCreate(ptap->Rd,ptap->AP_loc,NULL,&ptap->C_loc);CHKERRQ(ierr); 4414222ddf1SHong Zhang ierr = MatProductSetType(ptap->C_loc,MATPRODUCT_AB);CHKERRQ(ierr); 4424222ddf1SHong Zhang ierr = MatProductSetAlgorithm(ptap->C_loc,"default");CHKERRQ(ierr); 4434222ddf1SHong Zhang ierr = MatProductSetFill(ptap->C_loc,fill);CHKERRQ(ierr); 4444222ddf1SHong Zhang 4454222ddf1SHong Zhang ierr = MatSetOptionsPrefix(ptap->C_loc,prefix);CHKERRQ(ierr); 4464222ddf1SHong Zhang ierr = MatAppendOptionsPrefix(ptap->C_loc,"inner_diag_");CHKERRQ(ierr); 4474222ddf1SHong Zhang 4484222ddf1SHong Zhang ierr = MatProductSetFromOptions(ptap->C_loc);CHKERRQ(ierr); 4494222ddf1SHong Zhang ierr = MatProductSymbolic(ptap->C_loc);CHKERRQ(ierr); 4504222ddf1SHong Zhang 451e953e47cSHong Zhang c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 452a3bb6f32SFande Kong ierr = ISLocalToGlobalMappingApply(ptap->ltog,c_loc->i[ptap->C_loc->rmap->n],c_loc->j,c_loc->j);CHKERRQ(ierr); 453e953e47cSHong Zhang 454e953e47cSHong Zhang /* receives coj are complete */ 455e953e47cSHong Zhang for (i=0; i<nrecv; i++) { 456ffc4695bSBarry Smith ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRMPI(ierr); 457e953e47cSHong Zhang } 458e953e47cSHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 459ffc4695bSBarry Smith if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRMPI(ierr);} 460e953e47cSHong Zhang 461e953e47cSHong Zhang /* add received column indices into ta to update Crmax */ 462e953e47cSHong Zhang for (k=0; k<nrecv; k++) {/* k-th received message */ 463e953e47cSHong Zhang Jptr = buf_rj[k]; 464e953e47cSHong Zhang for (j=0; j<len_r[k]; j++) { 465e953e47cSHong Zhang ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 466e953e47cSHong Zhang } 467e953e47cSHong Zhang } 468e953e47cSHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 469e953e47cSHong Zhang ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 470e953e47cSHong Zhang 471e953e47cSHong Zhang /* (4) send and recv coi */ 472e953e47cSHong Zhang /*-----------------------*/ 473e953e47cSHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 474e953e47cSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 475e953e47cSHong Zhang ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 476e953e47cSHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 477e953e47cSHong Zhang for (proc=0,k=0; proc<size; proc++) { 478e953e47cSHong Zhang if (!len_s[proc]) continue; 479e953e47cSHong Zhang /* form outgoing message for i-structure: 480e953e47cSHong Zhang buf_si[0]: nrows to be sent 481e953e47cSHong Zhang [1:nrows]: row index (global) 482e953e47cSHong Zhang [nrows+1:2*nrows+1]: i-structure index 483e953e47cSHong Zhang */ 484e953e47cSHong Zhang /*-------------------------------------------*/ 485e953e47cSHong Zhang nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 486e953e47cSHong Zhang buf_si_i = buf_si + nrows+1; 487e953e47cSHong Zhang buf_si[0] = nrows; 488e953e47cSHong Zhang buf_si_i[0] = 0; 489e953e47cSHong Zhang nrows = 0; 490e953e47cSHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 491e953e47cSHong Zhang nzi = coi[i+1] - coi[i]; 492e953e47cSHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 493e953e47cSHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 494e953e47cSHong Zhang nrows++; 495e953e47cSHong Zhang } 496ffc4695bSBarry Smith ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRMPI(ierr); 497e953e47cSHong Zhang k++; 498e953e47cSHong Zhang buf_si += len_si[proc]; 499e953e47cSHong Zhang } 500e953e47cSHong Zhang for (i=0; i<nrecv; i++) { 501ffc4695bSBarry Smith ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRMPI(ierr); 502e953e47cSHong Zhang } 503e953e47cSHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 504ffc4695bSBarry Smith if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRMPI(ierr);} 505e953e47cSHong Zhang 506e953e47cSHong Zhang ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 507e953e47cSHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 508e953e47cSHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 509e953e47cSHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 510b4e8d1b6SHong Zhang 511e953e47cSHong Zhang /* (5) compute the local portion of Cmpi */ 512e953e47cSHong Zhang /* ------------------------------------------ */ 513e953e47cSHong Zhang /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 514e953e47cSHong Zhang ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 515e953e47cSHong Zhang current_space = free_space; 516e953e47cSHong Zhang 517e953e47cSHong Zhang ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 518e953e47cSHong Zhang for (k=0; k<nrecv; k++) { 519e953e47cSHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 520e953e47cSHong Zhang nrows = *buf_ri_k[k]; 521e953e47cSHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 522e953e47cSHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 523e953e47cSHong Zhang } 524e953e47cSHong Zhang 525e953e47cSHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 526cf742fccSHong Zhang ierr = PetscLLCondensedCreate_Scalable(Crmax,&lnk);CHKERRQ(ierr); 527e953e47cSHong Zhang for (i=0; i<pn; i++) { 528e953e47cSHong Zhang /* add C_loc into Cmpi */ 529e953e47cSHong Zhang nzi = c_loc->i[i+1] - c_loc->i[i]; 530e953e47cSHong Zhang Jptr = c_loc->j + c_loc->i[i]; 531cf742fccSHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 532e953e47cSHong Zhang 533e953e47cSHong Zhang /* add received col data into lnk */ 534e953e47cSHong Zhang for (k=0; k<nrecv; k++) { /* k-th received message */ 535e953e47cSHong Zhang if (i == *nextrow[k]) { /* i-th row */ 536e953e47cSHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 537e953e47cSHong Zhang Jptr = buf_rj[k] + *nextci[k]; 538cf742fccSHong Zhang ierr = PetscLLCondensedAddSorted_Scalable(nzi,Jptr,lnk);CHKERRQ(ierr); 539e953e47cSHong Zhang nextrow[k]++; nextci[k]++; 540e953e47cSHong Zhang } 541e953e47cSHong Zhang } 542e953e47cSHong Zhang nzi = lnk[0]; 543e953e47cSHong Zhang 544e953e47cSHong Zhang /* copy data into free space, then initialize lnk */ 545cf742fccSHong Zhang ierr = PetscLLCondensedClean_Scalable(nzi,current_space->array,lnk);CHKERRQ(ierr); 546e953e47cSHong Zhang ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 547e953e47cSHong Zhang } 548e953e47cSHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 549cf742fccSHong Zhang ierr = PetscLLCondensedDestroy_Scalable(lnk);CHKERRQ(ierr); 550e953e47cSHong Zhang ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 551e953e47cSHong Zhang 552e953e47cSHong Zhang /* local sizes and preallocation */ 553e953e47cSHong Zhang ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 554ac94c67aSTristan Konolige if (P->cmap->bs > 0) { 555ac94c67aSTristan Konolige ierr = PetscLayoutSetBlockSize(Cmpi->rmap,P->cmap->bs);CHKERRQ(ierr); 556ac94c67aSTristan Konolige ierr = PetscLayoutSetBlockSize(Cmpi->cmap,P->cmap->bs);CHKERRQ(ierr); 557ac94c67aSTristan Konolige } 558e953e47cSHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 559e953e47cSHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 560e953e47cSHong Zhang 561e953e47cSHong Zhang /* members in merge */ 562e953e47cSHong Zhang ierr = PetscFree(id_r);CHKERRQ(ierr); 563e953e47cSHong Zhang ierr = PetscFree(len_r);CHKERRQ(ierr); 564e953e47cSHong Zhang ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 565e953e47cSHong Zhang ierr = PetscFree(buf_ri);CHKERRQ(ierr); 566e953e47cSHong Zhang ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 567e953e47cSHong Zhang ierr = PetscFree(buf_rj);CHKERRQ(ierr); 568e953e47cSHong Zhang ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 569e953e47cSHong Zhang 570a3bb6f32SFande Kong nout = 0; 571a3bb6f32SFande 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); 572a3bb6f32SFande 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); 573a3bb6f32SFande 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); 574a3bb6f32SFande 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); 575a3bb6f32SFande Kong 5766718818eSStefano Zampini /* attach the supporting struct to Cmpi for reuse */ 5776718818eSStefano Zampini Cmpi->product->data = ptap; 5786718818eSStefano Zampini Cmpi->product->view = MatView_MPIAIJ_PtAP; 5796718818eSStefano Zampini Cmpi->product->destroy = MatDestroy_MPIAIJ_PtAP; 5806718818eSStefano Zampini 5816718818eSStefano Zampini /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 5826718818eSStefano Zampini Cmpi->assembled = PETSC_FALSE; 5836718818eSStefano Zampini Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_scalable; 584e953e47cSHong Zhang PetscFunctionReturn(0); 585e953e47cSHong Zhang } 586e953e47cSHong Zhang 587bc8e477aSFande Kong PETSC_STATIC_INLINE PetscErrorCode MatPtAPSymbolicComputeOneRowOfAP_private(Mat A,Mat P,Mat P_oth,const PetscInt *map,PetscInt dof,PetscInt i,PetscHSetI dht,PetscHSetI oht) 5884a56b808SFande Kong { 5894a56b808SFande Kong Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 5904a56b808SFande 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; 5914a56b808SFande Kong PetscInt *ai,nzi,j,*aj,row,col,*pi,*pj,pnz,nzpi,*p_othcols,k; 592bc8e477aSFande Kong PetscInt pcstart,pcend,column,offset; 5934a56b808SFande Kong PetscErrorCode ierr; 5944a56b808SFande Kong 5954a56b808SFande Kong PetscFunctionBegin; 5964a56b808SFande Kong pcstart = P->cmap->rstart; 597bc8e477aSFande Kong pcstart *= dof; 5984a56b808SFande Kong pcend = P->cmap->rend; 599bc8e477aSFande Kong pcend *= dof; 6004a56b808SFande Kong /* diagonal portion: Ad[i,:]*P */ 6014a56b808SFande Kong ai = ad->i; 6024a56b808SFande Kong nzi = ai[i+1] - ai[i]; 6034a56b808SFande Kong aj = ad->j + ai[i]; 6044a56b808SFande Kong for (j=0; j<nzi; j++) { 6054a56b808SFande Kong row = aj[j]; 606bc8e477aSFande Kong offset = row%dof; 607bc8e477aSFande Kong row /= dof; 6084a56b808SFande Kong nzpi = pd->i[row+1] - pd->i[row]; 6094a56b808SFande Kong pj = pd->j + pd->i[row]; 6104a56b808SFande Kong for (k=0; k<nzpi; k++) { 611bc8e477aSFande Kong ierr = PetscHSetIAdd(dht,pj[k]*dof+offset+pcstart);CHKERRQ(ierr); 6124a56b808SFande Kong } 6134a56b808SFande Kong } 614bc8e477aSFande Kong /* off diag P */ 6154a56b808SFande Kong for (j=0; j<nzi; j++) { 6164a56b808SFande Kong row = aj[j]; 617bc8e477aSFande Kong offset = row%dof; 618bc8e477aSFande Kong row /= dof; 6194a56b808SFande Kong nzpi = po->i[row+1] - po->i[row]; 6204a56b808SFande Kong pj = po->j + po->i[row]; 6214a56b808SFande Kong for (k=0; k<nzpi; k++) { 622bc8e477aSFande Kong ierr = PetscHSetIAdd(oht,p->garray[pj[k]]*dof+offset);CHKERRQ(ierr); 6234a56b808SFande Kong } 6244a56b808SFande Kong } 6254a56b808SFande Kong 6264a56b808SFande Kong /* off diagonal part: Ao[i, :]*P_oth */ 6274a56b808SFande Kong if (ao) { 6284a56b808SFande Kong ai = ao->i; 6294a56b808SFande Kong pi = p_oth->i; 6304a56b808SFande Kong nzi = ai[i+1] - ai[i]; 6314a56b808SFande Kong aj = ao->j + ai[i]; 6324a56b808SFande Kong for (j=0; j<nzi; j++) { 6334a56b808SFande Kong row = aj[j]; 634bc8e477aSFande Kong offset = a->garray[row]%dof; 635bc8e477aSFande Kong row = map[row]; 6364a56b808SFande Kong pnz = pi[row+1] - pi[row]; 6374a56b808SFande Kong p_othcols = p_oth->j + pi[row]; 6384a56b808SFande Kong for (col=0; col<pnz; col++) { 639bc8e477aSFande Kong column = p_othcols[col] * dof + offset; 6404a56b808SFande Kong if (column>=pcstart && column<pcend) { 6414a56b808SFande Kong ierr = PetscHSetIAdd(dht,column);CHKERRQ(ierr); 6424a56b808SFande Kong } else { 6434a56b808SFande Kong ierr = PetscHSetIAdd(oht,column);CHKERRQ(ierr); 6444a56b808SFande Kong } 6454a56b808SFande Kong } 6464a56b808SFande Kong } 6474a56b808SFande Kong } /* end if (ao) */ 6484a56b808SFande Kong PetscFunctionReturn(0); 6494a56b808SFande Kong } 6504a56b808SFande Kong 651bc8e477aSFande Kong PETSC_STATIC_INLINE PetscErrorCode MatPtAPNumericComputeOneRowOfAP_private(Mat A,Mat P,Mat P_oth,const PetscInt *map,PetscInt dof,PetscInt i,PetscHMapIV hmap) 6524a56b808SFande Kong { 6534a56b808SFande Kong Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 6544a56b808SFande 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; 655bc8e477aSFande Kong PetscInt *ai,nzi,j,*aj,row,col,*pi,pnz,*p_othcols,pcstart,*pj,k,nzpi,offset; 6564a56b808SFande Kong PetscScalar ra,*aa,*pa; 6574a56b808SFande Kong PetscErrorCode ierr; 6584a56b808SFande Kong 6594a56b808SFande Kong PetscFunctionBegin; 6604a56b808SFande Kong pcstart = P->cmap->rstart; 661bc8e477aSFande Kong pcstart *= dof; 662bc8e477aSFande Kong 6634a56b808SFande Kong /* diagonal portion: Ad[i,:]*P */ 6644a56b808SFande Kong ai = ad->i; 6654a56b808SFande Kong nzi = ai[i+1] - ai[i]; 6664a56b808SFande Kong aj = ad->j + ai[i]; 6674a56b808SFande Kong aa = ad->a + ai[i]; 6684a56b808SFande Kong for (j=0; j<nzi; j++) { 6694a56b808SFande Kong ra = aa[j]; 6704a56b808SFande Kong row = aj[j]; 671bc8e477aSFande Kong offset = row%dof; 672bc8e477aSFande Kong row /= dof; 6734a56b808SFande Kong nzpi = pd->i[row+1] - pd->i[row]; 6744a56b808SFande Kong pj = pd->j + pd->i[row]; 6754a56b808SFande Kong pa = pd->a + pd->i[row]; 6764a56b808SFande Kong for (k=0; k<nzpi; k++) { 677bc8e477aSFande Kong ierr = PetscHMapIVAddValue(hmap,pj[k]*dof+offset+pcstart,ra*pa[k]);CHKERRQ(ierr); 6784a56b808SFande Kong } 67949c8f2b8SFande Kong ierr = PetscLogFlops(2.0*nzpi);CHKERRQ(ierr); 6804a56b808SFande Kong } 6814a56b808SFande Kong for (j=0; j<nzi; j++) { 6824a56b808SFande Kong ra = aa[j]; 6834a56b808SFande Kong row = aj[j]; 684bc8e477aSFande Kong offset = row%dof; 685bc8e477aSFande Kong row /= dof; 6864a56b808SFande Kong nzpi = po->i[row+1] - po->i[row]; 6874a56b808SFande Kong pj = po->j + po->i[row]; 6884a56b808SFande Kong pa = po->a + po->i[row]; 6894a56b808SFande Kong for (k=0; k<nzpi; k++) { 690bc8e477aSFande Kong ierr = PetscHMapIVAddValue(hmap,p->garray[pj[k]]*dof+offset,ra*pa[k]);CHKERRQ(ierr); 6914a56b808SFande Kong } 69249c8f2b8SFande Kong ierr = PetscLogFlops(2.0*nzpi);CHKERRQ(ierr); 6934a56b808SFande Kong } 6944a56b808SFande Kong 6954a56b808SFande Kong /* off diagonal part: Ao[i, :]*P_oth */ 6964a56b808SFande Kong if (ao) { 6974a56b808SFande Kong ai = ao->i; 6984a56b808SFande Kong pi = p_oth->i; 6994a56b808SFande Kong nzi = ai[i+1] - ai[i]; 7004a56b808SFande Kong aj = ao->j + ai[i]; 7014a56b808SFande Kong aa = ao->a + ai[i]; 7024a56b808SFande Kong for (j=0; j<nzi; j++) { 7034a56b808SFande Kong row = aj[j]; 704bc8e477aSFande Kong offset = a->garray[row]%dof; 705bc8e477aSFande Kong row = map[row]; 7064a56b808SFande Kong ra = aa[j]; 7074a56b808SFande Kong pnz = pi[row+1] - pi[row]; 7084a56b808SFande Kong p_othcols = p_oth->j + pi[row]; 7094a56b808SFande Kong pa = p_oth->a + pi[row]; 7104a56b808SFande Kong for (col=0; col<pnz; col++) { 711bc8e477aSFande Kong ierr = PetscHMapIVAddValue(hmap,p_othcols[col]*dof+offset,ra*pa[col]);CHKERRQ(ierr); 7124a56b808SFande Kong } 71349c8f2b8SFande Kong ierr = PetscLogFlops(2.0*pnz);CHKERRQ(ierr); 7144a56b808SFande Kong } 7154a56b808SFande Kong } /* end if (ao) */ 716bb5ddf68SFande Kong 7174a56b808SFande Kong PetscFunctionReturn(0); 7184a56b808SFande Kong } 7194a56b808SFande Kong 720bc8e477aSFande Kong PetscErrorCode MatGetBrowsOfAcols_MPIXAIJ(Mat,Mat,PetscInt dof,MatReuse,Mat*); 7215c65b9ecSFande Kong 722bc8e477aSFande Kong PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIXAIJ_allatonce(Mat A,Mat P,PetscInt dof,Mat C) 7234a56b808SFande Kong { 7244a56b808SFande Kong PetscErrorCode ierr; 7254a56b808SFande Kong Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 726bc8e477aSFande Kong Mat_SeqAIJ *cd,*co,*po=(Mat_SeqAIJ*)p->B->data,*pd=(Mat_SeqAIJ*)p->A->data; 7276718818eSStefano Zampini Mat_APMPI *ptap; 7284a56b808SFande Kong PetscHMapIV hmap; 7294a56b808SFande 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; 7304a56b808SFande Kong PetscScalar *c_rmta,*c_otha,*poa,*pda,*apvalues,*apvaluestmp,*c_rmtaa; 731bc8e477aSFande Kong PetscInt offset,ii,pocol; 732bc8e477aSFande Kong const PetscInt *mappingindices; 733bc8e477aSFande Kong IS map; 7344a56b808SFande Kong 7354a56b808SFande Kong PetscFunctionBegin; 7366718818eSStefano Zampini MatCheckProduct(C,4); 7376718818eSStefano Zampini ptap = (Mat_APMPI*)C->product->data; 7386718818eSStefano Zampini if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data"); 7396718818eSStefano Zampini if (!ptap->P_oth) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatProductClear()"); 7404a56b808SFande Kong 7414a56b808SFande Kong ierr = MatZeroEntries(C);CHKERRQ(ierr); 7424a56b808SFande Kong 7434a56b808SFande Kong /* Get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 7444a56b808SFande Kong /*-----------------------------------------------------*/ 7454a56b808SFande Kong if (ptap->reuse == MAT_REUSE_MATRIX) { 7464a56b808SFande Kong /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 747bc8e477aSFande Kong ierr = MatGetBrowsOfAcols_MPIXAIJ(A,P,dof,MAT_REUSE_MATRIX,&ptap->P_oth);CHKERRQ(ierr); 7484a56b808SFande Kong } 749bc8e477aSFande Kong ierr = PetscObjectQuery((PetscObject)ptap->P_oth,"aoffdiagtopothmapping",(PetscObject*)&map);CHKERRQ(ierr); 7504a56b808SFande Kong 7514a56b808SFande Kong ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr); 752bc8e477aSFande Kong pon *= dof; 7534a56b808SFande Kong ierr = PetscCalloc2(ptap->c_rmti[pon],&c_rmtj,ptap->c_rmti[pon],&c_rmta);CHKERRQ(ierr); 7544a56b808SFande Kong ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr); 7554a56b808SFande Kong cmaxr = 0; 7564a56b808SFande Kong for (i=0; i<pon; i++) { 7574a56b808SFande Kong cmaxr = PetscMax(cmaxr,ptap->c_rmti[i+1]-ptap->c_rmti[i]); 7584a56b808SFande Kong } 7594a56b808SFande Kong ierr = PetscCalloc4(cmaxr,&apindices,cmaxr,&apvalues,cmaxr,&apvaluestmp,pon,&c_rmtc);CHKERRQ(ierr); 7604a56b808SFande Kong ierr = PetscHMapIVCreate(&hmap);CHKERRQ(ierr); 7614a56b808SFande Kong ierr = PetscHMapIVResize(hmap,cmaxr);CHKERRQ(ierr); 762bc8e477aSFande Kong ierr = ISGetIndices(map,&mappingindices);CHKERRQ(ierr); 7634a56b808SFande Kong for (i=0; i<am && pon; i++) { 7644a56b808SFande Kong ierr = PetscHMapIVClear(hmap);CHKERRQ(ierr); 765bc8e477aSFande Kong offset = i%dof; 766bc8e477aSFande Kong ii = i/dof; 767bc8e477aSFande Kong nzi = po->i[ii+1] - po->i[ii]; 7684a56b808SFande Kong if (!nzi) continue; 769bc8e477aSFande Kong ierr = MatPtAPNumericComputeOneRowOfAP_private(A,P,ptap->P_oth,mappingindices,dof,i,hmap);CHKERRQ(ierr); 7704a56b808SFande Kong voff = 0; 7714a56b808SFande Kong ierr = PetscHMapIVGetPairs(hmap,&voff,apindices,apvalues);CHKERRQ(ierr); 7724a56b808SFande Kong if (!voff) continue; 7734a56b808SFande Kong 7744a56b808SFande Kong /* Form C(ii, :) */ 775bc8e477aSFande Kong poj = po->j + po->i[ii]; 776bc8e477aSFande Kong poa = po->a + po->i[ii]; 7774a56b808SFande Kong for (j=0; j<nzi; j++) { 778bc8e477aSFande Kong pocol = poj[j]*dof+offset; 779bc8e477aSFande Kong c_rmtjj = c_rmtj + ptap->c_rmti[pocol]; 780bc8e477aSFande Kong c_rmtaa = c_rmta + ptap->c_rmti[pocol]; 7814a56b808SFande Kong for (jj=0; jj<voff; jj++) { 7824a56b808SFande Kong apvaluestmp[jj] = apvalues[jj]*poa[j]; 7834a56b808SFande Kong /*If the row is empty */ 784bc8e477aSFande Kong if (!c_rmtc[pocol]) { 7854a56b808SFande Kong c_rmtjj[jj] = apindices[jj]; 7864a56b808SFande Kong c_rmtaa[jj] = apvaluestmp[jj]; 787bc8e477aSFande Kong c_rmtc[pocol]++; 7884a56b808SFande Kong } else { 789bc8e477aSFande Kong ierr = PetscFindInt(apindices[jj],c_rmtc[pocol],c_rmtjj,&loc);CHKERRQ(ierr); 7904a56b808SFande Kong if (loc>=0){ /* hit */ 7914a56b808SFande Kong c_rmtaa[loc] += apvaluestmp[jj]; 79249c8f2b8SFande Kong ierr = PetscLogFlops(1.0);CHKERRQ(ierr); 7934a56b808SFande Kong } else { /* new element */ 7944a56b808SFande Kong loc = -(loc+1); 7954a56b808SFande Kong /* Move data backward */ 796bc8e477aSFande Kong for (kk=c_rmtc[pocol]; kk>loc; kk--) { 7974a56b808SFande Kong c_rmtjj[kk] = c_rmtjj[kk-1]; 7984a56b808SFande Kong c_rmtaa[kk] = c_rmtaa[kk-1]; 7994a56b808SFande Kong }/* End kk */ 8004a56b808SFande Kong c_rmtjj[loc] = apindices[jj]; 8014a56b808SFande Kong c_rmtaa[loc] = apvaluestmp[jj]; 802bc8e477aSFande Kong c_rmtc[pocol]++; 8034a56b808SFande Kong } 8044a56b808SFande Kong } 80549c8f2b8SFande Kong ierr = PetscLogFlops(voff);CHKERRQ(ierr); 8064a56b808SFande Kong } /* End jj */ 8074a56b808SFande Kong } /* End j */ 8084a56b808SFande Kong } /* End i */ 8094a56b808SFande Kong 8104a56b808SFande Kong ierr = PetscFree4(apindices,apvalues,apvaluestmp,c_rmtc);CHKERRQ(ierr); 811d4e5d74dSLawrence Mitchell ierr = PetscHMapIVDestroy(&hmap);CHKERRQ(ierr); 8124a56b808SFande Kong 8134a56b808SFande Kong ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr); 814bc8e477aSFande Kong pn *= dof; 8154a56b808SFande Kong ierr = PetscCalloc2(ptap->c_othi[pn],&c_othj,ptap->c_othi[pn],&c_otha);CHKERRQ(ierr); 8164a56b808SFande Kong 81783df288dSJunchao Zhang ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPI_REPLACE);CHKERRQ(ierr); 81883df288dSJunchao Zhang ierr = PetscSFReduceBegin(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPI_REPLACE);CHKERRQ(ierr); 8194a56b808SFande Kong ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr); 820bc8e477aSFande Kong pcstart = pcstart*dof; 821bc8e477aSFande Kong pcend = pcend*dof; 8224a56b808SFande Kong cd = (Mat_SeqAIJ*)(c->A)->data; 8234a56b808SFande Kong co = (Mat_SeqAIJ*)(c->B)->data; 8244a56b808SFande Kong 8254a56b808SFande Kong cmaxr = 0; 8264a56b808SFande Kong for (i=0; i<pn; i++) { 8274a56b808SFande Kong cmaxr = PetscMax(cmaxr,(cd->i[i+1]-cd->i[i])+(co->i[i+1]-co->i[i])); 8284a56b808SFande Kong } 8294a56b808SFande Kong ierr = PetscCalloc5(cmaxr,&apindices,cmaxr,&apvalues,cmaxr,&apvaluestmp,pn,&dcc,pn,&occ);CHKERRQ(ierr); 8304a56b808SFande Kong ierr = PetscHMapIVCreate(&hmap);CHKERRQ(ierr); 8314a56b808SFande Kong ierr = PetscHMapIVResize(hmap,cmaxr);CHKERRQ(ierr); 8324a56b808SFande Kong for (i=0; i<am && pn; i++) { 8334a56b808SFande Kong ierr = PetscHMapIVClear(hmap);CHKERRQ(ierr); 834bc8e477aSFande Kong offset = i%dof; 835bc8e477aSFande Kong ii = i/dof; 836bc8e477aSFande Kong nzi = pd->i[ii+1] - pd->i[ii]; 8374a56b808SFande Kong if (!nzi) continue; 838bc8e477aSFande Kong ierr = MatPtAPNumericComputeOneRowOfAP_private(A,P,ptap->P_oth,mappingindices,dof,i,hmap);CHKERRQ(ierr); 8394a56b808SFande Kong voff = 0; 8404a56b808SFande Kong ierr = PetscHMapIVGetPairs(hmap,&voff,apindices,apvalues);CHKERRQ(ierr); 8414a56b808SFande Kong if (!voff) continue; 8424a56b808SFande Kong /* Form C(ii, :) */ 843bc8e477aSFande Kong pdj = pd->j + pd->i[ii]; 844bc8e477aSFande Kong pda = pd->a + pd->i[ii]; 8454a56b808SFande Kong for (j=0; j<nzi; j++) { 846bc8e477aSFande Kong row = pcstart + pdj[j] * dof + offset; 8474a56b808SFande Kong for (jj=0; jj<voff; jj++) { 8484a56b808SFande Kong apvaluestmp[jj] = apvalues[jj]*pda[j]; 8494a56b808SFande Kong } 85049c8f2b8SFande Kong ierr = PetscLogFlops(voff);CHKERRQ(ierr); 8514a56b808SFande Kong ierr = MatSetValues(C,1,&row,voff,apindices,apvaluestmp,ADD_VALUES);CHKERRQ(ierr); 8524a56b808SFande Kong } 8534a56b808SFande Kong } 854bc8e477aSFande Kong ierr = ISRestoreIndices(map,&mappingindices);CHKERRQ(ierr); 8554a56b808SFande Kong ierr = MatGetOwnershipRangeColumn(C,&ccstart,&ccend);CHKERRQ(ierr); 8564a56b808SFande Kong ierr = PetscFree5(apindices,apvalues,apvaluestmp,dcc,occ);CHKERRQ(ierr); 8574a56b808SFande Kong ierr = PetscHMapIVDestroy(&hmap);CHKERRQ(ierr); 85883df288dSJunchao Zhang ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPI_REPLACE);CHKERRQ(ierr); 85983df288dSJunchao Zhang ierr = PetscSFReduceEnd(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPI_REPLACE);CHKERRQ(ierr); 860bc8e477aSFande Kong ierr = PetscFree2(c_rmtj,c_rmta);CHKERRQ(ierr); 861bc8e477aSFande Kong 862bc8e477aSFande Kong /* Add contributions from remote */ 863bc8e477aSFande Kong for (i = 0; i < pn; i++) { 864bc8e477aSFande Kong row = i + pcstart; 865bc8e477aSFande 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); 866bc8e477aSFande Kong } 867bc8e477aSFande Kong ierr = PetscFree2(c_othj,c_otha);CHKERRQ(ierr); 868bc8e477aSFande Kong 869bc8e477aSFande Kong ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 870bc8e477aSFande Kong ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 871bc8e477aSFande Kong 872bc8e477aSFande Kong ptap->reuse = MAT_REUSE_MATRIX; 873bc8e477aSFande Kong PetscFunctionReturn(0); 874bc8e477aSFande Kong } 875bc8e477aSFande Kong 876bc8e477aSFande Kong PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce(Mat A,Mat P,Mat C) 877bc8e477aSFande Kong { 878bc8e477aSFande Kong PetscErrorCode ierr; 879bc8e477aSFande Kong 880bc8e477aSFande Kong PetscFunctionBegin; 88134bcad68SFande Kong 882bc8e477aSFande Kong ierr = MatPtAPNumeric_MPIAIJ_MPIXAIJ_allatonce(A,P,1,C);CHKERRQ(ierr); 883bc8e477aSFande Kong PetscFunctionReturn(0); 884bc8e477aSFande Kong } 885bc8e477aSFande Kong 886bc8e477aSFande Kong PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIXAIJ_allatonce_merged(Mat A,Mat P,PetscInt dof,Mat C) 887bc8e477aSFande Kong { 888bc8e477aSFande Kong PetscErrorCode ierr; 889bc8e477aSFande Kong Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data,*c=(Mat_MPIAIJ*)C->data; 890bc8e477aSFande Kong Mat_SeqAIJ *cd,*co,*po=(Mat_SeqAIJ*)p->B->data,*pd=(Mat_SeqAIJ*)p->A->data; 8916718818eSStefano Zampini Mat_APMPI *ptap; 892bc8e477aSFande Kong PetscHMapIV hmap; 893bc8e477aSFande 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; 894bc8e477aSFande Kong PetscScalar *c_rmta,*c_otha,*poa,*pda,*apvalues,*apvaluestmp,*c_rmtaa; 895bc8e477aSFande Kong PetscInt offset,ii,pocol; 896bc8e477aSFande Kong const PetscInt *mappingindices; 897bc8e477aSFande Kong IS map; 898bc8e477aSFande Kong 899bc8e477aSFande Kong PetscFunctionBegin; 9006718818eSStefano Zampini MatCheckProduct(C,4); 9016718818eSStefano Zampini ptap = (Mat_APMPI*)C->product->data; 9026718818eSStefano Zampini if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data"); 9036718818eSStefano Zampini if (!ptap->P_oth) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatProductClear()"); 904bc8e477aSFande Kong 905bc8e477aSFande Kong ierr = MatZeroEntries(C);CHKERRQ(ierr); 906bc8e477aSFande Kong 907bc8e477aSFande Kong /* Get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 908bc8e477aSFande Kong /*-----------------------------------------------------*/ 909bc8e477aSFande Kong if (ptap->reuse == MAT_REUSE_MATRIX) { 910bc8e477aSFande Kong /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 911bc8e477aSFande Kong ierr = MatGetBrowsOfAcols_MPIXAIJ(A,P,dof,MAT_REUSE_MATRIX,&ptap->P_oth);CHKERRQ(ierr); 912bc8e477aSFande Kong } 913bc8e477aSFande Kong ierr = PetscObjectQuery((PetscObject)ptap->P_oth,"aoffdiagtopothmapping",(PetscObject*)&map);CHKERRQ(ierr); 914bc8e477aSFande Kong ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr); 915bc8e477aSFande Kong pon *= dof; 916bc8e477aSFande Kong ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr); 917bc8e477aSFande Kong pn *= dof; 918bc8e477aSFande Kong 919bc8e477aSFande Kong ierr = PetscCalloc2(ptap->c_rmti[pon],&c_rmtj,ptap->c_rmti[pon],&c_rmta);CHKERRQ(ierr); 920bc8e477aSFande Kong ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr); 921bc8e477aSFande Kong ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr); 922bc8e477aSFande Kong pcstart *= dof; 923bc8e477aSFande Kong pcend *= dof; 924bc8e477aSFande Kong cmaxr = 0; 925bc8e477aSFande Kong for (i=0; i<pon; i++) { 926bc8e477aSFande Kong cmaxr = PetscMax(cmaxr,ptap->c_rmti[i+1]-ptap->c_rmti[i]); 927bc8e477aSFande Kong } 928bc8e477aSFande Kong cd = (Mat_SeqAIJ*)(c->A)->data; 929bc8e477aSFande Kong co = (Mat_SeqAIJ*)(c->B)->data; 930bc8e477aSFande Kong for (i=0; i<pn; i++) { 931bc8e477aSFande Kong cmaxr = PetscMax(cmaxr,(cd->i[i+1]-cd->i[i])+(co->i[i+1]-co->i[i])); 932bc8e477aSFande Kong } 933bc8e477aSFande Kong ierr = PetscCalloc4(cmaxr,&apindices,cmaxr,&apvalues,cmaxr,&apvaluestmp,pon,&c_rmtc);CHKERRQ(ierr); 934bc8e477aSFande Kong ierr = PetscHMapIVCreate(&hmap);CHKERRQ(ierr); 935bc8e477aSFande Kong ierr = PetscHMapIVResize(hmap,cmaxr);CHKERRQ(ierr); 936bc8e477aSFande Kong ierr = ISGetIndices(map,&mappingindices);CHKERRQ(ierr); 937bc8e477aSFande Kong for (i=0; i<am && (pon || pn); i++) { 938bc8e477aSFande Kong ierr = PetscHMapIVClear(hmap);CHKERRQ(ierr); 939bc8e477aSFande Kong offset = i%dof; 940bc8e477aSFande Kong ii = i/dof; 941bc8e477aSFande Kong nzi = po->i[ii+1] - po->i[ii]; 942bc8e477aSFande Kong dnzi = pd->i[ii+1] - pd->i[ii]; 943bc8e477aSFande Kong if (!nzi && !dnzi) continue; 944bc8e477aSFande Kong ierr = MatPtAPNumericComputeOneRowOfAP_private(A,P,ptap->P_oth,mappingindices,dof,i,hmap);CHKERRQ(ierr); 945bc8e477aSFande Kong voff = 0; 946bc8e477aSFande Kong ierr = PetscHMapIVGetPairs(hmap,&voff,apindices,apvalues);CHKERRQ(ierr); 947bc8e477aSFande Kong if (!voff) continue; 948bc8e477aSFande Kong 949bc8e477aSFande Kong /* Form remote C(ii, :) */ 950bc8e477aSFande Kong poj = po->j + po->i[ii]; 951bc8e477aSFande Kong poa = po->a + po->i[ii]; 952bc8e477aSFande Kong for (j=0; j<nzi; j++) { 953bc8e477aSFande Kong pocol = poj[j]*dof+offset; 954bc8e477aSFande Kong c_rmtjj = c_rmtj + ptap->c_rmti[pocol]; 955bc8e477aSFande Kong c_rmtaa = c_rmta + ptap->c_rmti[pocol]; 956bc8e477aSFande Kong for (jj=0; jj<voff; jj++) { 957bc8e477aSFande Kong apvaluestmp[jj] = apvalues[jj]*poa[j]; 958bc8e477aSFande Kong /*If the row is empty */ 959bc8e477aSFande Kong if (!c_rmtc[pocol]) { 960bc8e477aSFande Kong c_rmtjj[jj] = apindices[jj]; 961bc8e477aSFande Kong c_rmtaa[jj] = apvaluestmp[jj]; 962bc8e477aSFande Kong c_rmtc[pocol]++; 963bc8e477aSFande Kong } else { 964bc8e477aSFande Kong ierr = PetscFindInt(apindices[jj],c_rmtc[pocol],c_rmtjj,&loc);CHKERRQ(ierr); 965bc8e477aSFande Kong if (loc>=0){ /* hit */ 966bc8e477aSFande Kong c_rmtaa[loc] += apvaluestmp[jj]; 967bc8e477aSFande Kong ierr = PetscLogFlops(1.0);CHKERRQ(ierr); 968bc8e477aSFande Kong } else { /* new element */ 969bc8e477aSFande Kong loc = -(loc+1); 970bc8e477aSFande Kong /* Move data backward */ 971bc8e477aSFande Kong for (kk=c_rmtc[pocol]; kk>loc; kk--) { 972bc8e477aSFande Kong c_rmtjj[kk] = c_rmtjj[kk-1]; 973bc8e477aSFande Kong c_rmtaa[kk] = c_rmtaa[kk-1]; 974bc8e477aSFande Kong }/* End kk */ 975bc8e477aSFande Kong c_rmtjj[loc] = apindices[jj]; 976bc8e477aSFande Kong c_rmtaa[loc] = apvaluestmp[jj]; 977bc8e477aSFande Kong c_rmtc[pocol]++; 978bc8e477aSFande Kong } 979bc8e477aSFande Kong } 980bc8e477aSFande Kong } /* End jj */ 981bc8e477aSFande Kong ierr = PetscLogFlops(voff);CHKERRQ(ierr); 982bc8e477aSFande Kong } /* End j */ 983bc8e477aSFande Kong 984bc8e477aSFande Kong /* Form local C(ii, :) */ 985bc8e477aSFande Kong pdj = pd->j + pd->i[ii]; 986bc8e477aSFande Kong pda = pd->a + pd->i[ii]; 987bc8e477aSFande Kong for (j=0; j<dnzi; j++) { 988bc8e477aSFande Kong row = pcstart + pdj[j] * dof + offset; 989bc8e477aSFande Kong for (jj=0; jj<voff; jj++) { 990bc8e477aSFande Kong apvaluestmp[jj] = apvalues[jj]*pda[j]; 991bc8e477aSFande Kong }/* End kk */ 992bc8e477aSFande Kong ierr = PetscLogFlops(voff);CHKERRQ(ierr); 993bc8e477aSFande Kong ierr = MatSetValues(C,1,&row,voff,apindices,apvaluestmp,ADD_VALUES);CHKERRQ(ierr); 994bc8e477aSFande Kong }/* End j */ 995bc8e477aSFande Kong } /* End i */ 996bc8e477aSFande Kong 997bc8e477aSFande Kong ierr = ISRestoreIndices(map,&mappingindices);CHKERRQ(ierr); 998bc8e477aSFande Kong ierr = PetscFree4(apindices,apvalues,apvaluestmp,c_rmtc);CHKERRQ(ierr); 999bc8e477aSFande Kong ierr = PetscHMapIVDestroy(&hmap);CHKERRQ(ierr); 1000bc8e477aSFande Kong ierr = PetscCalloc2(ptap->c_othi[pn],&c_othj,ptap->c_othi[pn],&c_otha);CHKERRQ(ierr); 1001bc8e477aSFande Kong 100283df288dSJunchao Zhang ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPI_REPLACE);CHKERRQ(ierr); 100383df288dSJunchao Zhang ierr = PetscSFReduceBegin(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPI_REPLACE);CHKERRQ(ierr); 100483df288dSJunchao Zhang ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPI_REPLACE);CHKERRQ(ierr); 100583df288dSJunchao Zhang ierr = PetscSFReduceEnd(ptap->sf,MPIU_SCALAR,c_rmta,c_otha,MPI_REPLACE);CHKERRQ(ierr); 10064a56b808SFande Kong ierr = PetscFree2(c_rmtj,c_rmta);CHKERRQ(ierr); 10074a56b808SFande Kong 10084a56b808SFande Kong /* Add contributions from remote */ 10094a56b808SFande Kong for (i = 0; i < pn; i++) { 10104a56b808SFande Kong row = i + pcstart; 10114a56b808SFande 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); 10124a56b808SFande Kong } 10134a56b808SFande Kong ierr = PetscFree2(c_othj,c_otha);CHKERRQ(ierr); 10144a56b808SFande Kong 10154a56b808SFande Kong ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10164a56b808SFande Kong ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 10174a56b808SFande Kong 10184a56b808SFande Kong ptap->reuse = MAT_REUSE_MATRIX; 10194a56b808SFande Kong PetscFunctionReturn(0); 10204a56b808SFande Kong } 10214a56b808SFande Kong 10224a56b808SFande Kong PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged(Mat A,Mat P,Mat C) 10234a56b808SFande Kong { 10244a56b808SFande Kong PetscErrorCode ierr; 10254a56b808SFande Kong 10264a56b808SFande Kong PetscFunctionBegin; 102734bcad68SFande Kong 1028bc8e477aSFande Kong ierr = MatPtAPNumeric_MPIAIJ_MPIXAIJ_allatonce_merged(A,P,1,C);CHKERRQ(ierr); 10294a56b808SFande Kong PetscFunctionReturn(0); 10304a56b808SFande Kong } 10314a56b808SFande Kong 10326718818eSStefano Zampini /* TODO: move algorithm selection to MatProductSetFromOptions */ 10334222ddf1SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIXAIJ_allatonce(Mat A,Mat P,PetscInt dof,PetscReal fill,Mat Cmpi) 10344a56b808SFande Kong { 10354a56b808SFande Kong Mat_APMPI *ptap; 10366718818eSStefano Zampini Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data; 10374a56b808SFande Kong MPI_Comm comm; 1038bc8e477aSFande Kong Mat_SeqAIJ *pd=(Mat_SeqAIJ*)p->A->data,*po=(Mat_SeqAIJ*)p->B->data; 10394a56b808SFande Kong MatType mtype; 10404a56b808SFande Kong PetscSF sf; 10414a56b808SFande Kong PetscSFNode *iremote; 10424a56b808SFande Kong PetscInt rootspacesize,*rootspace,*rootspaceoffsets,nleaves; 10434a56b808SFande Kong const PetscInt *rootdegrees; 10444a56b808SFande Kong PetscHSetI ht,oht,*hta,*hto; 10454a56b808SFande Kong PetscInt pn,pon,*c_rmtc,i,j,nzi,htsize,htosize,*c_rmtj,off,*c_othj,rcvncols,sendncols,*c_rmtoffsets; 1046131c27b5Sprj- PetscInt lidx,*rdj,col,pcstart,pcend,*dnz,*onz,am,arstart,arend,*poj,*pdj; 1047bc8e477aSFande Kong PetscInt nalg=2,alg=0,offset,ii; 1048131c27b5Sprj- PetscMPIInt owner; 1049bc8e477aSFande Kong const PetscInt *mappingindices; 10504a56b808SFande Kong PetscBool flg; 10514a56b808SFande Kong const char *algTypes[2] = {"overlapping","merged"}; 1052bc8e477aSFande Kong IS map; 10534a56b808SFande Kong PetscErrorCode ierr; 10544a56b808SFande Kong 10554a56b808SFande Kong PetscFunctionBegin; 10566718818eSStefano Zampini MatCheckProduct(Cmpi,5); 10576718818eSStefano Zampini if (Cmpi->product->data) SETERRQ(PetscObjectComm((PetscObject)Cmpi),PETSC_ERR_PLIB,"Product data not empty"); 10584a56b808SFande Kong ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 105934bcad68SFande Kong 10604a56b808SFande Kong /* Create symbolic parallel matrix Cmpi */ 10614a56b808SFande Kong ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr); 1062bc8e477aSFande Kong pn *= dof; 10634a56b808SFande Kong ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 10644a56b808SFande Kong ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 10654a56b808SFande Kong ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 10664a56b808SFande Kong 10674a56b808SFande Kong ierr = PetscNew(&ptap);CHKERRQ(ierr); 10684a56b808SFande Kong ptap->reuse = MAT_INITIAL_MATRIX; 10694a56b808SFande Kong ptap->algType = 2; 10704a56b808SFande Kong 10714a56b808SFande Kong /* Get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 1072bc8e477aSFande Kong ierr = MatGetBrowsOfAcols_MPIXAIJ(A,P,dof,MAT_INITIAL_MATRIX,&ptap->P_oth);CHKERRQ(ierr); 1073bc8e477aSFande Kong ierr = PetscObjectQuery((PetscObject)ptap->P_oth,"aoffdiagtopothmapping",(PetscObject*)&map);CHKERRQ(ierr); 10744a56b808SFande Kong /* This equals to the number of offdiag columns in P */ 10754a56b808SFande Kong ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr); 1076bc8e477aSFande Kong pon *= dof; 10774a56b808SFande Kong /* offsets */ 10784a56b808SFande Kong ierr = PetscMalloc1(pon+1,&ptap->c_rmti);CHKERRQ(ierr); 10794a56b808SFande Kong /* The number of columns we will send to remote ranks */ 10804a56b808SFande Kong ierr = PetscMalloc1(pon,&c_rmtc);CHKERRQ(ierr); 10814a56b808SFande Kong ierr = PetscMalloc1(pon,&hta);CHKERRQ(ierr); 10824a56b808SFande Kong for (i=0; i<pon; i++) { 10834a56b808SFande Kong ierr = PetscHSetICreate(&hta[i]);CHKERRQ(ierr); 10844a56b808SFande Kong } 10854a56b808SFande Kong ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr); 10864a56b808SFande Kong ierr = MatGetOwnershipRange(A,&arstart,&arend);CHKERRQ(ierr); 10874a56b808SFande Kong /* Create hash table to merge all columns for C(i, :) */ 10884a56b808SFande Kong ierr = PetscHSetICreate(&ht);CHKERRQ(ierr); 10894a56b808SFande Kong 1090bc8e477aSFande Kong ierr = ISGetIndices(map,&mappingindices);CHKERRQ(ierr); 10914a56b808SFande Kong ptap->c_rmti[0] = 0; 10924a56b808SFande Kong /* 2) Pass 1: calculate the size for C_rmt (a matrix need to be sent to other processors) */ 10934a56b808SFande Kong for (i=0; i<am && pon; i++) { 10944a56b808SFande Kong /* Form one row of AP */ 10954a56b808SFande Kong ierr = PetscHSetIClear(ht);CHKERRQ(ierr); 1096bc8e477aSFande Kong offset = i%dof; 1097bc8e477aSFande Kong ii = i/dof; 10984a56b808SFande Kong /* If the off diag is empty, we should not do any calculation */ 1099bc8e477aSFande Kong nzi = po->i[ii+1] - po->i[ii]; 11004a56b808SFande Kong if (!nzi) continue; 11014a56b808SFande Kong 1102bc8e477aSFande Kong ierr = MatPtAPSymbolicComputeOneRowOfAP_private(A,P,ptap->P_oth,mappingindices,dof,i,ht,ht);CHKERRQ(ierr); 11034a56b808SFande Kong ierr = PetscHSetIGetSize(ht,&htsize);CHKERRQ(ierr); 11044a56b808SFande Kong /* If AP is empty, just continue */ 11054a56b808SFande Kong if (!htsize) continue; 11064a56b808SFande Kong /* Form C(ii, :) */ 1107bc8e477aSFande Kong poj = po->j + po->i[ii]; 11084a56b808SFande Kong for (j=0; j<nzi; j++) { 1109bc8e477aSFande Kong ierr = PetscHSetIUpdate(hta[poj[j]*dof+offset],ht);CHKERRQ(ierr); 11104a56b808SFande Kong } 11114a56b808SFande Kong } 11124a56b808SFande Kong 11134a56b808SFande Kong for (i=0; i<pon; i++) { 11144a56b808SFande Kong ierr = PetscHSetIGetSize(hta[i],&htsize);CHKERRQ(ierr); 11154a56b808SFande Kong ptap->c_rmti[i+1] = ptap->c_rmti[i] + htsize; 11164a56b808SFande Kong c_rmtc[i] = htsize; 11174a56b808SFande Kong } 11184a56b808SFande Kong 11194a56b808SFande Kong ierr = PetscMalloc1(ptap->c_rmti[pon],&c_rmtj);CHKERRQ(ierr); 11204a56b808SFande Kong 11214a56b808SFande Kong for (i=0; i<pon; i++) { 11224a56b808SFande Kong off = 0; 11234a56b808SFande Kong ierr = PetscHSetIGetElems(hta[i],&off,c_rmtj+ptap->c_rmti[i]);CHKERRQ(ierr); 11244a56b808SFande Kong ierr = PetscHSetIDestroy(&hta[i]);CHKERRQ(ierr); 11254a56b808SFande Kong } 11264a56b808SFande Kong ierr = PetscFree(hta);CHKERRQ(ierr); 11274a56b808SFande Kong 1128071fcb05SBarry Smith ierr = PetscMalloc1(pon,&iremote);CHKERRQ(ierr); 11294a56b808SFande Kong for (i=0; i<pon; i++) { 1130ef7a94f2SFande Kong owner = 0; lidx = 0; 1131bc8e477aSFande Kong offset = i%dof; 1132bc8e477aSFande Kong ii = i/dof; 1133bc8e477aSFande Kong ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[ii],&owner,&lidx);CHKERRQ(ierr); 1134bc8e477aSFande Kong iremote[i].index = lidx*dof + offset; 11354a56b808SFande Kong iremote[i].rank = owner; 11364a56b808SFande Kong } 11374a56b808SFande Kong 11384a56b808SFande Kong ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr); 11394a56b808SFande Kong ierr = PetscSFSetGraph(sf,pn,pon,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 11404a56b808SFande Kong /* Reorder ranks properly so that the data handled by gather and scatter have the same order */ 11414a56b808SFande Kong ierr = PetscSFSetRankOrder(sf,PETSC_TRUE);CHKERRQ(ierr); 11424a56b808SFande Kong ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 11434a56b808SFande Kong ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 11444a56b808SFande Kong /* How many neighbors have contributions to my rows? */ 11454a56b808SFande Kong ierr = PetscSFComputeDegreeBegin(sf,&rootdegrees);CHKERRQ(ierr); 11464a56b808SFande Kong ierr = PetscSFComputeDegreeEnd(sf,&rootdegrees);CHKERRQ(ierr); 11474a56b808SFande Kong rootspacesize = 0; 11484a56b808SFande Kong for (i = 0; i < pn; i++) { 11494a56b808SFande Kong rootspacesize += rootdegrees[i]; 11504a56b808SFande Kong } 1151071fcb05SBarry Smith ierr = PetscMalloc1(rootspacesize,&rootspace);CHKERRQ(ierr); 1152071fcb05SBarry Smith ierr = PetscMalloc1(rootspacesize+1,&rootspaceoffsets);CHKERRQ(ierr); 11534a56b808SFande Kong /* Get information from leaves 11544a56b808SFande Kong * Number of columns other people contribute to my rows 11554a56b808SFande Kong * */ 11564a56b808SFande Kong ierr = PetscSFGatherBegin(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr); 11574a56b808SFande Kong ierr = PetscSFGatherEnd(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr); 11584a56b808SFande Kong ierr = PetscFree(c_rmtc);CHKERRQ(ierr); 11594a56b808SFande Kong ierr = PetscCalloc1(pn+1,&ptap->c_othi);CHKERRQ(ierr); 11604a56b808SFande Kong /* The number of columns is received for each row */ 11614a56b808SFande Kong ptap->c_othi[0] = 0; 11624a56b808SFande Kong rootspacesize = 0; 11634a56b808SFande Kong rootspaceoffsets[0] = 0; 11644a56b808SFande Kong for (i = 0; i < pn; i++) { 11654a56b808SFande Kong rcvncols = 0; 11664a56b808SFande Kong for (j = 0; j<rootdegrees[i]; j++) { 11674a56b808SFande Kong rcvncols += rootspace[rootspacesize]; 11684a56b808SFande Kong rootspaceoffsets[rootspacesize+1] = rootspaceoffsets[rootspacesize] + rootspace[rootspacesize]; 11694a56b808SFande Kong rootspacesize++; 11704a56b808SFande Kong } 11714a56b808SFande Kong ptap->c_othi[i+1] = ptap->c_othi[i] + rcvncols; 11724a56b808SFande Kong } 11734a56b808SFande Kong ierr = PetscFree(rootspace);CHKERRQ(ierr); 11744a56b808SFande Kong 1175071fcb05SBarry Smith ierr = PetscMalloc1(pon,&c_rmtoffsets);CHKERRQ(ierr); 11764a56b808SFande Kong ierr = PetscSFScatterBegin(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr); 11774a56b808SFande Kong ierr = PetscSFScatterEnd(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr); 11784a56b808SFande Kong ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 11794a56b808SFande Kong ierr = PetscFree(rootspaceoffsets);CHKERRQ(ierr); 11804a56b808SFande Kong 11814a56b808SFande Kong ierr = PetscCalloc1(ptap->c_rmti[pon],&iremote);CHKERRQ(ierr); 11824a56b808SFande Kong nleaves = 0; 11834a56b808SFande Kong for (i = 0; i<pon; i++) { 1184bc8e477aSFande Kong owner = 0; 1185bc8e477aSFande Kong ii = i/dof; 1186bc8e477aSFande Kong ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[ii],&owner,NULL);CHKERRQ(ierr); 11874a56b808SFande Kong sendncols = ptap->c_rmti[i+1] - ptap->c_rmti[i]; 11884a56b808SFande Kong for (j=0; j<sendncols; j++) { 11894a56b808SFande Kong iremote[nleaves].rank = owner; 11904a56b808SFande Kong iremote[nleaves++].index = c_rmtoffsets[i] + j; 11914a56b808SFande Kong } 11924a56b808SFande Kong } 11934a56b808SFande Kong ierr = PetscFree(c_rmtoffsets);CHKERRQ(ierr); 11944a56b808SFande Kong ierr = PetscCalloc1(ptap->c_othi[pn],&c_othj);CHKERRQ(ierr); 11954a56b808SFande Kong 11964a56b808SFande Kong ierr = PetscSFCreate(comm,&ptap->sf);CHKERRQ(ierr); 11974a56b808SFande Kong ierr = PetscSFSetGraph(ptap->sf,ptap->c_othi[pn],nleaves,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 11984a56b808SFande Kong ierr = PetscSFSetFromOptions(ptap->sf);CHKERRQ(ierr); 11994a56b808SFande Kong /* One to one map */ 120083df288dSJunchao Zhang ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPI_REPLACE);CHKERRQ(ierr); 12014a56b808SFande Kong 1202071fcb05SBarry Smith ierr = PetscMalloc2(pn,&dnz,pn,&onz);CHKERRQ(ierr); 12034a56b808SFande Kong ierr = PetscHSetICreate(&oht);CHKERRQ(ierr); 12044a56b808SFande Kong ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr); 1205bc8e477aSFande Kong pcstart *= dof; 1206bc8e477aSFande Kong pcend *= dof; 12074a56b808SFande Kong ierr = PetscMalloc2(pn,&hta,pn,&hto);CHKERRQ(ierr); 12084a56b808SFande Kong for (i=0; i<pn; i++) { 12094a56b808SFande Kong ierr = PetscHSetICreate(&hta[i]);CHKERRQ(ierr); 12104a56b808SFande Kong ierr = PetscHSetICreate(&hto[i]);CHKERRQ(ierr); 12114a56b808SFande Kong } 12124a56b808SFande Kong /* Work on local part */ 12134a56b808SFande Kong /* 4) Pass 1: Estimate memory for C_loc */ 12144a56b808SFande Kong for (i=0; i<am && pn; i++) { 12154a56b808SFande Kong ierr = PetscHSetIClear(ht);CHKERRQ(ierr); 12164a56b808SFande Kong ierr = PetscHSetIClear(oht);CHKERRQ(ierr); 1217bc8e477aSFande Kong offset = i%dof; 1218bc8e477aSFande Kong ii = i/dof; 1219bc8e477aSFande Kong nzi = pd->i[ii+1] - pd->i[ii]; 12204a56b808SFande Kong if (!nzi) continue; 12214a56b808SFande Kong 1222bc8e477aSFande Kong ierr = MatPtAPSymbolicComputeOneRowOfAP_private(A,P,ptap->P_oth,mappingindices,dof,i,ht,oht);CHKERRQ(ierr); 12234a56b808SFande Kong ierr = PetscHSetIGetSize(ht,&htsize);CHKERRQ(ierr); 12244a56b808SFande Kong ierr = PetscHSetIGetSize(oht,&htosize);CHKERRQ(ierr); 12254a56b808SFande Kong if (!(htsize+htosize)) continue; 12264a56b808SFande Kong /* Form C(ii, :) */ 1227bc8e477aSFande Kong pdj = pd->j + pd->i[ii]; 12284a56b808SFande Kong for (j=0; j<nzi; j++) { 1229bc8e477aSFande Kong ierr = PetscHSetIUpdate(hta[pdj[j]*dof+offset],ht);CHKERRQ(ierr); 1230bc8e477aSFande Kong ierr = PetscHSetIUpdate(hto[pdj[j]*dof+offset],oht);CHKERRQ(ierr); 12314a56b808SFande Kong } 12324a56b808SFande Kong } 12334a56b808SFande Kong 1234bc8e477aSFande Kong ierr = ISRestoreIndices(map,&mappingindices);CHKERRQ(ierr); 1235bc8e477aSFande Kong 12364a56b808SFande Kong ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr); 12374a56b808SFande Kong ierr = PetscHSetIDestroy(&oht);CHKERRQ(ierr); 12384a56b808SFande Kong 12394a56b808SFande Kong /* Get remote data */ 124083df288dSJunchao Zhang ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPI_REPLACE);CHKERRQ(ierr); 12414a56b808SFande Kong ierr = PetscFree(c_rmtj);CHKERRQ(ierr); 12424a56b808SFande Kong 12434a56b808SFande Kong for (i = 0; i < pn; i++) { 12444a56b808SFande Kong nzi = ptap->c_othi[i+1] - ptap->c_othi[i]; 12454a56b808SFande Kong rdj = c_othj + ptap->c_othi[i]; 12464a56b808SFande Kong for (j = 0; j < nzi; j++) { 12474a56b808SFande Kong col = rdj[j]; 12484a56b808SFande Kong /* diag part */ 12494a56b808SFande Kong if (col>=pcstart && col<pcend) { 12504a56b808SFande Kong ierr = PetscHSetIAdd(hta[i],col);CHKERRQ(ierr); 12514a56b808SFande Kong } else { /* off diag */ 12524a56b808SFande Kong ierr = PetscHSetIAdd(hto[i],col);CHKERRQ(ierr); 12534a56b808SFande Kong } 12544a56b808SFande Kong } 12554a56b808SFande Kong ierr = PetscHSetIGetSize(hta[i],&htsize);CHKERRQ(ierr); 12564a56b808SFande Kong dnz[i] = htsize; 12574a56b808SFande Kong ierr = PetscHSetIDestroy(&hta[i]);CHKERRQ(ierr); 12584a56b808SFande Kong ierr = PetscHSetIGetSize(hto[i],&htsize);CHKERRQ(ierr); 12594a56b808SFande Kong onz[i] = htsize; 12604a56b808SFande Kong ierr = PetscHSetIDestroy(&hto[i]);CHKERRQ(ierr); 12614a56b808SFande Kong } 12624a56b808SFande Kong 12634a56b808SFande Kong ierr = PetscFree2(hta,hto);CHKERRQ(ierr); 12644a56b808SFande Kong ierr = PetscFree(c_othj);CHKERRQ(ierr); 12654a56b808SFande Kong 12664a56b808SFande Kong /* local sizes and preallocation */ 12674a56b808SFande Kong ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 126834bcad68SFande Kong ierr = MatSetBlockSizes(Cmpi,dof>1? dof: P->cmap->bs,dof>1? dof: P->cmap->bs);CHKERRQ(ierr); 12694a56b808SFande Kong ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 1270bc8e477aSFande Kong ierr = MatSetUp(Cmpi);CHKERRQ(ierr); 12714a56b808SFande Kong ierr = PetscFree2(dnz,onz);CHKERRQ(ierr); 12724a56b808SFande Kong 12734a56b808SFande Kong /* attach the supporting struct to Cmpi for reuse */ 12746718818eSStefano Zampini Cmpi->product->data = ptap; 12756718818eSStefano Zampini Cmpi->product->destroy = MatDestroy_MPIAIJ_PtAP; 12766718818eSStefano Zampini Cmpi->product->view = MatView_MPIAIJ_PtAP; 12774a56b808SFande Kong 12784a56b808SFande Kong /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 12794a56b808SFande Kong Cmpi->assembled = PETSC_FALSE; 12804a56b808SFande Kong /* pick an algorithm */ 12814a56b808SFande Kong ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatPtAP","Mat");CHKERRQ(ierr); 12824a56b808SFande Kong alg = 0; 12834a56b808SFande Kong ierr = PetscOptionsEList("-matptap_allatonce_via","PtAP allatonce numeric approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 12844a56b808SFande Kong ierr = PetscOptionsEnd();CHKERRQ(ierr); 12854a56b808SFande Kong switch (alg) { 12864a56b808SFande Kong case 0: 12874a56b808SFande Kong Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce; 12884a56b808SFande Kong break; 12894a56b808SFande Kong case 1: 12904a56b808SFande Kong Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged; 12914a56b808SFande Kong break; 12924a56b808SFande Kong default: 12934a56b808SFande Kong SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG," Unsupported allatonce numerical algorithm \n"); 12944a56b808SFande Kong } 12954a56b808SFande Kong PetscFunctionReturn(0); 12964a56b808SFande Kong } 12974a56b808SFande Kong 12984222ddf1SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce(Mat A,Mat P,PetscReal fill,Mat C) 1299bc8e477aSFande Kong { 1300bc8e477aSFande Kong PetscErrorCode ierr; 1301bc8e477aSFande Kong 1302bc8e477aSFande Kong PetscFunctionBegin; 1303bc8e477aSFande Kong ierr = MatPtAPSymbolic_MPIAIJ_MPIXAIJ_allatonce(A,P,1,fill,C);CHKERRQ(ierr); 1304bc8e477aSFande Kong PetscFunctionReturn(0); 1305bc8e477aSFande Kong } 1306bc8e477aSFande Kong 13074222ddf1SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIXAIJ_allatonce_merged(Mat A,Mat P,PetscInt dof,PetscReal fill,Mat Cmpi) 13084a56b808SFande Kong { 13094a56b808SFande Kong Mat_APMPI *ptap; 13106718818eSStefano Zampini Mat_MPIAIJ *p=(Mat_MPIAIJ*)P->data; 13114a56b808SFande Kong MPI_Comm comm; 1312bc8e477aSFande Kong Mat_SeqAIJ *pd=(Mat_SeqAIJ*)p->A->data,*po=(Mat_SeqAIJ*)p->B->data; 13134a56b808SFande Kong MatType mtype; 13144a56b808SFande Kong PetscSF sf; 13154a56b808SFande Kong PetscSFNode *iremote; 13164a56b808SFande Kong PetscInt rootspacesize,*rootspace,*rootspaceoffsets,nleaves; 13174a56b808SFande Kong const PetscInt *rootdegrees; 13184a56b808SFande Kong PetscHSetI ht,oht,*hta,*hto,*htd; 13194a56b808SFande Kong PetscInt pn,pon,*c_rmtc,i,j,nzi,dnzi,htsize,htosize,*c_rmtj,off,*c_othj,rcvncols,sendncols,*c_rmtoffsets; 1320131c27b5Sprj- PetscInt lidx,*rdj,col,pcstart,pcend,*dnz,*onz,am,arstart,arend,*poj,*pdj; 1321bc8e477aSFande Kong PetscInt nalg=2,alg=0,offset,ii; 1322131c27b5Sprj- PetscMPIInt owner; 13234a56b808SFande Kong PetscBool flg; 13244a56b808SFande Kong const char *algTypes[2] = {"merged","overlapping"}; 1325bc8e477aSFande Kong const PetscInt *mappingindices; 1326bc8e477aSFande Kong IS map; 13274a56b808SFande Kong PetscErrorCode ierr; 13284a56b808SFande Kong 13294a56b808SFande Kong PetscFunctionBegin; 13306718818eSStefano Zampini MatCheckProduct(Cmpi,5); 13316718818eSStefano Zampini if (Cmpi->product->data) SETERRQ(PetscObjectComm((PetscObject)Cmpi),PETSC_ERR_PLIB,"Product data not empty"); 13324a56b808SFande Kong ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 133334bcad68SFande Kong 13344a56b808SFande Kong /* Create symbolic parallel matrix Cmpi */ 13354a56b808SFande Kong ierr = MatGetLocalSize(P,NULL,&pn);CHKERRQ(ierr); 1336bc8e477aSFande Kong pn *= dof; 13374a56b808SFande Kong ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 13384a56b808SFande Kong ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 13394a56b808SFande Kong ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 13404a56b808SFande Kong 13414a56b808SFande Kong ierr = PetscNew(&ptap);CHKERRQ(ierr); 13424a56b808SFande Kong ptap->reuse = MAT_INITIAL_MATRIX; 13434a56b808SFande Kong ptap->algType = 3; 13444a56b808SFande Kong 13454a56b808SFande Kong /* 0) Get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 1346bc8e477aSFande Kong ierr = MatGetBrowsOfAcols_MPIXAIJ(A,P,dof,MAT_INITIAL_MATRIX,&ptap->P_oth);CHKERRQ(ierr); 1347bc8e477aSFande Kong ierr = PetscObjectQuery((PetscObject)ptap->P_oth,"aoffdiagtopothmapping",(PetscObject*)&map);CHKERRQ(ierr); 13484a56b808SFande Kong 13494a56b808SFande Kong /* This equals to the number of offdiag columns in P */ 13504a56b808SFande Kong ierr = MatGetLocalSize(p->B,NULL,&pon);CHKERRQ(ierr); 1351bc8e477aSFande Kong pon *= dof; 13524a56b808SFande Kong /* offsets */ 13534a56b808SFande Kong ierr = PetscMalloc1(pon+1,&ptap->c_rmti);CHKERRQ(ierr); 13544a56b808SFande Kong /* The number of columns we will send to remote ranks */ 13554a56b808SFande Kong ierr = PetscMalloc1(pon,&c_rmtc);CHKERRQ(ierr); 13564a56b808SFande Kong ierr = PetscMalloc1(pon,&hta);CHKERRQ(ierr); 13574a56b808SFande Kong for (i=0; i<pon; i++) { 13584a56b808SFande Kong ierr = PetscHSetICreate(&hta[i]);CHKERRQ(ierr); 13594a56b808SFande Kong } 13604a56b808SFande Kong ierr = MatGetLocalSize(A,&am,NULL);CHKERRQ(ierr); 13614a56b808SFande Kong ierr = MatGetOwnershipRange(A,&arstart,&arend);CHKERRQ(ierr); 13624a56b808SFande Kong /* Create hash table to merge all columns for C(i, :) */ 13634a56b808SFande Kong ierr = PetscHSetICreate(&ht);CHKERRQ(ierr); 13644a56b808SFande Kong ierr = PetscHSetICreate(&oht);CHKERRQ(ierr); 13654a56b808SFande Kong ierr = PetscMalloc2(pn,&htd,pn,&hto);CHKERRQ(ierr); 13664a56b808SFande Kong for (i=0; i<pn; i++) { 13674a56b808SFande Kong ierr = PetscHSetICreate(&htd[i]);CHKERRQ(ierr); 13684a56b808SFande Kong ierr = PetscHSetICreate(&hto[i]);CHKERRQ(ierr); 13694a56b808SFande Kong } 1370bc8e477aSFande Kong 1371bc8e477aSFande Kong ierr = ISGetIndices(map,&mappingindices);CHKERRQ(ierr); 13724a56b808SFande Kong ptap->c_rmti[0] = 0; 13734a56b808SFande Kong /* 2) Pass 1: calculate the size for C_rmt (a matrix need to be sent to other processors) */ 13744a56b808SFande Kong for (i=0; i<am && (pon || pn); i++) { 13754a56b808SFande Kong /* Form one row of AP */ 13764a56b808SFande Kong ierr = PetscHSetIClear(ht);CHKERRQ(ierr); 13774a56b808SFande Kong ierr = PetscHSetIClear(oht);CHKERRQ(ierr); 1378bc8e477aSFande Kong offset = i%dof; 1379bc8e477aSFande Kong ii = i/dof; 13804a56b808SFande Kong /* If the off diag is empty, we should not do any calculation */ 1381bc8e477aSFande Kong nzi = po->i[ii+1] - po->i[ii]; 1382bc8e477aSFande Kong dnzi = pd->i[ii+1] - pd->i[ii]; 13834a56b808SFande Kong if (!nzi && !dnzi) continue; 13844a56b808SFande Kong 1385bc8e477aSFande Kong ierr = MatPtAPSymbolicComputeOneRowOfAP_private(A,P,ptap->P_oth,mappingindices,dof,i,ht,oht);CHKERRQ(ierr); 13864a56b808SFande Kong ierr = PetscHSetIGetSize(ht,&htsize);CHKERRQ(ierr); 13874a56b808SFande Kong ierr = PetscHSetIGetSize(oht,&htosize);CHKERRQ(ierr); 13884a56b808SFande Kong /* If AP is empty, just continue */ 13894a56b808SFande Kong if (!(htsize+htosize)) continue; 13904a56b808SFande Kong 13914a56b808SFande Kong /* Form remote C(ii, :) */ 1392bc8e477aSFande Kong poj = po->j + po->i[ii]; 13934a56b808SFande Kong for (j=0; j<nzi; j++) { 1394bc8e477aSFande Kong ierr = PetscHSetIUpdate(hta[poj[j]*dof+offset],ht);CHKERRQ(ierr); 1395bc8e477aSFande Kong ierr = PetscHSetIUpdate(hta[poj[j]*dof+offset],oht);CHKERRQ(ierr); 13964a56b808SFande Kong } 13974a56b808SFande Kong 13984a56b808SFande Kong /* Form local C(ii, :) */ 1399bc8e477aSFande Kong pdj = pd->j + pd->i[ii]; 14004a56b808SFande Kong for (j=0; j<dnzi; j++) { 1401bc8e477aSFande Kong ierr = PetscHSetIUpdate(htd[pdj[j]*dof+offset],ht);CHKERRQ(ierr); 1402bc8e477aSFande Kong ierr = PetscHSetIUpdate(hto[pdj[j]*dof+offset],oht);CHKERRQ(ierr); 14034a56b808SFande Kong } 14044a56b808SFande Kong } 14054a56b808SFande Kong 1406bc8e477aSFande Kong ierr = ISRestoreIndices(map,&mappingindices);CHKERRQ(ierr); 1407bc8e477aSFande Kong 14084a56b808SFande Kong ierr = PetscHSetIDestroy(&ht);CHKERRQ(ierr); 14094a56b808SFande Kong ierr = PetscHSetIDestroy(&oht);CHKERRQ(ierr); 14104a56b808SFande Kong 14114a56b808SFande Kong for (i=0; i<pon; i++) { 14124a56b808SFande Kong ierr = PetscHSetIGetSize(hta[i],&htsize);CHKERRQ(ierr); 14134a56b808SFande Kong ptap->c_rmti[i+1] = ptap->c_rmti[i] + htsize; 14144a56b808SFande Kong c_rmtc[i] = htsize; 14154a56b808SFande Kong } 14164a56b808SFande Kong 14174a56b808SFande Kong ierr = PetscMalloc1(ptap->c_rmti[pon],&c_rmtj);CHKERRQ(ierr); 14184a56b808SFande Kong 14194a56b808SFande Kong for (i=0; i<pon; i++) { 14204a56b808SFande Kong off = 0; 14214a56b808SFande Kong ierr = PetscHSetIGetElems(hta[i],&off,c_rmtj+ptap->c_rmti[i]);CHKERRQ(ierr); 14224a56b808SFande Kong ierr = PetscHSetIDestroy(&hta[i]);CHKERRQ(ierr); 14234a56b808SFande Kong } 14244a56b808SFande Kong ierr = PetscFree(hta);CHKERRQ(ierr); 14254a56b808SFande Kong 1426071fcb05SBarry Smith ierr = PetscMalloc1(pon,&iremote);CHKERRQ(ierr); 14274a56b808SFande Kong for (i=0; i<pon; i++) { 1428ef7a94f2SFande Kong owner = 0; lidx = 0; 1429bc8e477aSFande Kong offset = i%dof; 1430bc8e477aSFande Kong ii = i/dof; 1431bc8e477aSFande Kong ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[ii],&owner,&lidx);CHKERRQ(ierr); 1432bc8e477aSFande Kong iremote[i].index = lidx*dof+offset; 14334a56b808SFande Kong iremote[i].rank = owner; 14344a56b808SFande Kong } 14354a56b808SFande Kong 14364a56b808SFande Kong ierr = PetscSFCreate(comm,&sf);CHKERRQ(ierr); 14374a56b808SFande Kong ierr = PetscSFSetGraph(sf,pn,pon,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 14384a56b808SFande Kong /* Reorder ranks properly so that the data handled by gather and scatter have the same order */ 14394a56b808SFande Kong ierr = PetscSFSetRankOrder(sf,PETSC_TRUE);CHKERRQ(ierr); 14404a56b808SFande Kong ierr = PetscSFSetFromOptions(sf);CHKERRQ(ierr); 14414a56b808SFande Kong ierr = PetscSFSetUp(sf);CHKERRQ(ierr); 14424a56b808SFande Kong /* How many neighbors have contributions to my rows? */ 14434a56b808SFande Kong ierr = PetscSFComputeDegreeBegin(sf,&rootdegrees);CHKERRQ(ierr); 14444a56b808SFande Kong ierr = PetscSFComputeDegreeEnd(sf,&rootdegrees);CHKERRQ(ierr); 14454a56b808SFande Kong rootspacesize = 0; 14464a56b808SFande Kong for (i = 0; i < pn; i++) { 14474a56b808SFande Kong rootspacesize += rootdegrees[i]; 14484a56b808SFande Kong } 1449071fcb05SBarry Smith ierr = PetscMalloc1(rootspacesize,&rootspace);CHKERRQ(ierr); 1450071fcb05SBarry Smith ierr = PetscMalloc1(rootspacesize+1,&rootspaceoffsets);CHKERRQ(ierr); 14514a56b808SFande Kong /* Get information from leaves 14524a56b808SFande Kong * Number of columns other people contribute to my rows 14534a56b808SFande Kong * */ 14544a56b808SFande Kong ierr = PetscSFGatherBegin(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr); 14554a56b808SFande Kong ierr = PetscSFGatherEnd(sf,MPIU_INT,c_rmtc,rootspace);CHKERRQ(ierr); 14564a56b808SFande Kong ierr = PetscFree(c_rmtc);CHKERRQ(ierr); 1457071fcb05SBarry Smith ierr = PetscMalloc1(pn+1,&ptap->c_othi);CHKERRQ(ierr); 14584a56b808SFande Kong /* The number of columns is received for each row */ 14594a56b808SFande Kong ptap->c_othi[0] = 0; 14604a56b808SFande Kong rootspacesize = 0; 14614a56b808SFande Kong rootspaceoffsets[0] = 0; 14624a56b808SFande Kong for (i = 0; i < pn; i++) { 14634a56b808SFande Kong rcvncols = 0; 14644a56b808SFande Kong for (j = 0; j<rootdegrees[i]; j++) { 14654a56b808SFande Kong rcvncols += rootspace[rootspacesize]; 14664a56b808SFande Kong rootspaceoffsets[rootspacesize+1] = rootspaceoffsets[rootspacesize] + rootspace[rootspacesize]; 14674a56b808SFande Kong rootspacesize++; 14684a56b808SFande Kong } 14694a56b808SFande Kong ptap->c_othi[i+1] = ptap->c_othi[i] + rcvncols; 14704a56b808SFande Kong } 14714a56b808SFande Kong ierr = PetscFree(rootspace);CHKERRQ(ierr); 14724a56b808SFande Kong 1473071fcb05SBarry Smith ierr = PetscMalloc1(pon,&c_rmtoffsets);CHKERRQ(ierr); 14744a56b808SFande Kong ierr = PetscSFScatterBegin(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr); 14754a56b808SFande Kong ierr = PetscSFScatterEnd(sf,MPIU_INT,rootspaceoffsets,c_rmtoffsets);CHKERRQ(ierr); 14764a56b808SFande Kong ierr = PetscSFDestroy(&sf);CHKERRQ(ierr); 14774a56b808SFande Kong ierr = PetscFree(rootspaceoffsets);CHKERRQ(ierr); 14784a56b808SFande Kong 14794a56b808SFande Kong ierr = PetscCalloc1(ptap->c_rmti[pon],&iremote);CHKERRQ(ierr); 14804a56b808SFande Kong nleaves = 0; 14814a56b808SFande Kong for (i = 0; i<pon; i++) { 1482071fcb05SBarry Smith owner = 0; 1483bc8e477aSFande Kong ii = i/dof; 1484bc8e477aSFande Kong ierr = PetscLayoutFindOwnerIndex(P->cmap,p->garray[ii],&owner,NULL);CHKERRQ(ierr); 14854a56b808SFande Kong sendncols = ptap->c_rmti[i+1] - ptap->c_rmti[i]; 14864a56b808SFande Kong for (j=0; j<sendncols; j++) { 14874a56b808SFande Kong iremote[nleaves].rank = owner; 14884a56b808SFande Kong iremote[nleaves++].index = c_rmtoffsets[i] + j; 14894a56b808SFande Kong } 14904a56b808SFande Kong } 14914a56b808SFande Kong ierr = PetscFree(c_rmtoffsets);CHKERRQ(ierr); 14924a56b808SFande Kong ierr = PetscCalloc1(ptap->c_othi[pn],&c_othj);CHKERRQ(ierr); 14934a56b808SFande Kong 14944a56b808SFande Kong ierr = PetscSFCreate(comm,&ptap->sf);CHKERRQ(ierr); 14954a56b808SFande Kong ierr = PetscSFSetGraph(ptap->sf,ptap->c_othi[pn],nleaves,NULL,PETSC_OWN_POINTER,iremote,PETSC_OWN_POINTER);CHKERRQ(ierr); 14964a56b808SFande Kong ierr = PetscSFSetFromOptions(ptap->sf);CHKERRQ(ierr); 14974a56b808SFande Kong /* One to one map */ 149883df288dSJunchao Zhang ierr = PetscSFReduceBegin(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPI_REPLACE);CHKERRQ(ierr); 14994a56b808SFande Kong /* Get remote data */ 150083df288dSJunchao Zhang ierr = PetscSFReduceEnd(ptap->sf,MPIU_INT,c_rmtj,c_othj,MPI_REPLACE);CHKERRQ(ierr); 15014a56b808SFande Kong ierr = PetscFree(c_rmtj);CHKERRQ(ierr); 1502071fcb05SBarry Smith ierr = PetscMalloc2(pn,&dnz,pn,&onz);CHKERRQ(ierr); 15034a56b808SFande Kong ierr = MatGetOwnershipRangeColumn(P,&pcstart,&pcend);CHKERRQ(ierr); 1504bc8e477aSFande Kong pcstart *= dof; 1505bc8e477aSFande Kong pcend *= dof; 15064a56b808SFande Kong for (i = 0; i < pn; i++) { 15074a56b808SFande Kong nzi = ptap->c_othi[i+1] - ptap->c_othi[i]; 15084a56b808SFande Kong rdj = c_othj + ptap->c_othi[i]; 15094a56b808SFande Kong for (j = 0; j < nzi; j++) { 15104a56b808SFande Kong col = rdj[j]; 15114a56b808SFande Kong /* diag part */ 15124a56b808SFande Kong if (col>=pcstart && col<pcend) { 15134a56b808SFande Kong ierr = PetscHSetIAdd(htd[i],col);CHKERRQ(ierr); 15144a56b808SFande Kong } else { /* off diag */ 15154a56b808SFande Kong ierr = PetscHSetIAdd(hto[i],col);CHKERRQ(ierr); 15164a56b808SFande Kong } 15174a56b808SFande Kong } 15184a56b808SFande Kong ierr = PetscHSetIGetSize(htd[i],&htsize);CHKERRQ(ierr); 15194a56b808SFande Kong dnz[i] = htsize; 15204a56b808SFande Kong ierr = PetscHSetIDestroy(&htd[i]);CHKERRQ(ierr); 15214a56b808SFande Kong ierr = PetscHSetIGetSize(hto[i],&htsize);CHKERRQ(ierr); 15224a56b808SFande Kong onz[i] = htsize; 15234a56b808SFande Kong ierr = PetscHSetIDestroy(&hto[i]);CHKERRQ(ierr); 15244a56b808SFande Kong } 15254a56b808SFande Kong 15264a56b808SFande Kong ierr = PetscFree2(htd,hto);CHKERRQ(ierr); 15274a56b808SFande Kong ierr = PetscFree(c_othj);CHKERRQ(ierr); 15284a56b808SFande Kong 15294a56b808SFande Kong /* local sizes and preallocation */ 15304a56b808SFande Kong ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 153134bcad68SFande Kong ierr = MatSetBlockSizes(Cmpi, dof>1? dof: P->cmap->bs,dof>1? dof: P->cmap->bs);CHKERRQ(ierr); 15324a56b808SFande Kong ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 15334a56b808SFande Kong ierr = PetscFree2(dnz,onz);CHKERRQ(ierr); 15344a56b808SFande Kong 15354a56b808SFande Kong /* attach the supporting struct to Cmpi for reuse */ 15366718818eSStefano Zampini Cmpi->product->data = ptap; 15376718818eSStefano Zampini Cmpi->product->destroy = MatDestroy_MPIAIJ_PtAP; 15386718818eSStefano Zampini Cmpi->product->view = MatView_MPIAIJ_PtAP; 15394a56b808SFande Kong 15404a56b808SFande Kong /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 15414a56b808SFande Kong Cmpi->assembled = PETSC_FALSE; 15424a56b808SFande Kong /* pick an algorithm */ 15434a56b808SFande Kong ierr = PetscOptionsBegin(PetscObjectComm((PetscObject)A),((PetscObject)A)->prefix,"MatPtAP","Mat");CHKERRQ(ierr); 15444a56b808SFande Kong alg = 0; 15454a56b808SFande Kong ierr = PetscOptionsEList("-matptap_allatonce_via","PtAP allatonce numeric approach","MatPtAP",algTypes,nalg,algTypes[alg],&alg,&flg);CHKERRQ(ierr); 15464a56b808SFande Kong ierr = PetscOptionsEnd();CHKERRQ(ierr); 15474a56b808SFande Kong switch (alg) { 15484a56b808SFande Kong case 0: 15494a56b808SFande Kong Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce_merged; 15504a56b808SFande Kong break; 15514a56b808SFande Kong case 1: 15524a56b808SFande Kong Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ_allatonce; 15534a56b808SFande Kong break; 15544a56b808SFande Kong default: 15554a56b808SFande Kong SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG," Unsupported allatonce numerical algorithm \n"); 15564a56b808SFande Kong } 15574a56b808SFande Kong PetscFunctionReturn(0); 15584a56b808SFande Kong } 15594a56b808SFande Kong 15604222ddf1SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce_merged(Mat A,Mat P,PetscReal fill,Mat C) 1561bc8e477aSFande Kong { 1562bc8e477aSFande Kong PetscErrorCode ierr; 1563bc8e477aSFande Kong 1564bc8e477aSFande Kong PetscFunctionBegin; 1565bc8e477aSFande Kong ierr = MatPtAPSymbolic_MPIAIJ_MPIXAIJ_allatonce_merged(A,P,1,fill,C);CHKERRQ(ierr); 1566bc8e477aSFande Kong PetscFunctionReturn(0); 1567bc8e477aSFande Kong } 1568bc8e477aSFande Kong 15694222ddf1SHong Zhang PetscErrorCode MatPtAPSymbolic_MPIAIJ_MPIAIJ(Mat A,Mat P,PetscReal fill,Mat Cmpi) 1570e953e47cSHong Zhang { 1571e953e47cSHong Zhang PetscErrorCode ierr; 15723cdca5ebSHong Zhang Mat_APMPI *ptap; 15736718818eSStefano Zampini Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 1574e953e47cSHong Zhang MPI_Comm comm; 1575e953e47cSHong Zhang PetscMPIInt size,rank; 1576e953e47cSHong Zhang PetscFreeSpaceList free_space=NULL,current_space=NULL; 1577e953e47cSHong Zhang PetscInt am=A->rmap->n,pm=P->rmap->n,pN=P->cmap->N,pn=P->cmap->n; 1578e953e47cSHong Zhang PetscInt *lnk,i,k,pnz,row,nsend; 1579e953e47cSHong Zhang PetscBT lnkbt; 1580ec4bef21SJose E. Roman PetscMPIInt tagi,tagj,*len_si,*len_s,*len_ri,nrecv; 1581ec4bef21SJose E. Roman PETSC_UNUSED PetscMPIInt icompleted=0; 1582e953e47cSHong Zhang PetscInt **buf_rj,**buf_ri,**buf_ri_k; 1583e953e47cSHong Zhang PetscInt len,proc,*dnz,*onz,*owners,nzi,nspacedouble; 1584e953e47cSHong Zhang PetscInt nrows,*buf_s,*buf_si,*buf_si_i,**nextrow,**nextci; 1585e953e47cSHong Zhang MPI_Request *swaits,*rwaits; 1586e953e47cSHong Zhang MPI_Status *sstatus,rstatus; 1587e953e47cSHong Zhang PetscLayout rowmap; 1588e953e47cSHong Zhang PetscInt *owners_co,*coi,*coj; /* i and j array of (p->B)^T*A*P - used in the communication */ 1589e953e47cSHong Zhang PetscMPIInt *len_r,*id_r; /* array of length of comm->size, store send/recv matrix values */ 1590e953e47cSHong Zhang PetscInt *api,*apj,*Jptr,apnz,*prmap=p->garray,con,j,ap_rmax=0,Crmax,*aj,*ai,*pi; 1591ec07b8f8SHong Zhang Mat_SeqAIJ *p_loc,*p_oth=NULL,*ad=(Mat_SeqAIJ*)(a->A)->data,*ao=NULL,*c_loc,*c_oth; 1592e953e47cSHong Zhang PetscScalar *apv; 1593e953e47cSHong Zhang PetscTable ta; 1594b92f168fSBarry Smith MatType mtype; 1595e83e5f86SFande Kong const char *prefix; 1596e953e47cSHong Zhang #if defined(PETSC_USE_INFO) 1597e953e47cSHong Zhang PetscReal apfill; 1598e953e47cSHong Zhang #endif 1599e953e47cSHong Zhang 1600e953e47cSHong Zhang PetscFunctionBegin; 16016718818eSStefano Zampini MatCheckProduct(Cmpi,4); 16026718818eSStefano Zampini if (Cmpi->product->data) SETERRQ(PetscObjectComm((PetscObject)Cmpi),PETSC_ERR_PLIB,"Product data not empty"); 1603b4e8d1b6SHong Zhang ierr = PetscObjectGetComm((PetscObject)A,&comm);CHKERRQ(ierr); 1604ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 1605ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 1606e953e47cSHong Zhang 160752f7967eSHong Zhang if (size > 1) ao = (Mat_SeqAIJ*)(a->B)->data; 1608ec07b8f8SHong Zhang 1609e953e47cSHong Zhang /* create symbolic parallel matrix Cmpi */ 1610b92f168fSBarry Smith ierr = MatGetType(A,&mtype);CHKERRQ(ierr); 1611b92f168fSBarry Smith ierr = MatSetType(Cmpi,mtype);CHKERRQ(ierr); 1612e953e47cSHong Zhang 1613e953e47cSHong Zhang /* Do dense axpy in MatPtAPNumeric_MPIAIJ_MPIAIJ() */ 1614e953e47cSHong Zhang Cmpi->ops->ptapnumeric = MatPtAPNumeric_MPIAIJ_MPIAIJ; 1615e953e47cSHong Zhang 16163cdca5ebSHong Zhang /* create struct Mat_APMPI and attached it to C later */ 161715a3b8e2SHong Zhang ierr = PetscNew(&ptap);CHKERRQ(ierr); 161815a3b8e2SHong Zhang ptap->reuse = MAT_INITIAL_MATRIX; 1619a4ffb656SHong Zhang ptap->algType = 1; 162015a3b8e2SHong Zhang 162115a3b8e2SHong Zhang /* get P_oth by taking rows of P (= non-zero cols of local A) from other processors */ 162215a3b8e2SHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_INITIAL_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 162315a3b8e2SHong Zhang /* get P_loc by taking all local rows of P */ 162415a3b8e2SHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_INITIAL_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 162515a3b8e2SHong Zhang 162667a12041SHong Zhang /* (0) compute Rd = Pd^T, Ro = Po^T */ 162767a12041SHong Zhang /* --------------------------------- */ 1628419ecdd9Sandi selinger ierr = MatTranspose(p->A,MAT_INITIAL_MATRIX,&ptap->Rd);CHKERRQ(ierr); 1629419ecdd9Sandi selinger ierr = MatTranspose(p->B,MAT_INITIAL_MATRIX,&ptap->Ro);CHKERRQ(ierr); 163015a3b8e2SHong Zhang 163167a12041SHong Zhang /* (1) compute symbolic AP = A_loc*P = Ad*P_loc + Ao*P_oth (api,apj) */ 163267a12041SHong Zhang /* ----------------------------------------------------------------- */ 163367a12041SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 163452f7967eSHong Zhang if (ptap->P_oth) p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 1635445158ffSHong Zhang 16369a6dcab7SHong Zhang /* create and initialize a linked list */ 163745d00d1dSHong Zhang ierr = PetscTableCreate(pn,pN,&ta);CHKERRQ(ierr); /* for compute AP_loc and Cmpi */ 16384b38b95cSHong Zhang MatRowMergeMax_SeqAIJ(p_loc,ptap->P_loc->rmap->N,ta); 16394b38b95cSHong Zhang MatRowMergeMax_SeqAIJ(p_oth,ptap->P_oth->rmap->N,ta); 164078d30b94SHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); /* Crmax = nnz(sum of Prows) */ 1641d0e9a2f7SHong 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); */ 164278d30b94SHong Zhang 164378d30b94SHong Zhang ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 1644445158ffSHong Zhang 16458cb82516SHong Zhang /* Initial FreeSpace size is fill*(nnz(A) + nnz(P)) */ 1646ec07b8f8SHong Zhang if (ao) { 1647f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],PetscIntSumTruncate(ao->i[am],p_loc->i[pm]))),&free_space);CHKERRQ(ierr); 1648ec07b8f8SHong Zhang } else { 1649ec07b8f8SHong Zhang ierr = PetscFreeSpaceGet(PetscRealIntMultTruncate(fill,PetscIntSumTruncate(ad->i[am],p_loc->i[pm])),&free_space);CHKERRQ(ierr); 1650ec07b8f8SHong Zhang } 1651445158ffSHong Zhang current_space = free_space; 165267a12041SHong Zhang nspacedouble = 0; 165367a12041SHong Zhang 165467a12041SHong Zhang ierr = PetscMalloc1(am+1,&api);CHKERRQ(ierr); 165567a12041SHong Zhang api[0] = 0; 1656445158ffSHong Zhang for (i=0; i<am; i++) { 165767a12041SHong Zhang /* diagonal portion: Ad[i,:]*P */ 165867a12041SHong Zhang ai = ad->i; pi = p_loc->i; 165967a12041SHong Zhang nzi = ai[i+1] - ai[i]; 166067a12041SHong Zhang aj = ad->j + ai[i]; 1661445158ffSHong Zhang for (j=0; j<nzi; j++) { 1662445158ffSHong Zhang row = aj[j]; 166367a12041SHong Zhang pnz = pi[row+1] - pi[row]; 166467a12041SHong Zhang Jptr = p_loc->j + pi[row]; 1665445158ffSHong Zhang /* add non-zero cols of P into the sorted linked list lnk */ 1666445158ffSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1667445158ffSHong Zhang } 166867a12041SHong Zhang /* off-diagonal portion: Ao[i,:]*P */ 1669ec07b8f8SHong Zhang if (ao) { 167067a12041SHong Zhang ai = ao->i; pi = p_oth->i; 167167a12041SHong Zhang nzi = ai[i+1] - ai[i]; 167267a12041SHong Zhang aj = ao->j + ai[i]; 1673445158ffSHong Zhang for (j=0; j<nzi; j++) { 1674445158ffSHong Zhang row = aj[j]; 167567a12041SHong Zhang pnz = pi[row+1] - pi[row]; 167667a12041SHong Zhang Jptr = p_oth->j + pi[row]; 1677445158ffSHong Zhang ierr = PetscLLCondensedAddSorted(pnz,Jptr,lnk,lnkbt);CHKERRQ(ierr); 1678445158ffSHong Zhang } 1679ec07b8f8SHong Zhang } 1680445158ffSHong Zhang apnz = lnk[0]; 1681445158ffSHong Zhang api[i+1] = api[i] + apnz; 1682445158ffSHong Zhang if (ap_rmax < apnz) ap_rmax = apnz; 1683445158ffSHong Zhang 1684445158ffSHong Zhang /* if free space is not available, double the total space in the list */ 1685445158ffSHong Zhang if (current_space->local_remaining<apnz) { 1686f91af8c7SBarry Smith ierr = PetscFreeSpaceGet(PetscIntSumTruncate(apnz,current_space->total_array_size),¤t_space);CHKERRQ(ierr); 1687445158ffSHong Zhang nspacedouble++; 1688445158ffSHong Zhang } 1689445158ffSHong Zhang 1690445158ffSHong Zhang /* Copy data into free space, then initialize lnk */ 1691445158ffSHong Zhang ierr = PetscLLCondensedClean(pN,apnz,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1692445158ffSHong Zhang 1693445158ffSHong Zhang current_space->array += apnz; 1694445158ffSHong Zhang current_space->local_used += apnz; 1695445158ffSHong Zhang current_space->local_remaining -= apnz; 1696445158ffSHong Zhang } 1697681d504bSHong Zhang /* Allocate space for apj and apv, initialize apj, and */ 1698445158ffSHong Zhang /* destroy list of free space and other temporary array(s) */ 1699445158ffSHong Zhang ierr = PetscMalloc2(api[am],&apj,api[am],&apv);CHKERRQ(ierr); 1700445158ffSHong Zhang ierr = PetscFreeSpaceContiguous(&free_space,apj);CHKERRQ(ierr); 17019a6dcab7SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 17029a6dcab7SHong Zhang 1703aa690a28SHong Zhang /* Create AP_loc for reuse */ 1704445158ffSHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,pN,api,apj,apv,&ptap->AP_loc);CHKERRQ(ierr); 1705cec0a6c6SStefano Zampini ierr = MatSetType(ptap->AP_loc,((PetscObject)p->A)->type_name);CHKERRQ(ierr); 1706aa690a28SHong Zhang #if defined(PETSC_USE_INFO) 1707ec07b8f8SHong Zhang if (ao) { 1708aa690a28SHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+ao->i[am]+p_loc->i[pm]+1); 1709ec07b8f8SHong Zhang } else { 1710ec07b8f8SHong Zhang apfill = (PetscReal)api[am]/(ad->i[am]+p_loc->i[pm]+1); 1711ec07b8f8SHong Zhang } 1712aa690a28SHong Zhang ptap->AP_loc->info.mallocs = nspacedouble; 1713aa690a28SHong Zhang ptap->AP_loc->info.fill_ratio_given = fill; 1714aa690a28SHong Zhang ptap->AP_loc->info.fill_ratio_needed = apfill; 1715aa690a28SHong Zhang 1716aa690a28SHong Zhang if (api[am]) { 1717b11c1ec8SHong 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); 1718aa690a28SHong Zhang ierr = PetscInfo1(ptap->AP_loc,"Use MatPtAP(A,B,MatReuse,%g,&C) for best AP_loc performance.;\n",(double)apfill);CHKERRQ(ierr); 1719aa690a28SHong Zhang } else { 1720b11c1ec8SHong Zhang ierr = PetscInfo(ptap->AP_loc,"Nonscalable algorithm, AP_loc is empty \n");CHKERRQ(ierr); 1721aa690a28SHong Zhang } 1722aa690a28SHong Zhang #endif 1723aa690a28SHong Zhang 1724681d504bSHong Zhang /* (2-1) compute symbolic Co = Ro*AP_loc */ 1725681d504bSHong Zhang /* ------------------------------------ */ 1726e83e5f86SFande Kong ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1727e83e5f86SFande Kong ierr = MatSetOptionsPrefix(ptap->Ro,prefix);CHKERRQ(ierr); 1728e83e5f86SFande Kong ierr = MatAppendOptionsPrefix(ptap->Ro,"inner_offdiag_");CHKERRQ(ierr); 17294222ddf1SHong Zhang ierr = MatProductCreate(ptap->Ro,ptap->AP_loc,NULL,&ptap->C_oth);CHKERRQ(ierr); 17304222ddf1SHong Zhang ierr = MatGetOptionsPrefix(Cmpi,&prefix);CHKERRQ(ierr); 17314222ddf1SHong Zhang ierr = MatSetOptionsPrefix(ptap->C_oth,prefix);CHKERRQ(ierr); 17324222ddf1SHong Zhang ierr = MatAppendOptionsPrefix(ptap->C_oth,"inner_C_oth_");CHKERRQ(ierr); 17334222ddf1SHong Zhang ierr = MatProductSetType(ptap->C_oth,MATPRODUCT_AB);CHKERRQ(ierr); 17344222ddf1SHong Zhang ierr = MatProductSetAlgorithm(ptap->C_oth,"default");CHKERRQ(ierr); 17354222ddf1SHong Zhang ierr = MatProductSetFill(ptap->C_oth,fill);CHKERRQ(ierr); 17364222ddf1SHong Zhang ierr = MatProductSetFromOptions(ptap->C_oth);CHKERRQ(ierr); 17374222ddf1SHong Zhang ierr = MatProductSymbolic(ptap->C_oth);CHKERRQ(ierr); 173815a3b8e2SHong Zhang 173967a12041SHong Zhang /* (3) send coj of C_oth to other processors */ 174067a12041SHong Zhang /* ------------------------------------------ */ 174115a3b8e2SHong Zhang /* determine row ownership */ 1742445158ffSHong Zhang ierr = PetscLayoutCreate(comm,&rowmap);CHKERRQ(ierr); 1743445158ffSHong Zhang rowmap->n = pn; 1744445158ffSHong Zhang rowmap->bs = 1; 1745445158ffSHong Zhang ierr = PetscLayoutSetUp(rowmap);CHKERRQ(ierr); 1746445158ffSHong Zhang owners = rowmap->range; 174715a3b8e2SHong Zhang 174815a3b8e2SHong Zhang /* determine the number of messages to send, their lengths */ 17498cb82516SHong Zhang ierr = PetscMalloc4(size,&len_s,size,&len_si,size,&sstatus,size+2,&owners_co);CHKERRQ(ierr); 1750580bdb30SBarry Smith ierr = PetscArrayzero(len_s,size);CHKERRQ(ierr); 1751580bdb30SBarry Smith ierr = PetscArrayzero(len_si,size);CHKERRQ(ierr); 175215a3b8e2SHong Zhang 175367a12041SHong Zhang c_oth = (Mat_SeqAIJ*)ptap->C_oth->data; 175467a12041SHong Zhang coi = c_oth->i; coj = c_oth->j; 175567a12041SHong Zhang con = ptap->C_oth->rmap->n; 175615a3b8e2SHong Zhang proc = 0; 175767a12041SHong Zhang for (i=0; i<con; i++) { 175815a3b8e2SHong Zhang while (prmap[i] >= owners[proc+1]) proc++; 175915a3b8e2SHong Zhang len_si[proc]++; /* num of rows in Co(=Pt*AP) to be sent to [proc] */ 176015a3b8e2SHong Zhang len_s[proc] += coi[i+1] - coi[i]; /* num of nonzeros in Co to be sent to [proc] */ 176115a3b8e2SHong Zhang } 176215a3b8e2SHong Zhang 176315a3b8e2SHong Zhang len = 0; /* max length of buf_si[], see (4) */ 176415a3b8e2SHong Zhang owners_co[0] = 0; 176567a12041SHong Zhang nsend = 0; 176615a3b8e2SHong Zhang for (proc=0; proc<size; proc++) { 176715a3b8e2SHong Zhang owners_co[proc+1] = owners_co[proc] + len_si[proc]; 176815a3b8e2SHong Zhang if (len_s[proc]) { 1769445158ffSHong Zhang nsend++; 177015a3b8e2SHong Zhang len_si[proc] = 2*(len_si[proc] + 1); /* length of buf_si to be sent to [proc] */ 177115a3b8e2SHong Zhang len += len_si[proc]; 177215a3b8e2SHong Zhang } 177315a3b8e2SHong Zhang } 177415a3b8e2SHong Zhang 177515a3b8e2SHong Zhang /* determine the number and length of messages to receive for coi and coj */ 1776445158ffSHong Zhang ierr = PetscGatherNumberOfMessages(comm,NULL,len_s,&nrecv);CHKERRQ(ierr); 1777445158ffSHong Zhang ierr = PetscGatherMessageLengths2(comm,nsend,nrecv,len_s,len_si,&id_r,&len_r,&len_ri);CHKERRQ(ierr); 177815a3b8e2SHong Zhang 177915a3b8e2SHong Zhang /* post the Irecv and Isend of coj */ 178015a3b8e2SHong Zhang ierr = PetscCommGetNewTag(comm,&tagj);CHKERRQ(ierr); 1781445158ffSHong Zhang ierr = PetscPostIrecvInt(comm,tagj,nrecv,id_r,len_r,&buf_rj,&rwaits);CHKERRQ(ierr); 1782445158ffSHong Zhang ierr = PetscMalloc1(nsend+1,&swaits);CHKERRQ(ierr); 178315a3b8e2SHong Zhang for (proc=0, k=0; proc<size; proc++) { 178415a3b8e2SHong Zhang if (!len_s[proc]) continue; 178515a3b8e2SHong Zhang i = owners_co[proc]; 1786ffc4695bSBarry Smith ierr = MPI_Isend(coj+coi[i],len_s[proc],MPIU_INT,proc,tagj,comm,swaits+k);CHKERRMPI(ierr); 178715a3b8e2SHong Zhang k++; 178815a3b8e2SHong Zhang } 178915a3b8e2SHong Zhang 1790681d504bSHong Zhang /* (2-2) compute symbolic C_loc = Rd*AP_loc */ 1791681d504bSHong Zhang /* ---------------------------------------- */ 1792e83e5f86SFande Kong ierr = MatSetOptionsPrefix(ptap->Rd,prefix);CHKERRQ(ierr); 1793e83e5f86SFande Kong ierr = MatAppendOptionsPrefix(ptap->Rd,"inner_diag_");CHKERRQ(ierr); 17944222ddf1SHong Zhang ierr = MatProductCreate(ptap->Rd,ptap->AP_loc,NULL,&ptap->C_loc);CHKERRQ(ierr); 17954222ddf1SHong Zhang ierr = MatGetOptionsPrefix(Cmpi,&prefix);CHKERRQ(ierr); 17964222ddf1SHong Zhang ierr = MatSetOptionsPrefix(ptap->C_loc,prefix);CHKERRQ(ierr); 17974222ddf1SHong Zhang ierr = MatAppendOptionsPrefix(ptap->C_loc,"inner_C_loc_");CHKERRQ(ierr); 17984222ddf1SHong Zhang ierr = MatProductSetType(ptap->C_loc,MATPRODUCT_AB);CHKERRQ(ierr); 17994222ddf1SHong Zhang ierr = MatProductSetAlgorithm(ptap->C_loc,"default");CHKERRQ(ierr); 18004222ddf1SHong Zhang ierr = MatProductSetFill(ptap->C_loc,fill);CHKERRQ(ierr); 18014222ddf1SHong Zhang ierr = MatProductSetFromOptions(ptap->C_loc);CHKERRQ(ierr); 18024222ddf1SHong Zhang ierr = MatProductSymbolic(ptap->C_loc);CHKERRQ(ierr); 18034222ddf1SHong Zhang 1804681d504bSHong Zhang c_loc = (Mat_SeqAIJ*)ptap->C_loc->data; 1805681d504bSHong Zhang 1806681d504bSHong Zhang /* receives coj are complete */ 1807445158ffSHong Zhang for (i=0; i<nrecv; i++) { 1808ffc4695bSBarry Smith ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRMPI(ierr); 180915a3b8e2SHong Zhang } 181015a3b8e2SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1811ffc4695bSBarry Smith if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRMPI(ierr);} 181215a3b8e2SHong Zhang 181378d30b94SHong Zhang /* add received column indices into ta to update Crmax */ 181478d30b94SHong Zhang for (k=0; k<nrecv; k++) {/* k-th received message */ 181578d30b94SHong Zhang Jptr = buf_rj[k]; 181678d30b94SHong Zhang for (j=0; j<len_r[k]; j++) { 181778d30b94SHong Zhang ierr = PetscTableAdd(ta,*(Jptr+j)+1,1,INSERT_VALUES);CHKERRQ(ierr); 181878d30b94SHong Zhang } 181978d30b94SHong Zhang } 182078d30b94SHong Zhang ierr = PetscTableGetCount(ta,&Crmax);CHKERRQ(ierr); 182178d30b94SHong Zhang ierr = PetscTableDestroy(&ta);CHKERRQ(ierr); 18229a6dcab7SHong Zhang 182315a3b8e2SHong Zhang /* (4) send and recv coi */ 182415a3b8e2SHong Zhang /*-----------------------*/ 182515a3b8e2SHong Zhang ierr = PetscCommGetNewTag(comm,&tagi);CHKERRQ(ierr); 1826445158ffSHong Zhang ierr = PetscPostIrecvInt(comm,tagi,nrecv,id_r,len_ri,&buf_ri,&rwaits);CHKERRQ(ierr); 182715a3b8e2SHong Zhang ierr = PetscMalloc1(len+1,&buf_s);CHKERRQ(ierr); 182815a3b8e2SHong Zhang buf_si = buf_s; /* points to the beginning of k-th msg to be sent */ 182915a3b8e2SHong Zhang for (proc=0,k=0; proc<size; proc++) { 183015a3b8e2SHong Zhang if (!len_s[proc]) continue; 183115a3b8e2SHong Zhang /* form outgoing message for i-structure: 183215a3b8e2SHong Zhang buf_si[0]: nrows to be sent 183315a3b8e2SHong Zhang [1:nrows]: row index (global) 183415a3b8e2SHong Zhang [nrows+1:2*nrows+1]: i-structure index 183515a3b8e2SHong Zhang */ 183615a3b8e2SHong Zhang /*-------------------------------------------*/ 183715a3b8e2SHong Zhang nrows = len_si[proc]/2 - 1; /* num of rows in Co to be sent to [proc] */ 183815a3b8e2SHong Zhang buf_si_i = buf_si + nrows+1; 183915a3b8e2SHong Zhang buf_si[0] = nrows; 184015a3b8e2SHong Zhang buf_si_i[0] = 0; 184115a3b8e2SHong Zhang nrows = 0; 184215a3b8e2SHong Zhang for (i=owners_co[proc]; i<owners_co[proc+1]; i++) { 184315a3b8e2SHong Zhang nzi = coi[i+1] - coi[i]; 184415a3b8e2SHong Zhang buf_si_i[nrows+1] = buf_si_i[nrows] + nzi; /* i-structure */ 184515a3b8e2SHong Zhang buf_si[nrows+1] = prmap[i] -owners[proc]; /* local row index */ 184615a3b8e2SHong Zhang nrows++; 184715a3b8e2SHong Zhang } 1848ffc4695bSBarry Smith ierr = MPI_Isend(buf_si,len_si[proc],MPIU_INT,proc,tagi,comm,swaits+k);CHKERRMPI(ierr); 184915a3b8e2SHong Zhang k++; 185015a3b8e2SHong Zhang buf_si += len_si[proc]; 185115a3b8e2SHong Zhang } 1852681d504bSHong Zhang for (i=0; i<nrecv; i++) { 1853ffc4695bSBarry Smith ierr = MPI_Waitany(nrecv,rwaits,&icompleted,&rstatus);CHKERRMPI(ierr); 185415a3b8e2SHong Zhang } 185515a3b8e2SHong Zhang ierr = PetscFree(rwaits);CHKERRQ(ierr); 1856ffc4695bSBarry Smith if (nsend) {ierr = MPI_Waitall(nsend,swaits,sstatus);CHKERRMPI(ierr);} 185715a3b8e2SHong Zhang 18588cb82516SHong Zhang ierr = PetscFree4(len_s,len_si,sstatus,owners_co);CHKERRQ(ierr); 185915a3b8e2SHong Zhang ierr = PetscFree(len_ri);CHKERRQ(ierr); 186015a3b8e2SHong Zhang ierr = PetscFree(swaits);CHKERRQ(ierr); 186115a3b8e2SHong Zhang ierr = PetscFree(buf_s);CHKERRQ(ierr); 1862a4ffb656SHong Zhang 186367a12041SHong Zhang /* (5) compute the local portion of Cmpi */ 186467a12041SHong Zhang /* ------------------------------------------ */ 186578d30b94SHong Zhang /* set initial free space to be Crmax, sufficient for holding nozeros in each row of Cmpi */ 186678d30b94SHong Zhang ierr = PetscFreeSpaceGet(Crmax,&free_space);CHKERRQ(ierr); 186715a3b8e2SHong Zhang current_space = free_space; 186815a3b8e2SHong Zhang 1869445158ffSHong Zhang ierr = PetscMalloc3(nrecv,&buf_ri_k,nrecv,&nextrow,nrecv,&nextci);CHKERRQ(ierr); 1870445158ffSHong Zhang for (k=0; k<nrecv; k++) { 187115a3b8e2SHong Zhang buf_ri_k[k] = buf_ri[k]; /* beginning of k-th recved i-structure */ 187215a3b8e2SHong Zhang nrows = *buf_ri_k[k]; 187315a3b8e2SHong Zhang nextrow[k] = buf_ri_k[k] + 1; /* next row number of k-th recved i-structure */ 187415a3b8e2SHong Zhang nextci[k] = buf_ri_k[k] + (nrows + 1); /* poins to the next i-structure of k-th recved i-structure */ 187515a3b8e2SHong Zhang } 1876445158ffSHong Zhang 187778d30b94SHong Zhang ierr = MatPreallocateInitialize(comm,pn,pn,dnz,onz);CHKERRQ(ierr); 187878d30b94SHong Zhang ierr = PetscLLCondensedCreate(Crmax,pN,&lnk,&lnkbt);CHKERRQ(ierr); 187915a3b8e2SHong Zhang for (i=0; i<pn; i++) { 188067a12041SHong Zhang /* add C_loc into Cmpi */ 188167a12041SHong Zhang nzi = c_loc->i[i+1] - c_loc->i[i]; 188267a12041SHong Zhang Jptr = c_loc->j + c_loc->i[i]; 188367a12041SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 188415a3b8e2SHong Zhang 188515a3b8e2SHong Zhang /* add received col data into lnk */ 1886445158ffSHong Zhang for (k=0; k<nrecv; k++) { /* k-th received message */ 188715a3b8e2SHong Zhang if (i == *nextrow[k]) { /* i-th row */ 188815a3b8e2SHong Zhang nzi = *(nextci[k]+1) - *nextci[k]; 188915a3b8e2SHong Zhang Jptr = buf_rj[k] + *nextci[k]; 189015a3b8e2SHong Zhang ierr = PetscLLCondensedAddSorted(nzi,Jptr,lnk,lnkbt);CHKERRQ(ierr); 189115a3b8e2SHong Zhang nextrow[k]++; nextci[k]++; 189215a3b8e2SHong Zhang } 189315a3b8e2SHong Zhang } 1894d0e9a2f7SHong Zhang nzi = lnk[0]; 18958cb82516SHong Zhang 189615a3b8e2SHong Zhang /* copy data into free space, then initialize lnk */ 1897d0e9a2f7SHong Zhang ierr = PetscLLCondensedClean(pN,nzi,current_space->array,lnk,lnkbt);CHKERRQ(ierr); 1898d0e9a2f7SHong Zhang ierr = MatPreallocateSet(i+owners[rank],nzi,current_space->array,dnz,onz);CHKERRQ(ierr); 189915a3b8e2SHong Zhang } 190015a3b8e2SHong Zhang ierr = PetscFree3(buf_ri_k,nextrow,nextci);CHKERRQ(ierr); 190115a3b8e2SHong Zhang ierr = PetscLLDestroy(lnk,lnkbt);CHKERRQ(ierr); 1902445158ffSHong Zhang ierr = PetscFreeSpaceDestroy(free_space);CHKERRQ(ierr); 190380bb4639SHong Zhang 1904ae5f0867Sstefano_zampini /* local sizes and preallocation */ 190515a3b8e2SHong Zhang ierr = MatSetSizes(Cmpi,pn,pn,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 1906ac94c67aSTristan Konolige if (P->cmap->bs > 0) { 1907ac94c67aSTristan Konolige ierr = PetscLayoutSetBlockSize(Cmpi->rmap,P->cmap->bs);CHKERRQ(ierr); 1908ac94c67aSTristan Konolige ierr = PetscLayoutSetBlockSize(Cmpi->cmap,P->cmap->bs);CHKERRQ(ierr); 1909ac94c67aSTristan Konolige } 191015a3b8e2SHong Zhang ierr = MatMPIAIJSetPreallocation(Cmpi,0,dnz,0,onz);CHKERRQ(ierr); 191115a3b8e2SHong Zhang ierr = MatPreallocateFinalize(dnz,onz);CHKERRQ(ierr); 191215a3b8e2SHong Zhang 1913445158ffSHong Zhang /* members in merge */ 1914445158ffSHong Zhang ierr = PetscFree(id_r);CHKERRQ(ierr); 1915445158ffSHong Zhang ierr = PetscFree(len_r);CHKERRQ(ierr); 1916445158ffSHong Zhang ierr = PetscFree(buf_ri[0]);CHKERRQ(ierr); 1917445158ffSHong Zhang ierr = PetscFree(buf_ri);CHKERRQ(ierr); 1918445158ffSHong Zhang ierr = PetscFree(buf_rj[0]);CHKERRQ(ierr); 1919445158ffSHong Zhang ierr = PetscFree(buf_rj);CHKERRQ(ierr); 1920445158ffSHong Zhang ierr = PetscLayoutDestroy(&rowmap);CHKERRQ(ierr); 192115a3b8e2SHong Zhang 19228cb82516SHong Zhang ierr = PetscCalloc1(pN,&ptap->apa);CHKERRQ(ierr); 19232259aa2eSHong Zhang 19246718818eSStefano Zampini /* attach the supporting struct to Cmpi for reuse */ 19256718818eSStefano Zampini Cmpi->product->data = ptap; 19266718818eSStefano Zampini Cmpi->product->destroy = MatDestroy_MPIAIJ_PtAP; 19276718818eSStefano Zampini Cmpi->product->view = MatView_MPIAIJ_PtAP; 19286718818eSStefano Zampini 19291a47ec54SHong Zhang /* Cmpi is not ready for use - assembly will be done by MatPtAPNumeric() */ 19301a47ec54SHong Zhang Cmpi->assembled = PETSC_FALSE; 19310d3441aeSHong Zhang PetscFunctionReturn(0); 19320d3441aeSHong Zhang } 19330d3441aeSHong Zhang 1934aa690a28SHong Zhang PetscErrorCode MatPtAPNumeric_MPIAIJ_MPIAIJ(Mat A,Mat P,Mat C) 19350d3441aeSHong Zhang { 19360d3441aeSHong Zhang PetscErrorCode ierr; 19376718818eSStefano Zampini Mat_MPIAIJ *a=(Mat_MPIAIJ*)A->data,*p=(Mat_MPIAIJ*)P->data; 19380d3441aeSHong Zhang Mat_SeqAIJ *ad=(Mat_SeqAIJ*)(a->A)->data,*ao=(Mat_SeqAIJ*)(a->B)->data; 19392dd9e643SHong Zhang Mat_SeqAIJ *ap,*p_loc,*p_oth=NULL,*c_seq; 19406718818eSStefano Zampini Mat_APMPI *ptap; 19419ce11a7cSHong Zhang Mat AP_loc,C_loc,C_oth; 19420d3441aeSHong Zhang PetscInt i,rstart,rend,cm,ncols,row; 19438cb82516SHong Zhang PetscInt *api,*apj,am = A->rmap->n,j,col,apnz; 19448cb82516SHong Zhang PetscScalar *apa; 19450d3441aeSHong Zhang const PetscInt *cols; 19460d3441aeSHong Zhang const PetscScalar *vals; 19470d3441aeSHong Zhang 19480d3441aeSHong Zhang PetscFunctionBegin; 19496718818eSStefano Zampini MatCheckProduct(C,3); 19506718818eSStefano Zampini ptap = (Mat_APMPI*)C->product->data; 19516718818eSStefano Zampini if (!ptap) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be computed. Missing data"); 19526718818eSStefano Zampini if (!ptap->AP_loc) SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_ARG_WRONGSTATE,"PtAP cannot be reused. Do not call MatProductClear()"); 1953a9899c97SHong Zhang 19540d3441aeSHong Zhang ierr = MatZeroEntries(C);CHKERRQ(ierr); 1955e497d3c8SHong Zhang /* 1) get R = Pd^T,Ro = Po^T */ 1956748c7196SHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 1957419ecdd9Sandi selinger ierr = MatTranspose(p->A,MAT_REUSE_MATRIX,&ptap->Rd);CHKERRQ(ierr); 1958419ecdd9Sandi selinger ierr = MatTranspose(p->B,MAT_REUSE_MATRIX,&ptap->Ro);CHKERRQ(ierr); 1959748c7196SHong Zhang } 19600d3441aeSHong Zhang 1961e497d3c8SHong Zhang /* 2) get AP_loc */ 19620d3441aeSHong Zhang AP_loc = ptap->AP_loc; 19638cb82516SHong Zhang ap = (Mat_SeqAIJ*)AP_loc->data; 19640d3441aeSHong Zhang 1965e497d3c8SHong Zhang /* 2-1) get P_oth = ptap->P_oth and P_loc = ptap->P_loc */ 19660d3441aeSHong Zhang /*-----------------------------------------------------*/ 1967748c7196SHong Zhang if (ptap->reuse == MAT_REUSE_MATRIX) { 1968748c7196SHong Zhang /* P_oth and P_loc are obtained in MatPtASymbolic() when reuse == MAT_INITIAL_MATRIX */ 19690d3441aeSHong Zhang ierr = MatGetBrowsOfAoCols_MPIAIJ(A,P,MAT_REUSE_MATRIX,&ptap->startsj_s,&ptap->startsj_r,&ptap->bufa,&ptap->P_oth);CHKERRQ(ierr); 19700d3441aeSHong Zhang ierr = MatMPIAIJGetLocalMat(P,MAT_REUSE_MATRIX,&ptap->P_loc);CHKERRQ(ierr); 1971e497d3c8SHong Zhang } 19720d3441aeSHong Zhang 19738cb82516SHong Zhang /* 2-2) compute numeric A_loc*P - dominating part */ 19748cb82516SHong Zhang /* ---------------------------------------------- */ 19750d3441aeSHong Zhang /* get data from symbolic products */ 19768cb82516SHong Zhang p_loc = (Mat_SeqAIJ*)(ptap->P_loc)->data; 19772dd9e643SHong Zhang if (ptap->P_oth) { 19788cb82516SHong Zhang p_oth = (Mat_SeqAIJ*)(ptap->P_oth)->data; 19792dd9e643SHong Zhang } 19808cb82516SHong Zhang apa = ptap->apa; 1981681d504bSHong Zhang api = ap->i; 1982681d504bSHong Zhang apj = ap->j; 1983e497d3c8SHong Zhang for (i=0; i<am; i++) { 1984445158ffSHong Zhang /* AP[i,:] = A[i,:]*P = Ad*P_loc Ao*P_oth */ 1985e497d3c8SHong Zhang AProw_nonscalable(i,ad,ao,p_loc,p_oth,apa); 1986e497d3c8SHong Zhang apnz = api[i+1] - api[i]; 1987e497d3c8SHong Zhang for (j=0; j<apnz; j++) { 1988e497d3c8SHong Zhang col = apj[j+api[i]]; 1989e497d3c8SHong Zhang ap->a[j+ap->i[i]] = apa[col]; 1990e497d3c8SHong Zhang apa[col] = 0.0; 19910d3441aeSHong Zhang } 19920d3441aeSHong Zhang } 1993976452c9SRichard Tran Mills /* We have modified the contents of local matrix AP_loc and must increase its ObjectState, since we are not doing AssemblyBegin/End on it. */ 1994976452c9SRichard Tran Mills ierr = PetscObjectStateIncrease((PetscObject)AP_loc);CHKERRQ(ierr); 19950d3441aeSHong Zhang 19968cb82516SHong Zhang /* 3) C_loc = Rd*AP_loc, C_oth = Ro*AP_loc */ 1997*65e4b4d4SStefano Zampini ierr = MatProductNumeric(ptap->C_loc);CHKERRQ(ierr); 1998*65e4b4d4SStefano Zampini ierr = MatProductNumeric(ptap->C_oth);CHKERRQ(ierr); 19999ce11a7cSHong Zhang C_loc = ptap->C_loc; 20009ce11a7cSHong Zhang C_oth = ptap->C_oth; 20010d3441aeSHong Zhang 20020d3441aeSHong Zhang /* add C_loc and Co to to C */ 20030d3441aeSHong Zhang ierr = MatGetOwnershipRange(C,&rstart,&rend);CHKERRQ(ierr); 20040d3441aeSHong Zhang 20050d3441aeSHong Zhang /* C_loc -> C */ 20060d3441aeSHong Zhang cm = C_loc->rmap->N; 20079ce11a7cSHong Zhang c_seq = (Mat_SeqAIJ*)C_loc->data; 20088cb82516SHong Zhang cols = c_seq->j; 20098cb82516SHong Zhang vals = c_seq->a; 2010904d1e70Sandi selinger 2011904d1e70Sandi selinger 2012e9ede7d0Sandi selinger /* The (fast) MatSetValues_MPIAIJ_CopyFromCSRFormat function can only be used when C->was_assembled is PETSC_FALSE and */ 2013904d1e70Sandi selinger /* when there are no off-processor parts. */ 20141de21080Sandi selinger /* If was_assembled is true, then the statement aj[rowstart_diag+dnz_row] = mat_j[col] - cstart; in MatSetValues_MPIAIJ_CopyFromCSRFormat */ 20151de21080Sandi selinger /* is no longer true. Then the more complex function MatSetValues_MPIAIJ() has to be used, where the column index is looked up from */ 20161de21080Sandi selinger /* a table, and other, more complex stuff has to be done. */ 2017904d1e70Sandi selinger if (C->assembled) { 2018904d1e70Sandi selinger C->was_assembled = PETSC_TRUE; 2019904d1e70Sandi selinger C->assembled = PETSC_FALSE; 2020904d1e70Sandi selinger } 2021904d1e70Sandi selinger if (C->was_assembled) { 20220d3441aeSHong Zhang for (i=0; i<cm; i++) { 20239ce11a7cSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 20240d3441aeSHong Zhang row = rstart + i; 2025904d1e70Sandi selinger ierr = MatSetValues_MPIAIJ(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 20268cb82516SHong Zhang cols += ncols; vals += ncols; 20270d3441aeSHong Zhang } 2028904d1e70Sandi selinger } else { 2029e9ede7d0Sandi selinger ierr = MatSetValues_MPIAIJ_CopyFromCSRFormat(C,c_seq->j,c_seq->i,c_seq->a);CHKERRQ(ierr); 2030904d1e70Sandi selinger } 20310d3441aeSHong Zhang 20320d3441aeSHong Zhang /* Co -> C, off-processor part */ 20339ce11a7cSHong Zhang cm = C_oth->rmap->N; 20349ce11a7cSHong Zhang c_seq = (Mat_SeqAIJ*)C_oth->data; 20358cb82516SHong Zhang cols = c_seq->j; 20368cb82516SHong Zhang vals = c_seq->a; 20379ce11a7cSHong Zhang for (i=0; i<cm; i++) { 20389ce11a7cSHong Zhang ncols = c_seq->i[i+1] - c_seq->i[i]; 20390d3441aeSHong Zhang row = p->garray[i]; 20400d3441aeSHong Zhang ierr = MatSetValues(C,1,&row,ncols,cols,vals,ADD_VALUES);CHKERRQ(ierr); 20418cb82516SHong Zhang cols += ncols; vals += ncols; 20420d3441aeSHong Zhang } 2043904d1e70Sandi selinger 20440d3441aeSHong Zhang ierr = MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20450d3441aeSHong Zhang ierr = MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20460d3441aeSHong Zhang 2047748c7196SHong Zhang ptap->reuse = MAT_REUSE_MATRIX; 20480d3441aeSHong Zhang PetscFunctionReturn(0); 20490d3441aeSHong Zhang } 20504222ddf1SHong Zhang 20514222ddf1SHong Zhang PETSC_INTERN PetscErrorCode MatProductSymbolic_PtAP_MPIAIJ_MPIAIJ(Mat C) 20524222ddf1SHong Zhang { 20534222ddf1SHong Zhang PetscErrorCode ierr; 20544222ddf1SHong Zhang Mat_Product *product = C->product; 20554222ddf1SHong Zhang Mat A=product->A,P=product->B; 20564222ddf1SHong Zhang MatProductAlgorithm alg=product->alg; 20574222ddf1SHong Zhang PetscReal fill=product->fill; 20584222ddf1SHong Zhang PetscBool flg; 20594222ddf1SHong Zhang 20604222ddf1SHong Zhang PetscFunctionBegin; 20614222ddf1SHong Zhang /* scalable: do R=P^T locally, then C=R*A*P */ 20624222ddf1SHong Zhang ierr = PetscStrcmp(alg,"scalable",&flg);CHKERRQ(ierr); 20634222ddf1SHong Zhang if (flg) { 20644222ddf1SHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_scalable(A,P,product->fill,C);CHKERRQ(ierr); 20654222ddf1SHong Zhang C->ops->productnumeric = MatProductNumeric_PtAP; 20664222ddf1SHong Zhang goto next; 20674222ddf1SHong Zhang } 20684222ddf1SHong Zhang 20694222ddf1SHong Zhang /* nonscalable: do R=P^T locally, then C=R*A*P */ 20704222ddf1SHong Zhang ierr = PetscStrcmp(alg,"nonscalable",&flg);CHKERRQ(ierr); 20714222ddf1SHong Zhang if (flg) { 20724222ddf1SHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ(A,P,fill,C);CHKERRQ(ierr); 20734222ddf1SHong Zhang goto next; 20744222ddf1SHong Zhang } 20754222ddf1SHong Zhang 20764222ddf1SHong Zhang /* allatonce */ 20774222ddf1SHong Zhang ierr = PetscStrcmp(alg,"allatonce",&flg);CHKERRQ(ierr); 20784222ddf1SHong Zhang if (flg) { 20794222ddf1SHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce(A,P,fill,C);CHKERRQ(ierr); 20804222ddf1SHong Zhang goto next; 20814222ddf1SHong Zhang } 20824222ddf1SHong Zhang 20834222ddf1SHong Zhang /* allatonce_merged */ 20844222ddf1SHong Zhang ierr = PetscStrcmp(alg,"allatonce_merged",&flg);CHKERRQ(ierr); 20854222ddf1SHong Zhang if (flg) { 20864222ddf1SHong Zhang ierr = MatPtAPSymbolic_MPIAIJ_MPIAIJ_allatonce_merged(A,P,fill,C);CHKERRQ(ierr); 20874222ddf1SHong Zhang goto next; 20884222ddf1SHong Zhang } 20894222ddf1SHong Zhang 20904e84afc0SStefano Zampini /* backend general code */ 20914e84afc0SStefano Zampini ierr = PetscStrcmp(alg,"backend",&flg);CHKERRQ(ierr); 20924e84afc0SStefano Zampini if (flg) { 20934e84afc0SStefano Zampini ierr = MatProductSymbolic_MPIAIJBACKEND(C);CHKERRQ(ierr); 20944e84afc0SStefano Zampini PetscFunctionReturn(0); 20954e84afc0SStefano Zampini } 20964e84afc0SStefano Zampini 20974222ddf1SHong Zhang /* hypre */ 20984222ddf1SHong Zhang #if defined(PETSC_HAVE_HYPRE) 20994222ddf1SHong Zhang ierr = PetscStrcmp(alg,"hypre",&flg);CHKERRQ(ierr); 21004222ddf1SHong Zhang if (flg) { 21014222ddf1SHong Zhang ierr = MatPtAPSymbolic_AIJ_AIJ_wHYPRE(A,P,fill,C);CHKERRQ(ierr); 21024222ddf1SHong Zhang C->ops->productnumeric = MatProductNumeric_PtAP; 21034222ddf1SHong Zhang PetscFunctionReturn(0); 21044222ddf1SHong Zhang } 21054222ddf1SHong Zhang #endif 21064222ddf1SHong Zhang SETERRQ(PetscObjectComm((PetscObject)C),PETSC_ERR_SUP,"Mat Product Algorithm is not supported"); 21074222ddf1SHong Zhang 21084222ddf1SHong Zhang next: 21094222ddf1SHong Zhang C->ops->productnumeric = MatProductNumeric_PtAP; 21104222ddf1SHong Zhang PetscFunctionReturn(0); 21114222ddf1SHong Zhang } 2112